geometry_base.f90
Go to the documentation of this file.
1 !#############################################################################
2 !# #
3 !# fosite - 3D hydrodynamical simulation program #
4 !# module: geometry_generic.f90 #
5 !# #
6 !# Copyright (C) 2007-2018 #
7 !# Tobias Illenseer <tillense@astrophysik.uni-kiel.de> #
8 !# Jannes Klee <jklee@astrophysik.uni-kiel.de> #
9 !# #
10 !# This program is free software; you can redistribute it and/or modify #
11 !# it under the terms of the GNU General Public License as published by #
12 !# the Free Software Foundation; either version 2 of the License, or (at #
13 !# your option) any later version. #
14 !# #
15 !# This program is distributed in the hope that it will be useful, but #
16 !# WITHOUT ANY WARRANTY; without even the implied warranty of #
17 !# MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, GOOD TITLE or #
18 !# NON INFRINGEMENT. See the GNU General Public License for more #
19 !# details. #
20 !# #
21 !# You should have received a copy of the GNU General Public License #
22 !# along with this program; if not, write to the Free Software #
23 !# Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA. #
24 !# #
25 !#############################################################################
26 
27 !----------------------------------------------------------------------------!
35 !----------------------------------------------------------------------------!
40  USE common_dict
41  !--------------------------------------------------------------------------!
42  PRIVATE
43  REAL, PARAMETER :: pi = 3.1415926535897932384626433832795028842
44 ! INTERFACE GetScale
45 ! MODULE PROCEDURE GetScale1, GetScale2
46 ! END INTERFACE
47 ! INTERFACE SetScale
48 ! MODULE PROCEDURE SetScale1, SetScale2, SetScale3
49 ! END INTERFACE
50 
51  TYPE, ABSTRACT, EXTENDS (logging_base) :: geometry_base
52  REAL,DIMENSION(3) :: geoparam
53  INTEGER, PRIVATE :: azimuthindex = 0
54  CONTAINS
55 
56 ! PRIVATE
57  PROCEDURE :: scalefactors_0
58  procedure(scalefactors_1), DEFERRED :: scalefactors_1
59  procedure(scalefactors_2), DEFERRED :: scalefactors_2
60  procedure(scalefactors_3), DEFERRED :: scalefactors_3
61  procedure(scalefactors_4), DEFERRED :: scalefactors_4
62  PROCEDURE :: radius_0
63  procedure(radius_1), DEFERRED :: radius_1
64  procedure(radius_2), DEFERRED :: radius_2
65  procedure(radius_3), DEFERRED :: radius_3
66  procedure(radius_4), DEFERRED :: radius_4
67  PROCEDURE :: positionvector_0
68  procedure(positionvector_1), DEFERRED :: positionvector_1
69  procedure(positionvector_2), DEFERRED :: positionvector_2
70  procedure(positionvector_3), DEFERRED :: positionvector_3
71  procedure(positionvector_4), DEFERRED :: positionvector_4
73  procedure(convert2cartesian_coords_1), DEFERRED :: convert2cartesian_coords_1
74  procedure(convert2cartesian_coords_2), DEFERRED :: convert2cartesian_coords_2
75  procedure(convert2cartesian_coords_3), DEFERRED :: convert2cartesian_coords_3
76  procedure(convert2cartesian_coords_4), DEFERRED :: convert2cartesian_coords_4
78  procedure(convert2curvilinear_coords_1), DEFERRED :: convert2curvilinear_coords_1
79  procedure(convert2curvilinear_coords_2), DEFERRED :: convert2curvilinear_coords_2
80  procedure(convert2curvilinear_coords_3), DEFERRED :: convert2curvilinear_coords_3
81  procedure(convert2curvilinear_coords_4), DEFERRED :: convert2curvilinear_coords_4
83  procedure(convert2cartesian_vectors_1), DEFERRED :: convert2cartesian_vectors_1
84  procedure(convert2cartesian_vectors_2), DEFERRED :: convert2cartesian_vectors_2
85  procedure(convert2cartesian_vectors_3), DEFERRED :: convert2cartesian_vectors_3
86  procedure(convert2cartesian_vectors_4), DEFERRED :: convert2cartesian_vectors_4
88  procedure(convert2curvilinear_vectors_1), DEFERRED :: convert2curvilinear_vectors_1
89  procedure(convert2curvilinear_vectors_2), DEFERRED :: convert2curvilinear_vectors_2
90  procedure(convert2curvilinear_vectors_3), DEFERRED :: convert2curvilinear_vectors_3
91  procedure(convert2curvilinear_vectors_4), DEFERRED :: convert2curvilinear_vectors_4
92  PROCEDURE :: setscale1
93  PROCEDURE :: setscale2
94  PROCEDURE :: setscale3
95  PROCEDURE :: getscale1
96  PROCEDURE :: getscale2
97  PROCEDURE :: setazimuthindex
98  PROCEDURE :: getazimuthindex
99 
100 ! PUBLIC
101  PROCEDURE, PUBLIC :: initgeometry
102  PROCEDURE, PUBLIC :: finalize_base
103  procedure(finalize), DEFERRED :: finalize
104 
106  generic, PUBLIC :: scalefactors => scalefactors_0, scalefactors_1, scalefactors_2, &
107  scalefactors_3, scalefactors_4
109  generic, PUBLIC :: radius => radius_0, radius_1, radius_2, radius_3, radius_4
111  generic, PUBLIC :: positionvector => positionvector_0, positionvector_1, positionvector_2, &
112  positionvector_3, positionvector_4
114  generic, PUBLIC :: convert2cartesian => convert2cartesian_coords, convert2cartesian_vectors, &
115  convert2cartesian_coords_1, convert2cartesian_coords_2,&
116  convert2cartesian_coords_3, convert2cartesian_coords_4, &
117  convert2cartesian_vectors_1, convert2cartesian_vectors_2,&
118  convert2cartesian_vectors_3, convert2cartesian_vectors_4
120  generic, PUBLIC :: convert2curvilinear => convert2curvilinear_coords, convert2curvilinear_vectors, &
121  convert2curvilinear_coords_1, convert2curvilinear_coords_2,&
122  convert2curvilinear_coords_3, convert2curvilinear_coords_4,&
123  convert2curvilinear_vectors_1, convert2curvilinear_vectors_2,&
124  convert2curvilinear_vectors_3, convert2curvilinear_vectors_4
125  generic :: setscale => setscale1, setscale2, setscale3
126  generic :: getscale => getscale1, getscale2
127  END TYPE geometry_base
128 
129  abstract INTERFACE
130  PURE SUBROUTINE scalefactors_1(this,coords,hx,hy,hz)
132  IMPLICIT NONE
133  CLASS(geometry_base), INTENT(IN) :: this
134  REAL, DIMENSION(:,:), INTENT(IN) :: coords
135  REAL, DIMENSION(:), INTENT(OUT) :: hx,hy,hz
136  END SUBROUTINE
137  PURE SUBROUTINE scalefactors_2(this,coords,hx,hy,hz)
139  IMPLICIT NONE
140  CLASS(geometry_base), INTENT(IN) :: this
141  REAL, DIMENSION(:,:,:), INTENT(IN) :: coords
142  REAL, DIMENSION(:,:), INTENT(OUT) :: hx,hy,hz
143  END SUBROUTINE
144  PURE SUBROUTINE scalefactors_3(this,coords,hx,hy,hz)
146  IMPLICIT NONE
147  CLASS(geometry_base), INTENT(IN) :: this
148  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: coords
149  REAL, DIMENSION(:,:,:), INTENT(OUT) :: hx,hy,hz
150  END SUBROUTINE
151  PURE SUBROUTINE scalefactors_4(this,coords,hx,hy,hz)
153  IMPLICIT NONE
154  CLASS(geometry_base), INTENT(IN) :: this
155  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: coords
156  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: hx,hy,hz
157  END SUBROUTINE
158  PURE SUBROUTINE radius_1(this,coords,r)
160  IMPLICIT NONE
161  CLASS(geometry_base), INTENT(IN) :: this
162  REAL, DIMENSION(:,:), INTENT(IN) :: coords
163  REAL, DIMENSION(:), INTENT(OUT) :: r
164  END SUBROUTINE
165  PURE SUBROUTINE radius_2(this,coords,r)
167  IMPLICIT NONE
168  CLASS(geometry_base), INTENT(IN) :: this
169  REAL, DIMENSION(:,:,:), INTENT(IN) :: coords
170  REAL, DIMENSION(:,:), INTENT(OUT) :: r
171  END SUBROUTINE
172  PURE SUBROUTINE radius_3(this,coords,r)
174  IMPLICIT NONE
175  CLASS(geometry_base), INTENT(IN) :: this
176  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: coords
177  REAL, DIMENSION(:,:,:), INTENT(OUT) :: r
178  END SUBROUTINE
179  PURE SUBROUTINE radius_4(this,coords,r)
181  IMPLICIT NONE
182  CLASS(geometry_base), INTENT(IN) :: this
183  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: coords
184  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: r
185  END SUBROUTINE
186  PURE SUBROUTINE positionvector_1(this,coords,posvec)
188  IMPLICIT NONE
189  CLASS(geometry_base), INTENT(IN) :: this
190  REAL, DIMENSION(:,:), INTENT(IN) :: coords
191  REAL, DIMENSION(:,:), INTENT(OUT) :: posvec
192  END SUBROUTINE
193  PURE SUBROUTINE positionvector_2(this,coords,posvec)
195  IMPLICIT NONE
196  CLASS(geometry_base), INTENT(IN) :: this
197  REAL, DIMENSION(:,:,:), INTENT(IN) :: coords
198  REAL, DIMENSION(:,:,:), INTENT(OUT) :: posvec
199  END SUBROUTINE
200  PURE SUBROUTINE positionvector_3(this,coords,posvec)
202  IMPLICIT NONE
203  CLASS(geometry_base), INTENT(IN) :: this
204  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: coords
205  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: posvec
206  END SUBROUTINE
207  PURE SUBROUTINE positionvector_4(this,coords,posvec)
209  IMPLICIT NONE
210  CLASS(geometry_base), INTENT(IN) :: this
211  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: coords
212  REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: posvec
213  END SUBROUTINE
214  PURE SUBROUTINE convert2cartesian_coords_1(this,curv,cart)
216  IMPLICIT NONE
217  CLASS(geometry_base), INTENT(IN) :: this
218  REAL, DIMENSION(:,:), INTENT(IN) :: curv
219  REAL, DIMENSION(:,:), INTENT(OUT) :: cart
220  END SUBROUTINE
221  PURE SUBROUTINE convert2cartesian_coords_2(this,curv,cart)
223  IMPLICIT NONE
224  CLASS(geometry_base), INTENT(IN) :: this
225  REAL, DIMENSION(:,:,:), INTENT(IN) :: curv
226  REAL, DIMENSION(:,:,:), INTENT(OUT) :: cart
227  END SUBROUTINE
228  PURE SUBROUTINE convert2cartesian_coords_3(this,curv,cart)
230  IMPLICIT NONE
231  CLASS(geometry_base), INTENT(IN) :: this
232  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: curv
233  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: cart
234  END SUBROUTINE
235  PURE SUBROUTINE convert2cartesian_coords_4(this,curv,cart)
237  IMPLICIT NONE
238  CLASS(geometry_base), INTENT(IN) :: this
239  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: curv
240  REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: cart
241  END SUBROUTINE
242  PURE SUBROUTINE convert2curvilinear_coords_1(this,cart,curv)
244  IMPLICIT NONE
245  CLASS(geometry_base), INTENT(IN) :: this
246  REAL, DIMENSION(:,:), INTENT(IN) :: cart
247  REAL, DIMENSION(:,:), INTENT(OUT) :: curv
248  END SUBROUTINE
249  PURE SUBROUTINE convert2curvilinear_coords_2(this,cart,curv)
251  IMPLICIT NONE
252  CLASS(geometry_base), INTENT(IN) :: this
253  REAL, DIMENSION(:,:,:), INTENT(IN) :: cart
254  REAL, DIMENSION(:,:,:), INTENT(OUT) :: curv
255  END SUBROUTINE
256  PURE SUBROUTINE convert2curvilinear_coords_3(this,cart,curv)
258  IMPLICIT NONE
259  CLASS(geometry_base), INTENT(IN) :: this
260  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: cart
261  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: curv
262  END SUBROUTINE
263  PURE SUBROUTINE convert2curvilinear_coords_4(this,cart,curv)
265  IMPLICIT NONE
266  CLASS(geometry_base), INTENT(IN) :: this
267  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: cart
268  REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: curv
269  END SUBROUTINE
270  PURE SUBROUTINE convert2cartesian_vectors_1(this,curv,v_curv,v_cart)
272  IMPLICIT NONE
273  CLASS(geometry_base), INTENT(IN) :: this
274  REAL, DIMENSION(:,:), INTENT(IN) :: curv
275  REAL, DIMENSION(:,:), INTENT(IN) :: v_curv
276  REAL, DIMENSION(:,:), INTENT(OUT) :: v_cart
277  END SUBROUTINE
278  PURE SUBROUTINE convert2cartesian_vectors_2(this,curv,v_curv,v_cart)
280  IMPLICIT NONE
281  CLASS(geometry_base), INTENT(IN) :: this
282  REAL, DIMENSION(:,:,:), INTENT(IN) :: curv
283  REAL, DIMENSION(:,:,:), INTENT(IN) :: v_curv
284  REAL, DIMENSION(:,:,:), INTENT(OUT) :: v_cart
285  END SUBROUTINE
286  PURE SUBROUTINE convert2cartesian_vectors_3(this,curv,v_curv,v_cart)
288  IMPLICIT NONE
289  CLASS(geometry_base), INTENT(IN) :: this
290  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: curv
291  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: v_curv
292  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: v_cart
293  END SUBROUTINE
294  PURE SUBROUTINE convert2cartesian_vectors_4(this,curv,v_curv,v_cart)
296  IMPLICIT NONE
297  CLASS(geometry_base), INTENT(IN) :: this
298  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: curv
299  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: v_curv
300  REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: v_cart
301  END SUBROUTINE
302  PURE SUBROUTINE convert2curvilinear_vectors_1(this,curv,v_cart,v_curv)
304  IMPLICIT NONE
305  CLASS(geometry_base), INTENT(IN) :: this
306  REAL, DIMENSION(:,:), INTENT(IN) :: curv
307  REAL, DIMENSION(:,:), INTENT(IN) :: v_cart
308  REAL, DIMENSION(:,:), INTENT(OUT) :: v_curv
309  END SUBROUTINE
310  PURE SUBROUTINE convert2curvilinear_vectors_2(this,curv,v_cart,v_curv)
312  IMPLICIT NONE
313  CLASS(geometry_base), INTENT(IN) :: this
314  REAL, DIMENSION(:,:,:), INTENT(IN) :: curv
315  REAL, DIMENSION(:,:,:), INTENT(IN) :: v_cart
316  REAL, DIMENSION(:,:,:), INTENT(OUT) :: v_curv
317  END SUBROUTINE
318  PURE SUBROUTINE convert2curvilinear_vectors_3(this,curv,v_cart,v_curv)
320  IMPLICIT NONE
321  CLASS(geometry_base), INTENT(IN) :: this
322  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: curv
323  REAL, DIMENSION(:,:,:,:), INTENT(IN) :: v_cart
324  REAL, DIMENSION(:,:,:,:), INTENT(OUT) :: v_curv
325  END SUBROUTINE
326  PURE SUBROUTINE convert2curvilinear_vectors_4(this,curv,v_cart,v_curv)
328  IMPLICIT NONE
329  CLASS(geometry_base), INTENT(IN) :: this
330  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: curv
331  REAL, DIMENSION(:,:,:,:,:), INTENT(IN) :: v_cart
332  REAL, DIMENSION(:,:,:,:,:), INTENT(OUT) :: v_curv
333  END SUBROUTINE
334  SUBROUTINE finalize(this)
336  IMPLICIT NONE
337  CLASS(geometry_base), INTENT(INOUT) :: this
338  END SUBROUTINE
339  END INTERFACE
340 
343  !--------------------------------------------------------------------------!
344  INTEGER, PARAMETER :: cartesian = 1
345 ! INTEGER, PARAMETER :: SINHCARTESIAN = 2
346 ! INTEGER, PARAMETER :: POLAR = 20
347 ! INTEGER, PARAMETER :: LOGPOLAR = 21
348 ! INTEGER, PARAMETER :: TANPOLAR = 22
349 ! INTEGER, PARAMETER :: SINHPOLAR = 23
350 ! INTEGER, PARAMETER :: SINHTANHPOLAR = 24
351 ! INTEGER, PARAMETER :: POLYPOLAR = 25
352 ! INTEGER, PARAMETER :: ELLIPTIC = 27
353 ! INTEGER, PARAMETER :: BIPOLAR = 29
354  INTEGER, PARAMETER :: cylindrical = 30
355  INTEGER, PARAMETER :: logcylindrical = 31
356 ! INTEGER, PARAMETER :: TANCYLINDRICAL = 32
357 ! INTEGER, PARAMETER :: LNCOSHCYLINDRICAL = 33
358  INTEGER, PARAMETER :: spherical = 40
359  INTEGER, PARAMETER :: logspherical = 41
360 ! INTEGER, PARAMETER :: SINHSPHERICAL = 42
361 ! INTEGER, PARAMETER :: BIANGLESPHERICAL = 43
362 ! INTEGER, PARAMETER :: OBLATE_SPHEROIDAL = 50
363 ! INTEGER, PARAMETER :: CHANNEL = 60
365  !--------------------------------------------------------------------------!
366  PUBLIC :: &
367  ! types
368  geometry_base, &
369  pi, &
370  cartesian, &
373  !--------------------------------------------------------------------------!
374 
375 CONTAINS
376 
378  SUBROUTINE initgeometry(this,gnum,gname,config)
379  IMPLICIT NONE
380  !------------------------------------------------------------------------!
381  CLASS(geometry_base), INTENT(INOUT) :: this
382  TYPE(DICT_TYP),POINTER :: config
383  INTEGER :: gnum, gt !\todo{same variable}
384  CHARACTER(LEN=*) :: gname
385  REAL :: gs_def,gs_def2,gs_def3
386  !------------------------------------------------------------------------!
387  CHARACTER(LEN=8) :: gs_str
388  !------------------------------------------------------------------------!
389  INTENT(IN) :: gnum,gname
390  !------------------------------------------------------------------------!
391  CALL this%logging_base%InitLogging(gnum,gname)
392 
393  CALL getattr(config, "geometry", gt)
394 
395  ! check if geometry parameters were given
396  ! and set to defaults if not
397  SELECT CASE(gnum)
399  ! do nothing (no parameters needed)
400  CASE DEFAULT
401  ! geometries with at least one parameter
402  ! gparam defaults to 1 except for CHANNEL
403 ! IF (gt.EQ.CHANNEL) THEN
404 ! CALL GetAttr(config, "gparam", gs_def, 0.5)
405 ! ELSE
406  CALL getattr(config, "gparam", gs_def, 1.0)
407 ! END IF
408  ! geometries with three parameters
409 ! SELECT CASE(gt)
410 ! CASE(SINHPOLAR)
411 ! ! gp2,gp3 defaults to 0
412 ! CALL GetAttr(config, "gparam2", gs_def2, 0.0)
413 ! CALL GetAttr(config, "gparam3", gs_def3, 0.0)
414 ! CASE(SINHTANHPOLAR,LNCOSHCYLINDRICAL)
415 ! ! gp2,gp3 defaults to 1
416 ! CALL GetAttr(config, "gparam2", gs_def2, 1.0)
417 ! CALL GetAttr(config, "gparam3", gs_def3, 1.0)
418 ! END SELECT
419  END SELECT
420 
421  ! print some information
422  CALL this%Info(" GEOMETRY-> coordinates: " // trim(this%GetName()))
423  SELECT CASE(gt)
425  WRITE (gs_str,'(ES8.1)') this%GetScale()
426  CALL this%Info( " geometry scale: " // trim(gs_str))
427 ! CASE(LNCOSHCYLINDRICAL,SINHTANHPOLAR)
428 ! WRITE (gs_str,'(ES8.1)') this%GetScale()
429 ! CALL this%Info( " geometry scale: " // TRIM(gs_str))
430 ! WRITE (gs_str,'(ES8.1)') this%GetScale(2)
431 ! CALL this%Info(" geometry scale: " // TRIM(gs_str))
432 ! WRITE (gs_str,'(ES8.1)') this%GetScale(3)
433 ! CALL this%Info(" geometry scale: " // TRIM(gs_str))
434 
435  END SELECT
436  END SUBROUTINE initgeometry
437 
438  PURE FUNCTION getscale1(this) RESULT(gp)
439  IMPLICIT NONE
440  !------------------------------------------------------------------------!
441  CLASS(geometry_base), INTENT(IN) :: this
442  REAL :: gp
443  !------------------------------------------------------------------------!
444  gp = this%geoparam(1)
445  END FUNCTION getscale1
446 
447  PURE FUNCTION getscale2(this,i) RESULT(gp)
448  IMPLICIT NONE
449  !------------------------------------------------------------------------!
450  CLASS(geometry_base), INTENT(IN) :: this
451  INTEGER, INTENT(IN) :: i
452  REAL :: gp
453  !------------------------------------------------------------------------!
454  gp = this%geoparam(i)
455  END FUNCTION getscale2
456 
457  PURE SUBROUTINE setscale1(this,gp)
458  IMPLICIT NONE
459  !------------------------------------------------------------------------!
460  CLASS(geometry_base), INTENT(INOUT) :: this
461  REAL, INTENT(IN) :: gp
462  !------------------------------------------------------------------------!
463  this%geoparam(1) = gp
464  END SUBROUTINE setscale1
465 
466  PURE SUBROUTINE setscale2(this,gp,gp2)
467  IMPLICIT NONE
468  !------------------------------------------------------------------------!
469  CLASS(geometry_base), INTENT(INOUT) :: this
470  REAL, INTENT(IN) :: gp,gp2
471  !------------------------------------------------------------------------!
472  this%geoparam(1) = gp
473  this%geoparam(2) = gp2
474  END SUBROUTINE setscale2
475 
476  PURE SUBROUTINE setscale3(this,gp,gp2,gp3)
477  IMPLICIT NONE
478  !------------------------------------------------------------------------!
479  CLASS(geometry_base), INTENT(INOUT) :: this
480  REAL, INTENT(IN) :: gp,gp2,gp3
481  !------------------------------------------------------------------------!
482  this%geoparam(1) = gp
483  this%geoparam(2) = gp2
484  this%geoparam(3) = gp3
485  END SUBROUTINE setscale3
486 
487 
498  PURE SUBROUTINE scalefactors_0(this,coords,hx,hy,hz)
499  IMPLICIT NONE
500  !------------------------------------------------------------------------!
501  CLASS(geometry_base), INTENT(IN) :: this
502  TYPE(marray_cellvector), INTENT(IN) :: coords
503  TYPE(marray_cellscalar), INTENT(INOUT) :: hx,hy,hz
504  !------------------------------------------------------------------------!
505  CALL this%ScaleFactors(coords%data2d(:,:),hx%data1d(:),hy%data1d(:),hz%data1d(:))
506  END SUBROUTINE scalefactors_0
507 
514  PURE SUBROUTINE radius_0(this,coords,radius)
515  IMPLICIT NONE
516  !------------------------------------------------------------------------!
517  CLASS(geometry_base), INTENT(IN) :: this
518  TYPE(marray_cellvector), INTENT(IN) :: coords
519  TYPE(marray_cellscalar), INTENT(INOUT) :: radius
520  !------------------------------------------------------------------------!
521  CALL this%Radius(coords%data2d(:,:),radius%data1d(:))
522  END SUBROUTINE radius_0
523 
530  PURE SUBROUTINE positionvector_0(this,coords,posvec)
531  IMPLICIT NONE
532  !------------------------------------------------------------------------!
533  CLASS(geometry_base), INTENT(IN) :: this
534  TYPE(marray_cellvector), INTENT(IN) :: coords
535  TYPE(marray_cellvector), INTENT(INOUT) :: posvec
536  !------------------------------------------------------------------------!
537  CALL this%PositionVector(coords%data2d(:,:),posvec%data2d(:,:))
538  END SUBROUTINE positionvector_0
539 
544  PURE SUBROUTINE convert2cartesian_coords(this,curv,cart)
545  IMPLICIT NONE
546  !------------------------------------------------------------------------!
547  CLASS(geometry_base), INTENT(IN) :: this
548  TYPE(marray_cellvector), INTENT(IN) :: curv
549  TYPE(marray_cellvector), INTENT(INOUT) :: cart
550  !------------------------------------------------------------------------!
551  CALL this%Convert2Cartesian(curv%data2d(:,:),cart%data2d(:,:))
552  END SUBROUTINE convert2cartesian_coords
553 
558  PURE SUBROUTINE convert2curvilinear_coords(this,cart,curv)
559  IMPLICIT NONE
560  !------------------------------------------------------------------------!
561  CLASS(geometry_base), INTENT(IN) :: this
562  TYPE(marray_cellvector), INTENT(IN) :: cart
563  TYPE(marray_cellvector), INTENT(INOUT) :: curv
564  !------------------------------------------------------------------------!
565  CALL this%Convert2Curvilinear(cart%data2d(:,:),curv%data2d(:,:))
566  END SUBROUTINE convert2curvilinear_coords
567 
572  PURE SUBROUTINE convert2cartesian_vectors(this,curv,v_curv,v_cart)
573  IMPLICIT NONE
574  !------------------------------------------------------------------------!
575  CLASS(geometry_base), INTENT(IN) :: this
576  TYPE(marray_cellvector), INTENT(IN) :: curv
577  TYPE(marray_cellvector), INTENT(IN) :: v_curv
578  TYPE(marray_cellvector), INTENT(INOUT) :: v_cart
579  !------------------------------------------------------------------------!
580  CALL this%Convert2Cartesian(curv%data2d(:,:),v_curv%data2d(:,:),v_cart%data2d(:,:))
581  END SUBROUTINE convert2cartesian_vectors
582 
583 
585  PURE SUBROUTINE convert2curvilinear_vectors(this,curv,v_cart,v_curv)
586  IMPLICIT NONE
587  !------------------------------------------------------------------------!
588  CLASS(geometry_base), INTENT(IN) :: this
589  TYPE(marray_cellvector), INTENT(IN) :: curv
590  TYPE(marray_cellvector), INTENT(IN) :: v_cart
591  TYPE(marray_cellvector), INTENT(INOUT) :: v_curv
592  !------------------------------------------------------------------------!
593  CALL this%Convert2Curvilinear(curv%data2d(:,:),v_cart%data2d(:,:),v_curv%data2d(:,:))
594  END SUBROUTINE convert2curvilinear_vectors
595 
599  PURE SUBROUTINE setazimuthindex(this,idx)
600  IMPLICIT NONE
601  !------------------------------------------------------------------------!
602  CLASS(geometry_base), INTENT(INOUT) :: this
603  INTEGER, INTENT(IN) :: idx
604  !------------------------------------------------------------------------!
605  this%azimuthIndex = idx
606  END SUBROUTINE setazimuthindex
607 
611  PURE FUNCTION getazimuthindex(this) RESULT(idx)
612  IMPLICIT NONE
613  !------------------------------------------------------------------------!
614  CLASS(geometry_base), INTENT(IN) :: this
615  INTEGER :: idx
616  !------------------------------------------------------------------------!
617  idx = this%azimuthIndex
618  END FUNCTION getazimuthindex
619 
621  SUBROUTINE finalize_base(this)
622  IMPLICIT NONE
623  !------------------------------------------------------------------------!
624  CLASS(geometry_base), INTENT(INOUT) :: this
625  !------------------------------------------------------------------------!
626  IF (.NOT.this%Initialized()) &
627  CALL this%Error("CloseGeometry","not initialized")
628  !CALL this%logging_base%Finalize()
629  END SUBROUTINE finalize_base
630 
631 END MODULE geometry_base_mod
subroutine finalize_base(this)
Destructor of generic geometry module.
subroutine finalize(this)
Destructor of common class.
pure subroutine scalefactors_0(this, coords, hx, hy, hz)
Compute scale factors.
pure subroutine setazimuthindex(this, idx)
sets the coordinate index of the azimuthal angle
type(logging_base), save this
subroutine initgeometry(this, gnum, gname, config)
Constructor of generic geometry module.
derived mesh array class for scalar cell data
pure subroutine convert2cartesian_vectors(this, curv, v_curv, v_cart)
Convert curvilinear vector components to cartesian vector components.
integer, parameter, public logspherical
Basic fosite module.
common data structure
integer, parameter, public cartesian
real, parameter, public pi
base class for geometrical properties
integer, parameter, public logcylindrical
pure subroutine convert2cartesian_coords(this, curv, cart)
Convert curvilinear to cartesian coordinates.
derived mesh array class for vector cell data
integer, parameter, public spherical
pure subroutine radius_0(this, coords, radius)
Compute radial distances to the origin.
pure subroutine setscale2(this, gp, gp2)
Dictionary for generic data types.
Definition: common_dict.f90:61
pure subroutine setscale1(this, gp)
pure subroutine convert2curvilinear_vectors(this, curv, v_cart, v_curv)
Convert cartesian vector components to curvilinear vector components.
pure integer function getazimuthindex(this)
returns the coordinate index of the azimuthal angle
pure real function getscale1(this)
pure subroutine setscale3(this, gp, gp2, gp3)
integer, parameter, public cylindrical
pure real function getscale2(this, i)
pure subroutine convert2curvilinear_coords(this, cart, curv)
Convert cartesian to curvilinear coordinates.
pure subroutine positionvector_0(this, coords, posvec)
compute position vector components for all cell positions