2D AMRCLAW
saveqc.f
Go to the documentation of this file.
1 c
12 c ================================================================
13  subroutine saveqc(level,nvar,naux)
14 c ================================================================
15 c
16  use amr_module
17  implicit double precision (a-h,o-z)
18 
19  logical sticksout, found
20 ! make fliparray largest possible grid size
21  dimension fliparray(2*max1d*nghost*(nvar+naux))
22 c
23 c ::::::::::::::::::::::::: SAVEQC :::::::::::::::::::::::::::::::::
24 c prepare new fine grids to save fluxes after each integration step
25 c for future conservative fix-up on coarse grids.
26 c save all boundary fluxes of fine grid (even if on a phys. bndry.) -
27 c but only save space for every intrat of them.
28 c:::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
29 c
30  levc = level - 1
31  hxc = hxposs(levc)
32  hyc = hyposs(levc)
33  ng = 0 ! no ghost cells on coarsened enlarged patch
34 
35  mkid = lstart(level)
36  10 if (mkid .eq. 0) go to 99
37  nx = node(ndihi,mkid)-node(ndilo,mkid) + 1
38  ny = node(ndjhi,mkid)-node(ndjlo,mkid) + 1
39  ikeep = nx/intratx(level-1)
40  jkeep = ny/intraty(level-1)
41  lenbc = 2*(ikeep+jkeep)
42  ist = node(ffluxptr,mkid)
43  time = rnode(timemult,mkid)
44 
45 c make coarsened enlarged patch for conservative fixup
46  ilo = node(ndilo,mkid)
47  jlo = node(ndjlo,mkid)
48  ihi = node(ndihi,mkid)
49  jhi = node(ndjhi,mkid)
50  iclo = ilo/intratx(level-1) - 1
51  jclo = jlo/intraty(level-1) - 1
52  ichi = (ihi+1)/intratx(level-1)
53  jchi = (jhi+1)/intraty(level-1)
54  nrow = ichi-iclo+1
55  ncol = jchi-jclo+1
56  xl = rnode(cornxlo,mkid) - hxc
57  yb = rnode(cornylo,mkid) - hyc
58  xr = rnode(cornxhi,mkid) + hxc
59  yt = rnode(cornyhi,mkid) + hyc
60  loctmp = igetsp(nrow*ncol*(nvar+naux))
61  loctx = loctmp + nrow*ncol*nvar
62  do i = 1, nrow*ncol*naux
63  alloc(loctx+i-1) = needs_to_be_set
64  end do
65  locaux = node(storeaux,mkid)
66 
67  if (iclo .lt. 0 .or. ichi .eq. iregsz(levc) .or.
68  & jclo .lt. 0 .or. jchi .eq. jregsz(levc)) then
69  sticksout = .true.
70  else
71  sticksout = .false.
72  endif
73 
74  if (sticksout .and. (xperdom.or.yperdom.or.spheredom)) then
75  !iperim = nrow+ncol
76  !locflip = igetsp(iperim*nghost*(nvar+naux))
77  call preicall(alloc(loctmp),alloc(loctx),nrow,ncol,nvar,
78  . naux,iclo,ichi,jclo,jchi,level-1,
79  . fliparray)
80 ! . alloc(locflip))
81 ! call reclam(locflip,iperim*nghost*(nvar+naux))
82  else
83  call icall(alloc(loctmp),alloc(loctx),nrow,ncol,nvar,naux,
84  . iclo,ichi,jclo,jchi,level-1,1,1)
85  endif
86 ! in case any part sticks out of domain still need to set remaining aux
87 ! cells
88  if (naux .gt. 0 .and. sticksout) then
89  call setaux(ng,nrow,ncol,xl,yb,hxc,hyc,naux,alloc(loctx))
90  endif
91 !-- found = .false.
92 !-- do i = 1, naux*nrow*ncol, naux
93 !-- if (alloc(loctx+i-1) .eq. NEEDS_TO_BE_SET) then
94 !-- found = .true.
95 !-- endif
96 !-- end do
97 !-- if (found) write(*,*) "still have unset aux vals in qad"
98  call bc2amr(alloc(loctmp),alloc(loctx),nrow,ncol,nvar,naux,
99  . hxc,hyc,level,time,
100  . xl,xr,yb,yt)
101  call cstore(alloc(loctmp),nrow,ncol,nvar,
102  . alloc(ist+nvar*lenbc),lenbc,naux,alloc(loctx),
103  . alloc(ist+2*nvar*lenbc))
104  call reclam(loctmp,nrow*ncol*(nvar+naux))
105 
106  mkid = node(levelptr,mkid)
107  go to 10
108  99 return
109  end
subroutine bc2amr(val, aux, nrow, ncol, meqn, naux,
Take a grid patch with mesh widths hx,hy, of dimensions nrow by ncol, and set the values of any piece...
Definition: bc2amr.f:88
integer, parameter timemult
current simulation time on this grid
Definition: amr_module.f90:151
subroutine saveqc(level, nvar, naux)
Prepare new fine grids to save fluxes after each integration step for future conservative fix-up on c...
Definition: saveqc.f:14
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
function igetsp(nwords)
Allocate contiguous space of length nword in main storage array alloc.
Definition: igetsp.f:9
integer, parameter max1d
Definition: amr_module.f90:181
subroutine reclam(index, nwords)
Definition: reclam.f:5
integer, dimension(maxlv) iregsz
Definition: amr_module.f90:198
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(maxlv) jregsz
Definition: amr_module.f90:198
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
subroutine cstore(qc, nrow, ncol, nvar, qc1d, lenbc, naux, auxc, auxc1d)
Definition: cstore.f:5
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
logical yperdom
Definition: amr_module.f90:230
subroutine icall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, iputst, jputst)
For a rectangle defined on level level and bound by ilo, ihi, jlo, jhi, find intersecting grids at th...
Definition: icall.f:10
integer, parameter ndilo
global i index of left border of this grid
Definition: amr_module.f90:108
real(kind=8), parameter needs_to_be_set
Definition: amr_module.f90:168
integer, parameter ndjlo
global j index of lower border of this grid
Definition: amr_module.f90:114
logical spheredom
Definition: amr_module.f90:230
integer, dimension(maxlv) lstart
Definition: amr_module.f90:198
integer, dimension(maxlv) intraty
Definition: amr_module.f90:198
subroutine preicall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, fliparray)
Definition: preicall.f:6
integer, parameter ndjhi
global j index of upper border of this grid
Definition: amr_module.f90:117
logical xperdom
Definition: amr_module.f90:230
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
subroutine setaux(mbc, mx, my, xlower, ylower, dx, dy, maux, aux)
Definition: setaux.f90:2
integer, parameter cornyhi
y-coordinate of the upper border of this grid
Definition: amr_module.f90:149
integer, parameter storeaux
pointer to the address of memory storing auxiliary data on this grid
Definition: amr_module.f90:120
real(kind=8), dimension(:), allocatable alloc
Definition: amr_module.f90:218
integer, parameter ffluxptr
pointer to the address of memory storing fluxes in a layer around the grid, to be used in conservatio...
Definition: amr_module.f90:97