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