2D AMRCLAW
filval.f90
Go to the documentation of this file.
1 !
2 ! :::::::::::::::::::::::::::::: FILVAL ::::::::::::::::::::::::::
3 !
4 !
5 ! ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
6 !
7 !
8 ! ------------------------------------------------------------------
9 !
10 subroutine filval(val, mitot, mjtot, dx, dy, level, time, mic, &
11  mjc, xleft, xright, ybot, ytop, nvar, mptr, ilo, ihi, &
12  jlo, jhi, aux, naux)
13 
14  use amr_module, only: xlower, ylower, intratx, intraty, nghost, xperdom
15  use amr_module, only: yperdom, spheredom, xupper, yupper, alloc
16  use amr_module, only: outunit, needs_to_be_set, mcapa
17  use amr_module, only: newstl, iregsz, jregsz
18 
19  implicit none
20 
21  ! Input
22  integer, intent(in) :: mitot, mjtot, level, mic, mjc, nvar, mptr, ilo, ihi
23  integer, intent(in) :: jlo, jhi, naux
24  real(kind=8), intent(in) :: dx, dy, time, xleft, xright, ybot, ytop
25 
26  ! Output
27  real(kind=8), intent(in out) :: val(nvar,mitot,mjtot), aux(naux,mitot,mjtot)
28 
29  ! Local storage
30  integer :: refinement_ratio_x, refinement_ratio_y, iclo, jclo, ichi, jchi, ng
31  integer :: ivar, i, j, ico, jco, ifine, jfine, nx, ny
32  real(kind=8) :: valc(nvar,mic,mjc), auxc(naux,mic,mjc)
33  real(kind=8) :: dx_coarse, dy_coarse, xl, xr, yb, yt, area
34  real(kind=8) :: s1m, s1p, slopex, slopey, xoff, yoff
35  real(kind=8) :: fliparray((mitot+mjtot)*nghost*(nvar+naux))
36  real(kind=8) :: setflags(mitot,mjtot),maxauxdif,aux2(naux,mitot,mjtot)
37  integer :: mjb
38  integer :: jm, im, nm
39  logical :: sticksoutxfine, sticksoutyfine,sticksoutxcrse,sticksoutycrse
40 
41  !for setaux timing
42  integer :: clock_start, clock_finish, clock_rate
43  real(kind=8) :: cpu_start, cpu_finish
44 
45 
46  ! External function definitions
47  real(kind=8) :: get_max_speed
48 
49 
50  refinement_ratio_x = intratx(level-1)
51  refinement_ratio_y = intraty(level-1)
52  dx_coarse = dx * refinement_ratio_x
53  dy_coarse = dy * refinement_ratio_y
54  xl = xleft - dx_coarse
55  xr = xright + dx_coarse
56  yb = ybot - dy_coarse
57  yt = ytop + dy_coarse
58 
59  ! if topo not yet final then aux is set outside filval (in gfixup)
60  ! and so aux has real data already, (ie dont overwrite here)
61 
62  ! set integer indices for coarser patch enlarged by 1 cell
63  ! (can stick out of domain). proper nesting will insure this one
64  ! call is sufficient.
65  iclo = ilo / refinement_ratio_x - 1
66  jclo = jlo / refinement_ratio_y - 1
67  ichi = (ihi + 1) / refinement_ratio_x - 1 + 1
68  jchi = (jhi + 1) / refinement_ratio_y - 1 + 1
69  ng = 0
70 
71  sticksoutxfine = ( (ilo .lt. 0) .or. (ihi .ge. iregsz(level)))
72  sticksoutyfine = ( (jlo .lt. 0) .or. (jhi .ge. jregsz(level)))
73  sticksoutxcrse = ((iclo .lt. 0) .or. (ichi .ge. iregsz(level-1)))
74  sticksoutycrse = ((jclo .lt. 0) .or. (jchi .ge. jregsz(level-1)))
75 
76  if (naux == 0) then
77  if ((xperdom .and. sticksoutxcrse).or. (yperdom.and.sticksoutycrse) .or. spheredom) then
78  call preintcopy(valc,mic,mjc,nvar,iclo,ichi,jclo,jchi,level-1,fliparray)
79  else
80  call intcopy(valc,mic,mjc,nvar,iclo,ichi,jclo,jchi,level-1,1,1)
81  endif
82  else
83  ! intersect grids and copy all (soln and aux)
84  auxc(1,:,:) = needs_to_be_set
85  if ((xperdom .and.sticksoutxcrse) .or. (yperdom.and.sticksoutycrse) .or. spheredom) then
86  call preicall(valc,auxc,mic,mjc,nvar,naux,iclo,ichi,jclo,jchi, &
87  level-1,fliparray)
88  else
89  call icall(valc,auxc,mic,mjc,nvar,naux,iclo,ichi,jclo,jchi,level-1,1,1)
90  endif
91 !!$ do i = 1, mic
92 !!$ do j = 1, mjc
93 !!$ if (auxc(1,:,:) == NEEDS_TO_BE_SET) then
94 !!$ write(*,*)" *** coarsenened new fine grid not completely set from previously" &
95 !!$ "existing coarse grids ***"
96 !!$ stop
97 !!$ endif
98 !!$ end do
99 !!$ end do
100  ! no ghost cells on coarse enlarged patch. set any remaining
101  ! vals. should only be bcs that stick out of domain.
102  call setaux(ng,mic,mjc,xl,yb,dx_coarse,dy_coarse,naux,auxc)
103  endif
104 
105  call bc2amr(valc,auxc,mic,mjc,nvar,naux,dx_coarse,dy_coarse,level-1,time,xl,xr,yb, &
106  yt,xlower,ylower,xupper,yupper,xperdom,yperdom,spheredom)
107 
108 
109 
110 ! NOTE change in order of code. Since the interp from coarse to fine needs the aux
111 ! arrays set already, the fine copy is done first, to set up the aux arrays.
112 ! we can do this since we have the flag array to test where to overwrite.
113 
114 ! SO this is no longer overwriting but setting for the first time.
115 ! overwrite interpolated values with fine grid values, if available.
116  nx = mitot - 2*nghost
117  ny = mjtot - 2*nghost
118 
119  if (naux .gt. 0) then
120 ! ## NEEDS_TO_BE_SET is signal that aux array not set.
121 ! ## after calling icall to copy aux from other grids
122 ! ## any remaining NEEDS_TO_BE_SET signals will be set in setaux.
123 ! ## it also signals where soln was copied, so it wont be
124 ! ## overwritten with coarse grid interpolation
125  aux(1,:,:) = needs_to_be_set ! indicates fine cells not yet set.
126 
127  if ((xperdom.and.sticksoutxfine) .or. (yperdom.and.sticksoutyfine) .or. spheredom) then
128  call preicall(val,aux,mitot,mjtot,nvar,naux,ilo-nghost,ihi+nghost, &
129  jlo-nghost,jhi+nghost,level,fliparray)
130  else
131  call icall(val,aux,mitot,mjtot,nvar,naux,ilo-nghost,ihi+nghost, &
132  jlo-nghost,jhi+nghost,level,1,1)
133  endif
134  setflags = aux(1,:,:) ! save since will overwrite in setaux when setting all aux vals
135  ! need this so we know where to use coarse grid to set fine solution w/o overwriting
136  ! set remaining aux vals not set by copying from prev existing grids
137  call setaux(nghost,nx,ny,xleft,ybot,dx,dy,naux,aux)
138  else ! either no aux exists, or cant reuse yet
139  ! so only call intcopy (which copies soln) and not icall.
140  ! in this case flag q(1,:) to NEEDS_TO_BE_SET flag so wont be overwritten
141  ! by coarse grid interp. this is needed due to reversing order of
142  ! work - first copy from fine grids, then interpolate from coarse grids
143  val(1,:,:) = needs_to_be_set
144  if ((xperdom.and.sticksoutxfine) .or. (yperdom.and.sticksoutyfine) .or. spheredom) then
145  call preintcopy(val,mitot,mjtot,nvar,ilo-nghost,ihi+nghost, &
146  jlo-nghost,jhi+nghost,level,fliparray)
147  else
148  call intcopy(val,mitot,mjtot,nvar,ilo-nghost,ihi+nghost, &
149  jlo-nghost,jhi+nghost,level,1,1)
150  endif
151  setflags = val(1,:,:) ! remaining flags signals need to set
152  endif
153 
154 
155  ! Prepare slopes - use min-mod limiters
156 
157  do j=2, mjc-1
158  do i=2, mic-1
159  do ivar = 1, nvar
160 
161  s1p = valc(ivar,i+1,j) - valc(ivar,i,j)
162  s1m = valc(ivar,i,j) - valc(ivar,i-1,j)
163  slopex = min(abs(s1p), abs(s1m)) &
164  * sign(1.d0,valc(ivar,i+1,j) - valc(ivar,i-1,j))
165  ! if there's a sign change, set slope to 0.
166  if ( s1m*s1p <= 0.d0) slopex = 0.d0
167 
168  s1p = valc(ivar,i,j+1) - valc(ivar,i,j)
169  s1m = valc(ivar,i,j) - valc(ivar,i,j-1)
170  slopey = min(abs(s1p), abs(s1m)) &
171  * sign(1.0d0, valc(ivar,i,j+1) - valc(ivar,i,j-1))
172  if ( s1m*s1p <= 0.d0) slopey = 0.d0
173 
174  ! Interpolate from coarse cells to fine grid to find depth
175  do jco = 1,refinement_ratio_y
176  yoff = (real(jco,kind=8) - 0.5d0) / refinement_ratio_y - 0.5d0
177  jfine = (j-2) * refinement_ratio_y + nghost + jco
178 
179  do ico = 1,refinement_ratio_x
180  xoff = (real(ico,kind=8) - 0.5d0) / refinement_ratio_x - 0.5d0
181  ifine = (i-2) * refinement_ratio_x + nghost + ico
182 
183  if (setflags(ifine,jfine) .eq. needs_to_be_set) then
184  val(ivar,ifine,jfine) = valc(ivar,i,j) + xoff*slopex + yoff*slopey
185  endif
186 
187  end do
188  end do
189 
190  enddo !end of ivar loop
191  enddo !end of coarse i loop
192  enddo !end of coarse j loop
193 
194  ! adjust to conserve kappa*q, but only where coarse grid was interpolated
195  ! so now need to pass setflags to this subr.
196  if (mcapa .ne. 0) then
197  call fixcapaq(val,aux,mitot,mjtot,valc,auxc,mic,mjc,nvar,naux,level-1,setflags)
198  endif
199 
200 !!$! CHECK BY CALLING SETAUX AND SETTING ALL, THEN DIFFING
201 !!$ mjb = 0
202 !!$ if (naux .gt. 0 .and. mjb .eq. 1) then
203 !!$ aux2(1,:,:) = NEEDS_TO_BE_SET ! indicates fine cells not yet set
204 !!$ call setaux(nghost,nx,ny,xleft,ybot,dx,dy,naux,aux2)
205 !!$ maxauxdif = 1.d-13
206 !!$ do i = 1, mitot
207 !!$ do j = 1, mjtot
208 !!$ if (abs(aux(1,i,j)-aux2(1,i,j)) .gt. maxauxdif) then
209 !!$ maxauxdif = abs(aux(1,i,j)-aux2(1,i,j))
210 !!$ write(*,444)i,j,aux(1,i,j),aux2(1,i,j),maxauxdif
211 !!$444 format("i,j = ",2i4," auxs ",2e15.7," maxauxdif ",e12.5)
212 !!$ endif
213 !!$ end do
214 !!$ end do
215 !!$ if (maxauxdif .gt. 2.d-13) then
216 !!$ write(*,*)" maxauxdif = ",maxauxdif," with mitot,mjtot ",mitot,mjtot, &
217 !!$ " on grid ",mptr," level ",level
218 !!$ endif
219 !!$ endif
220 
221 end subroutine filval
222 
223 ! ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
224 
225 subroutine dumpaux(aux,naux,mitot,mjtot)
226  implicit none
227  real(kind=8) :: aux(naux,mitot,mjtot)
228  integer :: naux,mitot,mjtot,i,j,iaux
229 
230  do j = 1, mjtot
231  do i = 1, mitot
232  write(*,444) i,j,(aux(iaux,i,j),iaux=1,naux)
233  444 format(2i4,5e12.5)
234  end do
235  end do
236 
237 end subroutine dumpaux
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:87
subroutine filval(val, mitot, mjtot, dx, dy, level, time, mic, mjc, xleft, xright, ybot, ytop, nvar, mptr, ilo, ihi, jlo, jhi, aux, naux)
Definition: filval.f90:10
subroutine setflags(iflags, isize, jsize, rctold, idim3, mitot, mjtot, mptr)
Definition: setflags.f:4
subroutine icall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, iputst, jputst)
Definition: icall.f:4
subroutine intcopy(val, mitot, mjtot, nvar, ilo, ihi, jlo, jhi, level, iputst, jputst)
Definition: intcopy.f:4
subroutine fixcapaq(val, aux, mitot, mjtot, valc, auxc, mic, mjc, nvar, naux, levc, setflags)
Definition: fixcapaq.f:4
subroutine preintcopy(val, mitot, mjtot, nvar, ilo, ihi, jlo, jhi, level, fliparray)
Definition: preintcopy.f:4
subroutine dumpaux(aux, naux, mitot, mjtot)
Definition: filval.f90:225
subroutine preicall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, fliparray)
Definition: preicall.f:4
subroutine setaux(mbc, mx, my, xlower, ylower, dx, dy, maux, aux)
Definition: setaux.f90:1