aboutsummaryrefslogtreecommitdiff
path: root/src/conversion.F90
blob: 7dc6f6038125f4401dbbd9e11512a84b53db2bcd (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
! $Header$

#include "cctk.h"

module conversion
  implicit none
  private
  
  public make_cylindrical2cartesian
  public make_spherical2cartesian
  public make_cartesian2spherical
  
contains
  
  subroutine make_cylindrical2cartesian (xx, tt)
    CCTK_REAL, intent(in)  :: xx(3)
    CCTK_REAL, intent(out) :: tt(3,3)
    CCTK_REAL, parameter :: eps = 1.0d-20
    CCTK_REAL :: x,y,z
    CCTK_REAL :: rho
    
    ! cylindrical coordinates: rho phi z
    ! Cartesian   coordinates: x   y   z
    !    rho^2    = x^2 + y^2
    !    tan(phi) = y/x
    
    ! A[cart](p)_i = TT(p)_i^j A[cyl](p)_j
    !    where p = p[cart]
    
    ! Transformation:
    !    x = rho cos(phi)
    !    y = rho sin(phi)
    !    z = z
    !    rho = sqrt(x^2+y^2)
    !    phi = atan(y/x)
    !    z   = z
    
    ! Jacobian:
    !    drho/dx = x/rho
    !    drho/dy = y/rho
    !    drho/dz = 0
    !    dphi/dx = -y/rho^2
    !    dphi/dy = x/rho^2
    !    dphi/dz = 0
    !    dz  /dx = 0
    !    dz  /dy = 0
    !    dz  /dz = 1
    
    x = xx(1)
    y = xx(2)
    z = xx(3)
    
    rho = sqrt(x**2+y**2)
    
    tt(1,1) = x/(rho+eps)
    tt(2,1) = y/(rho+eps)
    tt(3,1) = 0
    tt(1,2) = -y/(rho**2+eps)
    tt(2,2) = x/(rho**2+eps)
    tt(3,2) = 0
    tt(1,3) = 0
    tt(2,3) = 0
    tt(3,3) = 1
    
  end subroutine make_cylindrical2cartesian
  
  subroutine make_spherical2cartesian (xx, tt)
    CCTK_REAL, intent(in)  :: xx(3)
    CCTK_REAL, intent(out) :: tt(3,3)
    CCTK_REAL, parameter :: eps = 1.0d-20
    CCTK_REAL :: x,y,z
    CCTK_REAL :: rho,r
    
    ! Cartesian coordinates: x y z
    ! spherical coordinates: r theta phi
    !    r^2        = x^2 + y^2 + z^2
    !    cos(theta) = z/r
    !    tan(phi)   = y/x
    
    ! A[cart](p)_i = TT(p)_i^j A[spher](p)_j
    !    where p = p[cart]
    
    ! Definition:
    !    rho   = r sin(theta)
    !    z/rho = cos(theta) / sin(theta)
    
    ! Transformation:
    !    x = r sin(theta) cos(phi)
    !    y = r sin(theta) sin(phi)
    !    z = r cos(theta)
    !    r     = sqrt(x^2+y^2+z^2)
    !    theta = acos(z/r)
    !    phi   = atan(y/x)
    
    ! Jacobian:
    !    dr    /dx =  x/r
    !    dr    /dy =  y/r
    !    dr    /dz =  z/r
    !    dtheta/dx =  xz/r^2rho           ! 0
    !    dtheta/dy =  yz/r^2rho           ! 0
    !    dtheta/dz = -rho/r^2             ! 1/rho
    !    dphi  /dx = -y/rho^2
    !    dphi  /dy =  x/rho^2
    !    dphi  /dz =  0
    
    x = xx(1)
    y = xx(2)
    z = xx(3)
    
    rho = sqrt(x**2+y**2)
    r = sqrt(x**2+y**2+z**2)
    
    tt(1,1) =  x/(r+eps)
    tt(2,1) =  y/(r+eps)
    tt(3,1) =  z/(r+eps)
    tt(1,2) =  x*z/(r**2*rho+eps)
    tt(2,2) =  y*z/(r**2*rho+eps)
    tt(3,2) = -rho**2/(r**2*rho+eps)
    tt(1,3) = -y/(rho**2+eps)
    tt(2,3) =  x/(rho**2+eps)
    tt(3,3) =  0
    
  end subroutine make_spherical2cartesian
  
  subroutine make_cartesian2spherical (xx, tt)
    CCTK_REAL, intent(in)  :: xx(3)
    CCTK_REAL, intent(out) :: tt(3,3)
    CCTK_REAL, parameter :: eps = 1.0d-20
    CCTK_REAL :: x,y,z
    CCTK_REAL :: rho,r
    
    ! Cartesian coordinates: x y z
    ! spherical coordinates: r theta phi
    !    r^2        = x^2 + y^2 + z^2
    !    cos(theta) = z/r
    !    tan(phi)   = y/x
    
    ! A[spher](p)_i = TT(p)_i^j A[cart](p)_j
    !    where p = p[cart]
    
    ! Definition:
    !    rho   = r sin(theta)
    !    z/rho = cos(theta) / sin(theta)
    
    ! Transformation:
    !    x = r sin(theta) cos(phi)
    !    y = r sin(theta) sin(phi)
    !    z = r cos(theta)
    !    r     = sqrt(x^2+y^2+z^2)
    !    theta = acos(z/r)
    !    phi   = atan(y/x)
    
    ! Jacobian:
    !    dx/dr     =     sin(theta) cos(phi) =  x/r
    !    dx/dtheta =   r cos(theta) cos(phi) =  x z/rho
    !    dx/dphi   = - r sin(theta) sin(phi) = -y
    !    dy/dr     =     sin(theta) sin(phi) =  y/r
    !    dy/dtheta =   r cos(theta) sin(phi) =  y z/rho
    !    dy/dphi   =   r sin(theta) cos(phi) =  x
    !    dz/dr     =     cos(theta)          =  z/r
    !    dz/dtheta = - r sin(theta)          = -rho
    !    dz/dphi   =   0                     =  0
    
    x = xx(1)
    y = xx(2)
    z = xx(3)
    
    rho = sqrt(x**2+y**2)
    r = sqrt(x**2+y**2+z**2)
    
    tt(1,1) =  x/(r+eps)
    tt(2,1) =  x * z/(rho+eps)
    tt(3,1) = -y
    tt(1,2) =  y/(r+eps)
    tt(2,2) =  y * z/(rho+eps)
    tt(3,2) =  x
    tt(1,3) =  z/(r+eps)
    tt(2,3) = -rho
    tt(3,3) =  0
    
  end subroutine make_cartesian2spherical
  
end module conversion