25 integer :: ig^
d, level, igrid, ipe
26 integer :: iside, i^
d, morton_no, isfc
53 {
do ig^db=1,
ng^db(1)\}
55 i^dd=
kr(^dd,^
d)*(2*iside-3);
62 call amr_morton_order()
71 integer,
intent(in) :: ig^D, level, igrid, ipe
72 logical,
intent(in) :: active
83 tree%node%active=active
85 nullify(tree%node%parent%node)
87 nullify(tree%node%child(ic^d)%node)
93 nullify({tree%node%neighbor(1,^d)%node},{tree%node%neighbor(2,^d)%node})
95 igrid_to_node(igrid,ipe)%node => tree%node
103 integer,
intent(in) :: igrid, ipe
104 integer,
dimension(2^D&),
intent(in) :: child_igrid, child_ipe
105 logical,
intent(out) :: active
108 integer :: level, ic^d, child_level, iside, iotherside, vote
111 tree%node =>
igrid_to_node(child_igrid(1^d&),child_ipe(1^d&))%node%parent%node
112 level=tree%node%level
120 child%node => tree%node%child(ic^d)%node
123 if(child%node%active) vote=vote+1
129 child_neighbor%node => child%node%neighbor(iside,^d)%node
130 if (
associated(child_neighbor%node))
then
131 if (child%node%ig^d==child_neighbor%node%ig^d)
then
132 nullify(child_neighbor%node%neighbor(iside,^d)%node)
135 nullify(child_neighbor%node%neighbor(iotherside,^d)%node)
139 nullify(tree%node%child(ic^d)%node)
140 deallocate(igrid_to_node(child_igrid(ic^d),child_ipe(ic^d))%node)
143 tree%node%leaf=.true.
144 tree%node%igrid=igrid
146 igrid_to_node(igrid,ipe)%node => tree%node
150 if (vote /= 2**^nd)
then
152 tree%node%active = .false.
153 nleafs_active = nleafs_active - vote
155 tree%node%active = .true.
156 nleafs_active = nleafs_active - vote + 1
158 active = tree%node%active
160 nleafs=nleafs-2**^nd+1
162 nleafs_level(child_level)=nleafs_level(child_level)-2**^nd
163 nleafs_level(level)=nleafs_level(level)+1
172 integer,
dimension(2^D&),
intent(in) :: child_igrid, child_ipe
173 integer,
intent(in) :: igrid, ipe
174 logical,
intent(out):: active
176 integer :: ig^
d, level, i^
d, ic^
d, child_ig^
d, child_level, iside
177 integer :: my_neighbor_type
178 logical,
dimension(ndim) :: pole
183 level=tree%node%level
184 active=tree%node%active
188 tree%node%leaf=.false.
189 tree%node%active=.true.
196 child_ig^
d=2*(ig^
d-1)+ic^
d;
198 child_igrid(ic^
d),child_ipe(ic^
d),active)
202 tree%node%child(ic^
d)%node => child%node
203 child%node%parent%node => tree%node
208 child%node => tree%node%child(ic^d)%node
210 i^dd=kr(^dd,^d)*(2*iside-3);
211 call find_neighbor(my_neighbor,my_neighbor_type,child,i^dd,pole)
212 select case (my_neighbor_type)
213 case (neighbor_sibling, neighbor_fine)
214 child%node%neighbor(iside,^d)%node => my_neighbor%node
216 my_neighbor%node%neighbor(iside,^d)%node => child%node
218 my_neighbor%node%neighbor(3-iside,^d)%node => child%node
221 nullify(child%node%neighbor(iside,^d)%node)
223 child%node%neighbor(3-ic^d,^d)%node=>tree%node%child(3-ic^d^d%ic^dd)%node\}
226 nleafs=nleafs+2**^nd-1
228 nleafs_level(child_level)=nleafs_level(child_level)+2**^nd
229 nleafs_level(level)=nleafs_level(level)-1
231 if (active) nleafs_active = nleafs_active + 2**^nd-1
239 integer,
intent(in) :: recv_igrid, recv_ipe, send_igrid, send_ipe
245 tree%node%igrid=recv_igrid
246 tree%node%ipe=recv_ipe
257 integer,
intent(in) :: level
260 nullify(tree%node%next%node)
268 nullify(tree%node%prev%node)
277 integer,
intent(in) :: level
282 prev%node => tree%node%prev%node
283 next%node => tree%node%next%node
284 if (
associated(next%node).and.
associated(prev%node))
then
285 prev%node%next%node => next%node
286 next%node%prev%node => prev%node
287 else if (
associated(prev%node))
then
289 nullify(prev%node%next%node)
290 else if (
associated(next%node))
then
292 nullify(next%node%prev%node)
304 integer,
intent(in) :: file_handle
306 integer,
dimension(MPI_STATUS_SIZE) :: status
323 call mpi_file_write(file_handle,tree%node%leaf,1,mpi_logical,status,
ierrmpi)
325 if (.not.tree%node%leaf)
then
343 integer,
intent(in) :: file_handle
345 integer,
dimension(MPI_STATUS_SIZE) :: status
347 integer :: ig^
d, level, morton_no, igrid, ipe, isfc
379 integer,
intent(in) :: ig^
d, level
381 integer :: ic^
d, child_ig^
d, child_level
385 call mpi_file_read(file_handle,leaf,1,mpi_logical, &
392 tree%node%level=level
393 tree%node%active=.true. .and. leaf
396 nullify(tree%node%child(ic^
d)%node)
398 nullify({tree%node%neighbor(1,^
d)%node},{tree%node%neighbor(2,^
d)%node})
399 nullify(tree%node%next%node,tree%node%prev%node)
407 morton_no=morton_no+1
410 tree%node%igrid=igrid
420 child_ig^
d=2*(ig^
d-1)+ic^
d;
421 allocate(tree%node%child(ic^
d)%node)
422 tree%node%child(ic^
d)%node%parent%node => tree%node
423 call read_node(tree%node%child(ic^
d),child_ig^
d,child_level)
recursive subroutine write_node(tree)
recursive subroutine read_node(tree, igD, level)
subroutine, public find_root_neighbor(tree_neighbor, tree, iD)
find neighors of level-one root blocks
subroutine, public asign_tree_neighbor(tree)
asign tree node neighor
subroutine, public find_neighbor(my_neighbor, my_neighbor_type, tree, iD, pole)
find neighors of all blocks
integer function, public getnode(ipe)
Get first available igrid on processor ipe.
Module with basic grid data structures.
integer, dimension(:), allocatable, save nleafs_level
How many leaves are present per refinement level.
integer, dimension(:), allocatable, save sfc_to_igrid
Go from a Morton number to an igrid index (for a single processor)
integer, save nleafs
Number of leaf block.
integer, dimension(:), allocatable, save igrid_to_sfc
Go from a grid index to Morton number (for a single processor)
integer, dimension(:,:), allocatable, save sfc_iglevel1
Space filling curve for level 1 grid. sfc_iglevel1(^D, MN) gives ig^D (the spatial index of the grid)
type(tree_node_ptr), dimension(:), allocatable, save level_tail
The tail pointer of the linked list per refinement level.
integer, dimension(:), allocatable, save morton_stop
Last Morton number per processor.
type(tree_node_ptr), dimension(:^d &), allocatable, save tree_root
Pointers to the coarse grid.
integer, save nparents
Number of parent blocks.
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 change_ipe_tree_leaf(recv_igrid, recv_ipe, send_igrid, send_ipe)
subroutine delete_from_linked_list(level, tree)
subroutine, public write_forest(file_handle)
subroutine, public read_forest(file_handle)
subroutine, public refine_tree_leaf(child_igrid, child_ipe, igrid, ipe, active)
subroutine add_to_linked_list(level, tree)
subroutine, public coarsen_tree_leaf(igrid, ipe, child_igrid, child_ipe, active)
subroutine init_tree_leaf(tree, igD, level, igrid, ipe, active)
subroutine, public init_forest_root
build root forest
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter nlevelshi
The maximum number of levels in the grid refinement.
integer icomm
The MPI communicator.
integer, dimension(:), allocatable ng
number of grid blocks in domain per dimension, in array over levels
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
subroutine get_morton_range
Set the Morton range for each processor.
subroutine amr_morton_order
Construct Morton-order as a global recursive lexicographic ordering.
subroutine level1_morton_order
build Morton space filling curve for level 1 grid