MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_amr_neighbors.t
Go to the documentation of this file.
2 
3  implicit none
4  private
5 
6  public :: find_neighbor
7  public :: find_root_neighbor
8  public :: asign_tree_neighbor
9 
10 contains
11 
12 
13  !> find neighors of level-one root blocks
14  subroutine find_root_neighbor(tree_neighbor,tree,i^D)
15  use mod_forest
17  use mod_geometry
18 
19  type(tree_node_ptr) :: tree_neighbor, tree
20  integer, intent(in) :: i^d
21 
22  integer :: jg^d
23 
24  jg^d=tree%node%ig^d+i^d;
25 
26  ! find the periodic grid indices, modulo(-1,10)=9
27  {if (periodb(^d)) jg^d=1+modulo(jg^d-1,ng^d(1))\}
28 
29  ! pi-periodicity at pole
30  select case (coordinate)
31  case (spherical) {^ifthreed
32  if (poleb(1,2).and.jg2==0) then ! northpole (theta=0)
33  jg2=1; jg3=1+modulo(jg3+ng3(1)/2-1,ng3(1))
34  end if
35  if (poleb(2,2).and.jg2==ng2(1)+1) then ! southpole (theta=pi)
36  jg2=ng2(1); jg3=1+modulo(jg3+ng3(1)/2-1,ng3(1))
37  end if}
38  case (cylindrical)
39  if (poleb(1,1).and.jg1==0) then ! cylindrical axis
40  jg1=1
41  {if (^d==phi_) jg^d=1+modulo(jg^d+ng^d(1)/2-1,ng^d(1))\}
42  end if
43  end select
44 
45  if (jg^d>=1.and.jg^d<=ng^d(1)|.and.) then
46  tree_neighbor%node => tree_root(jg^d)%node
47  else
48  nullify(tree_neighbor%node)
49  end if
50 
51  end subroutine find_root_neighbor
52 
53  !> find neighors of all blocks
54  subroutine find_neighbor(my_neighbor,my_neighbor_type,tree,i^D,pole)
55  use mod_forest
57 
58  type(tree_node_ptr) :: tree, my_neighbor
59  integer, intent(in) :: i^d
60  integer, intent(out) :: my_neighbor_type
61  logical, dimension(ndim), intent(out) :: pole
62 
63  integer :: level, ig^d, ic^d, n_ic^d, inp^d
64 
65  pole=.false.
66  level=tree%node%level
67  if (level==1) then
68  call find_root_neighbor(my_neighbor,tree,i^d)
69  if (associated(my_neighbor%node)) then
70  if (phi_ > 0) then
71  ig^d=tree%node%ig^d;
72  {if ((poleb(2,^d).and.ig^d==ng^d(1).and.i^d==1) .or. &
73  (poleb(1,^d).and.ig^d==1.and.i^d==-1)) pole(^d)=.true.\}
74  end if
75  if (my_neighbor%node%leaf) then
76  my_neighbor_type=3
77  else
78  my_neighbor_type=4
79  end if
80  else
81  my_neighbor_type=1
82  return
83  end if
84  else
85  ig^d=tree%node%ig^d;
86 
87  if (phi_ > 0) then
88  {if ((poleb(2,^d).and.ig^d==ng^d(level).and.i^d==1) .or. &
89  (poleb(1,^d).and.ig^d==1.and.i^d==-1)) pole(^d)=.true.\}
90  end if
91 
92  ! ic^D is 1 when ig^D is odd, is 2 when ig^D is even
93  ic^d=1+modulo(ig^d-1,2);
94  inp^d=int((ic^d+i^d+1)/2)-1;
95  my_neighbor%node => tree%node%parent%node
96  {if (inp^d/=0) then
97  my_neighbor%node => my_neighbor%node%neighbor(ic^d,^d)%node
98  if (.not.associated(my_neighbor%node)) then
99  my_neighbor_type=1
100  return
101  end if
102  end if\}
103  if (my_neighbor%node%leaf) then
104  my_neighbor_type=2
105  else
106  {if (i^d==0 .or. pole(^d)) then
107  n_ic^d=ic^d
108  else
109  n_ic^d=3-ic^d ! switch 1 <--> 2
110  end if\}
111  my_neighbor%node => my_neighbor%node%child(n_ic^d)%node
112  if (associated(my_neighbor%node)) then
113  if (my_neighbor%node%leaf) then
114  my_neighbor_type=3
115  else
116  my_neighbor_type=4
117  end if
118  else
119  my_neighbor_type=0
120  end if
121  end if
122  end if
123 
124  end subroutine find_neighbor
125 
126  !> asign tree node neighor
127  subroutine asign_tree_neighbor(tree)
128  use mod_forest
130 
131  type(tree_node_ptr) :: tree
132 
133  logical, dimension(ndim) :: pole
134  integer :: my_neighbor_type, i^d, iside
135  type(tree_node_ptr) :: my_neighbor
136 
137  {do iside=1,2
138  i^dd=kr(^dd,^d)*(2*iside-3);
139  call find_neighbor(my_neighbor,my_neighbor_type,tree,i^dd,pole)
140  select case (my_neighbor_type)
141  case (neighbor_sibling, neighbor_fine)
142  tree%node%neighbor(iside,^d)%node => my_neighbor%node
143  if (associated(my_neighbor%node)) then
144  if (pole(^d)) then
145  my_neighbor%node%neighbor(iside,^d)%node => tree%node
146  else
147  my_neighbor%node%neighbor(3-iside,^d)%node => tree%node
148  end if
149  end if
150  case default
151  nullify(tree%node%neighbor(iside,^d)%node)
152  end select
153  end do\}
154 
155  end subroutine asign_tree_neighbor
156 
157 end module mod_amr_neighbors
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
Module with basic grid data structures.
Definition: mod_forest.t:2
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:7
integer, parameter spherical
Definition: mod_geometry.t:11
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, dimension(:), allocatable ng
number of grid blocks in domain per dimension, in array over levels
integer, dimension(:), allocatable, parameter d
logical, dimension(ndim) periodb
True for dimensions with periodic boundaries.
logical, dimension(2, ndim) poleb
Indicates whether there is a pole at a boundary.
Pointer to a tree_node.
Definition: mod_forest.t:6