MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
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 
19 contains
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 
171  case (cartesian_expansion)
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 
317  end subroutine get_surface_area
318 
319  !> Calculate gradient of a scalar q within ixL in direction idir
320  subroutine gradient(q,ixI^L,ixO^L,idir,gradq)
322 
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)
326 
327  integer :: jxO^L, hxO^L
328 
329  hxo^l=ixo^l-kr(idir,^d);
330  jxo^l=ixo^l+kr(idir,^d);
331  select case(coordinate)
332  case(cartesian)
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))
336  case(spherical)
337  select case(idir)
338  case(1)
339  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((block%x(jxo^s,1)-block%x(hxo^s,1)))
340  {^nooned
341  case(2)
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))
343  }
344  {^ifthreed
345  case(3)
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)))
347  }
348  end select
349  case(cylindrical)
350  if(idir==phi_) then
351  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((block%x(jxo^s,phi_)-block%x(hxo^s,phi_))*block%x(ixo^s,r_))
352  else
353  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(block%x(jxo^s,idir)-block%x(hxo^s,idir))
354  end if
355  case default
356  call mpistop('Unknown geometry')
357  end select
358 
359  end subroutine gradient
360 
361  !> Calculate gradient of a scalar q in direction idir at cell interfaces
362  subroutine gradientx(q,x,ixI^L,ixO^L,idir,gradq,fourth_order)
364 
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
370 
371  hxo^l=ixo^l;
372  jxo^l=ixo^l+kr(idir,^d);
373  select case(coordinate)
374  case(cartesian)
375  if(fourth_order) then
376  ! Fourth order, stencil width is two
377  kxo^l=ixo^l^ladd2;
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)
384  else
385  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/dxlevel(idir)
386  end if
388  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
389  case(spherical)
390  select case(idir)
391  case(1)
392  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,1)-x(hxo^s,1)))
393  {^nooned
394  case(2)
395  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,2)-x(hxo^s,2))*x(ixo^s,1))
396  }
397  {^ifthreed
398  case(3)
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)))
400  }
401  end select
402  case(cylindrical)
403  if(idir==phi_) then
404  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,phi_)-x(hxo^s,phi_))*x(ixo^s,r_))
405  else
406  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
407  end if
408  case default
409  call mpistop('Unknown geometry')
410  end select
411 
412  end subroutine gradientx
413 
414  !> Calculate gradient of a scalar q in direction idir at cell interfaces
415  subroutine gradientq(q,x,ixI^L,ixO^L,idir,gradq)
417 
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
422 
423  jxo^l=ixo^l;
424  hxo^l=ixo^l-kr(idir,^d);
425  select case(coordinate)
426  case(cartesian)
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))
430  case(spherical)
431  select case(idir)
432  case(1)
433  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,1)-x(hxo^s,1)))
434  {^nooned
435  case(2)
436  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,2)-x(hxo^s,2))*x(ixo^s,1))
437  }
438  {^ifthreed
439  case(3)
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)))
441  }
442  end select
443  case(cylindrical)
444  if(idir==phi_) then
445  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/((x(jxo^s,phi_)-x(hxo^s,phi_))*x(ixo^s,r_))
446  else
447  gradq(ixo^s)=(q(jxo^s)-q(hxo^s))/(x(jxo^s,idir)-x(hxo^s,idir))
448  end if
449  case default
450  call mpistop('Unknown geometry')
451  end select
452 
453  end subroutine gradientq
454 
455  !> Calculate gradient of a scalar q within ixL in direction idir
456  !> first use limiter to go from cell center to edge
457  subroutine gradients(q,ixI^L,ixO^L,idir,gradq)
459  use mod_limiter
460  use mod_ppm
461 
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
466 
467  double precision :: x(ixI^S,1:ndim)
468  double precision :: invdx
469  integer :: hxO^L,ixC^L,jxC^L,gxC^L,hxC^L
470 
471  x(ixi^s,1:ndim)=block%x(ixi^s,1:ndim)
472 
473  invdx=1.d0/dxlevel(idir)
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);
479 
480  ! set the gradient limiter here
481  qr(gxc^s) = q(hxc^s)
482  ql(gxc^s) = q(gxc^s)
483  if (type_gradient_limiter(block%level)/=limiter_ppm) then
484  dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
485  call dwlimiter2(dqc,ixi^l,gxc^l,idir,type_gradient_limiter(block%level),ldw=ldq,rdw=rdq)
486  ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
487  qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
488  else
489  call ppmlimitervar(ixi^l,ixo^l,idir,q,q,ql,qr)
490  endif
491 
492  select case(coordinate)
493  case(cartesian)
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)
497  case(spherical)
498  gradq(ixo^s)=(qr(ixo^s)-ql(hxo^s))/block%dx(ixo^s,idir)
499  select case(idir)
500  case(2)
501  gradq(ixo^s)=gradq(ixo^s)/x(ixo^s,1)
502  {^ifthreed
503  case(3)
504  gradq(ixo^s)=gradq(ixo^s)/(x(ixo^s,1)*dsin(x(ixo^s,2)))
505  }
506  end select
507  case(cylindrical)
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)
510  end select
511 
512  end subroutine gradients
513 
514  !> Calculate divergence of a vector qvec within ixL
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 !< Default: false
521  logical, intent(in), optional :: sixthorder !< Default: false
522 
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
529 
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")
536 
537  if(use_4th_order) then
538  if (.not. slab_uniform) &
539  call mpistop("divvector: 4th order only supported for slab geometry")
540  ! Fourth order, stencil width is two
541  ix^l=ixo^l^ladd2;
542  else if(use_6th_order) then
543  ! Sixth order, stencil width is three
544  if (.not. slab_uniform) &
545  call mpistop("divvector: 6th order only supported for slab geometry")
546  ix^l=ixo^l^ladd3;
547  else
548  ! Second order, stencil width is one
549  ix^l=ixo^l^ladd1;
550  end if
551 
552  if (iximin^d>ixmin^d.or.iximax^d<ixmax^d|.or.) &
553  call mpistop("Error in divvector: Non-conforming input limits")
554 
555  invdx=1.d0/dxlevel
556  divq(ixo^s)=0.0d0
557 
558  if (slab_uniform) then
559  do idims=1,ndim
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))&
570  /(60.d0*dxlevel(idims))
571  else if(use_4th_order) then
572  ! Use fourth order scheme
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))
580  else
581  ! Use second order scheme
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)
586  end if
587  end do
588  else
589  do idims=1,ndim
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);
593  if(stretched_dim(idims) .and. stretch_uncentered) then
594  ! linear interpolation at cell interface along stretched dimension
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)))
597  else
598  qc(ixc^s)=block%surfaceC(ixc^s,idims)*half*(qvec(ixc^s,idims)+qvec(jxc^s,idims))
599  end if
600  divq(ixo^s)=divq(ixo^s)+qc(ixo^s)-qc(hxo^s)
601  end do
602  divq(ixo^s)=divq(ixo^s)/block%dvolume(ixo^s)
603  end if
604  end subroutine divvector
605 
606  !> Calculate divergence of a vector qvec within ixL
607  !> using limited extrapolation to cell edges
608  subroutine divvectors(qvec,ixI^L,ixO^L,divq)
610  use mod_limiter
611  use mod_ppm
612 
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
617 
618  double precision :: invdx(1:ndim)
619  integer :: hxO^L,ixC^L,jxC^L,idims,ix^L,gxC^L,hxC^L
620 
621  ix^l=ixo^l^ladd2;
622 
623  if (iximin^d>ixmin^d.or.iximax^d<ixmax^d|.or.) &
624  call mpistop("Error in divvectorS: Non-conforming input limits")
625 
626  invdx=1.d0/dxlevel
627  divq(ixo^s)=zero
628  do idims=1,ndim
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)
636  if(type_gradient_limiter(block%level)/=limiter_ppm) then
637  dqc(gxc^s)= qr(gxc^s)-ql(gxc^s)
638  call dwlimiter2(dqc,ixi^l,gxc^l,idims,type_gradient_limiter(block%level),ldw=ldq,rdw=rdq)
639  ql(ixc^s) = ql(ixc^s) + half*ldq(ixc^s)
640  qr(ixc^s) = qr(ixc^s) - half*rdq(jxc^s)
641  else
642  dqc(ixi^s)=qvec(ixi^s,idims)
643  call ppmlimitervar(ixi^l,ixo^l,idims,dqc,dqc,ql,qr)
644  endif
645 
646  if (slab_uniform) then
647  divq(ixo^s)=divq(ixo^s)+half*(qr(ixo^s)-ql(hxo^s))*invdx(idims)
648  else
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)
652  end if
653  end do
654  if(.not.slab_uniform) divq(ixo^s)=divq(ixo^s)/block%dvolume(ixo^s)
655 
656  end subroutine divvectors
657 
658  !> Calculate curl of a vector qvec within ixL
659  !> Options to
660  !> employ standard second order CD evaluations
661  !> use Gauss theorem for non-Cartesian grids
662  !> use Stokes theorem for non-Cartesian grids
663  subroutine curlvector(qvec,ixI^L,ixO^L,curlvec,idirmin,idirmin0,ndir0,fourthorder)
665 
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 !< Default: false
672 
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
677 
678  ! Calculate curl within ixL: CurlV_i=eps_ijk*d_j V_k
679  ! Curl can have components (idirmin:3)
680  ! Determine exact value of idirmin while doing the loop.
681 
682  use_4th_order = .false.
683  if (present(fourthorder)) use_4th_order = fourthorder
684 
685  if (use_4th_order) then
686  if (.not. slab_uniform) &
687  call mpistop("curlvector: 4th order only supported for slab geometry")
688  ! Fourth order, stencil width is two
689  ixa^l=ixo^l^ladd2;
690  else
691  ! Second order, stencil width is one
692  ixa^l=ixo^l^ladd1;
693  end if
694 
695  if (iximin^d>ixamin^d.or.iximax^d<ixamax^d|.or.) &
696  call mpistop("Error in curlvector: Non-conforming input limits")
697 
698  idirmin=4
699  curlvec(ixo^s,idirmin0:3)=zero
700 
701  ! all non-Cartesian cases
702  select case(coordinate)
703  case(cartesian) ! Cartesian grids
704  invdx=1.d0/dxlevel
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
708  ! Use second order scheme
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)
713  else
714  ! Use fourth order scheme
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))
721  end if
722  if(lvc(idir,jdir,kdir)==1)then
723  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp(ixo^s)
724  else
725  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp(ixo^s)
726  endif
727  if(idir<idirmin)idirmin=idir
728  endif
729  enddo; enddo; enddo;
730  case(cartesian_stretched) ! stretched Cartesian grids
731  do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
732  if(lvc(idir,jdir,kdir)/=0)then
733  select case(type_curl)
734  case(central)
735  tmp(ixa^s)=qvec(ixa^s,kdir)
736  hxo^l=ixo^l-kr(jdir,^d);
737  jxo^l=ixo^l+kr(jdir,^d);
738  ! second order centered differencing
739  tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,jdir)-block%x(hxo^s,jdir))
740  case(gaussbased)
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)
747  case(stokesbased)
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);
751  if(kdir<=ndim)then
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)))
754  else
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)))
757  endif
758  if(idir<=ndim)then
759  tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%surface(ixo^s,idir)
760  else ! essentially for 2.5D case, idir=3 and jdir,kdir<=2
761  tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/(block%ds(ixo^s,jdir)*block%ds(ixo^s,kdir))
762  endif
763  case default
764  call mpistop('no such curl evaluator')
765  end select
766  if(lvc(idir,jdir,kdir)==1)then
767  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
768  else
769  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
770  endif
771  if(idir<idirmin)idirmin=idir
772  endif
773  enddo; enddo; enddo;
774  case(spherical) ! possibly stretched spherical grids
775  select case(type_curl)
776  case(central) ! ok for any dimensionality
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);
782  select case(jdir)
783  case(1)
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))
786  {^nooned case(2)
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))
790  }
791  {^ifthreed case(3)
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)))
793  }
794  end select
795  if(lvc(idir,jdir,kdir)==1)then
796  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
797  else
798  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
799  endif
800  if(idir<idirmin)idirmin=idir
801  endif
802  enddo; enddo; enddo;
803  case(gaussbased)
804  do idir=idirmin0,3;
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)
815  else
816  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
817  endif
818  if(idir<idirmin)idirmin=idir
819  endif
820  enddo; enddo;
821  ! geometric terms
822  if(idir==2.and.phi_>0) curlvec(ixo^s,2)=curlvec(ixo^s,2)+qvec(ixo^s,phi_)/block%x(ixo^s,r_)
823  {^nooned
824  if(idir==phi_) curlvec(ixo^s,phi_)=curlvec(ixo^s,phi_)-qvec(ixo^s,2)/block%x(ixo^s,r_) &
825  +qvec(ixo^s,r_)*dcos(block%x(ixo^s,2))/(block%x(ixo^s,r_)*dsin(block%x(ixo^s,2)))
826  }
827  enddo;
828  case(stokesbased)
829  !if(ndim<3) call mpistop("Stokesbased for 3D spherical only")
830  do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
831  if(lvc(idir,jdir,kdir)/=0)then
832  select case(idir)
833  case(1)
834  if(jdir<kdir) then
835  ! idir=1,jdir=2,kdir=3
836  !! integral along 3rd dimension
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);
840  ! qvec(3) at cell interface along 2nd dimension
841  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
844  else
845  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
846  end if
847  ! 2nd coordinate at cell interface along 2nd dimension
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)
851  !! integral along 2nd dimension
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);
855  ! qvec(2) at cell interface along 3rd dimension
856  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
859  else
860  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
861  end if
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)
864  end if
865  case(2)
866  if(jdir<kdir) then
867  ! idir=2,jdir=1,kdir=3
868  !! integral along 1st dimension
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);
872  ! qvec(1) at cell interface along 3rd dimension
873  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
876  else
877  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
878  end if
879  curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,1)
880  !! integral along 3rd dimension
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);
884  ! qvec(3) at cell interface along 1st dimension
885  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
888  else
889  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
890  end if
891  ! 1st coordinate at cell interface along 1st dimension
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)
895  end if
896  case(3)
897  if(jdir<kdir) then
898  ! idir=3,jdir=1,kdir=2
899  !! integral along 1st dimension
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);
903  ! qvec(1) at cell interface along 2nd dimension
904  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
907  else
908  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
909  end if
910  curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
911  !! integral along 2nd dimension
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);
915  ! qvec(2) at cell interface along 1st dimension
916  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
919  else
920  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
921  end if
922  ! 1st coordinate at cell interface along 1st dimension
923  xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
924  if(ndim==3) then
925  surface(ixo^s)=block%surface(ixo^s,idir)
926  else
927  surface(ixo^s)=block%x(ixo^s,jdir)*block%dx(ixo^s,kdir)*block%dx(ixo^s,jdir)
928  end if
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))&
930  /surface(ixo^s)
931  end if
932  end select
933  if(idir<idirmin)idirmin=idir
934  endif
935  enddo; enddo; enddo;
936  case default
937  call mpistop('no such curl evaluator')
938  end select
939  case(cylindrical) ! possibly stretched cylindrical grids
940  select case(type_curl)
941  case(central) ! works for any dimensionality, polar/cylindrical
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
948  ! Case Polar_2D, Polar_2.5D or Polar_3D, i.e. R,phi,Z
949  select case(jdir)
950  case(1)
951  if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
952  ! computes d(R V_phi)/dR or d V_Z/dR
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) ! (1/R)*d(R V_phi)/dR
955  {^nooned case(2)
956  ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
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))
958  }
959  {^ifthreed case(3)
960  ! handles d V_phi/dZ or d V_R/dZ
961  tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,3)-block%x(hxo^s,3))
962  }
963  end select
964  end if
965  if(phi_==3.or.phi_==-1) then
966  ! Case Cylindrical_2D, Cylindrical_2.5D or Cylindrical_3D, i.e. R,Z,phi
967  select case(jdir)
968  case(1)
969  if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
970  ! computes d(R V_phi)/dR or d V_Z/dR
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) ! (1/R)*d(R V_phi)/dR
973  {^nooned case(2)
974  ! handles d V_phi/dZ or d V_R/dZ
975  tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,2)-block%x(hxo^s,2))
976  }
977  {^ifthreed case(3)
978  ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
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))
980  }
981  end select
982  end if
983  if(lvc(idir,jdir,kdir)==1)then
984  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
985  else
986  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
987  endif
988  if(idir<idirmin)idirmin=idir
989  endif
990  enddo; enddo; enddo;
991  case(gaussbased) ! works for any dimensionality, polar/cylindrical
992  if(ndim<2) call mpistop("Gaussbased for 2D, 2.5D or 3D polar or cylindrical only")
993  do idir=idirmin0,3;
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)
1004  else
1005  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1006  endif
1007  if(idir<idirmin)idirmin=idir
1008  endif
1009  enddo; enddo;
1010  ! geometric term from d e_R/d phi= e_phi for unit vectors e_R, e_phi
1011  ! but minus sign appears due to R,Z,phi ordering (?)
1012  ! note that in cylindrical 2D (RZ), phi_ is -1
1013  ! note that in polar 2D (Rphi), z_ is -1
1014  if((idir==phi_.or.(phi_==-1.and.idir==3)).and.z_>0) then
1015  ! cylindrical
1016  if( z_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-qvec(ixo^s,z_)/block%x(ixo^s,r_)
1017  ! polar
1018  if(phi_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+qvec(ixo^s,z_)/block%x(ixo^s,r_)
1019  endif
1020  enddo;
1021  case(stokesbased)
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
1025  if(idir==r_) then
1026  if(jdir==phi_) then
1027  ! idir=r,jdir=phi,kdir=z
1028  !! integral along z dimension
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);
1032  ! qvec(z) at cell interface along phi dimension
1033  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1036  else
1037  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1038  end if
1039  curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,kdir)
1040  !! integral along phi dimension
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);
1044  ! qvec(phi) at cell interface along z dimension
1045  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1048  else
1049  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1050  end if
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)
1053  end if
1054  else if(idir==phi_) then
1055  if(jdir<kdir) then
1056  ! idir=phi,jdir=r,kdir=z
1057  !! integral along r dimension
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);
1061  ! qvec(r) at cell interface along z dimension
1062  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1065  else
1066  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1067  end if
1068  curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,jdir)
1069  !! integral along z dimension
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);
1073  ! qvec(z) at cell interface along r dimension
1074  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1077  else
1078  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1079  end if
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)
1082  end if
1083  else ! idir==z_
1084  if(jdir<kdir) then
1085  ! idir=z,jdir=r,kdir=phi
1086  !! integral along r dimension
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);
1090  ! qvec(r) at cell interface along phi dimension
1091  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1094  else
1095  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1096  end if
1097  curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
1098  !! integral along phi dimension
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);
1102  ! qvec(phi) at cell interface along r dimension
1103  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1106  else
1107  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1108  end if
1109  ! r coordinate at cell interface along r dimension
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)
1113  end if
1114  end if
1115  if(idir<idirmin)idirmin=idir
1116  endif
1117  enddo; enddo; enddo;
1118  case default
1119  call mpistop('no such curl evaluator')
1120  end select
1121  case default
1122  call mpistop('not possible to calculate curl')
1123  end select
1124 
1125  end subroutine curlvector
1126 
1127  !> Calculate idim transverse components of curl of a vector qvec within ixL
1128  !> Options to
1129  !> employ standard second order CD evaluations
1130  !> use Gauss theorem for non-Cartesian grids
1131  !> use Stokes theorem for non-Cartesian grids
1132  subroutine curlvector_trans(qvec,qvecc,ixI^L,ixO^L,curlvec,idim,idirmin,idirmin0,ndir0)
1134 
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)
1140 
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
1144 
1145  ! Calculate curl within ixL: CurlV_i=eps_ijk*d_j V_k
1146  ! Curl can have components (idirmin:3)
1147  ! Determine exact value of idirmin while doing the loop.
1148 
1149  idirmin=4
1150  curlvec(ixo^s,idirmin0:3)=zero
1151  ! Second order, stencil width is one
1152  ixa^l=ixo^l^ladd1;
1153 
1154  ! all non-Cartesian cases
1155  select case(coordinate)
1156  case(cartesian) ! Cartesian grids
1157  invdx=1.d0/dxlevel
1158  do idir=idirmin0,3
1159  do jdir=1,ndim; do kdir=1,ndir0
1160  if(lvc(idir,jdir,kdir)/=0)then
1161  if(jdir/=idim) 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);
1165  else
1166  ! use two cell-center values for gradient at interface of the two cells
1167  ! because left and right neighbor interface values is unavailable at block boundary faces
1168  tmp(ixi^s)=qvecc(ixi^s,kdir)
1169  hxo^l=ixo^l;
1170  jxo^l=ixo^l+kr(jdir,^d);
1171  end if
1172  ! second order centered differencing
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)
1176  else
1177  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp(ixo^s)
1178  endif
1179  if(idir<idirmin)idirmin=idir
1180  endif
1181  enddo; enddo;
1182  end do
1183  case(cartesian_stretched) ! stretched Cartesian grids
1184  do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1185  if(lvc(idir,jdir,kdir)/=0)then
1186  select case(type_curl)
1187  case(central)
1188  tmp(ixi^s)=qvec(ixi^s,kdir)
1189  hxo^l=ixo^l-kr(jdir,^d);
1190  jxo^l=ixo^l+kr(jdir,^d);
1191  ! second order centered differencing
1192  tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,jdir)-block%x(hxo^s,jdir))
1193  case(gaussbased)
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)
1200  case(stokesbased)
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);
1204  if(kdir<=ndim)then
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)))
1207  else
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)))
1210  endif
1211  if(idir<=ndim)then
1212  tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/block%surface(ixo^s,idir)
1213  else ! essentially for 2.5D case, idir=3 and jdir,kdir<=2
1214  tmp2(ixo^s)=(tmp(ixo^s)-tmp(hxo^s))/(block%ds(ixo^s,jdir)*block%ds(ixo^s,kdir))
1215  endif
1216  case default
1217  call mpistop('no such curl evaluator')
1218  end select
1219  if(lvc(idir,jdir,kdir)==1)then
1220  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1221  else
1222  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1223  endif
1224  if(idir<idirmin)idirmin=idir
1225  endif
1226  enddo; enddo; enddo;
1227  case(spherical) ! possibly stretched spherical grids
1228  select case(type_curl)
1229  case(central) ! ok for any dimensionality
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);
1235  select case(jdir)
1236  case(1)
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))
1239  {^nooned case(2)
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))
1243  }
1244  {^ifthreed case(3)
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)))
1246  }
1247  end select
1248  if(lvc(idir,jdir,kdir)==1)then
1249  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1250  else
1251  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1252  endif
1253  if(idir<idirmin)idirmin=idir
1254  endif
1255  enddo; enddo; enddo;
1256  case(gaussbased)
1257  do idir=idirmin0,3;
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)
1268  else
1269  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1270  endif
1271  if(idir<idirmin)idirmin=idir
1272  endif
1273  enddo; enddo;
1274  ! geometric terms
1275  if(idir==2.and.phi_>0) curlvec(ixo^s,2)=curlvec(ixo^s,2)+qvec(ixo^s,phi_)/block%x(ixo^s,r_)
1276  {^nooned
1277  if(idir==phi_) curlvec(ixo^s,phi_)=curlvec(ixo^s,phi_)-qvec(ixo^s,2)/block%x(ixo^s,r_) &
1278  +qvec(ixo^s,r_)*dcos(block%x(ixo^s,2))/(block%x(ixo^s,r_)*dsin(block%x(ixo^s,2)))
1279  }
1280  enddo;
1281  case(stokesbased)
1282  !if(ndim<3) call mpistop("Stokesbased for 3D spherical only")
1283  do idir=idirmin0,3; do jdir=1,ndim; do kdir=1,ndir0
1284  if(lvc(idir,jdir,kdir)/=0)then
1285  select case(idir)
1286  case(1)
1287  if(jdir<kdir) then
1288  ! idir=1,jdir=2,kdir=3
1289  !! integral along 3rd dimension
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);
1293  ! qvec(3) at cell interface along 2nd dimension
1294  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1297  else
1298  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1299  end if
1300  ! 2nd coordinate at cell interface along 2nd dimension
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)
1304  !! integral along 2nd dimension
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);
1308  ! qvec(2) at cell interface along 3rd dimension
1309  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1312  else
1313  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1314  end if
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)
1317  end if
1318  case(2)
1319  if(jdir<kdir) then
1320  ! idir=2,jdir=1,kdir=3
1321  !! integral along 1st dimension
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);
1325  ! qvec(1) at cell interface along 3rd dimension
1326  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1329  else
1330  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1331  end if
1332  curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,1)
1333  !! integral along 3rd dimension
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);
1337  ! qvec(3) at cell interface along 1st dimension
1338  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1341  else
1342  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1343  end if
1344  ! 1st coordinate at cell interface along 1st dimension
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)
1348  end if
1349  case(3)
1350  if(jdir<kdir) then
1351  ! idir=3,jdir=1,kdir=2
1352  !! integral along 1st dimension
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);
1356  ! qvec(1) at cell interface along 2nd dimension
1357  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1360  else
1361  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1362  end if
1363  curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
1364  !! integral along 2nd dimension
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);
1368  ! qvec(2) at cell interface along 1st dimension
1369  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1372  else
1373  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1374  end if
1375  ! 1st coordinate at cell interface along 1st dimension
1376  xc(ixc^s)=block%x(ixc^s,jdir)+0.5d0*block%dx(ixc^s,jdir)
1377  if(ndim==3) then
1378  surface(ixo^s)=block%surface(ixo^s,idir)
1379  else
1380  surface(ixo^s)=block%x(ixo^s,jdir)*block%dx(ixo^s,kdir)*block%dx(ixo^s,jdir)
1381  end if
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))&
1383  /surface(ixo^s)
1384  end if
1385  end select
1386  if(idir<idirmin)idirmin=idir
1387  endif
1388  enddo; enddo; enddo;
1389  case default
1390  call mpistop('no such curl evaluator')
1391  end select
1392  case(cylindrical) ! possibly stretched cylindrical grids
1393  select case(type_curl)
1394  case(central) ! works for any dimensionality, polar/cylindrical
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
1401  ! Case Polar_2D, Polar_2.5D or Polar_3D, i.e. R,phi,Z
1402  select case(jdir)
1403  case(1)
1404  if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
1405  ! computes d(R V_phi)/dR or d V_Z/dR
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) ! (1/R)*d(R V_phi)/dR
1408  {^nooned case(2)
1409  ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
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))
1411  }
1412  {^ifthreed case(3)
1413  ! handles d V_phi/dZ or d V_R/dZ
1414  tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,3)-block%x(hxo^s,3))
1415  }
1416  end select
1417  end if
1418  if(phi_==3.or.phi_==-1) then
1419  ! Case Cylindrical_2D, Cylindrical_2.5D or Cylindrical_3D, i.e. R,Z,phi
1420  select case(jdir)
1421  case(1)
1422  if(idir==z_) tmp(ixa^s)=tmp(ixa^s)*block%x(ixa^s,1) ! R V_phi
1423  ! computes d(R V_phi)/dR or d V_Z/dR
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) ! (1/R)*d(R V_phi)/dR
1426  {^nooned case(2)
1427  ! handles d V_phi/dZ or d V_R/dZ
1428  tmp2(ixo^s)=(tmp(jxo^s)-tmp(hxo^s))/(block%x(jxo^s,2)-block%x(hxo^s,2))
1429  }
1430  {^ifthreed case(3)
1431  ! handles (1/R)d V_Z/dphi or (1/R)d V_R/dphi
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))
1433  }
1434  end select
1435  end if
1436  if(lvc(idir,jdir,kdir)==1)then
1437  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+tmp2(ixo^s)
1438  else
1439  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1440  endif
1441  if(idir<idirmin)idirmin=idir
1442  endif
1443  enddo; enddo; enddo;
1444  case(gaussbased) ! works for any dimensionality, polar/cylindrical
1445  if(ndim<2) call mpistop("Gaussbased for 2D, 2.5D or 3D polar or cylindrical only")
1446  do idir=idirmin0,3;
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)
1457  else
1458  curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-tmp2(ixo^s)
1459  endif
1460  if(idir<idirmin)idirmin=idir
1461  endif
1462  enddo; enddo;
1463  ! geometric term from d e_R/d phi= e_phi for unit vectors e_R, e_phi
1464  ! but minus sign appears due to R,Z,phi ordering (?)
1465  ! note that in cylindrical 2D (RZ), phi_ is -1
1466  ! note that in polar 2D (Rphi), z_ is -1
1467  if((idir==phi_.or.(phi_==-1.and.idir==3)).and.z_>0) then
1468  ! cylindrical
1469  if( z_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)-qvec(ixo^s,z_)/block%x(ixo^s,r_)
1470  ! polar
1471  if(phi_==2) curlvec(ixo^s,idir)=curlvec(ixo^s,idir)+qvec(ixo^s,z_)/block%x(ixo^s,r_)
1472  endif
1473  enddo;
1474  case(stokesbased)
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
1478  if(idir==r_) then
1479  if(jdir==phi_) then
1480  ! idir=r,jdir=phi,kdir=z
1481  !! integral along z dimension
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);
1485  ! qvec(z) at cell interface along phi dimension
1486  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1489  else
1490  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1491  end if
1492  curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,kdir)
1493  !! integral along phi dimension
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);
1497  ! qvec(phi) at cell interface along z dimension
1498  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1501  else
1502  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1503  end if
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)
1506  end if
1507  else if(idir==phi_) then
1508  if(jdir<kdir) then
1509  ! idir=phi,jdir=r,kdir=z
1510  !! integral along r dimension
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);
1514  ! qvec(r) at cell interface along z dimension
1515  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1518  else
1519  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1520  end if
1521  curlvec(ixo^s,idir)=(tmp(ixo^s)-tmp(hxo^s))*block%dx(ixo^s,jdir)
1522  !! integral along z dimension
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);
1526  ! qvec(z) at cell interface along r dimension
1527  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1530  else
1531  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1532  end if
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)
1535  end if
1536  else ! idir==z_
1537  if(jdir<kdir) then
1538  ! idir=z,jdir=r,kdir=phi
1539  !! integral along r dimension
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);
1543  ! qvec(r) at cell interface along phi dimension
1544  if(stretched_dim(kdir) .and. stretch_uncentered) then
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))
1547  else
1548  tmp(ixc^s)=0.5d0*(qvec(ixc^s,jdir)+qvec(jxc^s,jdir))
1549  end if
1550  curlvec(ixo^s,idir)=(tmp(hxo^s)-tmp(ixo^s))*block%dx(ixo^s,jdir)
1551  !! integral along phi dimension
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);
1555  ! qvec(phi) at cell interface along r dimension
1556  if(stretched_dim(jdir) .and. stretch_uncentered) then
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))
1559  else
1560  tmp(ixc^s)=0.5d0*(qvec(ixc^s,kdir)+qvec(jxc^s,kdir))
1561  end if
1562  ! r coordinate at cell interface along r dimension
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)
1566  end if
1567  end if
1568  if(idir<idirmin)idirmin=idir
1569  endif
1570  enddo; enddo; enddo;
1571  case default
1572  call mpistop('no such curl evaluator')
1573  end select
1574  case default
1575  call mpistop('not possible to calculate curl')
1576  end select
1577 
1578  end subroutine curlvector_trans
1579 
1580 end module mod_geometry
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
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.
Definition: mod_geometry.t:23
subroutine gradientq(q, x, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q in direction idir at cell interfaces.
Definition: mod_geometry.t:416
subroutine divvectors(qvec, ixIL, ixOL, divq)
Calculate divergence of a vector qvec within ixL using limited extrapolation to cell edges.
Definition: mod_geometry.t:609
subroutine get_surface_area(s, ixGL)
calculate area of surfaces of cells
Definition: mod_geometry.t:159
integer coordinate
Definition: mod_geometry.t:7
subroutine set_pole
Definition: mod_geometry.t:101
integer, parameter spherical
Definition: mod_geometry.t:11
integer, parameter cartesian
Definition: mod_geometry.t:8
integer, parameter stokesbased
Definition: mod_geometry.t:17
integer type_curl
Definition: mod_geometry.t:14
integer, parameter cylindrical
Definition: mod_geometry.t:10
integer, parameter gaussbased
Definition: mod_geometry.t:16
integer, parameter cartesian_expansion
Definition: mod_geometry.t:12
integer, parameter cartesian_stretched
Definition: mod_geometry.t:9
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
Definition: mod_geometry.t:321
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...
Definition: mod_geometry.t:664
integer, parameter central
Definition: mod_geometry.t:15
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...
Definition: mod_geometry.t:458
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.
Definition: mod_geometry.t:516
subroutine putgridgeo(igrid)
Deallocate geometry-related variables.
Definition: mod_geometry.t:149
subroutine gradientx(q, x, ixIL, ixOL, idir, gradq, fourth_order)
Calculate gradient of a scalar q in direction idir at cell interfaces.
Definition: mod_geometry.t:363
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.
Definition: mod_limiter.t:2
integer, parameter limiter_ppm
Definition: mod_limiter.t:27
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...
Definition: mod_limiter.t:129
double precision d
Definition: mod_limiter.t:14
Definition: mod_ppm.t:1
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