MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_amr_solution_node.t
Go to the documentation of this file.
2 use mod_comm_lib, only: mpistop
3
4 implicit none
5 private
6
7 public :: getnode, putnode
8 public :: alloc_node, alloc_state
9 public :: dealloc_node
10
11contains
12
13
14 !> Get first available igrid on processor ipe
15 integer function getnode(ipe)
16 use mod_forest, only: igrid_inuse
18
19 integer, intent(in) :: ipe
20 integer :: igrid, igrid_available
21
22 igrid_available=0
23
24 do igrid=1,max_blocks
25 if (igrid_inuse(igrid,ipe)) cycle
26
27 igrid_available=igrid
28 exit
29 end do
30
31 if (igrid_available == 0) then
32 getnode = -1
33 print *, "Current maximum number of grid blocks:", max_blocks
34 call mpistop("Insufficient grid blocks; increase max_blocks in meshlist")
35 else
36 getnode=igrid_available
37 igrid_inuse(igrid,ipe)=.true.
38 end if
39
40 if (ipe==mype) then
41 ! initialize nodal block
42 node(1:nodehi,getnode) = 0
43 rnode(1:rnodehi,getnode) = zero
44 end if
45
46 end function getnode
47
48 ! put igrid on processor ipe to be not in use
49 subroutine putnode(igrid,ipe)
50 use mod_forest
51 implicit none
52
53 integer, intent(in) :: igrid, ipe
54
55 igrid_inuse(igrid,ipe)=.false.
56
57 end subroutine putnode
58
59 !> allocate arrays on igrid node
60 subroutine alloc_node(igrid)
61 use mod_forest
63 use mod_geometry
66 use mod_b0, only: set_b0_grid
67
68 integer, intent(in) :: igrid
69
70 double precision :: dx^d, summeddx, sizeuniformpart^d
71 double precision :: xext(ixglo^d-1:ixghi^d+1,1:ndim)
72 double precision :: delx_ext(ixglo1-1:ixghi1+1)
73 double precision :: exp_factor_ext(ixglo1-1:ixghi1+1),del_exp_factor_ext(ixglo1-1:ixghi1+1),exp_factor_primitive_ext(ixglo1-1:ixghi1+1)
74 double precision :: xc(ixglo1:ixghi1),delxc(ixglo1:ixghi1)
75 double precision :: exp_factor_coarse(ixglo1:ixghi1),del_exp_factor_coarse(ixglo1:ixghi1),exp_factor_primitive_coarse(ixglo1:ixghi1)
76 integer :: level, ig^d, ign^d, ixcog^l, ix, i^d
77 integer :: imin, imax, index, igco^d, ixshift, offset, ifirst
78 integer :: icase, ixgext^l
79
80 ixcogmin^d=1;
81 ixcogmax^d=(ixghi^d-2*nghostcells)/2+2*nghostcells;
82
83 icase=mod(nghostcells,2)
84 if(stagger_grid) icase=1
85 select case(icase)
86 case(0)
87 ixgext^l=ixg^ll;
88 case(1)
89 ! for ghost cell related prolongations, we need
90 ! an extra layer with known volumes and dx-intervals
91 ! in case the number of ghost cells is odd
92 ixgext^l=ixg^ll^ladd1;
93 case default
94 call mpistop("no such case")
95 end select
96
97 ! set level information
98 level=igrid_to_node(igrid,mype)%node%level
99
100 if(.not. allocated(ps(igrid)%w)) then
101
102 ! allocate arrays for solution and space
103 call alloc_state(igrid, ps(igrid), ixg^ll, ixgext^l, .true.)
104 ! allocate arrays for one level coarser solution
105 call alloc_state_coarse(igrid, psc(igrid), ixcog^l, ixcog^l)
106 if(.not.convert) then
107 ! allocate arrays for temp solution 1
108 call alloc_state(igrid, ps1(igrid), ixg^ll, ixgext^l, .false.)
109
110 ! allocate temporary solution space
111 select case (t_integrator)
113 call alloc_state(igrid, ps2(igrid), ixg^ll, ixgext^l, .false.)
115 call alloc_state(igrid, ps2(igrid), ixg^ll, ixgext^l, .false.)
116 call alloc_state(igrid, ps3(igrid), ixg^ll, ixgext^l, .false.)
117 case(imex_ars3,imex_232)
118 call alloc_state(igrid, ps2(igrid), ixg^ll, ixgext^l, .false.)
119 call alloc_state(igrid, ps3(igrid), ixg^ll, ixgext^l, .false.)
120 call alloc_state(igrid, ps4(igrid), ixg^ll, ixgext^l, .false.)
121 end select
122 end if
123
124 end if
125
126 ps(igrid)%level=level
127 psc(igrid)%level=level-1
128 if(phys_trac) ps(igrid)%special_values=0.d0
129 if(.not.convert) then
130 ps1(igrid)%level=level
131 select case (t_integrator)
133 ps2(igrid)%level=level
135 ps2(igrid)%level=level
136 ps3(igrid)%level=level
137 case(imex_ars3,imex_232)
138 ps2(igrid)%level=level
139 ps3(igrid)%level=level
140 ps4(igrid)%level=level
141 end select
142 end if
143
144 ! block pointer to current block
145 block=>ps(igrid)
146
147 ig^d=igrid_to_node(igrid,mype)%node%ig^d;
148
149 node(plevel_,igrid)=level
150 ^d&node(pig^d_,igrid)=ig^d\
151
152 ! set dx information
153 ^d&rnode(rpdx^d_,igrid)=dx(^d,level)\
154 dxlevel(:)=dx(:,level)
155
156 ! uniform cartesian case as well as all unstretched coordinates
157 ! determine the minimal and maximal corners
158 ^d&rnode(rpxmin^d_,igrid)=xprobmin^d+dble(ig^d-1)*dg^d(level)\
159 ^d&rnode(rpxmax^d_,igrid)=xprobmin^d+dble(ig^d)*dg^d(level)\
160 {if(rnode(rpxmax^d_,igrid)>xprobmax^d) rnode(rpxmax^d_,igrid)=xprobmax^d\}
161
162 ^d&dx^d=rnode(rpdx^d_,igrid)\
163 {do ix=ixglo^d,ixmhi^d-nghostcells
164 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmin^d_,igrid)+(dble(ix-nghostcells)-half)*dx^d
165 end do\}
166 ! update overlap cells of neighboring blocks in the same way to get the same values
167 {do ix=ixmhi^d-nghostcells+1,ixghi^d
168 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmax^d_,igrid)+(dble(ix-ixmhi^d)-half)*dx^d
169 end do\}
170
171 ^d&dx^d=2.0d0*rnode(rpdx^d_,igrid)\
172 {do ix=ixcogmin^d,ixcogmax^d
173 psc(igrid)%x(ix^d%ixCoG^s,^d)=rnode(rpxmin^d_,igrid)+(dble(ix-nghostcells)-half)*dx^d
174 end do\}
175
176 ^d&ps(igrid)%dx(ixgext^s,^d)=rnode(rpdx^d_,igrid);
177 ^d&psc(igrid)%dx(ixcog^s,^d)=2.0d0*rnode(rpdx^d_,igrid);
178 ^d&dx^d=rnode(rpdx^d_,igrid)\
179 {do ix=ixgextmin^d,ixgextmax^d
180 xext(ix^d%ixGext^s,^d)=rnode(rpxmin^d_,igrid)+(dble(ix-nghostcells)-half)*dx^d
181 end do\}
182
183 if(any(stretched_dim)) then
184 {if(stretch_type(^d) == stretch_uni)then
185 imin=(ig^d-1)*block_nx^d
186 imax=ig^d*block_nx^d
187 rnode(rpxmin^d_,igrid)=xprobmin^d+dxfirst_1mq(level,^d) &
188 *(1.0d0-qstretch(level,^d)**imin)
189 rnode(rpxmax^d_,igrid)=xprobmin^d+dxfirst_1mq(level,^d) &
190 *(1.0d0-qstretch(level,^d)**imax)
191 ! fix possible out of bound due to precision
192 if(rnode(rpxmax^d_,igrid)>xprobmax^d) rnode(rpxmax^d_,igrid)=xprobmax^d
193 ixshift=(ig^d-1)*block_nx^d-nghostcells
194 do ix=ixgextmin^d,ixgextmax^d
195 index=ixshift+ix
196 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxfirst(level,^d)*qstretch(level,^d)**(index-1)
197 enddo
198 igco^d=(ig^d-1)/2
199 ixshift=igco^d*block_nx^d+(1-modulo(ig^d,2))*block_nx^d/2-nghostcells
200 do ix=ixcogmin^d,ixcogmax^d
201 index=ixshift+ix
202 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxfirst(level-1,^d)*qstretch(level-1,^d)**(index-1)
203 psc(igrid)%x(ix^d%ixCoG^s,^d)=xprobmin^d+dxfirst_1mq(level-1,^d)&
204 *(1.0d0-qstretch(level-1,^d)**(index-1))&
205 + 0.5d0*dxfirst(level-1,^d)*qstretch(level-1,^d)**(index-1)
206 end do
207 ! now that dx and grid boundaries are known: fill cell centers
208 ifirst=nghostcells+1
209 ! first fill the mesh
210 summeddx=0.0d0
211 do ix=ixmlo^d,ixmhi^d
212 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmin^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixG^t,^d)
213 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
214 enddo
215 ! then ghost cells to left
216 summeddx=0.0d0
217 do ix=nghostcells,1,-1
218 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmin^d_,igrid)-summeddx-0.5d0*ps(igrid)%dx(ix^d%ixG^t,^d)
219 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
220 enddo
221 ! then ghost cells to right
222 summeddx=0.0d0
223 do ix=ixghi^d-nghostcells+1,ixghi^d
224 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmax^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixG^t,^d)
225 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
226 enddo
227 select case(icase)
228 case(0)
229 ! if even number of ghost cells: xext is just copy of local x
230 xext(ixgext^s,^d)=ps(igrid)%x(ixgext^s,^d)
231 case(1)
232 ! if uneven number of ghost cells: extra layer left/right
233 summeddx=0.0d0
234 do ix=ixmlo^d,ixmhi^d
235 xext(ix^d%ixGext^s,^d)=rnode(rpxmin^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixGext^s,^d)
236 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
237 enddo
238 ! then ghost cells to left
239 summeddx=0.0d0
240 do ix=nghostcells,ixgextmin^d,-1
241 xext(ix^d%ixGext^s,^d)=rnode(rpxmin^d_,igrid)-summeddx-0.5d0*ps(igrid)%dx(ix^d%ixGext^s,^d)
242 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
243 enddo
244 ! then ghost cells to right
245 summeddx=0.0d0
246 do ix=ixghi^d-nghostcells+1,ixgextmax^d
247 xext(ix^d%ixGext^s,^d)=rnode(rpxmax^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixGext^s,^d)
248 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
249 enddo
250 case default
251 call mpistop("no such case")
252 end select
253 endif\}
254 {if(stretch_type(^d) == stretch_symm)then
255 ! here we distinguish three kinds of grid blocks
256 ! depending on their ig-index, set per level
257 ! the first n_stretchedblocks/2 will stretch to the left
258 ! the middle ntotal-n_stretchedblocks will be uniform
259 ! the last n_stretchedblocks/2 will stretch to the right
260 if(ig^d<=nstretchedblocks(level,^d)/2)then
261 ! stretch to the left
262 offset=block_nx^d*nstretchedblocks(level,^d)/2
263 imin=(ig^d-1)*block_nx^d
264 imax=ig^d*block_nx^d
265 rnode(rpxmin^d_,igrid)=xprobmin^d+xstretch^d-dxfirst_1mq(level,^d) &
266 *(1.0d0-qstretch(level,^d)**(offset-imin))
267 rnode(rpxmax^d_,igrid)=xprobmin^d+xstretch^d-dxfirst_1mq(level,^d) &
268 *(1.0d0-qstretch(level,^d)**(offset-imax))
269 ! fix possible out of bound due to precision
270 if(rnode(rpxmin^d_,igrid)<xprobmin^d) rnode(rpxmin^d_,igrid)=xprobmin^d
271 ixshift=(ig^d-1)*block_nx^d-nghostcells
272 do ix=ixgextmin^d,ixgextmax^d
273 index=ixshift+ix
274 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxfirst(level,^d)*qstretch(level,^d)**(offset-index)
275 enddo
276 ixshift=(nstretchedblocks(level,^d)/2-ig^d)*(block_nx^d/2)+block_nx^d/2+nghostcells
277 do ix=ixcogmin^d,ixcogmax^d
278 index=ixshift-ix
279 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxfirst(level-1,^d)*qstretch(level-1,^d)**index
280 enddo
281 ! last block: to modify ghost cells!!!
282 if(ig^d==nstretchedblocks(level,^d)/2)then
283 if(ng^d(level)==nstretchedblocks(level,^d))then
284 ! if middle blocks do not exist then use symmetry
285 do ix=ixghi^d-nghostcells+1,ixgextmax^d
286 ps(igrid)%dx(ix^d%ixGext^s,^d)= &
287 ps(igrid)%dx(2*(ixghi^d-nghostcells)+1-ix^d%ixGext^s,^d)
288 enddo
289 do ix=ixcogmax^d-nghostcells+1,ixcogmax^d
290 psc(igrid)%dx(ix^d%ixCoG^s,^d)= &
291 psc(igrid)%dx(2*(ixcogmax^d-nghostcells)+1-ix^d%ixCoG^s,^d)
292 enddo
293 else
294 ! if middle blocks exist then use same as middle blocks:
295 do ix=ixghi^d-nghostcells+1,ixgextmax^d
296 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxmid(level,^d)
297 enddo
298 do ix=ixcogmax^d-nghostcells+1,ixcogmax^d
299 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxmid(level-1,^d)
300 enddo
301 endif
302 endif
303 ! first block: make ghost cells symmetric (to allow periodicity)
304 if(ig^d==1)then
305 do ix=ixgextmin^d,nghostcells
306 ps(igrid)%dx(ix^d%ixGext^s,^d)=ps(igrid)%dx(2*nghostcells+1-ix^d%ixGext^s,^d)
307 enddo
308 do ix=1,nghostcells
309 psc(igrid)%dx(ix^d%ixCoG^s,^d)=psc(igrid)%dx(2*nghostcells+1-ix^d%ixCoG^s,^d)
310 enddo
311 endif
312 else
313 if(ig^d<=ng^d(level)-nstretchedblocks(level,^d)/2) then
314 ! keep uniform
315 ps(igrid)%dx(ixgext^s,^d)=dxmid(level,^d)
316 psc(igrid)%dx(ixcog^s,^d)=dxmid(level-1,^d)
317 rnode(rpxmin^d_,igrid)=xprobmin^d+xstretch^d+(ig^d-nstretchedblocks(level,^d)/2-1)*block_nx^d*dxmid(level,^d)
318 rnode(rpxmax^d_,igrid)=xprobmin^d+xstretch^d+(ig^d-nstretchedblocks(level,^d)/2) *block_nx^d*dxmid(level,^d)
319 ! first and last block: to modify the ghost cells!!!
320 if(ig^d==nstretchedblocks(level,^d)/2+1)then
321 do ix=ixgextmin^d,nghostcells
322 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxfirst(level,^d)*qstretch(level,^d)**(nghostcells-ix)
323 enddo
324 do ix=1,nghostcells
325 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxfirst(level-1,^d)*qstretch(level-1,^d)**(nghostcells-ix)
326 enddo
327 endif
328 if(ig^d==ng^d(level)-nstretchedblocks(level,^d))then
329 do ix=ixghi^d-nghostcells+1,ixgextmax^d
330 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxfirst(level,^d)*qstretch(level,^d)**(ix-block_nx^d-nghostcells-1)
331 enddo
332 do ix=ixcogmax^d-nghostcells+1,ixcogmax^d
333 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxfirst(level-1,^d)*qstretch(level-1,^d)**(ix-ixcogmax^d+nghostcells-1)
334 enddo
335 endif
336 else
337 ! stretch to the right
338 offset=block_nx^d*(ng^d(level)-nstretchedblocks(level,^d)/2)
339 sizeuniformpart^d=dxmid(1,^d)*(domain_nx^d-nstretchedblocks_baselevel(^d)*block_nx^d)
340 imin=(ig^d-1)*block_nx^d-offset
341 imax=ig^d*block_nx^d-offset
342 rnode(rpxmin^d_,igrid)=xprobmin^d+xstretch^d+sizeuniformpart^d+dxfirst_1mq(level,^d) &
343 *(1.0d0-qstretch(level,^d)**imin)
344 rnode(rpxmax^d_,igrid)=xprobmin^d+xstretch^d+sizeuniformpart^d+dxfirst_1mq(level,^d) &
345 *(1.0d0-qstretch(level,^d)**imax)
346 ! fix possible out of bound due to precision
347 if(rnode(rpxmax^d_,igrid)>xprobmax^d) rnode(rpxmax^d_,igrid)=xprobmax^d
348 ixshift=(ig^d-1)*block_nx^d-nghostcells-offset
349 do ix=ixgextmin^d,ixgextmax^d
350 index=ixshift+ix
351 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxfirst(level,^d)*qstretch(level,^d)**(index-1)
352 enddo
353 ixshift=(ig^d+nstretchedblocks(level,^d)/2-ng^d(level)-1)*(block_nx^d/2)-nghostcells
354 do ix=ixcogmin^d,ixcogmax^d
355 index=ixshift+ix
356 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxfirst(level-1,^d)*qstretch(level-1,^d)**(index-1)
357 enddo
358 ! first block: modify ghost cells!!!
359 if(ig^d==ng^d(level)-nstretchedblocks(level,^d)+1)then
360 if(ng^d(level)==nstretchedblocks(level,^d))then
361 ! if middle blocks do not exist then use symmetry
362 do ix=ixgextmin^d,nghostcells
363 ps(igrid)%dx(ix^d%ixGext^s,^d)=ps(igrid)%dx(2*nghostcells+1-ix^d%ixGext^s,^d)
364 enddo
365 do ix=1,nghostcells
366 psc(igrid)%dx(ix^d%ixCoG^s,^d)=psc(igrid)%dx(2*nghostcells+1-ix^d%ixCoG^s,^d)
367 enddo
368 else
369 ! if middle blocks exist then use same as middle blocks:
370 do ix=ixgextmin^d,nghostcells
371 ps(igrid)%dx(ix^d%ixGext^s,^d)=dxmid(level,^d)
372 enddo
373 do ix=1,nghostcells
374 psc(igrid)%dx(ix^d%ixCoG^s,^d)=dxmid(level-1,^d)
375 enddo
376 endif
377 endif
378 ! last block: make ghost cells symmetric (to allow periodicity)
379 if(ig^d==ng^d(level))then
380 do ix=ixghi^d-nghostcells+1,ixgextmax^d
381 ps(igrid)%dx(ix^d%ixGext^s,^d)=ps(igrid)%dx(2*(ixghi^d-nghostcells)+1-ix^d%ixGext^s,^d)
382 enddo
383 do ix=ixcogmax^d-nghostcells+1,ixcogmax^d
384 psc(igrid)%dx(ix^d%ixCoG^s,^d)=psc(igrid)%dx(2*(ixcogmax^d-nghostcells)+1-ix^d%ixCoG^s,^d)
385 enddo
386 endif
387 endif
388 endif
389 ! now that dx and grid boundaries are known: fill cell centers
390 ifirst=nghostcells+1
391 ! first fill the mesh
392 summeddx=0.0d0
393 do ix=ixmlo^d,ixmhi^d
394 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmin^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixG^t,^d)
395 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
396 enddo
397 ! then ghost cells to left
398 summeddx=0.0d0
399 do ix=nghostcells,1,-1
400 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmin^d_,igrid)-summeddx-0.5d0*ps(igrid)%dx(ix^d%ixG^t,^d)
401 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
402 enddo
403 ! then ghost cells to right
404 summeddx=0.0d0
405 do ix=ixghi^d-nghostcells+1,ixghi^d
406 ps(igrid)%x(ix^d%ixG^t,^d)=rnode(rpxmax^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixG^t,^d)
407 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
408 enddo
409 ! and next for the coarse representation
410 ! first fill the mesh
411 summeddx=0.0d0
412 do ix=nghostcells+1,ixcogmax^d-nghostcells
413 psc(igrid)%x(ix^d%ixCoG^s,^d)=rnode(rpxmin^d_,igrid)+summeddx+0.5d0*psc(igrid)%dx(ix^d%ixCoG^s,^d)
414 summeddx=summeddx+psc(igrid)%dx(ix^d%ifirst,^d)
415 enddo
416 ! then ghost cells to left
417 summeddx=0.0d0
418 do ix=nghostcells,1,-1
419 psc(igrid)%x(ix^d%ixCoG^s,^d)=rnode(rpxmin^d_,igrid)-summeddx-0.5d0*psc(igrid)%dx(ix^d%ixCoG^s,^d)
420 summeddx=summeddx+psc(igrid)%dx(ix^d%ifirst,^d)
421 enddo
422 ! then ghost cells to right
423 summeddx=0.0d0
424 do ix=ixcogmax^d-nghostcells+1,ixcogmax^d
425 psc(igrid)%x(ix^d%ixCoG^s,^d)=rnode(rpxmax^d_,igrid)+summeddx+0.5d0*psc(igrid)%dx(ix^d%ixCoG^s,^d)
426 summeddx=summeddx+psc(igrid)%dx(ix^d%ifirst,^d)
427 enddo
428 select case(icase)
429 case(0)
430 ! if even number of ghost cells: xext is just copy of local x
431 xext(ixgext^s,^d)=ps(igrid)%x(ixgext^s,^d)
432 case(1)
433 ! if uneven number of ghost cells: extra layer left/right
434 summeddx=0.0d0
435 do ix=ixmlo^d,ixmhi^d
436 xext(ix^d%ixGext^s,^d)=rnode(rpxmin^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixGext^s,^d)
437 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
438 enddo
439 ! then ghost cells to left
440 summeddx=0.0d0
441 do ix=nghostcells,ixgextmin^d,-1
442 xext(ix^d%ixGext^s,^d)=rnode(rpxmin^d_,igrid)-summeddx-0.5d0*ps(igrid)%dx(ix^d%ixGext^s,^d)
443 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
444 enddo
445 ! then ghost cells to right
446 summeddx=0.0d0
447 do ix=ixghi^d-nghostcells+1,ixgextmax^d
448 xext(ix^d%ixGext^s,^d)=rnode(rpxmax^d_,igrid)+summeddx+0.5d0*ps(igrid)%dx(ix^d%ixGext^s,^d)
449 summeddx=summeddx+ps(igrid)%dx(ix^d%ifirst,^d)
450 enddo
451 case default
452 call mpistop("no such case")
453 end select
454 endif\}
455 endif
456
457 ! calculate area of cell surfaces for standard block
458 call get_surface_area(ps(igrid),ixg^ll)
459 ! calculate area of cell surfaces for coarser representative block
460 call get_surface_area(psc(igrid),ixcog^l)
461 ! calculate volume and distance of cells
462 ps(igrid)%dsC=1.d0
463 select case (coordinate)
464 case (cartesian)
465 ps(igrid)%dvolume(ixgext^s)= {^d&rnode(rpdx^d_,igrid)|*}
466 ps(igrid)%ds(ixgext^s,1:ndim)=ps(igrid)%dx(ixgext^s,1:ndim)
467 ps(igrid)%dsC(ixgext^s,1:ndim)=ps(igrid)%dx(ixgext^s,1:ndim)
468 psc(igrid)%dvolume(ixcog^s)= {^d&2.d0*rnode(rpdx^d_,igrid)|*}
469 psc(igrid)%ds(ixcog^s,1:ndim)=psc(igrid)%dx(ixcog^s,1:ndim)
470 case (cartesian_stretched)
471 ps(igrid)%dvolume(ixgext^s)= {^d&ps(igrid)%dx(ixgext^s,^d)|*}
472 ps(igrid)%ds(ixgext^s,1:ndim)=ps(igrid)%dx(ixgext^s,1:ndim)
473 ps(igrid)%dsC(ixgext^s,1:ndim)=ps(igrid)%dx(ixgext^s,1:ndim)
474 psc(igrid)%dvolume(ixcog^s)= {^d&psc(igrid)%dx(ixcog^s,^d)|*}
475 psc(igrid)%ds(ixcog^s,1:ndim)=psc(igrid)%dx(ixcog^s,1:ndim)
476 case (cartesian_expansion)
477 {^ifoned
478 delx_ext(ixgext^s) = ps(igrid)%dx(ixgext^s,1)
479 if(associated(usr_set_surface)) call usr_set_surface(ixgext^l,xext(ixgext^s,1),delx_ext(ixgext^s),exp_factor_ext(ixgext^s),del_exp_factor_ext(ixgext^s),exp_factor_primitive_ext(ixgext^s))
480 ps(igrid)%dvolume(ixgext^s)= exp_factor_primitive_ext(ixgext^s)
481 ps(igrid)%ds(ixgext^s,1)=ps(igrid)%dx(ixgext^s,1)
482 ps(igrid)%dsC(ixgext^s,1)=ps(igrid)%dx(ixgext^s,1)
483 xc(ixcog^s) = psc(igrid)%x(ixcog^s,1)
484 delxc(ixcog^s) = psc(igrid)%dx(ixcog^s,1)
485 if(associated(usr_set_surface)) call usr_set_surface(ixcog^l,xc(ixcog^s),delxc(ixcog^s),exp_factor_coarse(ixcog^s),del_exp_factor_coarse(ixcog^s),exp_factor_primitive_coarse(ixcog^s))
486 psc(igrid)%dvolume(ixcog^s)= exp_factor_primitive_coarse(ixcog^s)
487 psc(igrid)%ds(ixcog^s,1)=psc(igrid)%dx(ixcog^s,1)
488 }
489 case (spherical)
490 ps(igrid)%dvolume(ixgext^s)=(xext(ixgext^s,1)**2 &
491 +ps(igrid)%dx(ixgext^s,1)**2/12.0d0)*&
492 ps(igrid)%dx(ixgext^s,1){^nooned &
493 *two*dabs(dsin(xext(ixgext^s,2))) &
494 *dsin(half*ps(igrid)%dx(ixgext^s,2))}{^ifthreed*ps(igrid)%dx(ixgext^s,3)}
495 psc(igrid)%dvolume(ixcog^s)=(psc(igrid)%x(ixcog^s,1)**2 &
496 +psc(igrid)%dx(ixcog^s,1)**2/12.0d0)*&
497 psc(igrid)%dx(ixcog^s,1){^nooned &
498 *two*dabs(dsin(psc(igrid)%x(ixcog^s,2))) &
499 *dsin(half*psc(igrid)%dx(ixcog^s,2))}{^ifthreed*psc(igrid)%dx(ixcog^s,3)}
500 ps(igrid)%ds(ixgext^s,1)=ps(igrid)%dx(ixgext^s,1)
501 {^nooned ps(igrid)%ds(ixgext^s,2)=xext(ixgext^s,1)*ps(igrid)%dx(ixgext^s,2)}
502 {^ifthreed ps(igrid)%ds(ixgext^s,3)=xext(ixgext^s,1)*dsin(xext(ixgext^s,2))*&
503 ps(igrid)%dx(ixgext^s,3)}
504 ps(igrid)%dsC(ixgext^s,1)=ps(igrid)%dx(ixgext^s,1)
505 {^nooned ps(igrid)%dsC(ixgext^s,2)=(xext(ixgext^s,1)+half*ps(igrid)%dx(ixgext^s,1))*&
506 ps(igrid)%dx(ixgext^s,2)
507 if(ndir>ndim) then
508 ps(igrid)%dsC(ixgext^s,3)=(xext(ixgext^s,1)+half*ps(igrid)%dx(ixgext^s,1))*&
509 dsin(xext(ixgext^s,2)+half*ps(igrid)%dx(ixgext^s,2))
510 end if
511 }
512 {^ifthreed ps(igrid)%dsC(ixgext^s,3)=(xext(ixgext^s,1)+half*ps(igrid)%dx(ixgext^s,1))*&
513 dsin(xext(ixgext^s,2)+half*ps(igrid)%dx(ixgext^s,2))*&
514 ps(igrid)%dx(ixgext^s,3)}
515 case (cylindrical)
516 ps(igrid)%dvolume(ixgext^s)=dabs(xext(ixgext^s,1)) &
517 *ps(igrid)%dx(ixgext^s,1){^de&*ps(igrid)%dx(ixgext^s,^de) }
518 psc(igrid)%dvolume(ixcog^s)=dabs(psc(igrid)%x(ixcog^s,1)) &
519 *psc(igrid)%dx(ixcog^s,1){^de&*psc(igrid)%dx(ixcog^s,^de) }
520 ps(igrid)%ds(ixgext^s,r_)=ps(igrid)%dx(ixgext^s,r_)
521 ps(igrid)%dsC(ixgext^s,r_)=ps(igrid)%dx(ixgext^s,r_)
522 if(z_>0.and.z_<=ndim) then
523 ps(igrid)%ds(ixgext^s,z_)=ps(igrid)%dx(ixgext^s,z_)
524 ps(igrid)%dsC(ixgext^s,z_)=ps(igrid)%dx(ixgext^s,z_)
525 if(phi_>z_.and.ndir>ndim) then
526 ps(igrid)%dsC(ixgext^s,phi_)=xext(ixgext^s,1)+half*ps(igrid)%dx(ixgext^s,1)
527 end if
528 end if
529 if(phi_>0.and.phi_<=ndim) then
530 ps(igrid)%ds(ixgext^s,phi_)=xext(ixgext^s,1)*ps(igrid)%dx(ixgext^s,phi_)
531 ps(igrid)%dsC(ixgext^s,phi_)=(xext(ixgext^s,1)+&
532 half*ps(igrid)%dx(ixgext^s,1))*ps(igrid)%dx(ixgext^s,phi_)
533 if(z_>phi_.and.ndir>ndim) ps(igrid)%dsC(ixgext^s,z_)=1.d0
534 end if
535 case default
536 call mpistop("Sorry, coordinate unknown")
537 end select
538
539 ! initialize background non-evolving solution
540 if (b0field) call set_b0_grid(igrid)
541 if (number_equi_vars>0) call phys_set_equi_vars(igrid)
542
543 ! find the blocks on the boundaries
544 ps(igrid)%is_physical_boundary=.false.
545 {
546 do i^d=-1,1
547 if(i^d==0) cycle
548 ign^d=ig^d+i^d
549 ! blocks at periodic boundary have neighbors in the physical domain
550 ! thus threated at internal blocks with no physical boundary
551 if (periodb(^d)) ign^d=1+modulo(ign^d-1,ng^d(level))
552 if (ign^d > ng^d(level)) then
553 if(phi_ > 0 .and. poleb(2,^d)) then
554 ! if at a pole, the boundary is not physical boundary
555 ps(igrid)%is_physical_boundary(2*^d)=.false.
556 else
557 ps(igrid)%is_physical_boundary(2*^d)=.true.
558 end if
559 else if (ign^d < 1) then
560 if(phi_ > 0 .and. poleb(1,^d)) then
561 ! if at a pole, the boundary is not physical boundary
562 ps(igrid)%is_physical_boundary(2*^d-1)=.false.
563 else
564 ps(igrid)%is_physical_boundary(2*^d-1)=.true.
565 end if
566 end if
567 end do
568 \}
569 if(any(ps(igrid)%is_physical_boundary)) then
570 phyboundblock(igrid)=.true.
571 else
572 phyboundblock(igrid)=.false.
573 end if
574
575 end subroutine alloc_node
576
577 !> allocate memory to physical state of igrid node
578 subroutine alloc_state(igrid, s, ixG^L, ixGext^L, alloc_once_for_ps)
580 type(state) :: s
581 integer, intent(in) :: igrid, ixg^l, ixgext^l
582 logical, intent(in) :: alloc_once_for_ps
583 integer :: ixgs^l
584
585 allocate(s%w(ixg^s,1:nw))
586 s%igrid=igrid
587 s%w=0.d0
588 s%ixG^l=ixg^l;
589 {^d& ixgsmin^d = ixgmin^d-1; ixgsmax^d = ixgmax^d|;}
590 if(stagger_grid) then
591 allocate(s%ws(ixgs^s,1:nws))
592 s%ws=0.d0
593 if(record_electric_field) then
594 allocate(s%we(ixgs^s,1:nws))
595 s%we=0.d0
596 end if
597 s%ixGs^l=ixgs^l;
598 end if
599 if(alloc_once_for_ps) then
600 ! allocate extra variables for ps state
601 if(nw_extra>0) allocate(s%wextra(ixg^s,1:nw_extra))
602 ! allocate coordinates
603 allocate(s%x(ixg^s,1:ndim))
604 allocate(s%dx(ixgext^s,1:ndim), &
605 s%ds(ixgext^s,1:ndim),s%dsC(ixgext^s,1:3))
606 allocate(s%dvolume(ixgext^s))
607 allocate(s%surfaceC(ixgs^s,1:ndim), &
608 s%surface(ixg^s,1:ndim))
609 ! allocate physical boundary flag
610 allocate(s%is_physical_boundary(2*ndim))
611 if(local_timestep) then
612 allocate(s%dt(ixg^s))
613 endif
614
615 if(b0field) then
616 allocate(s%B0(ixg^s,1:ndir,0:ndim))
617 allocate(s%J0(ixg^s,7-2*ndir:3))
618 end if
619 if(number_equi_vars > 0) then
620 allocate(s%equi_vars(ixg^s,1:number_equi_vars,0:ndim))
621 endif
622
623 ! allocate space for special values for each block state
624 if(phys_trac) then
625 ! special_values(1) Tcoff local
626 ! special_values(2) Tmax local
627 ! special_values(3:2+ndim) Bdir local
628 allocate(s%special_values(ndim+2))
629 end if
630 else
631 ! share common info from ps states to save memory
632 if(nw_extra>0) s%wextra=>ps(igrid)%wextra
633 s%x=>ps(igrid)%x
634 s%dx=>ps(igrid)%dx
635 s%ds=>ps(igrid)%ds
636 s%dsC=>ps(igrid)%dsC
637 s%dvolume=>ps(igrid)%dvolume
638 s%surfaceC=>ps(igrid)%surfaceC
639 s%surface=>ps(igrid)%surface
640 s%is_physical_boundary=>ps(igrid)%is_physical_boundary
641 if(b0field) then
642 s%B0=>ps(igrid)%B0
643 s%J0=>ps(igrid)%J0
644 end if
645 if(number_equi_vars > 0) then
646 s%equi_vars=>ps(igrid)%equi_vars
647 endif
648 if(phys_trac) s%special_values=>ps(igrid)%special_values
649 end if
650 end subroutine alloc_state
651
652 !> allocate memory to one-level coarser physical state of igrid node
653 subroutine alloc_state_coarse(igrid, s, ixG^L, ixGext^L)
655 type(state) :: s
656 integer, intent(in) :: igrid, ixg^l, ixgext^l
657 integer :: ixgs^l
658
659 allocate(s%w(ixg^s,1:nw))
660 s%igrid=igrid
661 s%w=0.d0
662 s%ixG^l=ixg^l;
663 {^d& ixgsmin^d = ixgmin^d-1; ixgsmax^d = ixgmax^d|;}
664 if(stagger_grid) then
665 allocate(s%ws(ixgs^s,1:nws))
666 s%ws=0.d0
667 s%ixGs^l=ixgs^l;
668 end if
669 if(b0fieldalloccoarse) then
670 allocate(s%B0(ixg^s,1:ndir,0:ndim))
671 end if
672 ! allocate coordinates
673 allocate(s%x(ixg^s,1:ndim))
674 allocate(s%dx(ixgext^s,1:ndim), &
675 s%ds(ixgext^s,1:ndim),s%dsC(ixgext^s,1:3))
676 allocate(s%dvolume(ixgext^s))
677 allocate(s%surfaceC(ixgs^s,1:ndim), &
678 s%surface(ixg^s,1:ndim))
679 ! allocate physical boundary flag
680 allocate(s%is_physical_boundary(2*ndim))
681 end subroutine alloc_state_coarse
682
683 subroutine dealloc_state(igrid, s,dealloc_x)
685 integer, intent(in) :: igrid
686 type(state) :: s
687 logical, intent(in) :: dealloc_x
688
689 deallocate(s%w)
690 if(stagger_grid) then
691 deallocate(s%ws)
692 if(record_electric_field) deallocate(s%we)
693 end if
694 if(dealloc_x) then
695 if(nw_extra>0) deallocate(s%wextra)
696 ! deallocate coordinates
697 deallocate(s%x)
698 deallocate(s%dx,s%dt,s%ds,s%dsC)
699 deallocate(s%dvolume)
700 deallocate(s%surfaceC,s%surface)
701 deallocate(s%is_physical_boundary)
702 if(b0field) then
703 deallocate(s%B0)
704 deallocate(s%J0)
705 end if
706 if(number_equi_vars > 0) then
707 deallocate(s%equi_vars)
708 end if
709 if(phys_trac) then
710 deallocate(s%special_values)
711 end if
712 end if
713 nullify(s%x,s%dx,s%dt,s%ds,s%dsC,s%dvolume,s%surfaceC,s%surface)
714 nullify(s%is_physical_boundary)
715 if(b0field) nullify(s%B0,s%J0)
716 if(number_equi_vars > 0) then
717 nullify(s%equi_vars)
718 end if
719 if(nw_extra>0) nullify(s%wextra)
720 if(phys_trac) nullify(s%special_values)
721 end subroutine dealloc_state
722
723 subroutine dealloc_state_coarse(igrid, s)
725 integer, intent(in) :: igrid
726 type(state) :: s
727
728 deallocate(s%w)
729 if(stagger_grid) then
730 deallocate(s%ws)
731 end if
732 if(b0fieldalloccoarse) then
733 deallocate(s%B0)
734 end if
735 ! deallocate coordinates
736 deallocate(s%x)
737 deallocate(s%dx,s%dt,s%ds,s%dsC)
738 deallocate(s%dvolume)
739 deallocate(s%surfaceC,s%surface)
740 nullify(s%x,s%dx,s%dt,s%ds,s%dsC,s%dvolume,s%surfaceC,s%surface)
741 nullify(s%is_physical_boundary)
742 end subroutine dealloc_state_coarse
743
744 subroutine dealloc_node(igrid)
746
747 integer, intent(in) :: igrid
748
749 if (igrid==0) then
750 call mpistop("trying to delete a non-existing grid in dealloc_node")
751 end if
752
753 call dealloc_state(igrid, ps(igrid),.true.)
754 call dealloc_state_coarse(igrid, psc(igrid))
755 if(.not.convert) then
756 call dealloc_state(igrid, ps1(igrid),.false.)
757 ! deallocate temporary solution space
758 select case (t_integrator)
760 call dealloc_state(igrid, ps2(igrid),.false.)
762 call dealloc_state(igrid, ps2(igrid),.false.)
763 call dealloc_state(igrid, ps3(igrid),.false.)
764 case(imex_ars3,imex_232)
765 call dealloc_state(igrid, ps2(igrid),.false.)
766 call dealloc_state(igrid, ps3(igrid),.false.)
767 call dealloc_state(igrid, ps4(igrid),.false.)
768 end select
769 end if
770
771 end subroutine dealloc_node
772
773end module mod_amr_solution_node
subroutine, public putnode(igrid, ipe)
subroutine, public alloc_state(igrid, s, ixgl, ixgextl, alloc_once_for_ps)
allocate memory to physical state of igrid node
subroutine, public dealloc_node(igrid)
subroutine, public alloc_node(igrid)
allocate arrays on igrid node
integer function, public getnode(ipe)
Get first available igrid on processor ipe.
Definition mod_b0.t:1
subroutine, public set_b0_grid(igrid)
Definition mod_b0.t:9
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with basic grid data structures.
Definition mod_forest.t:2
logical, dimension(:,:), allocatable, save igrid_inuse
Definition mod_forest.t:70
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Definition mod_forest.t:32
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
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.
integer ixghi
Upper index of grid block arrays.
double precision, dimension(:), allocatable dg
extent of grid blocks in domain per dimension, in array over levels
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
integer, parameter imex_trapezoidal
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
logical local_timestep
each cell has its own timestep or not
integer, parameter nodehi
grid hierarchy info (level and grid indices)
integer ndir
Number of spatial dimensions (components) for vector variables.
logical b0field
split magnetic field as background B0 field
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
double precision, dimension(:,:), allocatable dx
integer nghostcells
Number of ghost cells surrounding a grid.
logical phys_trac
Use TRAC for MHD or 1D HD.
logical convert
If true and restart_from_file is given, convert snapshots to other file formats.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
integer, parameter rnodehi
grid location info (corner coordinates and grid spacing)
integer max_blocks
The maximum number of grid blocks in a processor.
logical record_electric_field
True for record electric field.
integer t_integrator
time integrator method
integer number_equi_vars
number of equilibrium set variables, besides the mag field
integer, dimension(:,:), allocatable node
integer, parameter imex_midpoint
integer, parameter ixglo
Lower index of grid block arrays (always 1)
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_set_equi_vars), pointer phys_set_equi_vars
Definition mod_physics.t:86
Module with all the methods that users can customize in AMRVAC.
procedure(set_surface), pointer usr_set_surface