24 character(len=*),
intent(in) :: geom
30 case (
"Cartesian",
"Cartesian_1D",
"Cartesian_2D",
"Cartesian_3D")
33 case (
"Cartesian_1D_expansion")
34 if (
ndim /= 1)
call mpistop(
"Geometry Cartesian_1D_expansion but ndim /= 1")
37 case (
"Cartesian_1.5D")
38 if (
ndim /= 1)
call mpistop(
"Geometry Cartesian_1.5D but ndim /= 1")
41 case (
"Cartesian_1.75D")
42 if (
ndim /= 1)
call mpistop(
"Geometry Cartesian_1.75D but ndim /= 1")
45 case (
"Cartesian_2.5D")
46 if (
ndim /= 2)
call mpistop(
"Geometry Cartesian_2.5D but ndim /= 2")
49 case (
"cylindrical",
"cylindrical_2D",
"cylindrical_3D")
55 case (
"cylindrical_2.5D")
56 if (
ndim /= 2)
call mpistop(
"Geometry cylindrical_2.5D but ndim /= 2")
62 case (
"polar",
"polar_2D",
"polar_3D")
69 if (
ndim /= 1)
call mpistop(
"Geometry polar_1.5D but ndim /= 1")
75 if (
ndim /= 2)
call mpistop(
"Geometry polar_2.5D but ndim /= 2")
81 case (
"spherical",
"spherical_2D",
"spherical_3D")
87 case (
"spherical_2.5D")
89 call mpistop(
"Geometry spherical_2.5D requires ndim == 2")
96 call mpistop(
"Unknown geometry specified")
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
112 "Will apply pi-periodic conditions at northpole!"
117 if(abs(xprobmax2-dpi)<smalldouble)
then
119 "Will apply pi-periodic conditions at southpole!"
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!")
132 if(abs(xprobmin1)<smalldouble)
then
134 write(unitterm,*)
"Will apply pi-periodic conditions at r=0"
139 write(unitterm,*)
"There is no cylindrical axis!"
150 integer,
intent(in) :: igrid
152 deallocate(ps(igrid)%surfaceC,ps(igrid)%surface,ps(igrid)%dvolume,ps(igrid)%dx,&
153 psc(igrid)%dx,ps(igrid)%ds,psc(igrid)%dvolume)
163 integer,
intent(in) :: ixG^L
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)
172 drs(ixg^s)=s%dx(ixg^s,1)
173 x(ixg^s,1)=s%x(ixg^s,1)
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!")
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)
185 s%surfaceC(ixgmin1-1:ixgmax1,1)=exp_factor_ext(ixgmin1-1:ixgmax1)
189 drs(ixg^s)=s%dx(ixg^s,1)
191 dx2(ixg^s)=s%dx(ixg^s,2)}
193 dx3(ixg^s)=s%dx(ixg^s,3)}
196 s%surfaceC(ixg^s,1)=1.d0
197 s%surface(ixg^s,1) =1.d0
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)
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)
213 {s%surfaceC(0^
d%ixG^s,^
d)=s%surfaceC(1^
d%ixG^s,^
d);\}
215 x(ixg^s,1)=s%x(ixg^s,1)
217 x(ixg^s,2)=s%x(ixg^s,2)
219 drs(ixg^s)=s%dx(ixg^s,1)
221 dx2(ixg^s)=s%dx(ixg^s,2)
224 dx3(ixg^s)=s%dx(ixg^s,3)
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)}
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)}
235 s%surfaceC(ixg^s,3)=x(ixg^s,1)*drs(ixg^s)*dx2(ixg^s)
239 s%surfaceC(0,1)=dabs(x(1,1)-half*drs(1))**2
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,&
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))
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,&
260 s%surfaceC(ixgmin1:ixgmax1,ixgmin2:ixgmax2,0,3)=&
261 s%surfaceC(ixgmin1:ixgmax1,ixgmin2:ixgmax2,1,3)
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)}
267 s%surface(ixg^s,2)=x(ixg^s,1)*drs(ixg^s)&
268 *dsin(x(ixg^s,2))}{^ifthreed*dx3(ixg^s)}
271 s%surface(ixg^s,3)=x(ixg^s,1)*drs(ixg^s)*dx2(ixg^s)}
274 x(ixg^s,1)=s%x(ixg^s,1)
275 drs(ixg^s)=s%dx(ixg^s,1)
277 dx2(ixg^s)=s%dx(ixg^s,2)}
279 dx3(ixg^s)=s%dx(ixg^s,3)}
281 s%surfaceC(ixg^s,1)=dabs(x(ixg^s,1)+half*drs(ixg^s)){^de&*
dx^de(ixg^s) }
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)}
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)
291 s%surfaceC(0,1)=dabs(x(1,1)-half*drs(1))
294 s%surfaceC(0,ixgmin2:ixgmax2,1)=dabs(x(1,ixgmin2:ixgmax2,1)-half*drs(1,&
295 ixgmin2:ixgmax2))*dx2(1,ixgmin2:ixgmax2)
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,&
303 {s%surfaceC(0^de%ixG^s,^de)=s%surfaceC(1^de%ixG^s,^de);\}
305 s%surface(ixg^s,1)=dabs(x(ixg^s,1)){^de&*
dx^de(ixg^s) }
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)}}
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)}
314 call mpistop(
"Sorry, coordinate unknown")
323 integer,
intent(in) :: ixI^L, ixO^L, idir
324 double precision,
intent(in) :: q(ixI^S)
325 double precision,
intent(inout) :: gradq(ixI^S)
327 integer :: jxO^L, hxO^L
329 hxo^l=ixo^l-
kr(idir,^
d);
330 jxo^l=ixo^l+
kr(idir,^
d);
333 gradq(ixo^s)=half*(q(jxo^s)-q(hxo^s))/
dxlevel(idir)
335 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(
block%x(jxo^s,idir)-
block%x(hxo^s,idir))
339 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((
block%x(jxo^s,1)-
block%x(hxo^s,1)))
342 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((
block%x(jxo^s,2)-
block%x(hxo^s,2))*
block%x(ixo^s,1))
346 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)))
353 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(
block%x(jxo^s,idir)-
block%x(hxo^s,idir))
356 call mpistop(
'Unknown geometry')
362 subroutine gradientx(q,x,ixI^L,ixO^L,idir,gradq,fourth_order)
365 integer,
intent(in) :: ixI^L, ixO^L, idir
366 double precision,
intent(in) :: q(ixI^S), x(ixI^S,1:ndim)
367 double precision,
intent(inout) :: gradq(ixI^S)
368 logical,
intent(in) :: fourth_order
369 integer :: jxO^L, hxO^L, kxO^L
372 jxo^l=ixo^l+
kr(idir,^
d);
375 if(fourth_order)
then
378 if(iximin^
d>kxomin^
d.or.iximax^
d<kxomax^
d|.or.) &
379 call mpistop(
"Error in gradientx: Non-conforming input limits")
380 hxo^l=ixo^l-
kr(idir,^
d);
381 jxo^l=ixo^l+
kr(idir,^
d);
382 kxo^l=ixo^l+
kr(idir,^
d)*2;
383 gradq(ixo^s)=(27.d0*(q(jxo^s)-q(ixo^s))-q(kxo^s)+q(hxo^s))/24.d0/
dxlevel(idir)
385 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/
dxlevel(idir)
388 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
392 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,1)-x(hxo^s,1)))
395 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,2)-x(hxo^s,2))*x(ixo^s,1))
399 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)))
404 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,
phi_)-x(hxo^s,
phi_))*x(ixo^s,
r_))
406 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
409 call mpistop(
'Unknown geometry')
418 integer,
intent(in) :: ixI^L, ixO^L, idir
419 double precision,
intent(in) :: q(ixI^S), x(ixI^S,1:ndim)
420 double precision,
intent(inout) :: gradq(ixI^S)
421 integer :: jxO^L, hxO^L, kxO^L
424 hxo^l=ixo^l-
kr(idir,^
d);
427 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/
dxlevel(idir)
429 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
433 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,1)-x(hxo^s,1)))
436 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,2)-x(hxo^s,2))*x(ixo^s,1))
440 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)))
445 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,
phi_)-x(hxo^s,
phi_))*x(ixo^s,
r_))
447 gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
450 call mpistop(
'Unknown geometry')
462 integer,
intent(in) :: ixI^L, ixO^L, idir
463 double precision,
intent(in) :: q(ixI^S)
464 double precision,
intent(inout) :: gradq(ixI^S)
465 double precision ,
dimension(ixI^S) :: qC,qL,qR,dqC,ldq,rdq
467 double precision :: x(ixI^S,1:ndim)
468 double precision :: invdx
469 integer :: hxO^L,ixC^L,jxC^L,gxC^L,hxC^L
471 x(ixi^s,1:ndim)=
block%x(ixi^s,1:ndim)
474 hxo^l=ixo^l-
kr(idir,^
d);
475 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
476 jxc^l=ixc^l+
kr(idir,^
d);
477 gxcmin^
d=ixcmin^
d-
kr(idir,^
d);gxcmax^
d=jxcmax^
d;
478 hxc^l=gxc^l+
kr(idir,^
d);
484 dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
486 ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
487 qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
494 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))*invdx
496 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/
block%dx(ixo^s,idir)
498 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/
block%dx(ixo^s,idir)
501 gradq(ixo^s)=gradq(ixo^s)/x(ixo^s,1)
504 gradq(ixo^s)=gradq(ixo^s)/(x(ixo^s,1)*dsin(x(ixo^s,2)))
508 gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/
block%dx(ixo^s,idir)
509 if(idir==
phi_) gradq(ixo^s)=gradq(ixo^s)/x(ixo^s,1)
515 subroutine divvector(qvec,ixI^L,ixO^L,divq,fourthorder,sixthorder)
517 integer,
intent(in) :: ixI^L,ixO^L
518 double precision,
intent(in) :: qvec(ixI^S,1:ndir)
519 double precision,
intent(inout) :: divq(ixI^S)
520 logical,
intent(in),
optional :: fourthorder
521 logical,
intent(in),
optional :: sixthorder
523 double precision :: qC(ixI^S), invdx(1:ndim)
524 integer :: jxO^L, hxO^L, ixC^L, jxC^L
525 integer :: idims, ix^L, gxO^L, kxO^L
526 integer :: lxO^L, fxO^L
527 logical :: use_4th_order
528 logical :: use_6th_order
530 use_4th_order = .false.
531 use_6th_order = .false.
532 if (
present(fourthorder)) use_4th_order = fourthorder
533 if (
present(sixthorder)) use_6th_order = sixthorder
534 if(use_4th_order .and. use_6th_order) &
535 call mpistop(
"divvector: using 4th and 6th order at the same time")
537 if(use_4th_order)
then
539 call mpistop(
"divvector: 4th order only supported for slab geometry")
542 else if(use_6th_order)
then
545 call mpistop(
"divvector: 6th order only supported for slab geometry")
552 if (iximin^
d>ixmin^
d.or.iximax^
d<ixmax^
d|.or.) &
553 call mpistop(
"Error in divvector: Non-conforming input limits")
560 if(use_6th_order)
then
561 lxo^l=ixo^l+3*
kr(idims,^
d);
562 kxo^l=ixo^l+2*
kr(idims,^
d);
563 jxo^l=ixo^l+
kr(idims,^
d);
564 hxo^l=ixo^l-
kr(idims,^
d);
565 gxo^l=ixo^l-2*
kr(idims,^
d);
566 fxo^l=ixo^l-3*
kr(idims,^
d);
567 divq(ixo^s)=divq(ixo^s)+&
568 (-qvec(fxo^s,idims)+9.d0*qvec(gxo^s,idims)-45.d0*qvec(hxo^s,idims)&
569 +qvec(lxo^s,idims)-9.d0*qvec(kxo^s,idims)+45.d0*qvec(jxo^s,idims))&
571 else if(use_4th_order)
then
573 kxo^l=ixo^l+2*
kr(idims,^
d);
574 jxo^l=ixo^l+
kr(idims,^
d);
575 hxo^l=ixo^l-
kr(idims,^
d);
576 gxo^l=ixo^l-2*
kr(idims,^
d);
577 divq(ixo^s)=divq(ixo^s)+&
578 (-qvec(kxo^s,idims)+8.d0*qvec(jxo^s,idims)-8.d0*&
579 qvec(hxo^s,idims)+qvec(gxo^s,idims))/(12.d0*
dxlevel(idims))
582 jxo^l=ixo^l+
kr(idims,^
d);
583 hxo^l=ixo^l-
kr(idims,^
d);
584 divq(ixo^s)=divq(ixo^s)+half*(qvec(jxo^s,idims) &
585 - qvec(hxo^s,idims))*invdx(idims)
590 hxo^l=ixo^l-
kr(idims,^
d);
591 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
592 jxc^l=ixc^l+
kr(idims,^
d);
595 qc(ixc^s)=
block%surfaceC(ixc^s,idims)*(qvec(ixc^s,idims)+0.5d0*
block%dx(ixc^s,idims)*&
596 (qvec(jxc^s,idims)-qvec(ixc^s,idims))/(
block%x(jxc^s,idims)-
block%x(ixc^s,idims)))
598 qc(ixc^s)=
block%surfaceC(ixc^s,idims)*half*(qvec(ixc^s,idims)+qvec(jxc^s,idims))
600 divq(ixo^s)=divq(ixo^s)+qc(ixo^s)-qc(hxo^s)
602 divq(ixo^s)=divq(ixo^s)/
block%dvolume(ixo^s)
613 integer,
intent(in) :: ixI^L,ixO^L
614 double precision,
intent(in) :: qvec(ixI^S,1:ndir)
615 double precision,
intent(inout) :: divq(ixI^S)
616 double precision,
dimension(ixI^S) :: qL,qR,dqC,ldq,rdq
618 double precision :: invdx(1:ndim)
619 integer :: hxO^L,ixC^L,jxC^L,idims,ix^L,gxC^L,hxC^L
623 if (iximin^
d>ixmin^
d.or.iximax^
d<ixmax^
d|.or.) &
624 call mpistop(
"Error in divvectorS: Non-conforming input limits")
629 hxo^l=ixo^l-
kr(idims,^
d);
630 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
631 jxc^l=ixc^l+
kr(idims,^
d);
632 gxcmin^
d=ixcmin^
d-
kr(idims,^
d);gxcmax^
d=jxcmax^
d;
633 hxc^l=gxc^l+
kr(idims,^
d);
634 qr(gxc^s) = qvec(hxc^s,idims)
635 ql(gxc^s) = qvec(gxc^s,idims)
637 dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
639 ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
640 qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
642 dqc(ixi^s)=qvec(ixi^s,idims)
647 divq(ixo^s)=divq(ixo^s)+half*(qr(ixo^s)-ql(hxo^s))*invdx(idims)
649 qr(ixc^s)=
block%surfaceC(ixc^s,idims)*qr(ixc^s)
650 ql(ixc^s)=
block%surfaceC(ixc^s,idims)*ql(ixc^s)
651 divq(ixo^s)=divq(ixo^s)+qr(ixo^s)-ql(hxo^s)
663 subroutine curlvector(qvec,ixI^L,ixO^L,curlvec,idirmin,idirmin0,ndir0,fourthorder)
666 integer,
intent(in) :: ixI^L,ixO^L
667 integer,
intent(in) :: ndir0, idirmin0
668 integer,
intent(inout) :: idirmin
669 double precision,
intent(in) :: qvec(ixI^S,1:ndir0)
670 double precision,
intent(inout) :: curlvec(ixI^S,idirmin0:3)
671 logical,
intent(in),
optional :: fourthorder
673 double precision :: invdx(1:ndim)
674 double precision :: tmp(ixI^S),tmp2(ixI^S),xC(ixI^S),surface(ixI^S)
675 integer :: ixA^L,ixC^L,jxC^L,idir,jdir,kdir,hxO^L,jxO^L,kxO^L,gxO^L
676 logical :: use_4th_order
682 use_4th_order = .false.
683 if (
present(fourthorder)) use_4th_order = fourthorder
685 if (use_4th_order)
then
687 call mpistop(
"curlvector: 4th order only supported for slab geometry")
695 if (iximin^
d>ixamin^
d.or.iximax^
d<ixamax^
d|.or.) &
696 call mpistop(
"Error in curlvector: Non-conforming input limits")
699 curlvec(ixo^s,idirmin0:3)=zero
705 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
706 if(
lvc(idir,jdir,kdir)/=0)
then
707 if (.not. use_4th_order)
then
709 jxo^l=ixo^l+
kr(jdir,^
d);
710 hxo^l=ixo^l-
kr(jdir,^
d);
711 tmp(ixo^s)=half*(qvec(jxo^s,kdir) &
712 - qvec(hxo^s,kdir))*invdx(jdir)
715 kxo^l=ixo^l+2*
kr(jdir,^
d);
716 jxo^l=ixo^l+
kr(jdir,^
d);
717 hxo^l=ixo^l-
kr(jdir,^
d);
718 gxo^l=ixo^l-2*
kr(jdir,^
d);
719 tmp(ixo^s)=(-qvec(kxo^s,kdir) + 8.0d0 * qvec(jxo^s,kdir) - 8.0d0 * &
720 qvec(hxo^s,kdir) + qvec(gxo^s,kdir))/(12.0d0 *
dxlevel(jdir))
722 if(
lvc(idir,jdir,kdir)==1)
then
723 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp(ixo^s)
725 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp(ixo^s)
727 if(idir<idirmin)idirmin=idir
731 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
732 if(
lvc(idir,jdir,kdir)/=0)
then
735 tmp(ixa^s)=qvec(ixa^s,kdir)
736 hxo^l=ixo^l-
kr(jdir,^
d);
737 jxo^l=ixo^l+
kr(jdir,^
d);
739 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,jdir)-
block%x(hxo^s,jdir))
741 hxo^l=ixo^l-
kr(jdir,^
d);
742 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
743 jxc^l=ixc^l+
kr(jdir,^
d);
744 tmp(ixc^s)=
block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
745 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
746 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%dvolume(ixo^s)
748 hxo^l=ixo^l-
kr(jdir,^
d);
749 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
750 jxc^l=ixc^l+
kr(jdir,^
d);
752 tmp(ixc^s)=
block%ds(ixo^s,kdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
753 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
755 tmp(ixc^s)=(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
756 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
759 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%surface(ixo^s,idir)
761 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/(
block%ds(ixo^s,jdir)*
block%ds(ixo^s,kdir))
764 call mpistop(
'no such curl evaluator')
766 if(
lvc(idir,jdir,kdir)==1)
then
767 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
769 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
771 if(idir<idirmin)idirmin=idir
777 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
778 if(
lvc(idir,jdir,kdir)/=0)
then
779 tmp(ixa^s)=qvec(ixa^s,kdir)
780 hxo^l=ixo^l-
kr(jdir,^
d);
781 jxo^l=ixo^l+
kr(jdir,^
d);
784 tmp(ixa^s)=tmp(ixa^s)*
block%x(ixa^s,1)
785 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,1)-
block%x(hxo^s,1))*
block%x(ixo^s,1))
787 if(idir==1) tmp(ixa^s)=tmp(ixa^s)*dsin(
block%x(ixa^s,2))
788 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,2)-
block%x(hxo^s,2))*
block%x(ixo^s,1))
789 if(idir==1) tmp2(ixo^s)=tmp2(ixo^s)/dsin(
block%x(ixo^s,2))
792 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)))
795 if(
lvc(idir,jdir,kdir)==1)
then
796 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
798 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
800 if(idir<idirmin)idirmin=idir
805 do jdir=1,ndim;
do kdir=1,ndir0
806 if(
lvc(idir,jdir,kdir)/=0)
then
807 hxo^l=ixo^l-
kr(jdir,^
d);
808 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
809 jxc^l=ixc^l+
kr(jdir,^
d);
810 tmp(ixc^s)=
block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
811 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
812 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%dvolume(ixo^s)
813 if(
lvc(idir,jdir,kdir)==1)
then
814 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
816 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
818 if(idir<idirmin)idirmin=idir
822 if(idir==2.and.
phi_>0) curlvec(ixo^s,2)=curlvec(ixo^s,2)+qvec(ixo^s,
phi_)/
block%x(ixo^s,
r_)
830 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
831 if(
lvc(idir,jdir,kdir)/=0)
then
837 hxo^l=ixo^l-
kr(jdir,^
d);
838 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
839 jxc^l=ixc^l+
kr(jdir,^
d);
842 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
843 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
845 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
848 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
849 curlvec(ixo^s,idir)=(dsin(xc(ixo^s))*tmp(ixo^s)-&
850 dsin(xc(hxo^s))*tmp(hxo^s))*
block%dx(ixo^s,kdir)
852 hxo^l=ixo^l-
kr(kdir,^
d);
853 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
854 jxc^l=ixc^l+
kr(kdir,^
d);
857 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
858 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
860 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
862 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,jdir))&
863 /
block%surface(ixo^s,idir)*
block%x(ixo^s,idir)
869 hxo^l=ixo^l-
kr(kdir,^
d);
870 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
871 jxc^l=ixc^l+
kr(kdir,^
d);
874 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
875 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
877 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
879 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*
block%dx(ixo^s,1)
881 hxo^l=ixo^l-
kr(jdir,^
d);
882 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
883 jxc^l=ixc^l+
kr(jdir,^
d);
886 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
887 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
889 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
892 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
893 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(hxo^s)*tmp(hxo^s)-xc(ixo^s)*tmp(ixo^s))*&
894 dsin(
block%x(ixo^s,idir))*
block%dx(ixo^s,kdir))/
block%surface(ixo^s,idir)
900 hxo^l=ixo^l-
kr(kdir,^
d);
901 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
902 jxc^l=ixc^l+
kr(kdir,^
d);
905 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
906 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
908 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
910 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,jdir)
912 hxo^l=ixo^l-
kr(jdir,^
d);
913 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
914 jxc^l=ixc^l+
kr(jdir,^
d);
917 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
918 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
920 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
923 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
925 surface(ixo^s)=
block%surface(ixo^s,idir)
927 surface(ixo^s)=
block%x(ixo^s,jdir)*
block%dx(ixo^s,kdir)*
block%dx(ixo^s,jdir)
929 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))&
933 if(idir<idirmin)idirmin=idir
937 call mpistop(
'no such curl evaluator')
942 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
943 if(
lvc(idir,jdir,kdir)/=0)
then
944 tmp(ixa^s)=qvec(ixa^s,kdir)
945 hxo^l=ixo^l-
kr(jdir,^
d);
946 jxo^l=ixo^l+
kr(jdir,^
d);
947 if(
z_==3.or.
z_==-1)
then
951 if(idir==
z_) tmp(ixa^s)=tmp(ixa^s)*
block%x(ixa^s,1)
953 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,1)-
block%x(hxo^s,1))
954 if(idir==
z_) tmp2(ixo^s)=tmp2(ixo^s)/
block%x(ixo^s,1)
957 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,2)-
block%x(hxo^s,2))*
block%x(ixo^s,1))
961 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,3)-
block%x(hxo^s,3))
969 if(idir==
z_) tmp(ixa^s)=tmp(ixa^s)*
block%x(ixa^s,1)
971 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,1)-
block%x(hxo^s,1))
972 if(idir==
z_) tmp2(ixo^s)=tmp2(ixo^s)/
block%x(ixo^s,1)
975 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,2)-
block%x(hxo^s,2))
979 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,3)-
block%x(hxo^s,3))*
block%x(ixo^s,1))
983 if(
lvc(idir,jdir,kdir)==1)
then
984 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
986 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
988 if(idir<idirmin)idirmin=idir
992 if(ndim<2)
call mpistop(
"Gaussbased for 2D, 2.5D or 3D polar or cylindrical only")
994 do jdir=1,ndim;
do kdir=1,ndir0
995 if(
lvc(idir,jdir,kdir)/=0)
then
996 hxo^l=ixo^l-
kr(jdir,^
d);
997 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
998 jxc^l=ixc^l+
kr(jdir,^
d);
999 tmp(ixc^s)=
block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1000 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
1001 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%dvolume(ixo^s)
1002 if(
lvc(idir,jdir,kdir)==1)
then
1003 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1005 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1007 if(idir<idirmin)idirmin=idir
1014 if((idir==
phi_.or.(
phi_==-1.and.idir==3)).and.
z_>0)
then
1016 if(
z_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-qvec(ixo^s,
z_)/
block%x(ixo^s,
r_)
1018 if(
phi_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+qvec(ixo^s,
z_)/
block%x(ixo^s,
r_)
1022 if(ndim<3)
call mpistop(
"Stokesbased for 3D cylindrical only")
1023 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
1024 if(
lvc(idir,jdir,kdir)/=0)
then
1029 hxo^l=ixo^l-
kr(jdir,^
d);
1030 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1031 jxc^l=ixc^l+
kr(jdir,^
d);
1034 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1035 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1037 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1039 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*
block%dx(ixo^s,kdir)
1041 hxo^l=ixo^l-
kr(kdir,^
d);
1042 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1043 jxc^l=ixc^l+
kr(kdir,^
d);
1046 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1047 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1049 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1051 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*
block%x(ixo^s,idir)*
block%dx(ixo^s,jdir))&
1052 /
block%surface(ixo^s,idir)
1054 else if(idir==
phi_)
then
1058 hxo^l=ixo^l-
kr(kdir,^
d);
1059 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1060 jxc^l=ixc^l+
kr(kdir,^
d);
1063 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1064 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1066 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1068 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*
block%dx(ixo^s,jdir)
1070 hxo^l=ixo^l-
kr(jdir,^
d);
1071 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1072 jxc^l=ixc^l+
kr(jdir,^
d);
1075 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1076 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1078 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1080 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,kdir))&
1081 /
block%surface(ixo^s,idir)
1087 hxo^l=ixo^l-
kr(kdir,^
d);
1088 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1089 jxc^l=ixc^l+
kr(kdir,^
d);
1092 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1093 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1095 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1097 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,jdir)
1099 hxo^l=ixo^l-
kr(jdir,^
d);
1100 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1101 jxc^l=ixc^l+
kr(jdir,^
d);
1104 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1105 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1107 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1110 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
1111 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))&
1112 /
block%surface(ixo^s,idir)
1115 if(idir<idirmin)idirmin=idir
1117 enddo; enddo; enddo;
1119 call mpistop(
'no such curl evaluator')
1122 call mpistop(
'not possible to calculate curl')
1135 integer,
intent(in) :: ixI^L,ixO^L
1136 integer,
intent(in) :: idim, ndir0, idirmin0
1137 integer,
intent(inout) :: idirmin
1138 double precision,
intent(in) :: qvec(ixI^S,1:ndir0),qvecc(ixI^S,1:ndir0)
1139 double precision,
intent(inout) :: curlvec(ixI^S,idirmin0:3)
1141 double precision :: invdx(1:ndim)
1142 double precision :: tmp(ixI^S),tmp2(ixI^S),xC(ixI^S),surface(ixI^S)
1143 integer :: ixA^L,ixC^L,jxC^L,idir,jdir,kdir,hxO^L,jxO^L
1150 curlvec(ixo^s,idirmin0:3)=zero
1159 do jdir=1,ndim;
do kdir=1,ndir0
1160 if(
lvc(idir,jdir,kdir)/=0)
then
1162 tmp(ixi^s)=qvec(ixi^s,kdir)
1163 hxo^l=ixo^l-
kr(jdir,^
d);
1164 jxo^l=ixo^l+
kr(jdir,^
d);
1168 tmp(ixi^s)=qvecc(ixi^s,kdir)
1170 jxo^l=ixo^l+
kr(jdir,^
d);
1173 tmp(ixo^s)=half*(tmp(jxo^s)-tmp(hxo^s))*invdx(jdir)
1174 if(
lvc(idir,jdir,kdir)==1)
then
1175 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp(ixo^s)
1177 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp(ixo^s)
1179 if(idir<idirmin)idirmin=idir
1184 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
1185 if(
lvc(idir,jdir,kdir)/=0)
then
1188 tmp(ixi^s)=qvec(ixi^s,kdir)
1189 hxo^l=ixo^l-
kr(jdir,^
d);
1190 jxo^l=ixo^l+
kr(jdir,^
d);
1192 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,jdir)-
block%x(hxo^s,jdir))
1194 hxo^l=ixo^l-
kr(jdir,^
d);
1195 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1196 jxc^l=ixc^l+
kr(jdir,^
d);
1197 tmp(ixc^s)=
block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1198 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
1199 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%dvolume(ixo^s)
1201 hxo^l=ixo^l-
kr(jdir,^
d);
1202 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1203 jxc^l=ixc^l+
kr(jdir,^
d);
1205 tmp(ixc^s)=
block%ds(ixo^s,kdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1206 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
1208 tmp(ixc^s)=(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1209 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
1212 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%surface(ixo^s,idir)
1214 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/(
block%ds(ixo^s,jdir)*
block%ds(ixo^s,kdir))
1217 call mpistop(
'no such curl evaluator')
1219 if(
lvc(idir,jdir,kdir)==1)
then
1220 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1222 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1224 if(idir<idirmin)idirmin=idir
1226 enddo; enddo; enddo;
1230 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
1231 if(
lvc(idir,jdir,kdir)/=0)
then
1232 tmp(ixi^s)=qvec(ixi^s,kdir)
1233 hxo^l=ixo^l-
kr(jdir,^
d);
1234 jxo^l=ixo^l+
kr(jdir,^
d);
1237 tmp(ixa^s)=tmp(ixa^s)*
block%x(ixa^s,1)
1238 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,1)-
block%x(hxo^s,1))*
block%x(ixo^s,1))
1240 if(idir==1) tmp(ixa^s)=tmp(ixa^s)*dsin(
block%x(ixa^s,2))
1241 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,2)-
block%x(hxo^s,2))*
block%x(ixo^s,1))
1242 if(idir==1) tmp2(ixo^s)=tmp2(ixo^s)/dsin(
block%x(ixo^s,2))
1245 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)))
1248 if(
lvc(idir,jdir,kdir)==1)
then
1249 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1251 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1253 if(idir<idirmin)idirmin=idir
1255 enddo; enddo; enddo;
1258 do jdir=1,ndim;
do kdir=1,ndir0
1259 if(
lvc(idir,jdir,kdir)/=0)
then
1260 hxo^l=ixo^l-
kr(jdir,^
d);
1261 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1262 jxc^l=ixc^l+
kr(jdir,^
d);
1263 tmp(ixc^s)=
block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1264 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
1265 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%dvolume(ixo^s)
1266 if(
lvc(idir,jdir,kdir)==1)
then
1267 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1269 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1271 if(idir<idirmin)idirmin=idir
1275 if(idir==2.and.
phi_>0) curlvec(ixo^s,2)=curlvec(ixo^s,2)+qvec(ixo^s,
phi_)/
block%x(ixo^s,
r_)
1283 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
1284 if(
lvc(idir,jdir,kdir)/=0)
then
1290 hxo^l=ixo^l-
kr(jdir,^
d);
1291 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1292 jxc^l=ixc^l+
kr(jdir,^
d);
1295 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1296 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1298 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1301 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
1302 curlvec(ixo^s,idir)=(dsin(xc(ixo^s))*tmp(ixo^s)-&
1303 dsin(xc(hxo^s))*tmp(hxo^s))*
block%dx(ixo^s,kdir)
1305 hxo^l=ixo^l-
kr(kdir,^
d);
1306 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1307 jxc^l=ixc^l+
kr(kdir,^
d);
1310 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1311 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1313 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1315 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,jdir))&
1316 /
block%surface(ixo^s,idir)*
block%x(ixo^s,idir)
1322 hxo^l=ixo^l-
kr(kdir,^
d);
1323 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1324 jxc^l=ixc^l+
kr(kdir,^
d);
1327 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1328 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1330 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1332 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*
block%dx(ixo^s,1)
1334 hxo^l=ixo^l-
kr(jdir,^
d);
1335 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1336 jxc^l=ixc^l+
kr(jdir,^
d);
1339 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1340 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1342 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1345 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
1346 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(xc(hxo^s)*tmp(hxo^s)-xc(ixo^s)*tmp(ixo^s))*&
1347 dsin(
block%x(ixo^s,idir))*
block%dx(ixo^s,kdir))/
block%surface(ixo^s,idir)
1353 hxo^l=ixo^l-
kr(kdir,^
d);
1354 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1355 jxc^l=ixc^l+
kr(kdir,^
d);
1358 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1359 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1361 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1363 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,jdir)
1365 hxo^l=ixo^l-
kr(jdir,^
d);
1366 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1367 jxc^l=ixc^l+
kr(jdir,^
d);
1370 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1371 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1373 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1376 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
1378 surface(ixo^s)=
block%surface(ixo^s,idir)
1380 surface(ixo^s)=
block%x(ixo^s,jdir)*
block%dx(ixo^s,kdir)*
block%dx(ixo^s,jdir)
1382 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))&
1386 if(idir<idirmin)idirmin=idir
1388 enddo; enddo; enddo;
1390 call mpistop(
'no such curl evaluator')
1395 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
1396 if(
lvc(idir,jdir,kdir)/=0)
then
1397 tmp(ixi^s)=qvec(ixi^s,kdir)
1398 hxo^l=ixo^l-
kr(jdir,^
d);
1399 jxo^l=ixo^l+
kr(jdir,^
d);
1400 if(
z_==3.or.
z_==-1)
then
1404 if(idir==
z_) tmp(ixa^s)=tmp(ixa^s)*
block%x(ixa^s,1)
1406 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,1)-
block%x(hxo^s,1))
1407 if(idir==
z_) tmp2(ixo^s)=tmp2(ixo^s)/
block%x(ixo^s,1)
1410 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,2)-
block%x(hxo^s,2))*
block%x(ixo^s,1))
1414 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,3)-
block%x(hxo^s,3))
1422 if(idir==
z_) tmp(ixa^s)=tmp(ixa^s)*
block%x(ixa^s,1)
1424 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,1)-
block%x(hxo^s,1))
1425 if(idir==
z_) tmp2(ixo^s)=tmp2(ixo^s)/
block%x(ixo^s,1)
1428 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(
block%x(jxo^s,2)-
block%x(hxo^s,2))
1432 tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/((
block%x(jxo^s,3)-
block%x(hxo^s,3))*
block%x(ixo^s,1))
1436 if(
lvc(idir,jdir,kdir)==1)
then
1437 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1439 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1441 if(idir<idirmin)idirmin=idir
1443 enddo; enddo; enddo;
1445 if(ndim<2)
call mpistop(
"Gaussbased for 2D, 2.5D or 3D polar or cylindrical only")
1447 do jdir=1,ndim;
do kdir=1,ndir0
1448 if(
lvc(idir,jdir,kdir)/=0)
then
1449 hxo^l=ixo^l-
kr(jdir,^
d);
1450 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1451 jxc^l=ixc^l+
kr(jdir,^
d);
1452 tmp(ixc^s)=
block%surfaceC(ixc^s,jdir)*(qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1453 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir)))
1454 tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/
block%dvolume(ixo^s)
1455 if(
lvc(idir,jdir,kdir)==1)
then
1456 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1458 curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1460 if(idir<idirmin)idirmin=idir
1467 if((idir==
phi_.or.(
phi_==-1.and.idir==3)).and.
z_>0)
then
1469 if(
z_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-qvec(ixo^s,
z_)/
block%x(ixo^s,
r_)
1471 if(
phi_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+qvec(ixo^s,
z_)/
block%x(ixo^s,
r_)
1475 if(ndim<3)
call mpistop(
"Stokesbased for 3D cylindrical only")
1476 do idir=idirmin0,3;
do jdir=1,ndim;
do kdir=1,ndir0
1477 if(
lvc(idir,jdir,kdir)/=0)
then
1482 hxo^l=ixo^l-
kr(jdir,^
d);
1483 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1484 jxc^l=ixc^l+
kr(jdir,^
d);
1487 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1488 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1490 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1492 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*
block%dx(ixo^s,kdir)
1494 hxo^l=ixo^l-
kr(kdir,^
d);
1495 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1496 jxc^l=ixc^l+
kr(kdir,^
d);
1499 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1500 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1502 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1504 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*
block%x(ixo^s,idir)*
block%dx(ixo^s,jdir))&
1505 /
block%surface(ixo^s,idir)
1507 else if(idir==
phi_)
then
1511 hxo^l=ixo^l-
kr(kdir,^
d);
1512 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1513 jxc^l=ixc^l+
kr(kdir,^
d);
1516 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1517 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1519 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1521 curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*
block%dx(ixo^s,jdir)
1523 hxo^l=ixo^l-
kr(jdir,^
d);
1524 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1525 jxc^l=ixc^l+
kr(jdir,^
d);
1528 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1529 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1531 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1533 curlvec(ixo^s,idir)=(curlvec(ixo^s,idir)+(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,kdir))&
1534 /
block%surface(ixo^s,idir)
1540 hxo^l=ixo^l-
kr(kdir,^
d);
1541 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1542 jxc^l=ixc^l+
kr(kdir,^
d);
1545 tmp(ixc^s)=qvec(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,kdir)*&
1546 (qvec(jxc^s,jdir)-qvec(ixc^s,jdir))/(
block%x(jxc^s,kdir)-
block%x(ixc^s,kdir))
1548 tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1550 curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*
block%dx(ixo^s,jdir)
1552 hxo^l=ixo^l-
kr(jdir,^
d);
1553 ixcmin^
d=hxomin^
d;ixcmax^
d=ixomax^
d;
1554 jxc^l=ixc^l+
kr(jdir,^
d);
1557 tmp(ixc^s)=qvec(ixc^s,kdir)+0.5d0*
block%dx(ixc^s,jdir)*&
1558 (qvec(jxc^s,kdir)-qvec(ixc^s,kdir))/(
block%x(jxc^s,jdir)-
block%x(ixc^s,jdir))
1560 tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1563 xc(ixc^s)=
block%x(ixc^s,jdir)+0.5d0*
block%dx(ixc^s,jdir)
1564 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))&
1565 /
block%surface(ixo^s,idir)
1568 if(idir<idirmin)idirmin=idir
1570 enddo; enddo; enddo;
1572 call mpistop(
'no such curl evaluator')
1575 call mpistop(
'not possible to calculate curl')
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with geometry-related routines (e.g., divergence, curl)
subroutine set_coordinate_system(geom)
Set the coordinate system to be used.
subroutine gradientq(q, x, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q in direction idir at cell interfaces.
subroutine divvectors(qvec, ixIL, ixOL, divq)
Calculate divergence of a vector qvec within ixL using limited extrapolation to cell edges.
subroutine get_surface_area(s, ixGL)
calculate area of surfaces of cells
integer, parameter spherical
integer, parameter cartesian
integer, parameter stokesbased
integer, parameter cylindrical
integer, parameter gaussbased
integer, parameter cartesian_expansion
integer, parameter cartesian_stretched
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
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...
integer, parameter central
subroutine gradients(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir first use limiter to go from cell cente...
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...
subroutine divvector(qvec, ixIL, ixOL, divq, fourthorder, sixthorder)
Calculate divergence of a vector qvec within ixL.
subroutine putgridgeo(igrid)
Deallocate geometry-related variables.
subroutine gradientx(q, x, ixIL, ixOL, idir, gradq, fourth_order)
Calculate gradient of a scalar q in direction idir at cell interfaces.
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
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.
integer, parameter limiter_ppm
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)
Module with all the methods that users can customize in AMRVAC.
procedure(set_surface), pointer usr_set_surface