2D AMRCLAW
Functions/Subroutines
birect.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine birect (mptr1)
 Check each grid, starting with mptr1 (either newstl or lstart) to see that it has no more than max1d points in either dimensions. More...
 

Function/Subroutine Documentation

◆ birect()

subroutine birect (   mptr1)

Check each grid, starting with mptr1 (either newstl or lstart) to see that it has no more than max1d points in either dimensions.

needed so that scratch array space in stepgrid not exceeded.

Also check for too small grids - but has never happened.

Definition at line 11 of file birect.f.

References amr_module::cornxhi, amr_module::cornxlo, amr_module::cornyhi, amr_module::cornylo, amr_module::hxposs, amr_module::hyposs, amr_module::intratx, amr_module::intraty, amr_module::levelptr, amr_module::max1d, amr_module::ndihi, amr_module::ndilo, amr_module::ndjhi, amr_module::ndjlo, amr_module::nestlevel, amr_module::nghost, amr_module::node, nodget(), amr_module::rnode, and amr_module::timemult.

Referenced by domain(), and grdfit().

11 c
12  use amr_module
13  implicit double precision (a-h,o-z)
14 
15 
16 c
17 c ::::::::::::: BIRECT :::::::::::::::::::::::::::::::::::::::
18 c check each grid, starting with mptr1 (either newstl or lstart)
19 c to see that it has no more than max1d points in either dimensions.
20 c needed so that scratch array space in stepgrid not exceeded.
21 c
22 c also check for too small grids - but has never happened.
23 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
24 c
25  mptr = mptr1
26  level = node(nestlevel,mptr)
27  hx = hxposs(level)
28  hy = hyposs(level)
29 c
30 10 continue
31  cxlo = rnode(cornxlo,mptr)
32  cxhi = rnode(cornxhi,mptr)
33  cylo = rnode(cornylo,mptr)
34  cyhi = rnode(cornyhi,mptr)
35  nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
36  ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
37  minsize = 2*nghost
38 c
39 c check number of rows first - if too many, bisect grid with vertical
40 c line down the middle. make sure new grid corners are anchored
41 c on coarse grid point. make sure if bisecting coarse grid that
42 c new grids have even number of points
43 c
44  if (nx + 2*nghost .gt. max1d) then
45 
46  nxl = nx/2
47  if (level .gt. 1) then
48  lratio = intratx(level-1)
49  else
50  lratio = 2
51  endif
52  nxl = (nxl/lratio)*lratio
53  nxr = nx - nxl
54  cxmid = cxlo + nxl*hx
55 
56  mptrnx = nodget()
57  node(levelptr,mptrnx) = node(levelptr,mptr)
58  node(levelptr,mptr) = mptrnx
59 
60  rnode(cornxhi,mptr) = cxmid
61  node(ndihi,mptrnx) = node(ndihi,mptr)
62  node(ndihi,mptr) = node(ndilo,mptr) + nxl - 1
63  node(ndilo,mptrnx) = node(ndihi,mptr) + 1
64  node(ndjhi,mptrnx) = node(ndjhi,mptr)
65  node(ndjlo,mptrnx) = node(ndjlo,mptr)
66 
67  rnode(cornxlo,mptrnx) = cxmid
68  rnode(cornylo,mptrnx) = cylo
69  rnode(cornyhi,mptrnx) = cyhi
70  rnode(cornxhi,mptrnx) = cxhi
71  rnode(timemult,mptrnx) = rnode(timemult,mptr)
72  node(nestlevel,mptrnx) = node(nestlevel,mptr)
73 
74  go to 10
75 c
76 c check number of columns next - if too many, bisect grid with horizontal
77 c line down the middle
78 c
79  else if (ny + 2*nghost .gt. max1d) then
80 
81  nyl = ny/2
82  if (level .gt. 1) then
83  lratio = intraty(level-1)
84  else
85  lratio = 2
86  endif
87  nyl = (nyl/lratio)*lratio
88  nyr = ny - nyl
89  cymid = cylo + nyl*hy
90 
91  mptrnx = nodget()
92  node(levelptr,mptrnx) = node(levelptr,mptr)
93  node(levelptr,mptr) = mptrnx
94 
95  rnode(cornyhi,mptr) = cymid
96 
97  node(ndjhi,mptrnx) = node(ndjhi,mptr)
98  node(ndjhi,mptr) = node(ndjlo,mptr) + nyl - 1
99  node(ndjlo,mptrnx) = node(ndjhi,mptr) + 1
100  node(ndihi,mptrnx) = node(ndihi,mptr)
101  node(ndilo,mptrnx) = node(ndilo,mptr)
102 
103  rnode(cornxlo,mptrnx) = cxlo
104  rnode(cornylo,mptrnx) = cymid
105  rnode(cornyhi,mptrnx) = cyhi
106  rnode(cornxhi,mptrnx) = cxhi
107  node(nestlevel,mptrnx) = node(nestlevel,mptr)
108  rnode(timemult,mptrnx) = rnode(timemult,mptr)
109  go to 10
110 c
111 c grid ok - check the next
112 c
113  else
114  mptr = node(levelptr,mptr)
115  if (mptr.ne.0) go to 10
116 
117  endif
118 
119  return
integer, parameter timemult
current simulation time on this grid
Definition: amr_module.f90:151
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
integer, parameter max1d
Definition: amr_module.f90:181
real(kind=8), dimension(maxlv) hyposs
Definition: amr_module.f90:193
real(kind=8), dimension(maxlv) hxposs
Definition: amr_module.f90:193
integer, parameter ndihi
global i index of right border of this grid
Definition: amr_module.f90:111
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, parameter nestlevel
AMR level of the grid.
Definition: amr_module.f90:44
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
integer, parameter ndilo
global i index of left border of this grid
Definition: amr_module.f90:108
integer, parameter ndjlo
global j index of lower border of this grid
Definition: amr_module.f90:114
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
integer function nodget()
Get first free node of the linked list kept in node array.
Definition: nodget.f:11
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
integer, parameter cornylo
y-coordinate of the lower border of this grid
Definition: amr_module.f90:145
integer, dimension(maxlv) intratx
Definition: amr_module.f90:198
integer, parameter levelptr
node number (index) of next grid on the same level
Definition: amr_module.f90:35
integer, parameter cornxhi
x-coordinate of the right border of this grid
Definition: amr_module.f90:147
integer nghost
Definition: amr_module.f90:232
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
integer, parameter cornyhi
y-coordinate of the upper border of this grid
Definition: amr_module.f90:149
Here is the call graph for this function:
Here is the caller graph for this function: