MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_amr_fct.t
Go to the documentation of this file.
1 module mod_amr_fct
2  implicit none
3  private
4 
5  type facealloc
6  double precision, dimension(:^D&), pointer :: face
7  end type facealloc
8 
9  type fake_neighbors
10  integer :: igrid
11  integer :: ipe
12  end type fake_neighbors
13 
14  type(facealloc), dimension(:,:,:), allocatable, public :: pface
15 
16  type(fake_neighbors), dimension(:^D&,:,:), allocatable, public :: fine_neighbors
17 
18  integer, dimension(:,:^D&,:), allocatable, public :: old_neighbor
19 
20  double precision, allocatable :: recvbuffer(:), sendbuffer(:)
21  integer :: itag, isend, irecv
22  integer :: nrecv, nsend, ibuf_recv, ibuf_send, ibuf_send_next
23  integer, dimension(^ND) :: isize
24  integer, dimension(:), allocatable :: recvrequest, sendrequest
25  integer, dimension(:,:), allocatable :: recvstatus, sendstatus
26 
27  public :: store_faces
28  public :: comm_faces
29  public :: end_comm_faces
30  public :: deallocatebfaces
31  public :: old_neighbors
32  public :: prolong_2nd_stg
33  public :: already_fine
34 
35 contains
36  !> This subroutine performs a 2nd order prolongation for a staggered field F,
37  !> preserving the divergence of the coarse cell.
38  !> This is useful for preserving DivF=0.
39  !> If DivF=f(x), a different algorithm must be used.
40  subroutine prolong_2nd_stg(sCo,sFi,ixCo^Lin,ixFi^Lin,dxCo^D,xComin^D,dxFi^D,xFimin^D,ghost,fine_^Lin)
42  use mod_physics
43 
44  logical, intent(in) :: ghost
45  integer, intent(in) :: ixco^lin, ixfi^lin
46  double precision, intent(in) :: dxco^d, xcomin^d, dxfi^d, xfimin^d
47  type(state), intent(in) :: sco
48  type(state), intent(inout) :: sfi
49  logical, optional :: fine_^lin
50 
51  double precision :: eta^d, invdxco^d
52  double precision :: bfluxco(sco%ixgs^s,nws),bfluxfi(sfi%ixgs^s,nws)
53  double precision :: slopes(sco%ixgs^s,ndim),b_energy_change(ixg^t)
54  {^ifthreed
55  ! Directional bias, see pdf
56  ! These arrays should have ranges ixO^S
57  ! Since they are still not known, we give ranges ixG^T,
58  ! even though it is wasteful
59  double precision :: sigmau(ixg^t),sigmad(ixg^t),sigma(ixg^t,1:ndim), alpha(ixg^t,1:ndim)
60  ! Auxiliary arrays for magnetic fluxes
61  double precision :: f1(ixg^t),f2(ixg^t),f3(ixg^t),f4(ixg^t)
62  }
63  integer :: ixco^l,ixfi^l
64  integer :: idim1,idim2,ix^de,idim3,ixfis^l,ixgs^l,ixcos^l,ixfisc^l
65  integer :: ixcosv^l(1:ndim),ixfisv^l(1:ndim)
66  integer :: hxcos^l,jxcos^l,ixcose^l,ixfise^l,hxfisc^l,jxfisc^l,ipxfisc^l,ixcosc^l,imxfisc^l,jpxfisc^l,jmxfisc^l,hpxfisc^l
67  integer :: hxfi^l,jxfi^l,hijxfi^l,hjixfi^l,hjjxfi^l
68  integer :: iihxfi^l,iijxfi^l,ijhxfi^l,ijjxfi^l,ihixfi^l,ijixfi^l,ihjxfi^l
69  integer :: jihxfi^l,jijxfi^l,jjhxfi^l,jjjxfi^l,jhixfi^l,jjixfi^l,jhjxfi^l
70  logical :: fine_^l
71 
72  {^ifoned
73  call mpistop("CT prolongation not implemented in 1D. But CT is not needed.")
74  }
75 
76  {^nooned
77  ! Note on the indices:
78  ! ixCo Cells where
79  ! divergence-preserving prolongation will be applied.
80  ! ixFi Fine cells which correspond to that extent.
81  ! ixCoE For 'expanded', faces that need to be used to calculate
82  ! the slopes for the prolongation
83  ! ixCos
84  ! ixFis
85  ! ixCosV And
86  ! ixFisV For 'variable', since the ranges depend on the direction of the faces
87  ! ixCosC For 'copy',
88  ! ixFisC For 'copy', to fill the information from the fine grid needed in the internal faces prolongation step.
89  ! At the end, ixFisC is used to copy the assign the magnetic
90  ! fields components to the state structure.
91 
92  fine_^l=.false.;
93 
94  {if(present(fine_min^din)) fine_min^d=fine_min^din;}
95  {if(present(fine_max^din)) fine_max^d=fine_max^din;}
96 
97  ! When filling ghost cells, ixFi^L are given.
98  ! When refining the block, ixCo^L are given.
99  if (ghost) then
100  ixfi^l=ixfi^lin;
101  {ixcomin^d=int((ixfimin^d+nghostcells+1)/2);}
102  {ixcomax^d=int((ixfimax^d+nghostcells+1)/2);}
103  else
104  ixco^l=ixco^lin;
105  ixfi^l=ixm^ll;
106  end if
107 
108  ! Expanded range for staggered variables
109  ixcosmin^d=ixcomin^d-1;
110  ixcosmax^d=ixcomax^d;
111 
112  ixfismin^d=ixfimin^d-1;
113  ixfismax^d=ixfimax^d;
114 
115 
116  associate(wcos=>sco%ws, wfis=>sfi%ws,wco=>sco%w, wfi=>sfi%w)
117  ! Assemble general indices
118  ixgsmin^d=sfi%ixGsmin^d;
119  ixgsmax^d=sfi%ixGsmax^d;
120 
121  do idim1=1,ndim
122  ixcosvmin^d(idim1)=ixcomin^d-kr(^d,idim1);
123  ixcosvmax^d(idim1)=ixcomax^d;
124  ixfisvmin^d(idim1)=ixfimin^d-kr(^d,idim1);
125  ixfisvmax^d(idim1)=ixfimax^d;
126  end do
127 
128  ! Initialize auxiliary arrays at zero
129  bfluxco = zero
130  bfluxfi = zero
131  slopes = zero
132 
133 
134  invdxco^d=1.d0/dxco^d;
135 
136  ! Fill coarse magnetic flux array
137  do idim1=1,ndim
138  ! Fill information from parts already at the fine level
139  ! First set up indices
140  if (ghost) then
141  ixfisemin^d=max(1-kr(^d,idim1),ixfisvmin^d(idim1)-2*(1-kr(^d,idim1)));
142  ixfisemax^d=min(ixgsmax^d,ixfisvmax^d(idim1)+2*(1-kr(^d,idim1)));
143  ixcose^l=int((ixfise^l+nghostcells+1)/2);
144  else
145  ixcose^l=ixcosv^l(idim1)^ladd(1-kr(idim1,^d));
146  ixfise^l=ixfisv^l(idim1)^ladd2*(1-kr(idim1,^d));
147  end if
148  ! Convert fine fields to fluxes
149  bfluxfi(ixfise^s,idim1)=wfis(ixfise^s,idim1)*sfi%surfaceC(ixfise^s,idim1)
150  {^iftwod
151  idim2=1+mod(idim1,2)
152  }
153  {^ifthreed
154  idim2=1+mod(idim1,3)
155  idim3=1+mod(idim1+1,3)
156  }
157  bfluxco(ixcose^s,idim1) = zero
158  ! Add fine fluxes sharing the same fine face
159  {do ix^de=0,1\}
160  ixfisc^l=ixfise^l+ix2*kr(idim2,^d);
161  {^ifthreed
162  ixfisc^l=ixfisc^l+ix3*kr(idim3,^d);
163  }
164  bfluxco(ixcose^s,idim1)=bfluxco(ixcose^s,idim1)+bfluxfi(ixfiscmin^d:ixfiscmax^d:2,idim1)
165  {end do^DE&\}
166  end do
167 
168  ! Omit indices for already refined face, if any
169  {
170  if (fine_min^d) then
171  ixcosvmin^d(^d)=ixcosvmin^d(^d)+1
172  ixfisvmin^d(^d)=ixfisvmin^d(^d)+2
173  end if
174  if (fine_max^d) then
175  ixcosvmax^d(^d)=ixcosvmax^d(^d)-1
176  ixfisvmax^d(^d)=ixfisvmax^d(^d)-2
177  end if
178  \}
179 
180  do idim1=1,ndim
181  ixcose^l=ixcosv^l(idim1);
182  ! Omit part already refined
183  {
184  if (^d/=idim1) then
185  if ((.not.fine_min^d).or.(.not.ghost)) then
186  ixcosemin^d=ixcosvmin^d(idim1)-1
187  end if
188  if ((.not.fine_max^d).or.(.not.ghost)) then
189  ixcosemax^d=ixcosvmax^d(idim1)+1
190  end if
191  end if
192  \}
193  ! Fill coarse flux array from coarse field
194  bfluxco(ixcose^s,idim1)=wcos(ixcose^s,idim1)*sco%surfaceC(ixcose^s,idim1)
195  end do
196  ! Finished filling coarse flux array
197 
198  ! Distribute coarse fluxes among fine fluxes
199  ! There are too many loops here, perhaps optimize later
200  do idim1=1,ndim
201  ixcos^l=ixcosv^l(idim1);
202  do idim2=1,ndim
203  if(idim1==idim2) cycle
204  ! Calculate slope in direction idim2
205  ! Set up indices
206  jxcos^l=ixcos^l+kr(idim2,^d);
207  hxcos^l=ixcos^l-kr(idim2,^d);
208  slopes(ixcos^s,idim2)=0.125d0*(bfluxco(jxcos^s,idim1)-bfluxco(hxcos^s,idim1))
209  {^iftwod
210  do ix2=0,1
211  ixfiscmin^d=ixfisvmin^d(idim1)+ix2*kr(^d,idim2);
212  ixfiscmax^d=ixfisvmax^d(idim1)+ix2*kr(^d,idim2);
213  bfluxfi(ixfiscmin^d:ixfiscmax^d:2,idim1)=half*(bfluxco(ixcos^s,idim1)&
214  +(2*ix2-1)*slopes(ixcos^s,idim2))
215  end do
216  }
217  {^ifthreed
218  do idim3=1,ndim
219  ! Calculate slope in direction idim3
220  ! Set up indices
221  jxcos^l=ixcos^l+kr(idim3,^d);
222  hxcos^l=ixcos^l-kr(idim3,^d);
223  slopes(ixcos^s,idim3)=0.125d0*(bfluxco(jxcos^s,idim1)-bfluxco(hxcos^s,idim1))
224  if(lvc(idim1,idim2,idim3)<1) cycle
225  do ix2=0,1
226  do ix3=0,1
227  {ixfiscmin^d=ixfisvmin^d(idim1)+ix2*kr(^d,idim2)+ix3*kr(^d,idim3);}
228  {ixfiscmax^d=ixfisvmax^d(idim1)+ix2*kr(^d,idim2)+ix3*kr(^d,idim3);}
229  bfluxfi(ixfiscmin^d:ixfiscmax^d:2,idim1)=quarter*bfluxco(ixcos^s,idim1)&
230  +quarter*(2*ix2-1)*slopes(ixcos^s,idim2)&
231  +quarter*(2*ix3-1)*slopes(ixcos^s,idim3)
232  end do
233  end do
234  end do
235  }
236  end do
237  end do
238 
239  ! Calculate interior fine fluxes
240  {^ifthreed
241  ! Directional bias for nonlinear prolongation
242  do idim1=1,ndim
243  do idim2=1,ndim
244  do idim3=1,ndim
245  if (lvc(idim1,idim2,idim3)<1) cycle
246  ! Set up indices
247  hxfi^l=ixfi^l-kr(idim1,^d);
248  jxfi^l=ixfi^l+kr(idim1,^d);
249 
250  hijxfi^l=hxfi^l+kr(idim3,^d);
251  hjixfi^l=hxfi^l+kr(idim2,^d);
252  hjjxfi^l=hijxfi^l+kr(idim2,^d);
253 
254  iihxfi^l=ixfi^l-kr(idim3,^d);
255  iijxfi^l=ixfi^l+kr(idim3,^d);
256  ihixfi^l=ixfi^l-kr(idim2,^d);
257  ihjxfi^l=ihixfi^l+kr(idim3,^d);
258  ijixfi^l=ixfi^l+kr(idim2,^d);
259  ijhxfi^l=ijixfi^l-kr(idim3,^d);
260  ijjxfi^l=ijixfi^l+kr(idim3,^d);
261 
262  jihxfi^l=jxfi^l-kr(idim3,^d);
263  jijxfi^l=jxfi^l+kr(idim3,^d);
264  jhixfi^l=jxfi^l-kr(idim2,^d);
265  jhjxfi^l=jhixfi^l+kr(idim3,^d);
266  jjixfi^l=jxfi^l+kr(idim2,^d);
267  jjhxfi^l=jjixfi^l-kr(idim3,^d);
268  jjjxfi^l=jjixfi^l+kr(idim3,^d);
269 
270  sigmau(ixco^s)=&
271  abs(bfluxfi(jhixfimin^d:jhixfimax^d:2,idim2))&
272  + abs(bfluxfi(jhjxfimin^d:jhjxfimax^d:2,idim2))&
273  + abs(bfluxfi(jjixfimin^d:jjixfimax^d:2,idim2))&
274  + abs(bfluxfi(jjjxfimin^d:jjjxfimax^d:2,idim2))&
275  + abs(bfluxfi(jihxfimin^d:jihxfimax^d:2,idim3))&
276  + abs(bfluxfi(jijxfimin^d:jijxfimax^d:2,idim3))&
277  + abs(bfluxfi(jjhxfimin^d:jjhxfimax^d:2,idim3))&
278  + abs(bfluxfi(jjjxfimin^d:jjjxfimax^d:2,idim3))
279 
280  sigmad(ixco^s)=&
281  abs(bfluxfi(ihixfimin^d:ihixfimax^d:2,idim2))&
282  + abs(bfluxfi(ihjxfimin^d:ihjxfimax^d:2,idim2))&
283  + abs(bfluxfi(ijixfimin^d:ijixfimax^d:2,idim2))&
284  + abs(bfluxfi(ijjxfimin^d:ijjxfimax^d:2,idim2))&
285  + abs(bfluxfi(iihxfimin^d:iihxfimax^d:2,idim3))&
286  + abs(bfluxfi(iijxfimin^d:iijxfimax^d:2,idim3))&
287  + abs(bfluxfi(ijhxfimin^d:ijhxfimax^d:2,idim3))&
288  + abs(bfluxfi(ijjxfimin^d:ijjxfimax^d:2,idim3))
289 
290  sigma(ixco^s,idim1)=sigmau(ixco^s)+sigmad(ixco^s)
291  where(sigma(ixco^s,idim1)/=zero)
292  sigma(ixco^s,idim1)=abs(sigmau(ixco^s)-sigmad(ixco^s))/sigma(ixco^s,idim1)
293  elsewhere
294  sigma(ixco^s,idim1)=zero
295  end where
296 
297  end do
298  end do
299  end do
300 
301  ! Directional bias for Toth-Roe prolongation
302  do idim1=1,ndim
303  do idim2=1,ndim
304  do idim3=1,ndim
305  if (lvc(idim1,idim2,idim3)<1) cycle
306  ! Nonlinear
307  ! alpha(ixCo^S,idim1)=sigma(ixCo^S,idim2)-sigma(ixCo^S,idim3)
308  ! Homogeneous
309  ! alpha(ixCo^S,idim1)=0.d0
310  ! Toth-Roe
311  alpha(ixco^s,idim1)=(sco%dx(ixco^s,idim2)-sco%dx(ixco^s,idim3))/(sco%dx(ixco^s,idim2)+sco%dx(ixco^s,idim3))
312  end do
313  end do
314  end do
315  }
316 
317  do idim1=1,ndim
318  do idim2=1,ndim
319  {^iftwod
320  if (idim1==idim2) cycle
321  ixfiscmin^d=ixfismin^d+1;
322  ixfiscmax^d=ixfismax^d-1;
323  jxfisc^l=ixfisc^l+kr(idim1,^d);
324  hxfisc^l=ixfisc^l-kr(idim1,^d);
325  ipxfisc^l=ixfisc^l+kr(idim2,^d);
326  imxfisc^l=ixfisc^l-kr(idim2,^d);
327  jpxfisc^l=jxfisc^l+kr(idim2,^d);
328  jmxfisc^l=jxfisc^l-kr(idim2,^d);
329  hpxfisc^l=hxfisc^l+kr(idim2,^d);
330 
331  bfluxfi(ixfiscmin^d:ixfiscmax^d:2,idim1)=&
332  half*(bfluxfi(jxfiscmin^d:jxfiscmax^d:2,idim1)&
333  +bfluxfi(hxfiscmin^d:hxfiscmax^d:2,idim1))&
334  -quarter*(bfluxfi(ipxfiscmin^d:ipxfiscmax^d:2,idim2)&
335  -bfluxfi(jpxfiscmin^d:jpxfiscmax^d:2,idim2)&
336  -bfluxfi(imxfiscmin^d:imxfiscmax^d:2,idim2)&
337  +bfluxfi(jmxfiscmin^d:jmxfiscmax^d:2,idim2))
338 
339  bfluxfi(ipxfiscmin^d:ipxfiscmax^d:2,idim1)=&
340  half*(bfluxfi(jpxfiscmin^d:jpxfiscmax^d:2,idim1)&
341  +bfluxfi(hpxfiscmin^d:hpxfiscmax^d:2,idim1))&
342  -quarter*(bfluxfi(ipxfiscmin^d:ipxfiscmax^d:2,idim2)&
343  -bfluxfi(jpxfiscmin^d:jpxfiscmax^d:2,idim2)&
344  -bfluxfi(imxfiscmin^d:imxfiscmax^d:2,idim2)&
345  +bfluxfi(jmxfiscmin^d:jmxfiscmax^d:2,idim2))
346  }
347  {^ifthreed
348  do idim3=1,ndim
349  if (lvc(idim1,idim2,idim3)<1) cycle
350  ! Set up indices
351  hxfi^l=ixfi^l-kr(idim1,^d);
352  jxfi^l=ixfi^l+kr(idim1,^d);
353 
354  hijxfi^l=hxfi^l+kr(idim3,^d);
355  hjixfi^l=hxfi^l+kr(idim2,^d);
356  hjjxfi^l=hijxfi^l+kr(idim2,^d);
357 
358  iihxfi^l=ixfi^l-kr(idim3,^d);
359  iijxfi^l=ixfi^l+kr(idim3,^d);
360  ihixfi^l=ixfi^l-kr(idim2,^d);
361  ihjxfi^l=ihixfi^l+kr(idim3,^d);
362  ijixfi^l=ixfi^l+kr(idim2,^d);
363  ijhxfi^l=ijixfi^l-kr(idim3,^d);
364  ijjxfi^l=ijixfi^l+kr(idim3,^d);
365 
366  jihxfi^l=jxfi^l-kr(idim3,^d);
367  jijxfi^l=jxfi^l+kr(idim3,^d);
368  jhixfi^l=jxfi^l-kr(idim2,^d);
369  jhjxfi^l=jhixfi^l+kr(idim3,^d);
370  jjixfi^l=jxfi^l+kr(idim2,^d);
371  jjhxfi^l=jjixfi^l-kr(idim3,^d);
372  jjjxfi^l=jjixfi^l+kr(idim3,^d);
373 
374  ! Prolongation formulas
375  f1(ixco^s)=bfluxfi(ihixfimin^d:ihixfimax^d:2,idim2)&
376  -bfluxfi(jhixfimin^d:jhixfimax^d:2,idim2)&
377  -bfluxfi(ijixfimin^d:ijixfimax^d:2,idim2)&
378  +bfluxfi(jjixfimin^d:jjixfimax^d:2,idim2)
379 
380  f2(ixco^s)=bfluxfi(ihjxfimin^d:ihjxfimax^d:2,idim2)&
381  -bfluxfi(jhjxfimin^d:jhjxfimax^d:2,idim2)&
382  -bfluxfi(ijjxfimin^d:ijjxfimax^d:2,idim2)&
383  +bfluxfi(jjjxfimin^d:jjjxfimax^d:2,idim2)
384 
385  f3(ixco^s)=bfluxfi(iihxfimin^d:iihxfimax^d:2,idim3)&
386  -bfluxfi(jihxfimin^d:jihxfimax^d:2,idim3)&
387  -bfluxfi(iijxfimin^d:iijxfimax^d:2,idim3)&
388  +bfluxfi(jijxfimin^d:jijxfimax^d:2,idim3)
389 
390  f4(ixco^s)=bfluxfi(ijhxfimin^d:ijhxfimax^d:2,idim3)&
391  -bfluxfi(jjhxfimin^d:jjhxfimax^d:2,idim3)&
392  -bfluxfi(ijjxfimin^d:ijjxfimax^d:2,idim3)&
393  +bfluxfi(jjjxfimin^d:jjjxfimax^d:2,idim3)
394 
395  bfluxfi(ixfimin^d:ixfimax^d:2,idim1)=&
396  half*(bfluxfi(hxfimin^d:hxfimax^d:2,idim1)+bfluxfi(jxfimin^d:jxfimax^d:2,idim1))&
397  +6.25d-2*((3.d0+alpha(ixco^s,idim2))*f1(ixco^s)&
398  +(1.d0-alpha(ixco^s,idim2))*f2(ixco^s)&
399  +(3.d0-alpha(ixco^s,idim3))*f3(ixco^s)&
400  +(1.d0+alpha(ixco^s,idim3))*f4(ixco^s))
401 
402  bfluxfi(ijixfimin^d:ijixfimax^d:2,idim1)=&
403  half*(bfluxfi(hjixfimin^d:hjixfimax^d:2,idim1)+bfluxfi(jjixfimin^d:jjixfimax^d:2,idim1))&
404  +6.25d-2*((3.d0+alpha(ixco^s,idim2))*f1(ixco^s)&
405  +(1.d0-alpha(ixco^s,idim2))*f2(ixco^s)&
406  +(1.d0+alpha(ixco^s,idim3))*f3(ixco^s)&
407  +(3.d0-alpha(ixco^s,idim3))*f4(ixco^s))
408 
409  bfluxfi(iijxfimin^d:iijxfimax^d:2,idim1)=&
410  half*(bfluxfi(hijxfimin^d:hijxfimax^d:2,idim1)+bfluxfi(jijxfimin^d:jijxfimax^d:2,idim1))&
411  +6.25d-2*((1.d0-alpha(ixco^s,idim2))*f1(ixco^s)&
412  +(3.d0+alpha(ixco^s,idim2))*f2(ixco^s)&
413  +(3.d0-alpha(ixco^s,idim3))*f3(ixco^s)&
414  +(1.d0+alpha(ixco^s,idim3))*f4(ixco^s))
415 
416  bfluxfi(ijjxfimin^d:ijjxfimax^d:2,idim1)=&
417  half*(bfluxfi(hjjxfimin^d:hjjxfimax^d:2,idim1)+bfluxfi(jjjxfimin^d:jjjxfimax^d:2,idim1))&
418  +6.25d-2*((1.d0-alpha(ixco^s,idim2))*f1(ixco^s)&
419  +(3.d0+alpha(ixco^s,idim2))*f2(ixco^s)&
420  +(1.d0+alpha(ixco^s,idim3))*f3(ixco^s)&
421  +(3.d0-alpha(ixco^s,idim3))*f4(ixco^s))
422  end do
423  }
424  end do
425  end do
426 
427  ! Go back to magnetic fields
428  do idim1=1,ndim
429  ixfiscmax^d=ixfimax^d;
430  ixfiscmin^d=ixfimin^d-kr(^d,idim1);
431  where(sfi%surfaceC(ixfisc^s,idim1)/=zero)
432  wfis(ixfisc^s,idim1)=bfluxfi(ixfisc^s,idim1)/sfi%surfaceC(ixfisc^s,idim1)
433  elsewhere
434  wfis(ixfisc^s,idim1)=zero
435  end where
436  end do
437 
438  if(phys_total_energy.and. .not.prolongprimitive) then
439  b_energy_change(ixfi^s)=0.5d0*sum(wfi(ixfi^s,iw_mag(:))**2,dim=ndim+1)
440  end if
441  call phys_face_to_center(ixfi^l,sfi)
442  if(phys_total_energy.and. .not.prolongprimitive) then
443  b_energy_change(ixfi^s)=0.5d0*sum(wfi(ixfi^s,iw_mag(:))**2,dim=ndim+1)-&
444  b_energy_change(ixfi^s)
445  wfi(ixfi^s,iw_e)=wfi(ixfi^s,iw_e)+b_energy_change(ixfi^s)
446  end if
447 
448  end associate
449 
450  } ! END NOONED
451  end subroutine prolong_2nd_stg
452 
453  !> To achive consistency and thus conservation of divergence,
454  !> when refining a block we take into account the faces of the
455  !> already fine neighbours, if any. This routine stores them.
456  subroutine store_faces
457  use mod_forest, only: refine
459  integer :: igrid, iigrid, idims, iside, ineighbor, ipe_neighbor
460  integer :: nx^d, i^d, ic^d, inc^d
461 
462  if (npe>1) then
463  nsend_fc=0
464  nrecv_fc=0
465  end if
466 
467  ! Size of the block face
468  nx^d=ixmhi^d-ixmlo^d+1;
469 
470  do iigrid=1,igridstail; igrid=igrids(iigrid);
471  ! Check whether it is necessary to store any block face, i.e.
472  ! if any coarser neighbour is going to be refined
473  {do iside=1,2
474  i^dd=kr(^dd,^d)*(2*iside-3);
475  if (neighbor_pole(i^dd,igrid)/=0) cycle
476  if (neighbor_type(i^dd,igrid)==neighbor_coarse) then
477  ineighbor =neighbor(1,i^dd,igrid)
478  ipe_neighbor=neighbor(2,i^dd,igrid)
479  if (refine(ineighbor,ipe_neighbor)) then
480  allocate(pface(iside,^d,igrid)%face(1^d%1:nx^dd))
481  !! Store the faces
482  if (iside==1) then !! left side
483  pface(iside,^d,igrid)%face(1^d%1:nx^dd)=&
484  ps(igrid)%ws(ixmlo^d-1^d%ixM^t,^d)
485  else !! right side
486  pface(iside,^d,igrid)%face(1^d%1:nx^dd)=&
487  ps(igrid)%ws(ixmhi^d^d%ixM^t,^d)
488  end if
489  if (ipe_neighbor/=mype) nsend_fc(^d)=nsend_fc(^d)+1
490  end if
491  end if
492  end do\}
493 
494  ! If a grid is going to be refined,
495  ! remember what are its neighbours.
496  if (refine(igrid,mype)) then
497  ! Initialize neighbour array
498  fine_neighbors(:^d&,:,igrid)%igrid=-1
499  fine_neighbors(:^d&,:,igrid)%ipe=-1
500  do idims=1,ndim
501  do iside=1,2
502  i^d=kr(^d,idims)*(2*iside-3);
503  if (neighbor_pole(i^d,igrid)/=0) cycle
504  if (neighbor_type(i^d,igrid)==neighbor_fine) then
505  {do ic^db=1+int((1+i^db)/2),2-int((1-i^db)/2)
506  inc^db=ic^db+i^db\}
507  ineighbor=neighbor_child(1,inc^d,igrid)
508  ipe_neighbor=neighbor_child(2,inc^d,igrid)
509 
510  fine_neighbors(ic^d,idims,igrid)%igrid= ineighbor
511  fine_neighbors(ic^d,idims,igrid)%ipe=ipe_neighbor
512 
513  if (ipe_neighbor/=mype) nrecv_fc(idims)=nrecv_fc(idims)+1
514  {end do\}
515  end if
516  end do
517  end do
518  end if
519 
520  end do
521 
522  end subroutine store_faces
523 
524  !> When refining a coarse block with fine neighbours, it is necessary
525  !> prolong consistently with the already fine faces.
526  !> This routine takes care of the communication of such faces.
527  subroutine comm_faces
528  use mod_forest, only: refine
530 
531  integer :: iigrid,igrid,ineighbor,ipe_neighbor
532  integer :: idims,iside,i^d,ic^d,inc^d,nx^d
533  integer :: recvsize, sendsize
534 
535  ! Communicate the block faces to achieve consistency when refining
536  ! Initialize communication structures
537 
538  nrecv=0
539  nsend=0
540  recvsize=0
541  sendsize=0
542 
543  do idims=1,ndim
544  select case (idims)
545  {case (^d)
546  nrecv=nrecv+nrecv_fc(^d)
547  nsend=nsend+nsend_fc(^d)
548  nx^d=1^d%nx^dd=ixmhi^dd-ixmlo^dd+1;
549  isize(^d)={nx^dd*}
550  recvsize=recvsize+nrecv_fc(^d)*isize(^d)
551  sendsize=sendsize+nsend_fc(^d)*isize(^d)
552  \}
553  end select
554  end do
555 
556  if (nrecv>0) then
557  ! Allocate receive buffer
558  allocate(recvbuffer(recvsize),recvstatus(mpi_status_size,nrecv), &
559  recvrequest(nrecv))
560  recvrequest=mpi_request_null
561  ibuf_recv=1
562  irecv=0
563 
564  ! Receive
565  do iigrid=1,igridstail; igrid=igrids(iigrid);
566  if (refine(igrid,mype)) then
567  {do ic^db=1,2\}
568  ! Only one of the sides will be necessary,
569  ! so we do the loop only over dimensions, instead of
570  ! over dimensions and sizes as in the routines
571  ! old_neighbors and already_fine.
572  do idims=1,ndim
573  ipe_neighbor=fine_neighbors(ic^d,idims,igrid)%ipe
574  ineighbor =fine_neighbors(ic^d,idims,igrid)%igrid
575  if (ineighbor>0.and.ipe_neighbor/=mype) then
576  {if (idims==^d) iside=ic^d\}
577  !!! Check indices
578  i^d=kr(^d,idims)*(2*iside-3);
579  if (neighbor_pole(i^d,igrid)/=0) cycle
580  inc^d=ic^d+i^d;
581  irecv=irecv+1
582  itag=4**^nd*(igrid-1)+{inc^d*4**(^d-1)+}
583 
584  call mpi_irecv(recvbuffer(ibuf_recv),isize(idims), &
585  mpi_double_precision,ipe_neighbor,itag, &
586  icomm,recvrequest(irecv),ierrmpi)
587  ibuf_recv=ibuf_recv+isize(idims)
588  end if
589  end do
590  {end do\}
591  end if
592  end do
593 
594  end if
595 
596  if (nsend>0) then
597  ! Allocate send buffer
598  allocate(sendbuffer(sendsize),sendstatus(mpi_status_size,nsend),sendrequest(nsend))
599  sendrequest=mpi_request_null
600  isend=0
601  ibuf_send=1
602  do iigrid=1,igridstail; igrid=igrids(iigrid);
603  ! Check whether it is necessary to store any block face, i.e.
604  ! if any coarser neighbour is going to be refined
605  {do iside=1,2
606  i^dd=kr(^dd,^d)*(2*iside-3);
607  ! When there is a pole, faces are always zero and this is not necessary
608  if (neighbor_pole(i^dd,igrid)/=0) cycle
609  if (neighbor_type(i^dd,igrid)==neighbor_coarse) then
610  ineighbor =neighbor(1,i^dd,igrid)
611  ipe_neighbor=neighbor(2,i^dd,igrid)
612  if (refine(ineighbor,ipe_neighbor)) then
613  if (ipe_neighbor/=mype) then
614  ic^dd=1+modulo(node(pig^dd_,igrid)-1,2);
615  inc^d=-2*i^d+ic^d^d%inc^dd=ic^dd;
616  itag=4**^nd*(ineighbor-1)+{inc^dd*4**(^dd-1)+}
617  isend=isend+1
618  ibuf_send_next=ibuf_send+isize(^d)
619  sendbuffer(ibuf_send:ibuf_send_next-1)=&
620  reshape(pface(iside,^d,igrid)%face,(/isize(^d)/))
621  call mpi_isend(sendbuffer(ibuf_send),isize(^d), &
622  mpi_double_precision,ipe_neighbor,itag, &
623  icomm,sendrequest(isend),ierrmpi)
624  ibuf_send=ibuf_send_next
625  end if
626  end if
627  end if
628  end do\}
629  end do
630  end if
631 
632  ! Waitalls
633  if (nrecv>0) then
634  call mpi_waitall(nrecv,recvrequest,recvstatus,ierrmpi)
635  deallocate(recvstatus,recvrequest)
636  ibuf_recv=1
637  end if
638 
639  if (nsend>0) then
640  call mpi_waitall(nsend,sendrequest,sendstatus,ierrmpi)
641  deallocate(sendbuffer,sendstatus,sendrequest)
642  end if
643 
644  end subroutine comm_faces
645 
646  subroutine end_comm_faces
648  ! Deallocate receive buffer
649  if (nrecv>0) deallocate(recvbuffer)
650  end subroutine end_comm_faces
651 
652  subroutine deallocatebfaces
654  integer :: igrid, iigrid, iside
655 
656  do iigrid=1,igridstail; igrid=igrids(iigrid);
657  {do iside=1,2
658  if (associated(pface(iside,^d,igrid)%face)) then
659  deallocate(pface(iside,^d,igrid)%face)
660  end if
661  end do\}
662  end do
663 
664  end subroutine deallocatebfaces
665 
666  subroutine old_neighbors(child_igrid,child_ipe,igrid,ipe)
668  integer, dimension(2^D&), intent(in) :: child_igrid, child_ipe
669  integer, intent(in) :: igrid, ipe
670  integer :: iside, i^d, ic^d
671 
672  {do ic^db=1,2\}
673  old_neighbor(:,:^d&,child_igrid(ic^d))=-1
674  {do iside=1,2
675  if (ic^d==iside) then
676  i^dd=kr(^dd,^d)*(2*iside-3);
677  old_neighbor(1,i^dd,child_igrid(ic^dd))=&
678  fine_neighbors(ic^dd,^d,igrid)%igrid
679  old_neighbor(2,i^dd,child_igrid(ic^dd))=&
680  fine_neighbors(ic^dd,^d,igrid)%ipe
681  end if
682  end do\}
683  {end do\}
684 
685  end subroutine old_neighbors
686 
687  !> This routine fills the fine faces before prolonging.
688  !> It is the face equivalent of fix_conserve
689  subroutine already_fine(sFi,ichild,fine_^L)
690  use mod_forest
692  type(tree_node_ptr) :: tree
693  type(state) :: sfi
694  integer, intent(in) :: ichild
695  logical :: fine_^l
696 
697  integer :: ineighbor,ipe_neighbor,ibufnext
698  integer :: iside,iotherside,i^d,nx^d
699 
700  ! Size of the block face
701  nx^d=ixmhi^d-ixmlo^d+1;
702 
703  ! Initialise everything to zero and false
704  fine_min^d=.false.;
705  fine_max^d=.false.;
706  sfi%ws=zero
707 
708  {^nooned
709  ! This face communication is not needed in 1D
710  {do iside=1,2
711  i^dd=kr(^dd,^d)*(2*iside-3);
712  ! This is not necessary at the pole.
713  ! We are inside a loop over the children, so the grid index is ichild
714  if (neighbor_pole(i^dd,ichild)/=0) cycle
715 
716  ! Get old ipe and igrid of neighbour from array fake_neighbor
717  ! Then plug it into the structure pfaces and get the faces
718  ineighbor =old_neighbor(1,i^dd,ichild)
719  ipe_neighbor=old_neighbor(2,i^dd,ichild)
720 
721  iotherside=3-iside
722  if (ineighbor>0) then
723  if (iside==1) then ! left side
724  if (ipe_neighbor==mype) then
725  sfi%ws(ixmlo^d-1^d%ixM^t,^d)=pface(iotherside,^d,ineighbor)%face(1^d%1:nx^dd)
726  else
727  ibufnext=ibuf_recv+isize(^d)
728  sfi%ws(ixmlo^d-1^d%ixM^t,^d)=reshape(&
729  source=recvbuffer(ibuf_recv:ibufnext-1),&
730  shape=shape(sfi%ws(ixmlo^d-1^d%ixM^t,^d)))
731  ibuf_recv=ibufnext
732  end if
733  fine_min^d=.true.
734  else ! right side
735  if (ipe_neighbor==mype) then
736  sfi%ws(ixmhi^d^d%ixM^t,^d)=pface(iotherside,^d,ineighbor)%face(1^d%1:nx^dd)
737  else
738  ibufnext=ibuf_recv+isize(^d)
739  sfi%ws(ixmhi^d^d%ixM^t,^d)=reshape(&
740  source=recvbuffer(ibuf_recv:ibufnext-1),&
741  shape=shape(sfi%ws(ixmhi^d^d%ixM^t,^d)))
742  ibuf_recv=ibufnext
743  end if
744  fine_max^d=.true.
745  end if
746  end if
747  end do\}
748  }
749  end subroutine already_fine
750 
751 end module mod_amr_fct
type(fake_neighbors), dimension(:^d &,:,:), allocatable, public fine_neighbors
Definition: mod_amr_fct.t:16
subroutine, public old_neighbors(child_igrid, child_ipe, igrid, ipe)
Definition: mod_amr_fct.t:667
integer, dimension(:,:^d &,:), allocatable, public old_neighbor
Definition: mod_amr_fct.t:18
subroutine, public store_faces
To achive consistency and thus conservation of divergence, when refining a block we take into account...
Definition: mod_amr_fct.t:457
subroutine, public comm_faces
When refining a coarse block with fine neighbours, it is necessary prolong consistently with the alre...
Definition: mod_amr_fct.t:528
subroutine, public already_fine(sFi, ichild, fine_L)
This routine fills the fine faces before prolonging. It is the face equivalent of fix_conserve.
Definition: mod_amr_fct.t:690
type(facealloc), dimension(:,:,:), allocatable, public pface
Definition: mod_amr_fct.t:14
subroutine, public end_comm_faces
Definition: mod_amr_fct.t:647
subroutine, public prolong_2nd_stg(sCo, sFi, ixCoLin, ixFiLin, dxCoD, xCominD, dxFiD, xFiminD, ghost, fine_Lin)
This subroutine performs a 2nd order prolongation for a staggered field F, preserving the divergence ...
Definition: mod_amr_fct.t:41
subroutine, public deallocatebfaces
Definition: mod_amr_fct.t:653
Module with basic grid data structures.
Definition: mod_forest.t:2
logical, dimension(:,:), allocatable, save refine
Definition: mod_forest.t:70
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
integer nghostcells
Number of ghost cells surrounding a grid.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
Pointer to a tree_node.
Definition: mod_forest.t:6