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