6 integer :: itag, irecv, isend
7 integer,
dimension(:),
allocatable :: recvrequest, sendrequest
8 integer,
dimension(:,:),
allocatable :: recvstatus, sendstatus
11 integer,
dimension(:),
allocatable :: recvrequest_stg, sendrequest_stg
12 integer,
dimension(:,:),
allocatable :: recvstatus_stg, sendstatus_stg
39 integer :: iigrid, igrid, ipe, igridco, ipeco, level, ic^
d
40 integer,
dimension(2^D&) :: igridfi, ipefi
41 integer :: n_coarsen, n_refine
52 n_coarsen = count(
coarsen(:, :))
53 n_refine = count(
refine(:, :))
60 recvrequest=mpi_request_null
61 sendrequest=mpi_request_null
66 recvrequest_stg=mpi_request_null
67 sendrequest_stg=mpi_request_null
77 sibling%node => tree%node%child(ic^
d)%node
78 ipefi(ic^
d)=sibling%node%ipe
79 igridfi(ic^
d)=sibling%node%igrid
83 igridco=getnode(ipeco)
85 call coarsen_tree_leaf(igridco,ipeco,igridfi,ipefi,active)
87 call coarsen_grid_siblings(igridco,ipeco,igridfi,ipefi,active)
91 if (ipefi(ic^d)==ipeco)
then
92 call putnode(igridfi(ic^d),ipefi(ic^d))
93 coarsen(igridfi(ic^d),ipefi(ic^d))=.false.
101 call mpi_waitall(irecv,recvrequest,recvstatus,ierrmpi)
102 if(stagger_grid)
call mpi_waitall(irecv,recvrequest_stg,recvstatus_stg,ierrmpi)
105 call mpi_waitall(isend,sendrequest,sendstatus,ierrmpi)
106 if(stagger_grid)
call mpi_waitall(isend,sendrequest_stg,sendstatus_stg,ierrmpi)
109 deallocate(recvstatus,recvrequest,sendstatus,sendrequest)
110 if(stagger_grid)
deallocate(recvstatus_stg,recvrequest_stg,sendstatus_stg,sendrequest_stg)
114 do igrid=1,max_blocks
115 if (coarsen(igrid,ipe))
then
121 call putnode(igrid,ipe)
122 coarsen(igrid,ipe)=.false.
128 do igrid=1,max_blocks
129 if (refine(igrid,ipe))
then
132 igridfi(ic^d)=getnode(ipe)
136 call refine_tree_leaf(igridfi,ipefi,igrid,ipe,active)
138 if (ipe==mype)
call refine_grids(igridfi,ipefi,igrid,ipe,active)
141 call putnode(igrid,ipe)
142 refine(igrid,ipe)=.false.
151 if(.not.time_advance)
call mpi_barrier(icomm,ierrmpi)
153 if(stagger_grid)
call end_comm_faces
158 call amr_morton_order()
164 call build_connectivity
172 if (time_advance)
then
173 call getbc(global_time+dt,0.d0,ps,iwstart,nwgc)
175 call getbc(global_time,0.d0,ps,iwstart,nwgc)
179 if (use_multigrid)
call mg_update_refinement(n_coarsen, n_refine)
182 if (
associated(usr_after_refine))
then
183 call usr_after_refine(n_coarsen, n_refine)
195 logical,
dimension(:,:),
allocatable :: refine2
196 integer :: iigrid, igrid, level, ic^D, inp^D, i^D, my_neighbor_type,ipe
197 logical :: coarsening, pole(ndim), sendbuf(max_blocks)
198 type(
tree_node_ptr) :: tree, p_neighbor, my_parent, sibling, my_neighbor, &
202 allocate(refine2(max_blocks,
npe))
203 call mpi_allreduce(
refine,refine2,max_blocks*
npe,mpi_logical,mpi_lor, &
208 call mpi_allgather(sendbuf,max_blocks,mpi_logical,
refine,max_blocks, &
215 if (.not.
associated(tree%node))
exit
217 if (
refine(tree%node%igrid,tree%node%ipe))
then
218 ic^d=1+modulo(tree%node%ig^d-1,2);
219 {
do inp^db=ic^db-2,ic^db-1\}
220 if (inp^d==0|.and.) cycle
221 p_neighbor%node => tree%node%parent%node
223 p_neighbor%node => p_neighbor%node%neighbor(ic^d,^d)%node
224 if (.not.
associated(p_neighbor%node)) cycle
226 if (p_neighbor%node%leaf)
then
227 refine(p_neighbor%node%igrid,p_neighbor%node%ipe)=.true.
232 tree%node => tree%node%next%node
239 do iigrid=1,igridstail; igrid=igrids(iigrid);
240 if (refine(igrid,mype).and.coarsen(igrid,mype)) coarsen(igrid,mype)=.false.
244 sendbuf(:)=coarsen(:,mype)
245 call mpi_allgather(sendbuf,max_blocks,mpi_logical,coarsen,max_blocks, &
246 mpi_logical,icomm,ierrmpi)
248 do level=levmax,max(2,levmin),-1
249 tree%node => level_head(level)%node
251 if (.not.
associated(tree%node))
exit
253 if (coarsen(tree%node%igrid,tree%node%ipe))
then
255 my_parent%node => tree%node%parent%node
258 check1: {
do ic^db=1,2\}
259 sibling%node => my_parent%node%child(ic^d)%node
260 if (sibling%node%leaf)
then
261 if (coarsen(sibling%node%igrid,sibling%node%ipe)) cycle
270 check2: {
do ic^db=1,2\}
271 sibling%node => my_parent%node%child(ic^d)%node
272 {
do i^db=ic^db-2,ic^db-1\}
273 if (i^d==0|.and.) cycle
274 call find_neighbor(my_neighbor,my_neighbor_type, &
276 select case (my_neighbor_type)
277 case (neighbor_sibling)
278 if (refine(my_neighbor%node%igrid, &
279 my_neighbor%node%ipe))
then
286 neighborchild%node=>my_neighbor%node%child(1^d&)%node
287 if (neighborchild%node%leaf)
then
288 if (coarsen(neighborchild%node%igrid, &
289 neighborchild%node%ipe))
then
302 tree%node => tree%node%next%node
311 type(tree_node_ptr) :: sibling
314 sibling%node => my_parent%node%child(ic^d)%node
315 if (sibling%node%leaf)
then
316 coarsen(sibling%node%igrid,sibling%node%ipe)=.false.
326 subroutine coarsen_grid_siblings(igrid,ipe,child_igrid,child_ipe,active)
332 integer,
intent(in) :: igrid, ipe
333 integer,
dimension(2^D&),
intent(in) :: child_igrid, child_ipe
334 logical,
intent(in) :: active
336 integer :: igridFi, ipeFi, ixCo^L, ixCoG^L, ixCoM^L, ic^D, idir
341 if (.not. active)
then
342 if (ipe ==
mype)
then
345 igridfi=child_igrid(ic^d)
346 ipefi=child_ipe(ic^d)
357 igridfi=child_igrid(ic^d)
358 ipefi=child_ipe(ic^d)
360 if (ipefi==mype)
then
361 ^d&dxlevel(^d)=rnode(rpdx^d_,igridfi);
363 ixcomin^d=ixmlo^d+(ic^d-1)*(ixmhi^d-ixmlo^d+1)/2;
364 ixcomax^d=ixmhi^d+(ic^d-2)*(ixmhi^d-ixmlo^d+1)/2;
366 call coarsen_grid(ps(igridfi),ixg^ll,ixm^ll,ps(igrid),ixg^ll,ixco^l)
371 ixcogmax^d=ixghi^d/2+nghostcells;
372 ixcom^l=ixcog^l^lsubnghostcells;
373 call coarsen_grid(ps(igridfi),ixg^ll,ixm^ll,psc(igridfi), &
379 call mpi_isend(psc(igridfi)%w,1,type_coarse_block,ipe,itag, &
380 icomm,sendrequest(isend),ierrmpi)
381 if(stagger_grid)
then
384 itag_stg=(npe+ipefi+1)+igridfi*(ndir-1+idir)
385 call mpi_isend(psc(igridfi)%ws,1,type_coarse_block_stg(idir,ic^d),ipe,itag_stg, &
386 icomm,sendrequest_stg(isend),ierrmpi)
395 call mpi_irecv(ps(igrid)%w,1,type_sub_block(ic^d),ipefi,itag, &
396 icomm,recvrequest(irecv),ierrmpi)
397 if(stagger_grid)
then
400 itag_stg=(npe+ipefi+1)+igridfi*(ndir-1+idir)
401 call mpi_irecv(ps(igrid)%ws,1,type_sub_block_stg(idir,ic^d),ipefi,itag_stg, &
402 icomm,recvrequest_stg(irecv),ierrmpi)
409 end subroutine coarsen_grid_siblings
subroutine unflag_coarsen_siblings
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 find_neighbor(my_neighbor, my_neighbor_type, tree, iD, pole)
find neighors of all blocks
subroutine, public putnode(igrid, ipe)
subroutine, public alloc_node(igrid)
allocate arrays on igrid node
integer function, public getnode(ipe)
Get first available igrid on processor ipe.
Module to coarsen and refine grids for AMR.
subroutine proper_nesting
For all grids on all processors, do a check on refinement flags. Make sure that neighbors will not di...
subroutine, public amr_coarsen_refine
coarsen and refine blocks to update AMR grid
subroutine, public coarsen_grid(sFi, ixFiGL, ixFiL, sCo, ixCoGL, ixCoL)
coarsen one grid to its coarser representative
Module with basic grid data structures.
logical, dimension(:,:), allocatable, save refine
logical, dimension(:,:), allocatable, save coarsen
AMR flags and grids-in-use identifier per processor (igrid,ipe)
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
type(tree_node_ptr), dimension(:), allocatable, save level_head
The head pointer of the linked list per refinement level.
subroutine, public get_level_range
subroutine, public build_connectivity
subroutine, public getigrids
subroutine, public refine_tree_leaf(child_igrid, child_ipe, igrid, ipe, active)
subroutine, public coarsen_tree_leaf(igrid, ipe, child_igrid, child_ipe, active)
update ghost cells of all blocks including physical boundaries
This module contains definitions of global parameters and variables and some generic functions/subrou...
logical stagger_grid
True for using stagger grid.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
integer nbufferx
Number of cells as buffer zone.
integer refine_max_level
Maximal number of AMR levels.
integer max_blocks
The maximum number of grid blocks in a processor.
subroutine, public initial_condition(igrid)
fill in initial condition
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
subroutine, public refine_grids(child_igrid, child_ipe, igrid, ipe, active)
refine one block to its children blocks
subroutine, public selectgrids
Module with all the methods that users can customize in AMRVAC.
procedure(after_refine), pointer usr_after_refine