geometry_spherical.f90
Go to the documentation of this file.
1!#############################################################################
2!# #
3!# fosite - 3D hydrodynamical simulation program #
4!# module: geometry_spherical.f90 #
5!# #
6!# Copyright (C) 2007-2018 #
7!# Tobias Illenseer <tillense@astrophysik.uni-kiel.de> #
8!# Manuel Jung #
9!# Lars Bösch #
10!# Jannes Klee <jklee@astrophysik.uni-kiel.de> #
11!# #
12!# This program is free software; you can redistribute it and/or modify #
13!# it under the terms of the GNU General Public License as published by #
14!# the Free Software Foundation; either version 2 of the License, or (at #
15!# your option) any later version. #
16!# #
17!# This program is distributed in the hope that it will be useful, but #
18!# WITHOUT ANY WARRANTY; without even the implied warranty of #
19!# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
20!# NON INFRINGEMENT. See the GNU General Public License for more #
21!# details. #
22!# #
23!# You should have received a copy of the GNU General Public License #
24!# along with this program; if not, write to the Free Software #
25!# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
26!# #
27!#############################################################################
28
29!----------------------------------------------------------------------------!
41!----------------------------------------------------------------------------!
46 USE common_dict
47 IMPLICIT NONE
48 !--------------------------------------------------------------------------!
50 CONTAINS
52 PROCEDURE :: scalefactors_1
53 PROCEDURE :: scalefactors_2
54 PROCEDURE :: scalefactors_3
55 PROCEDURE :: scalefactors_4
56 PROCEDURE :: radius_1
57 PROCEDURE :: radius_2
58 PROCEDURE :: radius_3
59 PROCEDURE :: radius_4
60 PROCEDURE :: positionvector_1
61 PROCEDURE :: positionvector_2
62 PROCEDURE :: positionvector_3
63 PROCEDURE :: positionvector_4
80 PROCEDURE :: finalize
81 END TYPE
82 PRIVATE
83 CHARACTER(LEN=32), PARAMETER :: geometry_name = "spherical"
84 !--------------------------------------------------------------------------!
85 PUBLIC :: geometry_spherical
86 !--------------------------------------------------------------------------!
87
88CONTAINS
89
90 SUBROUTINE initgeometry_spherical(this,config)
91 IMPLICIT NONE
92 !------------------------------------------------------------------------!
93 CLASS(geometry_spherical), INTENT(INOUT) :: this
94 TYPE(dict_typ),POINTER :: config
95 !------------------------------------------------------------------------!
96 CALL this%InitGeometry(spherical,geometry_name,config)
97 CALL this%SetAzimuthIndex(3)
98 END SUBROUTINE initgeometry_spherical
99
100 PURE SUBROUTINE scalefactors_1(this,coords,hx,hy,hz)
101 IMPLICIT NONE
102 !------------------------------------------------------------------------!
103 CLASS(geometry_spherical), INTENT(IN) :: this
104 REAL, INTENT(IN), DIMENSION(:,:) :: coords
105 REAL, INTENT(OUT), DIMENSION(:) :: hx,hy,hz
106 !------------------------------------------------------------------------!
107 CALL scalefactors(coords(:,1),coords(:,2),hx(:),hy(:),hz(:))
108 END SUBROUTINE scalefactors_1
109
110 PURE SUBROUTINE scalefactors_2(this,coords,hx,hy,hz)
111 IMPLICIT NONE
112 !------------------------------------------------------------------------!
113 CLASS(geometry_spherical), INTENT(IN) :: this
114 REAL, INTENT(IN), DIMENSION(:,:,:) :: coords
115 REAL, INTENT(OUT), DIMENSION(:,:) :: hx,hy,hz
116 !------------------------------------------------------------------------!
117 CALL scalefactors(coords(:,:,1),coords(:,:,2),hx(:,:),hy(:,:),hz(:,:))
118 END SUBROUTINE scalefactors_2
119
120 PURE SUBROUTINE scalefactors_3(this,coords,hx,hy,hz)
121 IMPLICIT NONE
122 !------------------------------------------------------------------------!
123 CLASS(geometry_spherical), INTENT(IN) :: this
124 REAL, INTENT(IN), DIMENSION(:,:,:,:) :: coords
125 REAL, INTENT(OUT), DIMENSION(:,:,:) :: hx,hy,hz
126 !------------------------------------------------------------------------!
127 CALL scalefactors(coords(:,:,:,1),coords(:,:,:,2),hx(:,:,:), &
128 hy(:,:,:),hz(:,:,:))
129 END SUBROUTINE scalefactors_3
130
131 PURE SUBROUTINE scalefactors_4(this,coords,hx,hy,hz)
132 IMPLICIT NONE
133 !------------------------------------------------------------------------!
134 CLASS(geometry_spherical), INTENT(IN) :: this
135 REAL, INTENT(IN), DIMENSION(:,:,:,:,:) :: coords
136 REAL, INTENT(OUT), DIMENSION(:,:,:,:) :: hx,hy,hz
137 !------------------------------------------------------------------------!
138 CALL scalefactors(coords(:,:,:,:,1),coords(:,:,:,:,2),hx(:,:,:,:), &
139 hy(:,:,:,:),hz(:,:,:,:))
140 END SUBROUTINE scalefactors_4
141
142 PURE SUBROUTINE radius_1(this,coords,r)
143 IMPLICIT NONE
144 !------------------------------------------------------------------------!
145 CLASS(geometry_spherical), INTENT(IN) :: this
146 REAL, DIMENSION(:,:), INTENT(IN) :: coords
147 REAL, DIMENSION(:), INTENT(OUT) :: r
148 !------------------------------------------------------------------------!
149 r = radius(coords(:,1))
150 END SUBROUTINE radius_1
151
152 PURE SUBROUTINE radius_2(this,coords,r)
153 IMPLICIT NONE
154 !------------------------------------------------------------------------!
155 CLASS(geometry_spherical), INTENT(IN) :: this
156 REAL, DIMENSION(:,:,:), INTENT(IN) :: coords
157 REAL, DIMENSION(:,:), INTENT(OUT) :: r
158 !------------------------------------------------------------------------!
159 r = radius(coords(:,:,1))
160 END SUBROUTINE radius_2
161
162 PURE SUBROUTINE radius_3(this,coords,r)
163 IMPLICIT NONE
164 !------------------------------------------------------------------------!
165 CLASS(geometry_spherical), INTENT(IN) :: this
166 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: coords
167 REAL, DIMENSION(:,:,:), INTENT(OUT) :: r
168 !------------------------------------------------------------------------!
169 r = radius(coords(:,:,:,1))
170 END SUBROUTINE radius_3
171
172 PURE SUBROUTINE radius_4(this,coords,r)
173 IMPLICIT NONE
174 !------------------------------------------------------------------------!
175 CLASS(geometry_spherical), INTENT(IN) :: this
176 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: coords
177 REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: r
178 !------------------------------------------------------------------------!
179 r = radius(coords(:,:,:,:,1))
180 END SUBROUTINE radius_4
181
182 PURE SUBROUTINE positionvector_1(this,coords,posvec)
183 IMPLICIT NONE
184 !------------------------------------------------------------------------!
185 CLASS(geometry_spherical), INTENT(IN) :: this
186 REAL, DIMENSION(:,:), INTENT(IN) :: coords
187 REAL, DIMENSION(:,:), INTENT(OUT) :: posvec
188 !------------------------------------------------------------------------!
189 CALL positionvector(coords(:,1),posvec(:,1),posvec(:,2),posvec(:,3))
190 END SUBROUTINE positionvector_1
191
192 PURE SUBROUTINE positionvector_2(this,coords,posvec)
193 IMPLICIT NONE
194 !------------------------------------------------------------------------!
195 CLASS(geometry_spherical), INTENT(IN) :: this
196 REAL, DIMENSION(:,:,:), INTENT(IN) :: coords
197 REAL, DIMENSION(:,:,:), INTENT(OUT) :: posvec
198 !------------------------------------------------------------------------!
199 CALL positionvector(coords(:,:,1),posvec(:,:,1),posvec(:,:,2),posvec(:,:,3))
200 END SUBROUTINE positionvector_2
201
202 PURE SUBROUTINE positionvector_3(this,coords,posvec)
203 IMPLICIT NONE
204 !------------------------------------------------------------------------!
205 CLASS(geometry_spherical), INTENT(IN) :: this
206 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: coords
207 REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: posvec
208 !------------------------------------------------------------------------!
209 CALL positionvector(coords(:,:,:,1),posvec(:,:,:,1), &
210 posvec(:,:,:,2),posvec(:,:,:,3))
211 END SUBROUTINE positionvector_3
212
213 PURE SUBROUTINE positionvector_4(this,coords,posvec)
214 IMPLICIT NONE
215 !------------------------------------------------------------------------!
216 CLASS(geometry_spherical), INTENT(IN) :: this
217 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: coords
218 REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: posvec
219 !------------------------------------------------------------------------!
220 CALL positionvector(coords(:,:,:,:,1),posvec(:,:,:,:,1), &
221 posvec(:,:,:,:,2),posvec(:,:,:,:,3))
222 END SUBROUTINE positionvector_4
223
224 PURE SUBROUTINE convert2cartesian_coords_1(this,curv,cart)
225 IMPLICIT NONE
226 !------------------------------------------------------------------------!
227 CLASS(geometry_spherical), INTENT(IN) :: this
228 REAL, DIMENSION(:,:), INTENT(IN) :: curv
229 REAL, DIMENSION(:,:), INTENT(OUT) :: cart
230 !------------------------------------------------------------------------!
231 CALL convert2cartesian_coords(curv(:,1),curv(:,2),curv(:,3), &
232 cart(:,1),cart(:,2),cart(:,3))
233 END SUBROUTINE convert2cartesian_coords_1
234
235 PURE SUBROUTINE convert2cartesian_coords_2(this,curv,cart)
236 IMPLICIT NONE
237 !------------------------------------------------------------------------!
238 CLASS(geometry_spherical), INTENT(IN) :: this
239 REAL, DIMENSION(:,:,:), INTENT(IN) :: curv
240 REAL, DIMENSION(:,:,:), INTENT(OUT) :: cart
241 !------------------------------------------------------------------------!
242 CALL convert2cartesian_coords(curv(:,:,1),curv(:,:,2),curv(:,:,3), &
243 cart(:,:,1),cart(:,:,2),cart(:,:,3))
244 END SUBROUTINE convert2cartesian_coords_2
245
246 PURE SUBROUTINE convert2cartesian_coords_3(this,curv,cart)
247 IMPLICIT NONE
248 !------------------------------------------------------------------------!
249 CLASS(geometry_spherical), INTENT(IN) :: this
250 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: curv
251 REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: cart
252 !------------------------------------------------------------------------!
253 CALL convert2cartesian_coords(curv(:,:,:,1),curv(:,:,:,2),curv(:,:,:,3), &
254 cart(:,:,:,1),cart(:,:,:,2),cart(:,:,:,3))
255 END SUBROUTINE convert2cartesian_coords_3
256
257 PURE SUBROUTINE convert2cartesian_coords_4(this,curv,cart)
258 IMPLICIT NONE
259 !------------------------------------------------------------------------!
260 CLASS(geometry_spherical), INTENT(IN) :: this
261 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: curv
262 REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: cart
263 !------------------------------------------------------------------------!
264 CALL convert2cartesian_coords(curv(:,:,:,:,1),curv(:,:,:,:,2),curv(:,:,:,:,3), &
265 cart(:,:,:,:,1),cart(:,:,:,:,2),cart(:,:,:,:,3))
266 END SUBROUTINE convert2cartesian_coords_4
267
268 PURE SUBROUTINE convert2curvilinear_coords_1(this,cart,curv)
269 IMPLICIT NONE
270 !------------------------------------------------------------------------!
271 CLASS(geometry_spherical), INTENT(IN) :: this
272 REAL, DIMENSION(:,:), INTENT(IN) :: cart
273 REAL, DIMENSION(:,:), INTENT(OUT) :: curv
274 !------------------------------------------------------------------------!
275 CALL convert2curvilinear_coords(cart(:,1),cart(:,2),cart(:,3), &
276 curv(:,1),curv(:,2),curv(:,3))
277 END SUBROUTINE convert2curvilinear_coords_1
278
279 PURE SUBROUTINE convert2curvilinear_coords_2(this,cart,curv)
280 IMPLICIT NONE
281 !------------------------------------------------------------------------!
282 CLASS(geometry_spherical), INTENT(IN) :: this
283 REAL, DIMENSION(:,:,:), INTENT(IN) :: cart
284 REAL, DIMENSION(:,:,:), INTENT(OUT) :: curv
285 !------------------------------------------------------------------------!
286 CALL convert2curvilinear_coords(cart(:,:,1),cart(:,:,2),cart(:,:,3), &
287 curv(:,:,1),curv(:,:,2),curv(:,:,3))
288 END SUBROUTINE convert2curvilinear_coords_2
289
290 PURE SUBROUTINE convert2curvilinear_coords_3(this,cart,curv)
291 IMPLICIT NONE
292 !------------------------------------------------------------------------!
293 CLASS(geometry_spherical), INTENT(IN) :: this
294 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: cart
295 REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: curv
296 !------------------------------------------------------------------------!
297 CALL convert2curvilinear_coords(cart(:,:,:,1),cart(:,:,:,2),cart(:,:,:,3), &
298 curv(:,:,:,1),curv(:,:,:,2),curv(:,:,:,3))
299 END SUBROUTINE convert2curvilinear_coords_3
300
301 PURE SUBROUTINE convert2curvilinear_coords_4(this,cart,curv)
302 IMPLICIT NONE
303 !------------------------------------------------------------------------!
304 CLASS(geometry_spherical), INTENT(IN) :: this
305 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: cart
306 REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: curv
307 !------------------------------------------------------------------------!
308 CALL convert2curvilinear_coords(cart(:,:,:,:,1),cart(:,:,:,:,2),cart(:,:,:,:,3), &
309 curv(:,:,:,:,1),curv(:,:,:,:,2),curv(:,:,:,:,3))
310 END SUBROUTINE convert2curvilinear_coords_4
311
312 PURE SUBROUTINE convert2cartesian_vectors_1(this,curv,v_curv,v_cart)
313 IMPLICIT NONE
314 !------------------------------------------------------------------------!
315 CLASS(geometry_spherical), INTENT(IN) :: this
316 REAL, DIMENSION(:,:), INTENT(IN) :: curv
317 REAL, DIMENSION(:,:), INTENT(IN) :: v_curv
318 REAL, DIMENSION(:,:), INTENT(OUT) :: v_cart
319 !------------------------------------------------------------------------!
320 CALL convert2cartesian_vectors(curv(:,1),curv(:,2),curv(:,3), &
321 v_curv(:,1),v_curv(:,2),v_curv(:,3), &
322 v_cart(:,1),v_cart(:,2),v_cart(:,3))
323 END SUBROUTINE convert2cartesian_vectors_1
324
325 PURE SUBROUTINE convert2cartesian_vectors_2(this,curv,v_curv,v_cart)
326 IMPLICIT NONE
327 !------------------------------------------------------------------------!
328 CLASS(geometry_spherical), INTENT(IN) :: this
329 REAL, DIMENSION(:,:,:), INTENT(IN) :: curv
330 REAL, DIMENSION(:,:,:), INTENT(IN) :: v_curv
331 REAL, DIMENSION(:,:,:), INTENT(OUT) :: v_cart
332 !------------------------------------------------------------------------!
333 CALL convert2cartesian_vectors(curv(:,:,1),curv(:,:,2),curv(:,:,3), &
334 v_curv(:,:,1),v_curv(:,:,2),v_curv(:,:,3), &
335 v_cart(:,:,1),v_cart(:,:,2),v_cart(:,:,3))
336 END SUBROUTINE convert2cartesian_vectors_2
337
338 PURE SUBROUTINE convert2cartesian_vectors_3(this,curv,v_curv,v_cart)
339 IMPLICIT NONE
340 !------------------------------------------------------------------------!
341 CLASS(geometry_spherical), INTENT(IN) :: this
342 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: curv
343 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: v_curv
344 REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: v_cart
345 !------------------------------------------------------------------------!
346 CALL convert2cartesian_vectors(curv(:,:,:,1),curv(:,:,:,2),curv(:,:,:,3), &
347 v_curv(:,:,:,1),v_curv(:,:,:,2),v_curv(:,:,:,3), &
348 v_cart(:,:,:,1),v_cart(:,:,:,2),v_cart(:,:,:,3))
349 END SUBROUTINE convert2cartesian_vectors_3
350
351 PURE SUBROUTINE convert2cartesian_vectors_4(this,curv,v_curv,v_cart)
352 IMPLICIT NONE
353 !------------------------------------------------------------------------!
354 CLASS(geometry_spherical), INTENT(IN) :: this
355 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: curv
356 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: v_curv
357 REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: v_cart
358 !------------------------------------------------------------------------!
359 CALL convert2cartesian_vectors(curv(:,:,:,:,1),curv(:,:,:,:,2),curv(:,:,:,:,3), &
360 v_curv(:,:,:,:,1),v_curv(:,:,:,:,2),v_curv(:,:,:,:,3), &
361 v_cart(:,:,:,:,1),v_cart(:,:,:,:,2),v_cart(:,:,:,:,3))
362 END SUBROUTINE convert2cartesian_vectors_4
363
364 PURE SUBROUTINE convert2curvilinear_vectors_1(this,curv,v_cart,v_curv)
365 IMPLICIT NONE
366 !------------------------------------------------------------------------!
367 CLASS(geometry_spherical), INTENT(IN) :: this
368 REAL, DIMENSION(:,:), INTENT(IN) :: curv
369 REAL, DIMENSION(:,:), INTENT(IN) :: v_cart
370 REAL, DIMENSION(:,:), INTENT(OUT) :: v_curv
371 !------------------------------------------------------------------------!
372 CALL convert2curvilinear_vectors(curv(:,1),curv(:,2),curv(:,3), &
373 v_cart(:,1),v_cart(:,2),v_cart(:,3), &
374 v_curv(:,1),v_curv(:,2),v_curv(:,3))
375 END SUBROUTINE convert2curvilinear_vectors_1
376
377 PURE SUBROUTINE convert2curvilinear_vectors_2(this,curv,v_cart,v_curv)
378 IMPLICIT NONE
379 !------------------------------------------------------------------------!
380 CLASS(geometry_spherical), INTENT(IN) :: this
381 REAL, DIMENSION(:,:,:), INTENT(IN) :: curv
382 REAL, DIMENSION(:,:,:), INTENT(IN) :: v_cart
383 REAL, DIMENSION(:,:,:), INTENT(OUT) :: v_curv
384 !------------------------------------------------------------------------!
385 CALL convert2curvilinear_vectors(curv(:,:,1),curv(:,:,2),curv(:,:,3), &
386 v_cart(:,:,1),v_cart(:,:,2),v_cart(:,:,3), &
387 v_curv(:,:,1),v_curv(:,:,2),v_curv(:,:,3))
388 END SUBROUTINE convert2curvilinear_vectors_2
389
390 PURE SUBROUTINE convert2curvilinear_vectors_3(this,curv,v_cart,v_curv)
391 IMPLICIT NONE
392 !------------------------------------------------------------------------!
393 CLASS(geometry_spherical), INTENT(IN) :: this
394 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: curv
395 REAL, DIMENSION(:,:,:,:), INTENT(IN) :: v_cart
396 REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: v_curv
397 !------------------------------------------------------------------------!
398 CALL convert2curvilinear_vectors(curv(:,:,:,1),curv(:,:,:,2),curv(:,:,:,3), &
399 v_cart(:,:,:,1),v_cart(:,:,:,2),v_cart(:,:,:,3), &
400 v_curv(:,:,:,1),v_curv(:,:,:,2),v_curv(:,:,:,3))
401 END SUBROUTINE convert2curvilinear_vectors_3
402
403 PURE SUBROUTINE convert2curvilinear_vectors_4(this,curv,v_cart,v_curv)
404 IMPLICIT NONE
405 !------------------------------------------------------------------------!
406 CLASS(geometry_spherical), INTENT(IN) :: this
407 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: curv
408 REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: v_cart
409 REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: v_curv
410 !------------------------------------------------------------------------!
411 CALL convert2curvilinear_vectors(curv(:,:,:,:,1),curv(:,:,:,:,2),curv(:,:,:,:,3), &
412 v_cart(:,:,:,:,1),v_cart(:,:,:,:,2),v_cart(:,:,:,:,3), &
413 v_curv(:,:,:,:,1),v_curv(:,:,:,:,2),v_curv(:,:,:,:,3))
414 END SUBROUTINE convert2curvilinear_vectors_4
415
416 SUBROUTINE finalize(this)
417 IMPLICIT NONE
418 !------------------------------------------------------------------------!
419 CLASS(geometry_spherical), INTENT(INOUT) :: this
420 !------------------------------------------------------------------------!
421 CALL this%Finalize_base()
422 END SUBROUTINE finalize
423
424 ELEMENTAL SUBROUTINE scalefactors(r,theta,hr,htheta,hphi)
425 IMPLICIT NONE
426 !------------------------------------------------------------------------!
427 REAL, INTENT(IN) :: r,theta
428 REAL, INTENT(OUT) :: hr,htheta,hphi
429 !------------------------------------------------------------------------!
430 hr = 1.
431 htheta = r
432 hphi = r * sin(theta)
433 END SUBROUTINE scalefactors
434
435 ELEMENTAL FUNCTION radius(r)
436 IMPLICIT NONE
437 !------------------------------------------------------------------------!
438 REAL, INTENT(IN) :: r
439 REAL :: radius
440 !------------------------------------------------------------------------!
441 radius = r
442 END FUNCTION radius
443
444 ELEMENTAL SUBROUTINE positionvector(r,rx,ry,rz)
445 IMPLICIT NONE
446 !------------------------------------------------------------------------!
447 REAL, INTENT(IN) :: r
448 REAL, INTENT(OUT) :: rx,ry,rz
449 !------------------------------------------------------------------------!
450 rx = r
451 ry = 0.0
452 rz = 0.0
453 END SUBROUTINE positionvector
454
455 ELEMENTAL SUBROUTINE convert2cartesian_coords(r,theta,phi,x,y,z)
456 IMPLICIT NONE
457 !------------------------------------------------------------------------!
458 REAL, INTENT(IN) :: r,theta,phi
459 REAL, INTENT(OUT) :: x,y,z
460 !------------------------------------------------------------------------!
461 x = r*sin(theta)*cos(phi)
462 y = r*sin(theta)*sin(phi)
463 z = r*cos(theta)
464 END SUBROUTINE convert2cartesian_coords
465
466 ELEMENTAL SUBROUTINE convert2curvilinear_coords(x,y,z,r,theta,phi)
467 IMPLICIT NONE
468 !------------------------------------------------------------------------!
469 REAL, INTENT(IN) :: x,y,z
470 REAL, INTENT(OUT) :: r,theta,phi
471 !------------------------------------------------------------------------!
472 r = sqrt(x*x+y*y+z*z)
473 theta = acos(z/r)
474 phi = atan2(y,x)
475 IF (phi.LT.0.0) THEN
476 phi = phi + 2.0*pi
477 END IF
478 END SUBROUTINE convert2curvilinear_coords
479
481 ELEMENTAL SUBROUTINE convert2cartesian_vectors(r,theta,phi,vr,vtheta,vphi,vx,vy,vz)
482 IMPLICIT NONE
483 !------------------------------------------------------------------------!
484 REAL, INTENT(IN) :: r,theta,phi,vr,vtheta,vphi
485 REAL, INTENT(OUT) :: vx,vy,vz
486 !------------------------------------------------------------------------!
487 vx = vr*sin(theta)*cos(phi) + vtheta*cos(theta)*cos(phi) - vphi*sin(phi)
488 vy = vr*sin(theta)*sin(phi) + vtheta*cos(theta)*sin(phi) + vphi*cos(phi)
489 vz = vr*cos(theta) - vtheta*sin(theta)
490 END SUBROUTINE convert2cartesian_vectors
491
493 ELEMENTAL SUBROUTINE convert2curvilinear_vectors(r,theta,phi,vx,vy,vz,vr,vtheta,vphi)
494 IMPLICIT NONE
495 !------------------------------------------------------------------------!
496 REAL, INTENT(IN) :: r,theta,phi,vx,vy,vz
497 REAL, INTENT(OUT) :: vr,vtheta,vphi
498 !------------------------------------------------------------------------!
499 vr = vx*sin(theta)*cos(phi) + vy*sin(theta)*sin(phi) + vz*cos(theta)
500 vtheta = vx*cos(theta)*cos(phi) + vy*cos(theta)*sin(phi) - vz*sin(theta)
501 vphi = -vx*sin(phi) + vy*cos(phi)
502 END SUBROUTINE convert2curvilinear_vectors
503
504END MODULE geometry_spherical_mod
Dictionary for generic data types.
Definition: common_dict.f90:61
type(logging_base), save this
base class for geometrical properties
real, parameter, public pi
integer, parameter, public spherical
pure subroutine convert2cartesian_vectors(this, curv, v_curv, v_cart)
Convert curvilinear vector components to cartesian vector components.
pure subroutine convert2curvilinear_coords(this, cart, curv)
Convert cartesian to curvilinear coordinates.
pure subroutine convert2curvilinear_vectors(this, curv, v_cart, v_curv)
Convert cartesian vector components to curvilinear vector components.
pure subroutine convert2cartesian_coords(this, curv, cart)
Convert curvilinear to cartesian coordinates.
defines properties of a 3D spherical mesh
pure subroutine scalefactors_4(this, coords, hx, hy, hz)
pure subroutine convert2curvilinear_coords_4(this, cart, curv)
elemental subroutine positionvector(r, rx, ry, rz)
pure subroutine positionvector_1(this, coords, posvec)
pure subroutine convert2cartesian_vectors_3(this, curv, v_curv, v_cart)
pure subroutine convert2curvilinear_coords_2(this, cart, curv)
pure subroutine convert2curvilinear_vectors_2(this, curv, v_cart, v_curv)
pure subroutine positionvector_3(this, coords, posvec)
pure subroutine convert2cartesian_coords_2(this, curv, cart)
pure subroutine convert2curvilinear_vectors_4(this, curv, v_cart, v_curv)
pure subroutine convert2curvilinear_vectors_1(this, curv, v_cart, v_curv)
pure subroutine positionvector_4(this, coords, posvec)
pure subroutine convert2cartesian_vectors_2(this, curv, v_curv, v_cart)
pure subroutine convert2curvilinear_vectors_3(this, curv, v_cart, v_curv)
pure subroutine radius_3(this, coords, r)
pure subroutine radius_4(this, coords, r)
subroutine initgeometry_spherical(this, config)
pure subroutine convert2cartesian_coords_4(this, curv, cart)
character(len=32), parameter geometry_name
pure subroutine convert2curvilinear_coords_1(this, cart, curv)
pure subroutine scalefactors_1(this, coords, hx, hy, hz)
pure subroutine convert2cartesian_coords_3(this, curv, cart)
pure subroutine scalefactors_3(this, coords, hx, hy, hz)
pure subroutine positionvector_2(this, coords, posvec)
elemental subroutine scalefactors(r, theta, hr, htheta, hphi)
pure subroutine scalefactors_2(this, coords, hx, hy, hz)
pure subroutine convert2cartesian_vectors_4(this, curv, v_curv, v_cart)
pure subroutine radius_1(this, coords, r)
pure subroutine convert2curvilinear_coords_3(this, cart, curv)
pure subroutine radius_2(this, coords, r)
elemental real function radius(r)
pure subroutine convert2cartesian_vectors_1(this, curv, v_curv, v_cart)
pure subroutine convert2cartesian_coords_1(this, curv, cart)
derived mesh array class for scalar cell data
derived mesh array class for vector cell data