MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_selectgrids.t
Go to the documentation of this file.
2 
3  implicit none
4  private
5 
6  public :: selectgrids
7 
8 contains
9 
10 
11  !=============================================================================
12  subroutine selectgrids
13 
14  use mod_forest
16  integer :: iigrid, igrid, jgrid, kgrid, isave, my_isafety
17  integer, allocatable, dimension(:,:) :: isafety
18 
19  ! Set the number of safety-blocks (additional blocks after
20  ! flag_grid_usr):
21  integer, parameter :: nsafety = 1
22  integer :: ixo^l, ipe
23  type(tree_node_ptr) :: tree
24  integer :: userflag
25  !-----------------------------------------------------------------------------
26  if (.not. allocated(isafety)) &
27  allocate(isafety(max_blocks,0:npe-1))
28 
29  ! reset all grids to active:
30  neighbor_active = .true.
31 
32  jgrid=0
33  kgrid=0
34  isafety = -1
35  userflag = -1
36 
37  ! Check the user flag:
38  do iigrid=1,igridstail; igrid=igrids(iigrid);
39  userflag = igrid_active(igrid)
40  if (userflag <= 0) then
41  jgrid=jgrid+1
42  igrids_active(jgrid)=igrid
43  else
44  kgrid=kgrid+1
45  igrids_passive(kgrid)=igrid
46  end if
47  isafety(igrid,mype) = userflag
48  end do
49 
50  igridstail_active = jgrid
51  igridstail_passive = kgrid
52 
53  ! Check if user wants to deactivate grids at all and return if not:
54  if (userflag == -1) return
55 
56  ! Got the passive grids.
57  ! Now, we re-activate a safety belt of radius nsafety blocks.
58 
59  ! First communicate the current isafety buffer:
60  call mpi_allgather(isafety(:,mype),max_blocks,mpi_integer,isafety,&
61  max_blocks,mpi_integer,icomm,ierrmpi)
62 
63  ! Now check the distance of neighbors to the active zone:
64  do isave = 1, nsafety
65  do iigrid=1,igridstail_passive; igrid=igrids_passive(iigrid);
66  if ( isafety(igrid,mype) /= isave) cycle
67  ! Get the minimum neighbor isafety:
68  if(min_isafety_neighbor(igrid) >= isafety(igrid,mype)) then
69  ! Increment the buffer:
70  isafety(igrid,mype)=isafety(igrid,mype)+1
71  end if
72  end do
73 
74  ! Communicate the incremented buffers:
75  call mpi_allgather(isafety(:,mype),max_blocks,mpi_integer,isafety,&
76  max_blocks,mpi_integer,icomm,ierrmpi)
77  end do
78 
79  ! Update the active and passive arrays:
80  jgrid=0
81  kgrid=0
82  do iigrid=1,igridstail; igrid=igrids(iigrid);
83  if (isafety(igrid,mype) <= nsafety) then
84  jgrid=jgrid+1
85  igrids_active(jgrid)=igrid
86  else
87  kgrid=kgrid+1
88  igrids_passive(kgrid)=igrid
89  end if
90  ! Create the neighbor flags:
91  call set_neighbor_state(igrid)
92  end do
93  igridstail_active = jgrid
94  igridstail_passive = kgrid
95 
96  ! Update the tree:
97  nleafs_active = 0
98  do ipe=0,npe-1
99  do igrid=1,max_blocks
100  if (isafety(igrid,ipe) == -1) cycle
101  if (.not.associated(igrid_to_node(igrid,ipe)%node)) cycle
102  tree%node => igrid_to_node(igrid,ipe)%node
103  if (isafety(igrid,ipe) > nsafety) then
104  tree%node%active=.false.
105  else
106  tree%node%active=.true.
108  end if
109  end do
110  end do
111 
112  ! Just for output and testing:
113  ! ixO^L=ixG^LL^LSUBnghostcells;
114  ! do iigrid=1,igridstail; igrid=igrids(iigrid);
115  ! ps(igrid)%w(ixO^S,flg_) = dble(isafety(igrid,mype))
116  ! ps(igrid)%w(ixO^S,cpu_) = dble(mype)
117  ! end do
118 
119  contains
120  !=============================================================================
121  subroutine set_neighbor_state(igrid)
122 
123  integer, intent(in) :: igrid
124  integer :: my_neighbor_type, i^D, isafety_neighbor
125  !-----------------------------------------------------------------------------
126 
127 
128  {do i^db=-1,1\}
129  if (i^d==0|.and.) then
130  if (isafety(igrid,mype) > nsafety) &
131  neighbor_active(i^d,igrid) = .false.
132  end if
133  my_neighbor_type=neighbor_type(i^d,igrid)
134 
135  select case (my_neighbor_type)
136  case (neighbor_boundary) ! boundary
137  isafety_neighbor = nsafety+1
138  case (neighbor_coarse) ! fine-coarse
139  isafety_neighbor = isafety_fc(i^d,igrid)
140  case (neighbor_sibling) ! same level
141  isafety_neighbor = isafety_srl(i^d,igrid)
142  case (neighbor_fine) ! coarse-fine
143  isafety_neighbor = isafety_cf_max(i^d,igrid)
144  end select
145 
146  if (isafety_neighbor > nsafety) &
147  neighbor_active(i^d,igrid) = .false.
148 
149  {end do\}
150 
151  end subroutine set_neighbor_state
152  !=============================================================================
153  integer function min_isafety_neighbor(igrid)
154 
155  integer, intent(in) :: igrid
156  integer :: my_neighbor_type, i^D
157  !-----------------------------------------------------------------------------
158 
159  min_isafety_neighbor = biginteger
160 
161  {do i^db=-1,1\}
162  if (i^d==0|.and.) cycle
163  my_neighbor_type=neighbor_type(i^d,igrid)
164 
165  select case (my_neighbor_type)
166  case (neighbor_coarse) ! fine-coarse
167  min_isafety_neighbor = min(isafety_fc(i^d,igrid),&
168  min_isafety_neighbor)
169  case (neighbor_sibling) ! same level
170  min_isafety_neighbor = min(isafety_srl(i^d,igrid),&
171  min_isafety_neighbor)
172  case (neighbor_fine) ! coarse-fine
173  min_isafety_neighbor = min(isafety_cf_min(i^d,igrid),&
174  min_isafety_neighbor)
175  end select
176 
177  {end do\}
178 
179  end function min_isafety_neighbor
180  !=============================================================================
181  integer function isafety_fc(i^D,igrid)
182 
183  integer, intent(in) :: i^D, igrid
184  integer :: ineighbor, ipe_neighbor
185  !-----------------------------------------------------------------------------
186  ineighbor=neighbor(1,i^d,igrid)
187  ipe_neighbor=neighbor(2,i^d,igrid)
188 
189  isafety_fc = isafety(ineighbor,ipe_neighbor)
190 
191  end function isafety_fc
192  !=============================================================================
193  integer function isafety_srl(i^D,igrid)
194 
195  integer, intent(in) :: i^D, igrid
196  integer :: ineighbor, ipe_neighbor
197  !-----------------------------------------------------------------------------
198  ineighbor=neighbor(1,i^d,igrid)
199  ipe_neighbor=neighbor(2,i^d,igrid)
200 
201  isafety_srl = isafety(ineighbor,ipe_neighbor)
202 
203  end function isafety_srl
204  !=============================================================================
205  integer function isafety_cf_min(i^D,igrid)
206 
207  integer, intent(in) :: i^D, igrid
208  integer :: ic^D, inc^D
209  integer :: ineighbor, ipe_neighbor
210  !-----------------------------------------------------------------------------
211 
212  isafety_cf_min = biginteger
213 
214  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
215  inc^db=2*i^db+ic^db
216  \}
217  ineighbor = neighbor_child(1,inc^d,igrid)
218  ipe_neighbor = neighbor_child(2,inc^d,igrid)
219 
220  isafety_cf_min = min(isafety_cf_min,isafety(ineighbor,ipe_neighbor))
221 
222  {end do\}
223 
224  end function isafety_cf_min
225  !=============================================================================
226  integer function isafety_cf_max(i^D,igrid)
227 
228  integer, intent(in) :: i^D, igrid
229  integer :: ic^D, inc^D
230  integer :: ineighbor, ipe_neighbor
231  !-----------------------------------------------------------------------------
232 
233  isafety_cf_max = - biginteger
234 
235  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
236  inc^db=2*i^db+ic^db
237  \}
238  ineighbor = neighbor_child(1,inc^d,igrid)
239  ipe_neighbor = neighbor_child(2,inc^d,igrid)
240 
241  isafety_cf_max = max(isafety_cf_max,isafety(ineighbor,ipe_neighbor))
242 
243  {end do\}
244 
245  end function isafety_cf_max
246  !=============================================================================
247  end subroutine selectgrids
248  !=============================================================================
249 
250  integer function igrid_active(igrid)
251  use mod_usr_methods, only: usr_flag_grid
253 
254  integer, intent(in) :: igrid
255  integer :: ixO^L, flag
256  !-----------------------------------------------------------------------------
257  ixo^l=ixg^ll^lsubnghostcells;
258 
259  igrid_active = -1
260 
261  if (associated(usr_flag_grid)) then
262  call usr_flag_grid(global_time,ixg^ll,ixo^l,ps(igrid)%w,ps(igrid)%x,igrid_active)
263  end if
264 
265  end function igrid_active
266  !=============================================================================
267 
268 end module mod_selectgrids
subroutine set_neighbor_state(igrid)
Module with basic grid data structures.
Definition: mod_forest.t:2
integer nleafs_active
Definition: mod_forest.t:78
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
This module contains definitions of global parameters and variables and some generic functions/subrou...
double precision global_time
The global simulation time.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
integer max_blocks
The maximum number of grid blocks in a processor.
subroutine, public selectgrids
Module with all the methods that users can customize in AMRVAC.
procedure(flag_grid), pointer usr_flag_grid
Pointer to a tree_node.
Definition: mod_forest.t:6