93 CLASS(geometry_spherical),
INTENT(INOUT) :: this
94 TYPE(DICT_TYP),
POINTER :: config
97 CALL this%SetAzimuthIndex(3)
104 REAL,
INTENT(IN),
DIMENSION(:,:) :: coords
105 REAL,
INTENT(OUT),
DIMENSION(:) :: hx,hy,hz
107 CALL scalefactors(coords(:,1),coords(:,2),hx(:),hy(:),hz(:))
114 REAL,
INTENT(IN),
DIMENSION(:,:,:) :: coords
115 REAL,
INTENT(OUT),
DIMENSION(:,:) :: hx,hy,hz
117 CALL scalefactors(coords(:,:,1),coords(:,:,2),hx(:,:),hy(:,:),hz(:,:))
124 REAL,
INTENT(IN),
DIMENSION(:,:,:,:) :: coords
125 REAL,
INTENT(OUT),
DIMENSION(:,:,:) :: hx,hy,hz
127 CALL scalefactors(coords(:,:,:,1),coords(:,:,:,2),hx(:,:,:), &
135 REAL,
INTENT(IN),
DIMENSION(:,:,:,:,:) :: coords
136 REAL,
INTENT(OUT),
DIMENSION(:,:,:,:) :: hx,hy,hz
138 CALL scalefactors(coords(:,:,:,:,1),coords(:,:,:,:,2),hx(:,:,:,:), &
139 hy(:,:,:,:),hz(:,:,:,:))
142 PURE SUBROUTINE radius_1(this,coords,r)
146 REAL,
DIMENSION(:,:),
INTENT(IN) :: coords
147 REAL,
DIMENSION(:),
INTENT(OUT) :: r
152 PURE SUBROUTINE radius_2(this,coords,r)
156 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: coords
157 REAL,
DIMENSION(:,:),
INTENT(OUT) :: r
162 PURE SUBROUTINE radius_3(this,coords,r)
166 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: coords
167 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: r
169 r =
radius(coords(:,:,:,1))
172 PURE SUBROUTINE radius_4(this,coords,r)
176 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: coords
177 REAL,
DIMENSION(:,:,:,:),
INTENT(OUT) :: r
179 r =
radius(coords(:,:,:,:,1))
186 REAL,
DIMENSION(:,:),
INTENT(IN) :: coords
187 REAL,
DIMENSION(:,:),
INTENT(OUT) :: posvec
189 CALL positionvector(coords(:,1),posvec(:,1),posvec(:,2),posvec(:,3))
196 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: coords
197 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: posvec
199 CALL positionvector(coords(:,:,1),posvec(:,:,1),posvec(:,:,2),posvec(:,:,3))
206 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: coords
207 REAL,
DIMENSION(:,:,:,:),
INTENT(OUT) :: posvec
210 posvec(:,:,:,2),posvec(:,:,:,3))
217 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: coords
218 REAL,
DIMENSION(:,:,:,:,:),
INTENT(OUT) :: posvec
220 CALL positionvector(coords(:,:,:,:,1),posvec(:,:,:,:,1), &
221 posvec(:,:,:,:,2),posvec(:,:,:,:,3))
228 REAL,
DIMENSION(:,:),
INTENT(IN) :: curv
229 REAL,
DIMENSION(:,:),
INTENT(OUT) :: cart
232 cart(:,1),cart(:,2),cart(:,3))
239 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: curv
240 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: cart
243 cart(:,:,1),cart(:,:,2),cart(:,:,3))
250 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: curv
251 REAL,
DIMENSION(:,:,:,:),
INTENT(OUT) :: cart
253 CALL convert2cartesian_coords(curv(:,:,:,1),curv(:,:,:,2),curv(:,:,:,3), &
254 cart(:,:,:,1),cart(:,:,:,2),cart(:,:,:,3))
261 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: curv
262 REAL,
DIMENSION(:,:,:,:,:),
INTENT(OUT) :: cart
264 CALL convert2cartesian_coords(curv(:,:,:,:,1),curv(:,:,:,:,2),curv(:,:,:,:,3), &
265 cart(:,:,:,:,1),cart(:,:,:,:,2),cart(:,:,:,:,3))
272 REAL,
DIMENSION(:,:),
INTENT(IN) :: cart
273 REAL,
DIMENSION(:,:),
INTENT(OUT) :: curv
276 curv(:,1),curv(:,2),curv(:,3))
283 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: cart
284 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: curv
287 curv(:,:,1),curv(:,:,2),curv(:,:,3))
294 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: cart
295 REAL,
DIMENSION(:,:,:,:),
INTENT(OUT) :: curv
297 CALL convert2curvilinear_coords(cart(:,:,:,1),cart(:,:,:,2),cart(:,:,:,3), &
298 curv(:,:,:,1),curv(:,:,:,2),curv(:,:,:,3))
305 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: cart
306 REAL,
DIMENSION(:,:,:,:,:),
INTENT(OUT) :: curv
308 CALL convert2curvilinear_coords(cart(:,:,:,:,1),cart(:,:,:,:,2),cart(:,:,:,:,3), &
309 curv(:,:,:,:,1),curv(:,:,:,:,2),curv(:,:,:,:,3))
316 REAL,
DIMENSION(:,:),
INTENT(IN) :: curv
317 REAL,
DIMENSION(:,:),
INTENT(IN) :: v_curv
318 REAL,
DIMENSION(:,:),
INTENT(OUT) :: v_cart
321 v_curv(:,1),v_curv(:,2),v_curv(:,3), &
322 v_cart(:,1),v_cart(:,2),v_cart(:,3))
329 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: curv
330 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: v_curv
331 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: v_cart
334 v_curv(:,:,1),v_curv(:,:,2),v_curv(:,:,3), &
335 v_cart(:,:,1),v_cart(:,:,2),v_cart(:,:,3))
342 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: curv
343 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: v_curv
344 REAL,
DIMENSION(:,:,:,:),
INTENT(OUT) :: v_cart
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))
355 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: curv
356 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: v_curv
357 REAL,
DIMENSION(:,:,:,:,:),
INTENT(OUT) :: v_cart
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))
368 REAL,
DIMENSION(:,:),
INTENT(IN) :: curv
369 REAL,
DIMENSION(:,:),
INTENT(IN) :: v_cart
370 REAL,
DIMENSION(:,:),
INTENT(OUT) :: v_curv
373 v_cart(:,1),v_cart(:,2),v_cart(:,3), &
374 v_curv(:,1),v_curv(:,2),v_curv(:,3))
381 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: curv
382 REAL,
DIMENSION(:,:,:),
INTENT(IN) :: v_cart
383 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: v_curv
386 v_cart(:,:,1),v_cart(:,:,2),v_cart(:,:,3), &
387 v_curv(:,:,1),v_curv(:,:,2),v_curv(:,:,3))
394 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: curv
395 REAL,
DIMENSION(:,:,:,:),
INTENT(IN) :: v_cart
396 REAL,
DIMENSION(:,:,:,:),
INTENT(OUT) :: v_curv
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))
407 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: curv
408 REAL,
DIMENSION(:,:,:,:,:),
INTENT(IN) :: v_cart
409 REAL,
DIMENSION(:,:,:,:,:),
INTENT(OUT) :: v_curv
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))
419 CLASS(geometry_spherical),
INTENT(INOUT) :: this
421 CALL this%Finalize_base()
424 ELEMENTAL SUBROUTINE scalefactors(r,theta,hr,htheta,hphi)
427 REAL,
INTENT(IN) :: r,theta
428 REAL,
INTENT(OUT) :: hr,htheta,hphi
432 hphi = r * sin(theta)
435 ELEMENTAL FUNCTION radius(r)
438 REAL,
INTENT(IN) :: r
447 REAL,
INTENT(IN) :: r
448 REAL,
INTENT(OUT) :: rx,ry,rz
458 REAL,
INTENT(IN) :: r,theta,phi
459 REAL,
INTENT(OUT) :: x,y,z
461 x = r*sin(theta)*cos(phi)
462 y = r*sin(theta)*sin(phi)
469 REAL,
INTENT(IN) :: x,y,z
470 REAL,
INTENT(OUT) :: r,theta,phi
472 r = sqrt(x*x+y*y+z*z)
484 REAL,
INTENT(IN) :: r,theta,phi,vr,vtheta,vphi
485 REAL,
INTENT(OUT) :: vx,vy,vz
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)
496 REAL,
INTENT(IN) :: r,theta,phi,vx,vy,vz
497 REAL,
INTENT(OUT) :: vr,vtheta,vphi
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)
subroutine finalize(this)
Destructor of common class.
pure subroutine scalefactors_1(this, coords, hx, hy, hz)
pure subroutine convert2cartesian_coords_4(this, curv, cart)
type(logging_base), save this
derived mesh array class for scalar cell data
pure subroutine convert2cartesian_vectors_1(this, curv, v_curv, v_cart)
pure subroutine convert2cartesian_vectors_4(this, curv, v_curv, v_cart)
pure subroutine convert2curvilinear_vectors_2(this, curv, v_cart, v_curv)
pure subroutine convert2cartesian_vectors(this, curv, v_curv, v_cart)
Convert curvilinear vector components to cartesian vector components.
pure subroutine convert2curvilinear_vectors_1(this, curv, v_cart, v_curv)
character(len=32), parameter geometry_name
pure subroutine scalefactors_2(this, coords, hx, hy, hz)
pure subroutine convert2cartesian_vectors_3(this, curv, v_curv, v_cart)
defines properties of a 3D spherical mesh
pure subroutine convert2curvilinear_coords_4(this, cart, curv)
pure subroutine scalefactors_3(this, coords, hx, hy, hz)
pure subroutine radius_3(this, coords, r)
real, parameter, public pi
pure subroutine convert2curvilinear_coords_2(this, cart, curv)
base class for geometrical properties
pure subroutine convert2curvilinear_coords_3(this, cart, curv)
pure subroutine positionvector_2(this, coords, posvec)
pure subroutine positionvector_4(this, coords, posvec)
elemental subroutine positionvector(r, rx, ry, rz)
pure subroutine radius_1(this, coords, r)
pure subroutine convert2cartesian_coords(this, curv, cart)
Convert curvilinear to cartesian coordinates.
derived mesh array class for vector cell data
pure subroutine convert2cartesian_vectors_2(this, curv, v_curv, v_cart)
pure subroutine positionvector_3(this, coords, posvec)
integer, parameter, public spherical
pure subroutine convert2curvilinear_vectors_3(this, curv, v_cart, v_curv)
pure subroutine convert2cartesian_coords_2(this, curv, cart)
Dictionary for generic data types.
pure subroutine convert2curvilinear_vectors_4(this, curv, v_cart, v_curv)
pure subroutine convert2curvilinear_vectors(this, curv, v_cart, v_curv)
Convert cartesian vector components to curvilinear vector components.
pure subroutine radius_2(this, coords, r)
pure subroutine positionvector_1(this, coords, posvec)
pure subroutine radius_4(this, coords, r)
elemental real function radius(r)
pure subroutine convert2curvilinear_coords(this, cart, curv)
Convert cartesian to curvilinear coordinates.
pure subroutine convert2cartesian_coords_1(this, curv, cart)
pure subroutine convert2cartesian_coords_3(this, curv, cart)
elemental subroutine scalefactors(r, theta, hr, htheta, hphi)
subroutine initgeometry_spherical(this, config)
pure subroutine scalefactors_4(this, coords, hx, hy, hz)
pure subroutine convert2curvilinear_coords_1(this, cart, curv)