MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_geometry.t
Go to the documentation of this file.
1!> Module with geometry-related routines (e.g., divergence, curl)
3 use mod_comm_lib, only: mpistop
4 implicit none
5 public
6
7 integer :: coordinate=-1
8 integer, parameter :: cartesian = 0
9 integer, parameter :: cartesian_stretched= 1
10 integer, parameter :: cylindrical = 2
11 integer, parameter :: spherical = 3
12 integer, parameter :: cartesian_expansion= 4
13
14 integer :: type_curl=0
15 integer, parameter :: central=1
16 integer, parameter :: gaussbased=2
17 integer, parameter :: stokesbased=3
18
19contains
20
21 !> Set the coordinate system to be used
22 subroutine set_coordinate_system(geom)
24 character(len=*), intent(in) :: geom !< Name of the coordinate system
25
26 ! Store the geometry name
27 geometry_name = geom
28
29 select case (geom)
30 case ("Cartesian","Cartesian_1D","Cartesian_2D","Cartesian_3D")
31 ndir = ndim
33 case ("Cartesian_1D_expansion")
34 if (ndim /= 1) call mpistop("Geometry Cartesian_1D_expansion but ndim /= 1")
35 ndir = ndim
37 case ("Cartesian_1.5D")
38 if (ndim /= 1) call mpistop("Geometry Cartesian_1.5D but ndim /= 1")
40 ndir = 2
41 case ("Cartesian_1.75D")
42 if (ndim /= 1) call mpistop("Geometry Cartesian_1.75D but ndim /= 1")
44 ndir = 3
45 case ("Cartesian_2.5D")
46 if (ndim /= 2) call mpistop("Geometry Cartesian_2.5D but ndim /= 2")
48 ndir = 3
49 case ("cylindrical","cylindrical_2D","cylindrical_3D")
50 ndir = ndim
51 r_ = 1
52 z_ = 2
53 if(ndir==3) phi_ = 3
55 case ("cylindrical_2.5D")
56 if (ndim /= 2) call mpistop("Geometry cylindrical_2.5D but ndim /= 2")
57 ndir = 3
58 r_ = 1
59 z_ = 2
60 phi_ = 3
62 case ("polar","polar_2D","polar_3D")
63 ndir = ndim
64 r_ = 1
65 phi_ = 2
66 if(ndir==3) z_ = 3
68 case ("polar_1.5D")
69 if (ndim /= 1) call mpistop("Geometry polar_1.5D but ndim /= 1")
70 ndir = 2
71 r_ = 1
72 phi_ = 2
74 case ("polar_2.5D")
75 if (ndim /= 2) call mpistop("Geometry polar_2.5D but ndim /= 2")
76 ndir = 3
77 r_ = 1
78 phi_ = 2
79 z_ = 3
81 case ("spherical","spherical_2D","spherical_3D")
82 ndir = ndim
83 r_ = 1
84 if(ndir==3) phi_ = 3
85 z_ = -1
87 case ("spherical_2.5D")
88 if (ndim /= 2) &
89 call mpistop("Geometry spherical_2.5D requires ndim == 2")
90 ndir = 3
91 r_ = 1
92 phi_ = 3
93 z_ = -1
95 case default
96 call mpistop("Unknown geometry specified")
97 end select
98 end subroutine set_coordinate_system
99
100 subroutine set_pole
102
103 select case (coordinate)
104 case (spherical) {^ifthreed
105 ! For spherical grid, check whether phi-direction is periodic
106 if(periodb(ndim)) then
107 if(phi_/=3) call mpistop("phi_ should be 3 in 3D spherical coord!")
108 if(mod(ng3(1),2)/=0) &
109 call mpistop("Number of meshes in phi-direction should be even!")
110 if(abs(xprobmin2)<smalldouble) then
111 if(mype==0) write(unitterm,*) &
112 "Will apply pi-periodic conditions at northpole!"
113 poleb(1,2)=.true.
114 else
115 if(mype==0) write(unitterm,*) "There is no northpole!"
116 end if
117 if(abs(xprobmax2-dpi)<smalldouble) then
118 if(mype==0) write(unitterm,*) &
119 "Will apply pi-periodic conditions at southpole!"
120 poleb(2,2)=.true.
121 else
122 if(mype==0) write(unitterm,*) "There is no southpole!"
123 end if
124 end if}
125 case (cylindrical)
126 {
127 if (^d == phi_ .and. periodb(^d)) then
128 if(mod(ng^d(1),2)/=0) then
129 call mpistop("Number of meshes in phi-direction should be even!")
130 end if
131
132 if(abs(xprobmin1)<smalldouble) then
133 if (mype==0) then
134 write(unitterm,*) "Will apply pi-periodic conditions at r=0"
135 end if
136 poleb(1,1)=.true.
137 else
138 if (mype==0) then
139 write(unitterm,*) "There is no cylindrical axis!"
140 end if
141 end if
142 end if\}
143 end select
144
145 end subroutine set_pole
146
147 !> Deallocate geometry-related variables
148 subroutine putgridgeo(igrid)
150 integer, intent(in) :: igrid
151
152 deallocate(ps(igrid)%surfaceC,ps(igrid)%surface,ps(igrid)%dvolume,ps(igrid)%dx,&
153 psc(igrid)%dx,ps(igrid)%ds,psc(igrid)%dvolume)
154
155 end subroutine putgridgeo
156
157 !> calculate area of surfaces of cells
158 subroutine get_surface_area(s,ixG^L)
161
162 type(state) :: s
163 integer, intent(in) :: ixG^L
164
165 double precision :: x(ixG^S,ndim), xext(ixGmin1-1:ixGmax1,1), drs(ixG^S), drs_ext(ixGmin1-1:ixGmax1), dx2(ixG^S), dx3(ixG^S)
166 double precision :: exp_factor_ext(ixGmin1-1:ixGmax1),del_exp_factor_ext(ixGmin1-1:ixGmax1),exp_factor_primitive_ext(ixGmin1-1:ixGmax1)
167 double precision :: exp_factor(ixGmin1:ixGmax1),del_exp_factor(ixGmin1:ixGmax1),exp_factor_primitive(ixGmin1:ixGmax1)
168
169 select case (coordinate)
170
172 drs(ixg^s)=s%dx(ixg^s,1)
173 x(ixg^s,1)=s%x(ixg^s,1)
174 {^ifoned
175 if(associated(usr_set_surface))then
176 call usr_set_surface(ixgmin1,ixgmax1,x,drs,exp_factor,del_exp_factor,exp_factor_primitive)
177 if (any(exp_factor <= zero)) call mpistop("The area must always be strictly positive!")
178 endif
179 s%surface(ixg^s,1)=exp_factor(ixg^s)
180 xext(ixgmin1-1,1)=x(1,1)-half*drs(1)
181 xext(ixg^s,1)=x(ixg^s,1)+half*drs(ixg^s)
182 drs_ext(ixgmin1-1)=drs(1)
183 drs_ext(ixg^s)=drs(ixg^s)
184 if(associated(usr_set_surface)) call usr_set_surface(ixgmin1-1,ixgmax1,xext,drs_ext,exp_factor_ext,del_exp_factor_ext,exp_factor_primitive_ext)
185 s%surfaceC(ixgmin1-1:ixgmax1,1)=exp_factor_ext(ixgmin1-1:ixgmax1)
186 }
187
189 drs(ixg^s)=s%dx(ixg^s,1)
190 {^nooned
191 dx2(ixg^s)=s%dx(ixg^s,2)}
192 {^ifthreed
193 dx3(ixg^s)=s%dx(ixg^s,3)}
194
195 {^ifoned
196 s%surfaceC(ixg^s,1)=1.d0
197 s%surface(ixg^s,1) =1.d0
198 }
199 {^iftwod
200 s%surfaceC(ixg^s,1)=dx2(ixg^s)
201 s%surfaceC(ixg^s,2)=drs(ixg^s)
202 s%surface(ixg^s,1) =dx2(ixg^s)
203 s%surface(ixg^s,2)=drs(ixg^s)
204 }
205 {^ifthreed
206 s%surfaceC(ixg^s,1)= dx2(ixg^s)*dx3(ixg^s)
207 s%surfaceC(ixg^s,2)= drs(ixg^s)*dx3(ixg^s)
208 s%surfaceC(ixg^s,3)= drs(ixg^s)*dx2(ixg^s)
209 s%surface(ixg^s,1)=s%surfaceC(ixg^s,1)
210 s%surface(ixg^s,2)=s%surfaceC(ixg^s,2)
211 s%surface(ixg^s,3)=s%surfaceC(ixg^s,3)
212 }
213 {s%surfaceC(0^d%ixG^s,^d)=s%surfaceC(1^d%ixG^s,^d);\}
214 case (spherical)
215 x(ixg^s,1)=s%x(ixg^s,1)
216 {^nooned
217 x(ixg^s,2)=s%x(ixg^s,2)
218 }
219 drs(ixg^s)=s%dx(ixg^s,1)
220 {^nooned
221 dx2(ixg^s)=s%dx(ixg^s,2)
222 }
223 {^ifthreed
224 dx3(ixg^s)=s%dx(ixg^s,3)
225 }
226
227 s%surfaceC(ixg^s,1)=(x(ixg^s,1)+half*drs(ixg^s))**2 {^nooned &
228 *two*dabs(dsin(x(ixg^s,2)))*dsin(half*dx2(ixg^s))}{^ifthreed*dx3(ixg^s)}
229 ! change negative theta to positive to preserve positive area near pole
230
231 {^nooned
232 s%surfaceC(ixg^s,2)=x(ixg^s,1)*drs(ixg^s)&
233 *dabs(dsin(x(ixg^s,2)+half*dx2(ixg^s)))}{^ifthreed*dx3(ixg^s)}
234
235 {^ifthreed
236 s%surfaceC(ixg^s,3)=x(ixg^s,1)*drs(ixg^s)*dx2(ixg^s)
237 }
238
239 {^ifoned
240 s%surfaceC(0,1)=(x(1,1)-half*drs(1))**2
241 }
242 {^iftwod
243 s%surfaceC(0,ixgmin2:ixgmax2,1)=(x(1,ixgmin2:ixgmax2,1)-half*drs(1,&
244 ixgmin2:ixgmax2))**2*two*dabs(dsin(x(1,ixgmin2:ixgmax2,2)))*dsin(half*dx2(1,&
245 ixgmin2:ixgmax2))
246 s%surfaceC(ixgmin1:ixgmax1,0,2)=x(ixgmin1:ixgmax1,1,&
247 1)*drs(ixgmin1:ixgmax1,1)*dabs(dsin(x(ixgmin1:ixgmax1,1,&
248 2)-half*dx2(ixgmin1:ixgmax1,1)))
249 }
250 {^ifthreed
251 s%surfaceC(0,ixgmin2:ixgmax2,ixgmin3:ixgmax3,1)=(x(1,ixgmin2:ixgmax2,&
252 ixgmin3:ixgmax3,1)-half*drs(1,ixgmin2:ixgmax2,&
253 ixgmin3:ixgmax3))**2*two*dabs(dsin(x(1,ixgmin2:ixgmax2,ixgmin3:ixgmax3,&
254 2)))*dsin(half*dx2(1,ixgmin2:ixgmax2,ixgmin3:ixgmax3))*dx3(1,&
255 ixgmin2:ixgmax2,ixgmin3:ixgmax3)
256 s%surfaceC(ixgmin1:ixgmax1,0,ixgmin3:ixgmax3,2)=x(ixgmin1:ixgmax1,1,&
257 ixgmin3:ixgmax3,1)*drs(ixgmin1:ixgmax1,1,&
258 ixgmin3:ixgmax3)*dabs(dsin(x(ixgmin1:ixgmax1,1,ixgmin3:ixgmax3,&
259 2)-half*dx2(ixgmin1:ixgmax1,1,ixgmin3:ixgmax3)))*dx3(ixgmin1:ixgmax1,1,&
260 ixgmin3:ixgmax3)
261 s%surfaceC(ixgmin1:ixgmax1,ixgmin2:ixgmax2,0,3)=&
262 s%surfaceC(ixgmin1:ixgmax1,ixgmin2:ixgmax2,1,3)
263 }
264
265 s%surface(ixg^s,1)=x(ixg^s,1)**2 {^nooned &
266 *two*dabs(dsin(x(ixg^s,2)))*dsin(half*dx2(ixg^s))}{^ifthreed*dx3(ixg^s)}
267 {^nooned
268 s%surface(ixg^s,2)=x(ixg^s,1)*drs(ixg^s)&
269 *dabs(dsin(x(ixg^s,2)))}{^ifthreed*dx3(ixg^s)}
270
271 {^ifthreed
272 s%surface(ixg^s,3)=x(ixg^s,1)*drs(ixg^s)*dx2(ixg^s)}
273
274 case (cylindrical)
275 x(ixg^s,1)=s%x(ixg^s,1)
276 drs(ixg^s)=s%dx(ixg^s,1)
277 {^nooned
278 dx2(ixg^s)=s%dx(ixg^s,2)}
279 {^ifthreed
280 dx3(ixg^s)=s%dx(ixg^s,3)}
281
282 s%surfaceC(ixg^s,1)=dabs(x(ixg^s,1)+half*drs(ixg^s)){^de&*dx^de(ixg^s) }
283 ! change negative r to positive to preserve positive area near pole
284 {^nooned
285 if (z_==2) s%surfaceC(ixg^s,2)=dabs(x(ixg^s,1))*drs(ixg^s){^ifthreed*dx3(ixg^s)}
286 if (phi_==2) s%surfaceC(ixg^s,2)=drs(ixg^s){^ifthreed*dx3(ixg^s)}
287 }
288 {^ifthreed
289 if (z_==3) s%surfaceC(ixg^s,3)=dabs(x(ixg^s,1))*drs(ixg^s)*dx2(ixg^s)
290 if (phi_==3) s%surfaceC(ixg^s,3)=drs(ixg^s)*dx2(ixg^s)
291 }
292 {^ifoned
293 s%surfaceC(0,1)=dabs(x(1,1)-half*drs(1))
294 }
295 {^iftwod
296 s%surfaceC(0,ixgmin2:ixgmax2,1)=dabs(x(1,ixgmin2:ixgmax2,1)-half*drs(1,&
297 ixgmin2:ixgmax2))*dx2(1,ixgmin2:ixgmax2)
298 }
299 {^ifthreed
300 s%surfaceC(0,ixgmin2:ixgmax2,ixgmin3:ixgmax3,1)=dabs(x(1,ixgmin2:ixgmax2,&
301 ixgmin3:ixgmax3,1)-half*drs(1,ixgmin2:ixgmax2,ixgmin3:ixgmax3))*dx2(1,&
302 ixgmin2:ixgmax2,ixgmin3:ixgmax3)*dx3(1,ixgmin2:ixgmax2,&
303 ixgmin3:ixgmax3)
304 }
305 {s%surfaceC(0^de%ixG^s,^de)=s%surfaceC(1^de%ixG^s,^de);\}
306
307 s%surface(ixg^s,1)=dabs(x(ixg^s,1)){^de&*dx^de(ixg^s) }
308 {^nooned
309 if (z_==2) s%surface(ixg^s,2)=dabs(x(ixg^s,1))*drs(ixg^s){^ifthreed*dx3(ixg^s)}
310 if (phi_==2) s%surface(ixg^s,2)=drs(ixg^s){^ifthreed*dx3(ixg^s)}}
311 {^ifthreed
312 if (z_==3) s%surface(ixg^s,3)=dabs(x(ixg^s,1))*drs(ixg^s)*dx2(ixg^s)
313 if (phi_==3) s%surface(ixg^s,3)=drs(ixg^s)*dx2(ixg^s)}
314
315 case default
316 call mpistop("Sorry, coordinate unknown")
317 end select
318 end subroutine get_surface_area
319
320 !**************************************************************************
321 ! Purpose: Computes the gradient of a scalar field q(ixI^S) at cell
322 ! centers in the idir direction and outputs it in gradq(ixO^S),
323 ! which is also defined at cell **centers**. The gradient is
324 ! computed using a 2nd-order central difference scheme,
325 ! utilizing cell center values on both sides.
326 !
327 ! For uniform Cartesian coordinates, an optional input
328 ! parameter nth_in allows increasing the order of the
329 ! central difference scheme to 2*nth_in.
330 !**************************************************************************
331 subroutine gradient(q,ixI^L,ixO^L,idir,gradq,nth_in)
333 integer, intent(in) :: ixI^L, ixO^L, idir
334 integer, intent(in), optional :: nth_in
335 double precision, intent(in) :: q(ixI^S)
336 double precision, intent(inout) :: gradq(ixI^S)
337 integer :: jxO^L, hxO^L, nth
338
339 if(present(nth_in)) then
340 nth = nth_in
341 else
342 nth = 1
343 endif
344 if(nth .gt. nghostcells) then
345 call mpistop("gradient stencil too wide")
346 endif
347
348 hxo^l=ixo^l-kr(idir,^d);
349 jxo^l=ixo^l+kr(idir,^d);
350 select case(coordinate)
351 case(cartesian)
352 select case(nth)
353 case(1)
354 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))*half
355 case(2)
356 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))*2.d0/3.d0
357 hxo^l=ixo^l-2*kr(idir,^d);
358 jxo^l=ixo^l+2*kr(idir,^d);
359 gradq(ixo^s)=gradq(ixo^s)-(q(jxo^s)-q(hxo^s))/12.d0
360 case(3)
361 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))*3.d0/4.d0
362 hxo^l=ixo^l-2*kr(idir,^d);
363 jxo^l=ixo^l+2*kr(idir,^d);
364 gradq(ixo^s)=gradq(ixo^s)-(q(jxo^s)-q(hxo^s))*3.d0/20.d0
365 hxo^l=ixo^l-3*kr(idir,^d);
366 jxo^l=ixo^l+3*kr(idir,^d);
367 gradq(ixo^s)=gradq(ixo^s)+(q(jxo^s)-q(hxo^s))/60.d0
368 case default
369 call mpistop("unknown stencil gradient")
370 end select
371 gradq(ixo^s)=gradq(ixo^s)/dxlevel(idir)
373 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(block%x(jxo^s,idir)-block%x(hxo^s,idir))
374 case(spherical)
375 select case(idir)
376 case(1)
377 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((block%x(jxo^s,1)-block%x(hxo^s,1)))
378 {^nooned
379 case(2)
380 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((block%x(jxo^s,2)-block%x(hxo^s,2))*block%x(ixo^s,1))
381 }
382 {^ifthreed
383 case(3)
384 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((block%x(jxo^s,3)-block%x(hxo^s,3))*block%x(ixo^s,1)*dsin(block%x(ixo^s,2)))
385 }
386 end select
387 case(cylindrical)
388 if(idir==phi_) then
389 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((block%x(jxo^s,phi_)-block%x(hxo^s,phi_))*block%x(ixo^s,r_))
390 else
391 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(block%x(jxo^s,idir)-block%x(hxo^s,idir))
392 end if
393 case default
394 call mpistop('Unknown geometry')
395 end select
396 end subroutine gradient
397
398 !**************************************************************************
399 ! Purpose: Originally named as gradientS.
400 ! Computes the gradient of a scalar field q(ixI^S) at cell
401 ! centers in the idir direction and outputs it in gradq(ixO^S),
402 ! which is also defined at cell **centers**. Unlike the standard
403 ! gradient computation, this routine first reconstructs cell
404 ! face values using a slope **Limiter** and then computes the
405 ! gradient from these face values.
406 !**************************************************************************
407 subroutine gradientl(q,ixI^L,ixO^L,idir,gradq)
409 use mod_limiter
410 use mod_ppm
411 integer, intent(in) :: ixI^L, ixO^L, idir
412 double precision, intent(in) :: q(ixI^S)
413 double precision, intent(inout) :: gradq(ixI^S)
414 double precision ,dimension(ixI^S) :: qC,qL,qR,dqC,ldq,rdq
415 double precision :: x(ixI^S,1:ndim)
416 double precision :: invdx
417 integer :: hxO^L,ixC^L,jxC^L,gxC^L,hxC^L
418
419 x(ixi^s,1:ndim)=block%x(ixi^s,1:ndim)
420 invdx=1.d0/dxlevel(idir)
421 hxo^l=ixo^l-kr(idir,^d);
422 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
423 jxc^l=ixc^l+kr(idir,^d);
424 gxcmin^d=ixcmin^d-kr(idir,^d);gxcmax^d=jxcmax^d;
425 hxc^l=gxc^l+kr(idir,^d);
426 ! set the gradient limiter here
427 qr(gxc^s) = q(hxc^s)
428 ql(gxc^s) = q(gxc^s)
430 dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
431 call dwlimiter2(dqc,ixi^l,gxc^l,idir,type_gradient_limiter(block%level),ldw=ldq,rdw=rdq)
432 ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
433 qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
434 else
435 call ppmlimitervar(ixi^l,ixo^l,idir,q,q,ql,qr)
436 endif
437
438 select case(coordinate)
439 case(cartesian)
440 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))*invdx
442 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/block%dx(ixo^s,idir)
443 case(spherical)
444 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/block%dx(ixo^s,idir)
445 select case(idir)
446 case(2)
447 gradq(ixo^s)=gradq(ixo^s)/x(ixo^s,1)
448 {^ifthreed
449 case(3)
450 gradq(ixo^s)=gradq(ixo^s)/(x(ixo^s,1)*dsin(x(ixo^s,2)))
451 }
452 end select
453 case(cylindrical)
454 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/block%dx(ixo^s,idir)
455 if(idir==phi_) gradq(ixo^s)=gradq(ixo^s)/x(ixo^s,1)
456 end select
457 end subroutine gradientl
458
459 !**************************************************************************
460 ! Purpose: Merged from gradientx, gradientq and gradientC.
461 ! Computes the gradient of a scalar field q(ixI^S) at cell
462 ! centers in the idir direction and outputs it in gradq(ixO^S),
463 ! which is defined at cell **Faces**. The gradient is computed
464 ! using a 2nd-order central difference scheme, utilizing
465 ! the cell center values on both sides of a given cell face.
466 !
467 ! For uniform Cartesian coordinates, an optional input
468 ! parameter nth_in allows increasing the order of the
469 ! central difference scheme to 2*nth_in.
470 !
471 ! By default, the index is shifted one position to the left.
472 ! This structure also allows the routine to be used for
473 ! one-sided difference calculations of the gradient at cell
474 ! **centers**. In such cases, the optional logical parameter pm_in
475 ! controls whether the index shifts to the left (default)
476 ! or to the right.
477 !**************************************************************************
478 subroutine gradientf(q,x,ixI^L,ixO^L,idir,gradq,nth_in,pm_in)
480 integer, intent(in) :: ixI^L, ixO^L, idir
481 double precision, intent(in) :: q(ixI^S),x(ixI^S,1:ndim)
482 double precision, intent(inout) :: gradq(ixI^S)
483 integer, intent(in), optional :: nth_in
484 logical, intent(in), optional :: pm_in
485 logical :: pm
486 integer :: nth,jxO^L,hxO^L
487
488 if(present(nth_in)) then
489 nth = nth_in
490 else
491 nth = 1
492 endif
493 if(nth .gt. nghostcells) then
494 call mpistop("gradient stencil too wide")
495 endif
496 if(present(pm_in)) then
497 pm = pm_in
498 else
499 pm = .true.
500 endif
501
502 if(pm) then
503 hxo^l=ixo^l;
504 jxo^l=ixo^l+kr(idir,^d);
505 else
506 hxo^l=ixo^l-kr(idir,^d);
507 jxo^l=ixo^l;
508 endif
509 select case(coordinate)
510 case(cartesian)
511 select case(nth)
512 case(1)
513 gradq(ixo^s)=q(jxo^s)-q(hxo^s)
514 case(2)
515 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))*27.d0/24.d0
516 if(pm) then
517 hxo^l=ixo^l-kr(idir,^d);
518 jxo^l=ixo^l+2*kr(idir,^d);
519 else
520 hxo^l=ixo^l-2*kr(idir,^d);
521 jxo^l=ixo^l+kr(idir,^d);
522 endif
523 gradq(ixo^s)=gradq(ixo^s)-(q(jxo^s)-q(hxo^s))/24.d0
524 case(3)
525 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))*225.d0/192.d0
526 if(pm) then
527 hxo^l=ixo^l-kr(idir,^d);
528 jxo^l=ixo^l+2*kr(idir,^d);
529 else
530 hxo^l=ixo^l-2*kr(idir,^d);
531 jxo^l=ixo^l+kr(idir,^d);
532 endif
533 gradq(ixo^s)=gradq(ixo^s)-(q(jxo^s)-q(hxo^s))*25.d0/384.d0
534 if(pm) then
535 hxo^l=ixo^l-2*kr(idir,^d);
536 jxo^l=ixo^l+3*kr(idir,^d);
537 else
538 hxo^l=ixo^l-3*kr(idir,^d);
539 jxo^l=ixo^l+2*kr(idir,^d);
540 endif
541 gradq(ixo^s)=gradq(ixo^s)+(q(jxo^s)-q(hxo^s))*3.0/640.d0
542 case default
543 call mpistop("unknown stencil gradient")
544 end select
545 gradq(ixo^s)=gradq(ixo^s)/dxlevel(idir)
547 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
548 case(spherical)
549 select case(idir)
550 case(1)
551 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,1)-x(hxo^s,1))
552 {^nooned
553 case(2)
554 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,2)-x(hxo^s,2))*x(ixo^s,1))
555 }
556 {^ifthreed
557 case(3)
558 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,3)-x(hxo^s,3))*x(ixo^s,1)*dsin(x(ixo^s,2)))
559 }
560 end select
561 case(cylindrical)
562 if(idir==phi_) then
563 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,phi_)-x(hxo^s,phi_))*x(ixo^s,r_))
564 else
565 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
566 end if
567 case default
568 call mpistop('Unknown geometry')
569 end select
570 end subroutine gradientf
571
572 !**************************************************************************
573 ! Purpose: Computes the divergence of a vector field qvec(ixI^S,1:ndim)
574 ! at cell centers and stores the result in q(ixO^S), which is
575 ! also defined at cell centers. The divergence is computed using
576 ! a 2nd-order central difference scheme, utilizing the vector
577 ! components at adjacent cell centers.
578 !
579 ! For uniform Cartesian coordinates, an optional input parameter
580 ! nth_in allows increasing the order of the central difference
581 ! scheme to 2*nth_in, improving the accuracy of the divergence
582 ! computation.
583 !**************************************************************************
584 subroutine divvector(qvec,ixI^L,ixO^L,divq,nth_in)
586 integer, intent(in) :: ixI^L,ixO^L
587 double precision, intent(in) :: qvec(ixI^S,1:ndir)
588 double precision, intent(inout) :: divq(ixI^S)
589 integer, intent(in), optional :: nth_in
590 double precision :: qC(ixI^S), invdx(1:ndim)
591 integer :: jxO^L, hxO^L, ixC^L, jxC^L
592 integer :: idims, gxO^L, kxO^L
593 integer :: lxO^L, fxO^L, nth
594
595 if(present(nth_in)) then
596 nth = nth_in
597 else
598 nth = 1
599 endif
600 if(nth .gt. nghostcells) then
601 call mpistop("divvector stencil too wide")
602 endif
603
604 invdx=1.d0/dxlevel
605 divq(ixo^s)=0.0d0
606 if (slab_uniform) then
607 do idims=1,ndim
608 select case(nth)
609 case(1)
610 jxo^l=ixo^l+kr(idims,^d);
611 hxo^l=ixo^l-kr(idims,^d);
612 divq(ixo^s)=divq(ixo^s)+half*(qvec(jxo^s,idims) &
613 - qvec(hxo^s,idims))*invdx(idims)
614 case(2)
615 kxo^l=ixo^l+2*kr(idims,^d);
616 jxo^l=ixo^l+kr(idims,^d);
617 hxo^l=ixo^l-kr(idims,^d);
618 gxo^l=ixo^l-2*kr(idims,^d);
619 divq(ixo^s)=divq(ixo^s)+&
620 (-qvec(kxo^s,idims)+8.d0*qvec(jxo^s,idims)-8.d0*&
621 qvec(hxo^s,idims)+qvec(gxo^s,idims))/(12.d0*dxlevel(idims))
622 case(3)
623 lxo^l=ixo^l+3*kr(idims,^d);
624 kxo^l=ixo^l+2*kr(idims,^d);
625 jxo^l=ixo^l+kr(idims,^d);
626 hxo^l=ixo^l-kr(idims,^d);
627 gxo^l=ixo^l-2*kr(idims,^d);
628 fxo^l=ixo^l-3*kr(idims,^d);
629 divq(ixo^s)=divq(ixo^s)+&
630 (-qvec(fxo^s,idims)+9.d0*qvec(gxo^s,idims)-45.d0*qvec(hxo^s,idims)&
631 +qvec(lxo^s,idims)-9.d0*qvec(kxo^s,idims)+45.d0*qvec(jxo^s,idims))&
632 /(60.d0*dxlevel(idims))
633 case default
634 call mpistop("unknown stencil in divvector")
635 end select
636 enddo
637 else
638 do idims=1,ndim
639 hxo^l=ixo^l-kr(idims,^d);
640 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
641 jxc^l=ixc^l+kr(idims,^d);
642 if(stretched_dim(idims) .and. stretch_uncentered) then
643 ! linear interpolation at cell interface along stretched dimension
644 qc(ixc^s)=block%surfaceC(ixc^s,idims)*(qvec(ixc^s,idims)+0.5d0*block%dx(ixc^s,idims)*&
645 (qvec(jxc^s,idims)-qvec(ixc^s,idims))/(block%x(jxc^s,idims)-block%x(ixc^s,idims)))
646 else
647 qc(ixc^s)=block%surfaceC(ixc^s,idims)*half*(qvec(ixc^s,idims)+qvec(jxc^s,idims))
648 end if
649 divq(ixo^s)=divq(ixo^s)+qc(ixo^s)-qc(hxo^s)
650 end do
651 divq(ixo^s)=divq(ixo^s)/block%dvolume(ixo^s)
652 end if
653 end subroutine divvector
654
655 !> Calculate divergence of a vector qvec within ixL
656 !> using limited extrapolation to cell edges
657 subroutine divvectors(qvec,ixI^L,ixO^L,divq)
659 use mod_limiter
660 use mod_ppm
661 integer, intent(in) :: ixI^L,ixO^L
662 double precision, intent(in) :: qvec(ixI^S,1:ndir)
663 double precision, intent(inout) :: divq(ixI^S)
664 double precision, dimension(ixI^S) :: qL,qR,dqC,ldq,rdq
665
666 double precision :: invdx(1:ndim)
667 integer :: hxO^L,ixC^L,jxC^L,idims,ix^L,gxC^L,hxC^L
668
669 ix^l=ixo^l^ladd2;
670
671 if (iximin^d>ixmin^d.or.iximax^d<ixmax^d|.or.) &
672 call mpistop("Error in divvectorS: Non-conforming input limits")
673
674 invdx=1.d0/dxlevel
675 divq(ixo^s)=zero
676 do idims=1,ndim
677 hxo^l=ixo^l-kr(idims,^d);
678 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
679 jxc^l=ixc^l+kr(idims,^d);
680 gxcmin^d=ixcmin^d-kr(idims,^d);gxcmax^d=jxcmax^d;
681 hxc^l=gxc^l+kr(idims,^d);
682 qr(gxc^s) = qvec(hxc^s,idims)
683 ql(gxc^s) = qvec(gxc^s,idims)
685 dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
686 call dwlimiter2(dqc,ixi^l,gxc^l,idims,type_gradient_limiter(block%level),ldw=ldq,rdw=rdq)
687 ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
688 qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
689 else
690 dqc(ixi^s)=qvec(ixi^s,idims)
691 call ppmlimitervar(ixi^l,ixo^l,idims,dqc,dqc,ql,qr)
692 endif
693
694 if (slab_uniform) then
695 divq(ixo^s)=divq(ixo^s)+half*(qr(ixo^s)-ql(hxo^s))*invdx(idims)
696 else
697 qr(ixc^s)=block%surfaceC(ixc^s,idims)*qr(ixc^s)
698 ql(ixc^s)=block%surfaceC(ixc^s,idims)*ql(ixc^s)
699 divq(ixo^s)=divq(ixo^s)+qr(ixo^s)-ql(hxo^s)
700 end if
701 end do
702 if(.not.slab_uniform) divq(ixo^s)=divq(ixo^s)/block%dvolume(ixo^s)
703 end subroutine divvectors
704
705 !> Calculate curl of a vector qvec within ixL
706 !> Options to
707 !> employ standard second order CD evaluations
708 !> use Gauss theorem for non-Cartesian grids
709 !> use Stokes theorem for non-Cartesian grids
710 subroutine curlvector(qvec,ixI^L,ixO^L,curlvec,idirmin,idirmin0,ndir0,fourthorder)
712
713 integer, intent(in) :: ixI^L,ixO^L
714 integer, intent(in) :: ndir0, idirmin0
715 integer, intent(inout) :: idirmin
716 double precision, intent(in) :: qvec(ixI^S,1:ndir0)
717 double precision, intent(inout) :: curlvec(ixI^S,idirmin0:3)
718 logical, intent(in), optional :: fourthorder !< Default: false
719
720 double precision :: invdx(1:ndim)
721 double precision :: tmp(ixI^S),tmp2(ixI^S),xC(ixI^S),surface(ixI^S)
722 integer :: ixA^L,ixC^L,jxC^L,idir,jdir,kdir,hxO^L,jxO^L,kxO^L,gxO^L
723 logical :: use_4th_order
724
725 ! Calculate curl within ixL: CurlV_i=eps_ijk*d_j V_k
726 ! Curl can have components (idirmin:3)
727 ! Determine exact value of idirmin while doing the loop.
728
729 use_4th_order = .false.
730 if (present(fourthorder)) use_4th_order = fourthorder
731
732 if (use_4th_order) then
733 if (.not. slab_uniform) &
734 call mpistop("curlvector: 4th order only supported for slab geometry")
735 ! Fourth order, stencil width is two
736 ixa^l=ixo^l^ladd2;
737 else
738 ! Second order, stencil width is one
739 ixa^l=ixo^l^ladd1;
740 end if
741
742 if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
743 call mpistop("Error in curlvector: Non-conforming input limits")
744
745 idirmin=4
746 curlvec(ixo^s,idirmin0:3)=zero
747
748 ! all non-Cartesian cases
749 select case(coordinate)
750 case(cartesian) ! Cartesian grids
751 invdx=1.d0/dxlevel
752 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
753 if(lvc(idir,jdir,kdir)/=0)then
754 if (.not. use_4th_order) then
755 ! Use second order scheme
756 jxo^l=ixo^l+kr(jdir,^d);
757 hxo^l=ixo^l-kr(jdir,^d);
758 tmp(ixo^s)=half*(qvec(jxo^s,kdir) &
759 - qvec(hxo^s,kdir))*invdx(jdir)
760 else
761 ! Use fourth order scheme
762 kxo^l=ixo^l+2*kr(jdir,^d);
763 jxo^l=ixo^l+kr(jdir,^d);
764 hxo^l=ixo^l-kr(jdir,^d);
765 gxo^l=ixo^l-2*kr(jdir,^d);
766 tmp(ixo^s)=(-qvec(kxo^s,kdir) + 8.0d0 * qvec(jxo^s,kdir) - 8.0d0 * &
767 qvec(hxo^s,kdir) + qvec(gxo^s,kdir))/(12.0d0 * dxlevel(jdir))
768 end if
769 if(lvc(idir,jdir,kdir)==1)then
770 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp(ixo^s)
771 else
772 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp(ixo^s)
773 endif
774 if(idir<idirmin)idirmin=idir
775 endif
776 enddo; enddo; enddo;
777 case(cartesian_stretched) ! stretched Cartesian grids
778 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
779 if(lvc(idir,jdir,kdir)/=0)then
780 select case(type_curl)
781 case(central)
782 tmp(ixa^s)=qvec(ixa^s,kdir)
783 hxo^l=ixo^l-kr(jdir,^d);
784 jxo^l=ixo^l+kr(jdir,^d);
785 ! second order centered differencing
786 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,jdir)-block%x(hxo^s,jdir))
787 case(gaussbased)
788 hxo^l=ixo^l-kr(jdir,^d);
789 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
790 jxc^l=ixc^l+kr(jdir,^d);
791 tmp(ixc^s)=block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
792 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
793 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%dvolume(ixo^s)
794 case(stokesbased)
795 hxo^l=ixo^l-kr(jdir,^d);
796 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
797 jxc^l=ixc^l+kr(jdir,^d);
798 if(kdir<=ndim)then
799 tmp(ixc^s)=block%ds(ixo^s,kdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
800 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
801 else
802 tmp(ixc^s)=(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
803 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
804 endif
805 if(idir<=ndim)then
806 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%surface(ixo^s,idir)
807 else ! essentially for 2.5D case, idir=3 and jdir,kdir<=2
808 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/(block%ds(ixo^s,jdir)*block%ds(ixo^s,kdir))
809 endif
810 case default
811 call mpistop('no such curl evaluator')
812 end select
813 if(lvc(idir,jdir,kdir)==1)then
814 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
815 else
816 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
817 endif
818 if(idir<idirmin)idirmin=idir
819 endif
820 enddo; enddo; enddo;
821 case(spherical) ! possibly stretched spherical grids
822 select case(type_curl)
823 case(central) ! ok for any dimensionality
824 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
825 if(lvc(idir,jdir,kdir)/=0)then
826 tmp(ixa^s)=qvec(ixa^s,kdir)
827 hxo^l=ixo^l-kr(jdir,^d);
828 jxo^l=ixo^l+kr(jdir,^d);
829 select case(jdir)
830 case(1)
831 tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1)
832 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,1)-block%x(hxo^s,1))*block%x(ixo^s,1))
833 {^nooned case(2)
834 if(idir==1) tmp(ixa^s)=tmp(ixa^s)*dsin(block%x(ixa^s,2))
835 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,2)-block%x(hxo^s,2))*block%x(ixo^s,1))
836 if(idir==1) tmp2(ixo^s)=tmp2(ixo^s)/dsin(block%x(ixo^s,2))
837 }
838 {^ifthreed case(3)
839 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,3)-block%x(hxo^s,3))*block%x(ixo^s,1)*dsin(block%x(ixo^s,2)))
840 }
841 end select
842 if(lvc(idir,jdir,kdir)==1)then
843 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
844 else
845 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
846 endif
847 if(idir<idirmin)idirmin=idir
848 endif
849 enddo; enddo; enddo;
850 case(gaussbased)
851 do idir=idirmin0,3;
852 do jdir=1,ndim; do kdir=1,ndir0
853 if(lvc(idir,jdir,kdir)/=0)then
854 hxo^l=ixo^l-kr(jdir,^d);
855 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
856 jxc^l=ixc^l+kr(jdir,^d);
857 tmp(ixc^s)=block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
858 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
859 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%dvolume(ixo^s)
860 if(lvc(idir,jdir,kdir)==1)then
861 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
862 else
863 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
864 endif
865 if(idir<idirmin)idirmin=idir
866 endif
867 enddo; enddo;
868 ! geometric terms
869 if(idir==2.and.phi_>0) curlvec(ixo^s,2)=curlvec(ixo^s,2)+qvec(ixo^s,phi_)/block%x(ixo^s,r_)
870 {^nooned
871 if(idir==phi_) curlvec(ixo^s,phi_)=curlvec(ixo^s,phi_)-qvec(ixo^s,2)/block%x(ixo^s,r_) &
872 +qvec(ixo^s,r_)*dcos(block%x(ixo^s,2))/(block%x(ixo^s,r_)*dsin(block%x(ixo^s,2)))
873 }
874 enddo;
875 case(stokesbased)
876 !if(ndim<3) call mpistop("Stokesbased for 3D spherical only")
877 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
878 if(lvc(idir,jdir,kdir)/=0)then
879 select case(idir)
880 case(1)
881 if(jdir<kdir) then
882 ! idir=1,jdir=2,kdir=3
883 !! integral along 3rd dimension
884 hxo^l=ixo^l-kr(jdir,^d);
885 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
886 jxc^l=ixc^l+kr(jdir,^d);
887 ! qvec(3) at cell interface along 2nd dimension
888 if(stretched_dim(jdir) .and. stretch_uncentered) then
889 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
890 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
891 else
892 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
893 end if
894 ! 2nd coordinate at cell interface along 2nd dimension
895 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
896 curlvec(ixo^s,idir)=(dsin(xc(ixo^s))*tmp(ixo^s)-&
897 dsin(xc(hxo^s))*tmp(hxo^s))*block%dx(ixo^s,kdir)
898 !! integral along 2nd dimension
899 hxo^l=ixo^l-kr(kdir,^d);
900 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
901 jxc^l=ixc^l+kr(kdir,^d);
902 ! qvec(2) at cell interface along 3rd dimension
903 if(stretched_dim(kdir) .and. stretch_uncentered) then
904 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
905 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
906 else
907 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
908 end if
909 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir))&
910 /block%surface(ixo^s,idir)*block%x(ixo^s,idir)
911 end if
912 case(2)
913 if(jdir<kdir) then
914 ! idir=2,jdir=1,kdir=3
915 !! integral along 1st dimension
916 hxo^l=ixo^l-kr(kdir,^d);
917 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
918 jxc^l=ixc^l+kr(kdir,^d);
919 ! qvec(1) at cell interface along 3rd dimension
920 if(stretched_dim(kdir) .and. stretch_uncentered) then
921 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
922 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
923 else
924 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
925 end if
926 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,1)
927 !! integral along 3rd dimension
928 hxo^l=ixo^l-kr(jdir,^d);
929 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
930 jxc^l=ixc^l+kr(jdir,^d);
931 ! qvec(3) at cell interface along 1st dimension
932 if(stretched_dim(jdir) .and. stretch_uncentered) then
933 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
934 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
935 else
936 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
937 end if
938 ! 1st coordinate at cell interface along 1st dimension
939 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
940 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(hxo^s)*tmp(hxo^s)-xc(ixo^s)*tmp(ixo^s))*&
941 dsin(block%x(ixo^s,idir))*block%dx(ixo^s,kdir))/block%surface(ixo^s,idir)
942 end if
943 case(3)
944 if(jdir<kdir) then
945 ! idir=3,jdir=1,kdir=2
946 !! integral along 1st dimension
947 hxo^l=ixo^l-kr(kdir,^d);
948 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
949 jxc^l=ixc^l+kr(kdir,^d);
950 ! qvec(1) at cell interface along 2nd dimension
951 if(stretched_dim(kdir) .and. stretch_uncentered) then
952 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
953 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
954 else
955 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
956 end if
957 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
958 !! integral along 2nd dimension
959 hxo^l=ixo^l-kr(jdir,^d);
960 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
961 jxc^l=ixc^l+kr(jdir,^d);
962 ! qvec(2) at cell interface along 1st dimension
963 if(stretched_dim(jdir) .and. stretch_uncentered) then
964 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
965 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
966 else
967 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
968 end if
969 ! 1st coordinate at cell interface along 1st dimension
970 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
971 if(ndim==3) then
972 surface(ixo^s)=block%surface(ixo^s,idir)
973 else
974 surface(ixo^s)=block%x(ixo^s,jdir)*block%dx(ixo^s,kdir)*block%dx(ixo^s,jdir)
975 end if
976 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(ixo^s)*tmp(ixo^s)-xc(hxo^s)*tmp(hxo^s))*block%dx(ixo^s,kdir))&
977 /surface(ixo^s)
978 end if
979 end select
980 if(idir<idirmin)idirmin=idir
981 endif
982 enddo; enddo; enddo;
983 case default
984 call mpistop('no such curl evaluator')
985 end select
986 case(cylindrical) ! possibly stretched cylindrical grids
987 select case(type_curl)
988 case(central) ! works for any dimensionality, polar/cylindrical
989 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
990 if(lvc(idir,jdir,kdir)/=0)then
991 tmp(ixa^s)=qvec(ixa^s,kdir)
992 hxo^l=ixo^l-kr(jdir,^d);
993 jxo^l=ixo^l+kr(jdir,^d);
994 if(z_==3.or.z_==-1) then
995 ! Case Polar_2D, Polar_2.5D or Polar_3D, i.e. R,phi,Z
996 select case(jdir)
997 case(1)
998 if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
999 ! computes d(R V_phi)/dR or d V_Z/dR
1000 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,1)-block%x(hxo^s,1))
1001 if(idir==z_) tmp2(ixo^s)=tmp2(ixo^s)/block%x(ixo^s,1) ! (1/R)*d(R V_phi)/dR
1002 {^nooned case(2)
1003 ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
1004 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,2)-block%x(hxo^s,2))*block%x(ixo^s,1))
1005 }
1006 {^ifthreed case(3)
1007 ! handles d V_phi/dZ or d V_R/dZ
1008 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,3)-block%x(hxo^s,3))
1009 }
1010 end select
1011 end if
1012 if(phi_==3.or.phi_==-1) then
1013 ! Case Cylindrical_2D, Cylindrical_2.5D or Cylindrical_3D, i.e. R,Z,phi
1014 select case(jdir)
1015 case(1)
1016 if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
1017 ! computes d(R V_phi)/dR or d V_Z/dR
1018 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,1)-block%x(hxo^s,1))
1019 if(idir==z_) tmp2(ixo^s)=tmp2(ixo^s)/block%x(ixo^s,1) ! (1/R)*d(R V_phi)/dR
1020 {^nooned case(2)
1021 ! handles d V_phi/dZ or d V_R/dZ
1022 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,2)-block%x(hxo^s,2))
1023 }
1024 {^ifthreed case(3)
1025 ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
1026 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,3)-block%x(hxo^s,3))*block%x(ixo^s,1))
1027 }
1028 end select
1029 end if
1030 if(lvc(idir,jdir,kdir)==1)then
1031 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1032 else
1033 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1034 endif
1035 if(idir<idirmin)idirmin=idir
1036 endif
1037 enddo; enddo; enddo;
1038 case(gaussbased) ! works for any dimensionality, polar/cylindrical
1039 if(ndim<2) call mpistop("Gaussbased for 2D, 2.5D or 3D polar or cylindrical only")
1040 do idir=idirmin0,3;
1041 do jdir=1,ndim; do kdir=1,ndir0
1042 if(lvc(idir,jdir,kdir)/=0)then
1043 hxo^l=ixo^l-kr(jdir,^d);
1044 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1045 jxc^l=ixc^l+kr(jdir,^d);
1046 tmp(ixc^s)=block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1047 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
1048 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%dvolume(ixo^s)
1049 if(lvc(idir,jdir,kdir)==1)then
1050 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1051 else
1052 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1053 endif
1054 if(idir<idirmin)idirmin=idir
1055 endif
1056 enddo; enddo;
1057 ! geometric term from d e_R/d phi= e_phi for unit vectors e_R, e_phi
1058 ! but minus sign appears due to R,Z,phi ordering (?)
1059 ! note that in cylindrical 2D (RZ), phi_ is -1
1060 ! note that in polar 2D (Rphi), z_ is -1
1061 if((idir==phi_.or.(phi_==-1.and.idir==3)).and.z_>0) then
1062 ! cylindrical
1063 if( z_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-qvec(ixo^s,z_)/block%x(ixo^s,r_)
1064 ! polar
1065 if(phi_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+qvec(ixo^s,z_)/block%x(ixo^s,r_)
1066 endif
1067 enddo;
1068 case(stokesbased)
1069 if(ndim<3) call mpistop("Stokesbased for 3D cylindrical only")
1070 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1071 if(lvc(idir,jdir,kdir)/=0)then
1072 if(idir==r_) then
1073 if(jdir==phi_) then
1074 ! idir=r,jdir=phi,kdir=z
1075 !! integral along z dimension
1076 hxo^l=ixo^l-kr(jdir,^d);
1077 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1078 jxc^l=ixc^l+kr(jdir,^d);
1079 ! qvec(z) at cell interface along phi dimension
1080 if(stretched_dim(jdir) .and. stretch_uncentered) then
1081 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1082 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1083 else
1084 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1085 end if
1086 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,kdir)
1087 !! integral along phi dimension
1088 hxo^l=ixo^l-kr(kdir,^d);
1089 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1090 jxc^l=ixc^l+kr(kdir,^d);
1091 ! qvec(phi) at cell interface along z dimension
1092 if(stretched_dim(kdir) .and. stretch_uncentered) then
1093 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1094 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1095 else
1096 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1097 end if
1098 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*block%x(ixo^s,idir)*block%dx(ixo^s,jdir))&
1099 /block%surface(ixo^s,idir)
1100 end if
1101 else if(idir==phi_) then
1102 if(jdir<kdir) then
1103 ! idir=phi,jdir=r,kdir=z
1104 !! integral along r dimension
1105 hxo^l=ixo^l-kr(kdir,^d);
1106 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1107 jxc^l=ixc^l+kr(kdir,^d);
1108 ! qvec(r) at cell interface along z dimension
1109 if(stretched_dim(kdir) .and. stretch_uncentered) then
1110 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1111 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1112 else
1113 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1114 end if
1115 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,jdir)
1116 !! integral along z dimension
1117 hxo^l=ixo^l-kr(jdir,^d);
1118 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1119 jxc^l=ixc^l+kr(jdir,^d);
1120 ! qvec(z) at cell interface along r dimension
1121 if(stretched_dim(jdir) .and. stretch_uncentered) then
1122 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1123 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1124 else
1125 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1126 end if
1127 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,kdir))&
1128 /block%surface(ixo^s,idir)
1129 end if
1130 else ! idir==z_
1131 if(jdir<kdir) then
1132 ! idir=z,jdir=r,kdir=phi
1133 !! integral along r dimension
1134 hxo^l=ixo^l-kr(kdir,^d);
1135 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1136 jxc^l=ixc^l+kr(kdir,^d);
1137 ! qvec(r) at cell interface along phi dimension
1138 if(stretched_dim(kdir) .and. stretch_uncentered) then
1139 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1140 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1141 else
1142 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1143 end if
1144 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
1145 !! integral along phi dimension
1146 hxo^l=ixo^l-kr(jdir,^d);
1147 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1148 jxc^l=ixc^l+kr(jdir,^d);
1149 ! qvec(phi) at cell interface along r dimension
1150 if(stretched_dim(jdir) .and. stretch_uncentered) then
1151 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1152 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1153 else
1154 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1155 end if
1156 ! r coordinate at cell interface along r dimension
1157 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
1158 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(ixo^s)*tmp(ixo^s)-xc(hxo^s)*tmp(hxo^s))*block%dx(ixo^s,kdir))&
1159 /block%surface(ixo^s,idir)
1160 end if
1161 end if
1162 if(idir<idirmin)idirmin=idir
1163 endif
1164 enddo; enddo; enddo;
1165 case default
1166 call mpistop('no such curl evaluator')
1167 end select
1168 case default
1169 call mpistop('not possible to calculate curl')
1170 end select
1171
1172 end subroutine curlvector
1173
1174 !> Calculate idim transverse components of curl of a vector qvec within ixL
1175 !> Options to
1176 !> employ standard second order CD evaluations
1177 !> use Gauss theorem for non-Cartesian grids
1178 !> use Stokes theorem for non-Cartesian grids
1179 subroutine curlvector_trans(qvec,qvecc,ixI^L,ixO^L,curlvec,idim,idirmin,idirmin0,ndir0)
1181
1182 integer, intent(in) :: ixI^L,ixO^L
1183 integer, intent(in) :: idim, ndir0, idirmin0
1184 integer, intent(inout) :: idirmin
1185 double precision, intent(in) :: qvec(ixI^S,1:ndir0),qvecc(ixI^S,1:ndir0)
1186 double precision, intent(inout) :: curlvec(ixI^S,idirmin0:3)
1187
1188 double precision :: invdx(1:ndim)
1189 double precision :: tmp(ixI^S),tmp2(ixI^S),xC(ixI^S),surface(ixI^S)
1190 integer :: ixA^L,ixC^L,jxC^L,idir,jdir,kdir,hxO^L,jxO^L
1191
1192 ! Calculate curl within ixL: CurlV_i=eps_ijk*d_j V_k
1193 ! Curl can have components (idirmin:3)
1194 ! Determine exact value of idirmin while doing the loop.
1195
1196 idirmin=4
1197 curlvec(ixo^s,idirmin0:3)=zero
1198 ! Second order, stencil width is one
1199 ixa^l=ixo^l^ladd1;
1200
1201 ! all non-Cartesian cases
1202 select case(coordinate)
1203 case(cartesian) ! Cartesian grids
1204 invdx=1.d0/dxlevel
1205 do idir=idirmin0,3
1206 do jdir=1,ndim; do kdir=1,ndir0
1207 if(lvc(idir,jdir,kdir)/=0)then
1208 if(jdir/=idim) then
1209 tmp(ixi^s)=qvec(ixi^s,kdir)
1210 hxo^l=ixo^l-kr(jdir,^d);
1211 jxo^l=ixo^l+kr(jdir,^d);
1212 else
1213 ! use two cell-center values for gradient at interface of the two cells
1214 ! because left and right neighbor interface values is unavailable at block boundary faces
1215 tmp(ixi^s)=qvecc(ixi^s,kdir)
1216 hxo^l=ixo^l;
1217 jxo^l=ixo^l+kr(jdir,^d);
1218 end if
1219 ! second order centered differencing
1220 tmp(ixo^s)=half*(tmp(jxo^s)-tmp(hxo^s))*invdx(jdir)
1221 if(lvc(idir,jdir,kdir)==1)then
1222 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp(ixo^s)
1223 else
1224 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp(ixo^s)
1225 endif
1226 if(idir<idirmin)idirmin=idir
1227 endif
1228 enddo; enddo;
1229 end do
1230 case(cartesian_stretched) ! stretched Cartesian grids
1231 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1232 if(lvc(idir,jdir,kdir)/=0)then
1233 select case(type_curl)
1234 case(central)
1235 tmp(ixi^s)=qvec(ixi^s,kdir)
1236 hxo^l=ixo^l-kr(jdir,^d);
1237 jxo^l=ixo^l+kr(jdir,^d);
1238 ! second order centered differencing
1239 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,jdir)-block%x(hxo^s,jdir))
1240 case(gaussbased)
1241 hxo^l=ixo^l-kr(jdir,^d);
1242 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1243 jxc^l=ixc^l+kr(jdir,^d);
1244 tmp(ixc^s)=block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1245 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
1246 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%dvolume(ixo^s)
1247 case(stokesbased)
1248 hxo^l=ixo^l-kr(jdir,^d);
1249 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1250 jxc^l=ixc^l+kr(jdir,^d);
1251 if(kdir<=ndim)then
1252 tmp(ixc^s)=block%ds(ixo^s,kdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1253 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
1254 else
1255 tmp(ixc^s)=(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1256 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
1257 endif
1258 if(idir<=ndim)then
1259 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%surface(ixo^s,idir)
1260 else ! essentially for 2.5D case, idir=3 and jdir,kdir<=2
1261 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/(block%ds(ixo^s,jdir)*block%ds(ixo^s,kdir))
1262 endif
1263 case default
1264 call mpistop('no such curl evaluator')
1265 end select
1266 if(lvc(idir,jdir,kdir)==1)then
1267 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1268 else
1269 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1270 endif
1271 if(idir<idirmin)idirmin=idir
1272 endif
1273 enddo; enddo; enddo;
1274 case(spherical) ! possibly stretched spherical grids
1275 select case(type_curl)
1276 case(central) ! ok for any dimensionality
1277 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1278 if(lvc(idir,jdir,kdir)/=0)then
1279 tmp(ixi^s)=qvec(ixi^s,kdir)
1280 hxo^l=ixo^l-kr(jdir,^d);
1281 jxo^l=ixo^l+kr(jdir,^d);
1282 select case(jdir)
1283 case(1)
1284 tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1)
1285 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,1)-block%x(hxo^s,1))*block%x(ixo^s,1))
1286 {^nooned case(2)
1287 if(idir==1) tmp(ixa^s)=tmp(ixa^s)*dsin(block%x(ixa^s,2))
1288 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,2)-block%x(hxo^s,2))*block%x(ixo^s,1))
1289 if(idir==1) tmp2(ixo^s)=tmp2(ixo^s)/dsin(block%x(ixo^s,2))
1290 }
1291 {^ifthreed case(3)
1292 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,3)-block%x(hxo^s,3))*block%x(ixo^s,1)*dsin(block%x(ixo^s,2)))
1293 }
1294 end select
1295 if(lvc(idir,jdir,kdir)==1)then
1296 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1297 else
1298 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1299 endif
1300 if(idir<idirmin)idirmin=idir
1301 endif
1302 enddo; enddo; enddo;
1303 case(gaussbased)
1304 do idir=idirmin0,3;
1305 do jdir=1,ndim; do kdir=1,ndir0
1306 if(lvc(idir,jdir,kdir)/=0)then
1307 hxo^l=ixo^l-kr(jdir,^d);
1308 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1309 jxc^l=ixc^l+kr(jdir,^d);
1310 tmp(ixc^s)=block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1311 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
1312 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%dvolume(ixo^s)
1313 if(lvc(idir,jdir,kdir)==1)then
1314 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1315 else
1316 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1317 endif
1318 if(idir<idirmin)idirmin=idir
1319 endif
1320 enddo; enddo;
1321 ! geometric terms
1322 if(idir==2.and.phi_>0) curlvec(ixo^s,2)=curlvec(ixo^s,2)+qvec(ixo^s,phi_)/block%x(ixo^s,r_)
1323 {^nooned
1324 if(idir==phi_) curlvec(ixo^s,phi_)=curlvec(ixo^s,phi_)-qvec(ixo^s,2)/block%x(ixo^s,r_) &
1325 +qvec(ixo^s,r_)*dcos(block%x(ixo^s,2))/(block%x(ixo^s,r_)*dsin(block%x(ixo^s,2)))
1326 }
1327 enddo;
1328 case(stokesbased)
1329 !if(ndim<3) call mpistop("Stokesbased for 3D spherical only")
1330 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1331 if(lvc(idir,jdir,kdir)/=0)then
1332 select case(idir)
1333 case(1)
1334 if(jdir<kdir) then
1335 ! idir=1,jdir=2,kdir=3
1336 !! integral along 3rd dimension
1337 hxo^l=ixo^l-kr(jdir,^d);
1338 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1339 jxc^l=ixc^l+kr(jdir,^d);
1340 ! qvec(3) at cell interface along 2nd dimension
1341 if(stretched_dim(jdir) .and. stretch_uncentered) then
1342 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1343 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1344 else
1345 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1346 end if
1347 ! 2nd coordinate at cell interface along 2nd dimension
1348 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
1349 curlvec(ixo^s,idir)=(dsin(xc(ixo^s))*tmp(ixo^s)-&
1350 dsin(xc(hxo^s))*tmp(hxo^s))*block%dx(ixo^s,kdir)
1351 !! integral along 2nd dimension
1352 hxo^l=ixo^l-kr(kdir,^d);
1353 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1354 jxc^l=ixc^l+kr(kdir,^d);
1355 ! qvec(2) at cell interface along 3rd dimension
1356 if(stretched_dim(kdir) .and. stretch_uncentered) then
1357 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1358 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1359 else
1360 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1361 end if
1362 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir))&
1363 /block%surface(ixo^s,idir)*block%x(ixo^s,idir)
1364 end if
1365 case(2)
1366 if(jdir<kdir) then
1367 ! idir=2,jdir=1,kdir=3
1368 !! integral along 1st dimension
1369 hxo^l=ixo^l-kr(kdir,^d);
1370 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1371 jxc^l=ixc^l+kr(kdir,^d);
1372 ! qvec(1) at cell interface along 3rd dimension
1373 if(stretched_dim(kdir) .and. stretch_uncentered) then
1374 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1375 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1376 else
1377 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1378 end if
1379 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,1)
1380 !! integral along 3rd dimension
1381 hxo^l=ixo^l-kr(jdir,^d);
1382 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1383 jxc^l=ixc^l+kr(jdir,^d);
1384 ! qvec(3) at cell interface along 1st dimension
1385 if(stretched_dim(jdir) .and. stretch_uncentered) then
1386 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1387 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1388 else
1389 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1390 end if
1391 ! 1st coordinate at cell interface along 1st dimension
1392 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
1393 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(hxo^s)*tmp(hxo^s)-xc(ixo^s)*tmp(ixo^s))*&
1394 dsin(block%x(ixo^s,idir))*block%dx(ixo^s,kdir))/block%surface(ixo^s,idir)
1395 end if
1396 case(3)
1397 if(jdir<kdir) then
1398 ! idir=3,jdir=1,kdir=2
1399 !! integral along 1st dimension
1400 hxo^l=ixo^l-kr(kdir,^d);
1401 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1402 jxc^l=ixc^l+kr(kdir,^d);
1403 ! qvec(1) at cell interface along 2nd dimension
1404 if(stretched_dim(kdir) .and. stretch_uncentered) then
1405 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1406 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1407 else
1408 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1409 end if
1410 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
1411 !! integral along 2nd dimension
1412 hxo^l=ixo^l-kr(jdir,^d);
1413 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1414 jxc^l=ixc^l+kr(jdir,^d);
1415 ! qvec(2) at cell interface along 1st dimension
1416 if(stretched_dim(jdir) .and. stretch_uncentered) then
1417 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1418 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1419 else
1420 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1421 end if
1422 ! 1st coordinate at cell interface along 1st dimension
1423 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
1424 if(ndim==3) then
1425 surface(ixo^s)=block%surface(ixo^s,idir)
1426 else
1427 surface(ixo^s)=block%x(ixo^s,jdir)*block%dx(ixo^s,kdir)*block%dx(ixo^s,jdir)
1428 end if
1429 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(ixo^s)*tmp(ixo^s)-xc(hxo^s)*tmp(hxo^s))*block%dx(ixo^s,kdir))&
1430 /surface(ixo^s)
1431 end if
1432 end select
1433 if(idir<idirmin)idirmin=idir
1434 endif
1435 enddo; enddo; enddo;
1436 case default
1437 call mpistop('no such curl evaluator')
1438 end select
1439 case(cylindrical) ! possibly stretched cylindrical grids
1440 select case(type_curl)
1441 case(central) ! works for any dimensionality, polar/cylindrical
1442 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1443 if(lvc(idir,jdir,kdir)/=0)then
1444 tmp(ixi^s)=qvec(ixi^s,kdir)
1445 hxo^l=ixo^l-kr(jdir,^d);
1446 jxo^l=ixo^l+kr(jdir,^d);
1447 if(z_==3.or.z_==-1) then
1448 ! Case Polar_2D, Polar_2.5D or Polar_3D, i.e. R,phi,Z
1449 select case(jdir)
1450 case(1)
1451 if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
1452 ! computes d(R V_phi)/dR or d V_Z/dR
1453 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,1)-block%x(hxo^s,1))
1454 if(idir==z_) tmp2(ixo^s)=tmp2(ixo^s)/block%x(ixo^s,1) ! (1/R)*d(R V_phi)/dR
1455 {^nooned case(2)
1456 ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
1457 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,2)-block%x(hxo^s,2))*block%x(ixo^s,1))
1458 }
1459 {^ifthreed case(3)
1460 ! handles d V_phi/dZ or d V_R/dZ
1461 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,3)-block%x(hxo^s,3))
1462 }
1463 end select
1464 end if
1465 if(phi_==3.or.phi_==-1) then
1466 ! Case Cylindrical_2D, Cylindrical_2.5D or Cylindrical_3D, i.e. R,Z,phi
1467 select case(jdir)
1468 case(1)
1469 if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
1470 ! computes d(R V_phi)/dR or d V_Z/dR
1471 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,1)-block%x(hxo^s,1))
1472 if(idir==z_) tmp2(ixo^s)=tmp2(ixo^s)/block%x(ixo^s,1) ! (1/R)*d(R V_phi)/dR
1473 {^nooned case(2)
1474 ! handles d V_phi/dZ or d V_R/dZ
1475 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,2)-block%x(hxo^s,2))
1476 }
1477 {^ifthreed case(3)
1478 ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
1479 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((block%x(jxo^s,3)-block%x(hxo^s,3))*block%x(ixo^s,1))
1480 }
1481 end select
1482 end if
1483 if(lvc(idir,jdir,kdir)==1)then
1484 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1485 else
1486 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1487 endif
1488 if(idir<idirmin)idirmin=idir
1489 endif
1490 enddo; enddo; enddo;
1491 case(gaussbased) ! works for any dimensionality, polar/cylindrical
1492 if(ndim<2) call mpistop("Gaussbased for 2D, 2.5D or 3D polar or cylindrical only")
1493 do idir=idirmin0,3;
1494 do jdir=1,ndim; do kdir=1,ndir0
1495 if(lvc(idir,jdir,kdir)/=0)then
1496 hxo^l=ixo^l-kr(jdir,^d);
1497 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1498 jxc^l=ixc^l+kr(jdir,^d);
1499 tmp(ixc^s)=block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1500 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir)))
1501 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%dvolume(ixo^s)
1502 if(lvc(idir,jdir,kdir)==1)then
1503 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1504 else
1505 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1506 endif
1507 if(idir<idirmin)idirmin=idir
1508 endif
1509 enddo; enddo;
1510 ! geometric term from d e_R/d phi= e_phi for unit vectors e_R, e_phi
1511 ! but minus sign appears due to R,Z,phi ordering (?)
1512 ! note that in cylindrical 2D (RZ), phi_ is -1
1513 ! note that in polar 2D (Rphi), z_ is -1
1514 if((idir==phi_.or.(phi_==-1.and.idir==3)).and.z_>0) then
1515 ! cylindrical
1516 if( z_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-qvec(ixo^s,z_)/block%x(ixo^s,r_)
1517 ! polar
1518 if(phi_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+qvec(ixo^s,z_)/block%x(ixo^s,r_)
1519 endif
1520 enddo;
1521 case(stokesbased)
1522 if(ndim<3) call mpistop("Stokesbased for 3D cylindrical only")
1523 do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1524 if(lvc(idir,jdir,kdir)/=0)then
1525 if(idir==r_) then
1526 if(jdir==phi_) then
1527 ! idir=r,jdir=phi,kdir=z
1528 !! integral along z dimension
1529 hxo^l=ixo^l-kr(jdir,^d);
1530 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1531 jxc^l=ixc^l+kr(jdir,^d);
1532 ! qvec(z) at cell interface along phi dimension
1533 if(stretched_dim(jdir) .and. stretch_uncentered) then
1534 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1535 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1536 else
1537 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1538 end if
1539 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,kdir)
1540 !! integral along phi dimension
1541 hxo^l=ixo^l-kr(kdir,^d);
1542 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1543 jxc^l=ixc^l+kr(kdir,^d);
1544 ! qvec(phi) at cell interface along z dimension
1545 if(stretched_dim(kdir) .and. stretch_uncentered) then
1546 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1547 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1548 else
1549 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1550 end if
1551 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*block%x(ixo^s,idir)*block%dx(ixo^s,jdir))&
1552 /block%surface(ixo^s,idir)
1553 end if
1554 else if(idir==phi_) then
1555 if(jdir<kdir) then
1556 ! idir=phi,jdir=r,kdir=z
1557 !! integral along r dimension
1558 hxo^l=ixo^l-kr(kdir,^d);
1559 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1560 jxc^l=ixc^l+kr(kdir,^d);
1561 ! qvec(r) at cell interface along z dimension
1562 if(stretched_dim(kdir) .and. stretch_uncentered) then
1563 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1564 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1565 else
1566 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1567 end if
1568 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,jdir)
1569 !! integral along z dimension
1570 hxo^l=ixo^l-kr(jdir,^d);
1571 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1572 jxc^l=ixc^l+kr(jdir,^d);
1573 ! qvec(z) at cell interface along r dimension
1574 if(stretched_dim(jdir) .and. stretch_uncentered) then
1575 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1576 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1577 else
1578 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1579 end if
1580 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,kdir))&
1581 /block%surface(ixo^s,idir)
1582 end if
1583 else ! idir==z_
1584 if(jdir<kdir) then
1585 ! idir=z,jdir=r,kdir=phi
1586 !! integral along r dimension
1587 hxo^l=ixo^l-kr(kdir,^d);
1588 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1589 jxc^l=ixc^l+kr(kdir,^d);
1590 ! qvec(r) at cell interface along phi dimension
1591 if(stretched_dim(kdir) .and. stretch_uncentered) then
1592 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*block%dx(ixc^s,kdir)*&
1593 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(block%x(jxc^s,kdir)-block%x(ixc^s,kdir))
1594 else
1595 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1596 end if
1597 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
1598 !! integral along phi dimension
1599 hxo^l=ixo^l-kr(jdir,^d);
1600 ixcmin^d=hxomin^d;ixcmax^d=ixomax^d;
1601 jxc^l=ixc^l+kr(jdir,^d);
1602 ! qvec(phi) at cell interface along r dimension
1603 if(stretched_dim(jdir) .and. stretch_uncentered) then
1604 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*block%dx(ixc^s,jdir)*&
1605 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(block%x(jxc^s,jdir)-block%x(ixc^s,jdir))
1606 else
1607 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1608 end if
1609 ! r coordinate at cell interface along r dimension
1610 xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
1611 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(ixo^s)*tmp(ixo^s)-xc(hxo^s)*tmp(hxo^s))*block%dx(ixo^s,kdir))&
1612 /block%surface(ixo^s,idir)
1613 end if
1614 end if
1615 if(idir<idirmin)idirmin=idir
1616 endif
1617 enddo; enddo; enddo;
1618 case default
1619 call mpistop('no such curl evaluator')
1620 end select
1621 case default
1622 call mpistop('not possible to calculate curl')
1623 end select
1624
1625 end subroutine curlvector_trans
1626
1627end module mod_geometry
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
subroutine set_coordinate_system(geom)
Set the coordinate system to be used.
subroutine divvector(qvec, ixil, ixol, divq, nth_in)
integer coordinate
Definition mod_geometry.t:7
subroutine set_pole
subroutine get_surface_area(s, ixgl)
calculate area of surfaces of cells
integer, parameter spherical
integer, parameter cartesian
Definition mod_geometry.t:8
integer, parameter stokesbased
integer type_curl
integer, parameter cylindrical
integer, parameter gaussbased
integer, parameter cartesian_expansion
integer, parameter cartesian_stretched
Definition mod_geometry.t:9
subroutine curlvector(qvec, ixil, ixol, curlvec, idirmin, idirmin0, ndir0, fourthorder)
Calculate curl of a vector qvec within ixL Options to employ standard second order CD evaluations use...
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
subroutine divvectors(qvec, ixil, ixol, divq)
Calculate divergence of a vector qvec within ixL using limited extrapolation to cell edges.
subroutine gradientf(q, x, ixil, ixol, idir, gradq, nth_in, pm_in)
integer, parameter central
subroutine gradientl(q, ixil, ixol, idir, gradq)
subroutine putgridgeo(igrid)
Deallocate geometry-related variables.
subroutine curlvector_trans(qvec, qvecc, ixil, ixol, curlvec, idim, idirmin, idirmin0, ndir0)
Calculate idim transverse components of curl of a vector qvec within ixL Options to employ standard s...
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
character(len=std_len) geometry_name
integer, dimension(3, 3, 3) lvc
Levi-Civita tensor.
integer, dimension(3, 3) kr
Kronecker delta tensor.
logical stretch_uncentered
If true, adjust mod_geometry routines to account for grid stretching (but the flux computation will n...
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
logical, dimension(ndim) stretched_dim
True if a dimension is stretched.
integer, parameter unitterm
Unit for standard output.
logical, dimension(ndim) periodb
True for dimensions with periodic boundaries.
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer nghostcells
Number of ghost cells surrounding a grid.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
logical, dimension(2, ndim) poleb
Indicates whether there is a pole at a boundary.
integer, dimension(:), allocatable type_gradient_limiter
Type of slope limiter used for computing gradients or divergences, when typegrad or typediv are set t...
Module with slope/flux limiters.
Definition mod_limiter.t:2
integer, parameter limiter_ppm
Definition mod_limiter.t:27
double precision d
Definition mod_limiter.t:14
subroutine dwlimiter2(dwc, ixil, ixcl, idims, typelim, ldw, rdw, a2max)
Limit the centered dwC differences within ixC for iw in direction idim. The limiter is chosen accordi...
subroutine, public ppmlimitervar(ixil, ixl, idims, q, qct, qlc, qrc)
Definition mod_ppm.t:12
Module with all the methods that users can customize in AMRVAC.
procedure(set_surface), pointer usr_set_surface