98 subroutine intfil(val,mi,mj,time,flaguse,nrowst,ncolst,ilo,ihi,jlo,jhi,level,nvar,naux,msrc)
109 integer,
intent(in) :: mi, mj, nrowst, ncolst, ilo, ihi, jlo, jhi, level, nvar, naux,msrc
110 real(kind=8),
intent(in) :: time
113 integer(kind=1),
intent(in out) :: flaguse(ilo:ihi, jlo:jhi)
114 real(kind=8),
intent(in out) :: val(nvar,mi,mj)
117 integer :: imlo, jmlo, imhi, jmhi, nx, ny, mitot, mjtot
118 integer :: ixlo, ixhi, jxlo, jxhi, locold, locnew, nextSpot
119 integer :: icount, bndNum, bndLoc, levSt
120 integer :: ivar, i, j, mptr, mstart, loc, numg
121 real(kind=8) :: dt, alphac, alphai
122 logical :: t_interpolate
124 integer :: patch_rect(4)
126 real(kind=8),
parameter :: t_epsilon = 1.0d-4
129 character(len=*),
parameter :: missing_error = &
130 "(' time wanted ',e15.7,' not available from grid ',i4,'level',i4)" 131 character(len=*),
parameter :: time_error = &
132 "(' trying to interpolate from previous time values ',/," // &
133 "' for a patch with corners ilo,ihi,jlo,jhi:'" // &
135 "' from source grid ',i4,' at level ',i4,/," // &
136 "' time wanted ',e24.16,' source time is ',e24.16,/," // &
137 "' alphai, t_epsilon ',2e24.16)" 139 patch_rect = [ilo,ihi,jlo,jhi]
153 if (msrc .eq. -1)
then 157 bndloc = node(bndlistst,msrc)
158 bndnum = node(bndlistnum,msrc)
159 nextspot = node(bndlistst, msrc)
165 if (msrc .eq. -1)
then 168 mptr =
bndlist(nextspot,gridnbor)
172 imlo = node(ndilo, mptr)
173 jmlo = node(ndjlo, mptr)
174 imhi = node(ndihi, mptr)
175 jmhi = node(ndjhi, mptr)
177 nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
178 ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
180 mitot = nx + 2 * nghost
181 mjtot = ny + 2 * nghost
183 ixlo = max(imlo,patch_rect(1))
184 ixhi = min(imhi,patch_rect(2))
185 jxlo = max(jmlo,patch_rect(3))
186 jxhi = min(jmhi,patch_rect(4))
190 if (ixlo <= ixhi .and. jxlo <= jxhi)
then 194 alphac = (rnode(timemult,mptr) - time)/dt
195 alphai = 1.d0 - alphac
197 if ((alphai < -t_epsilon) .or. (alphai > 1.d0 + t_epsilon))
then 198 write(outunit,missing_error) time, mptr, level
199 print missing_error, time, mptr, level
200 write(outunit,
'(A,E24.16)')
'Line 80', dt
201 write(outunit,time_error) patch_rect,mptr,level,time,rnode(timemult,mptr),alphai,t_epsilon
202 print time_error, patch_rect,mptr,level,time,rnode(timemult,mptr),alphai,t_epsilon
203 call outtre(mstart,.false.,nvar,naux)
208 t_interpolate = .false.
209 if (abs(alphai - 1.d0) < t_epsilon)
then 211 loc = node(store1,mptr)
212 else if (dabs(alphai) .lt. t_epsilon)
then 213 loc = node(store2,mptr)
214 if (level == mxnest)
then 215 write(outunit,
'(A,E24.16)')
'Line 95', dt
216 write(outunit,time_error) patch_rect,mptr,level,time, &
217 rnode(timemult,mptr),alphai,t_epsilon
218 write(*,time_error) patch_rect,mptr,level,time, &
219 rnode(timemult,mptr),alphai,t_epsilon
223 locold = node(store2,mptr)
224 locnew = node(store1,mptr)
225 t_interpolate = .true.
228 if (level == mxnest)
then 229 write(outunit,
'(A,E24.16)')
'Line 107',dt
230 write(outunit,time_error) patch_rect,mptr,level,time,rnode(timemult,mptr),alphai,t_epsilon
231 print time_error, patch_rect,mptr,level,time,rnode(timemult,mptr),alphai,t_epsilon
237 if (.not. t_interpolate)
then 242 val(ivar,i-patch_rect(1)+nrowst,j-jlo+ncolst) = &
243 alloc(
iadd(ivar,i-imlo+nghost+1,j-jmlo+nghost+1))
253 val(ivar,i-patch_rect(1)+nrowst,j-jlo+ncolst) = &
254 alloc(
iadnew(ivar,i-imlo+nghost+1,j-jmlo+nghost+1))*alphai + &
255 alloc(
iadold(ivar,i-imlo+nghost+1,j-jmlo+nghost+1))*alphac
266 if (msrc .ne. -1)
then 267 nextspot =
bndlist(nextspot,nextfree)
275 if (jhi >= jregsz(level))
then 276 flaguse(patch_rect(1):ihi, max(jregsz(level),jlo):jhi) = 1
280 flaguse(patch_rect(1):ihi, jlo:min(-1,ncolst + jhi - jlo )) = 1
283 if (patch_rect(1) < 0)
then 284 flaguse(patch_rect(1):min(-1,nrowst + ihi - patch_rect(1)), jlo:jhi) = 1
287 if (ihi >= iregsz(level))
then 288 flaguse(max(iregsz(level),patch_rect(1)):ihi, jlo:jhi) = 1
293 integer pure function iadd(ivar,i,j)
295 integer,
intent(in) :: ivar, i, j
296 iadd = loc + ivar-1 + nvar*((j-1)*mitot+i-1)
299 integer pure function iadnew(ivar,i,j)
301 integer,
intent(in) :: ivar, i, j
302 iadnew = locnew + ivar-1 + nvar*((j-1)*mitot+i-1)
305 integer pure function iadold(ivar,i,j)
307 integer,
intent(in) :: ivar, i, j
308 iadold = locold + ivar-1 + nvar*((j-1)*mitot+i-1)
integer, parameter timemult
current simulation time on this grid
integer, parameter bndlistnum
number of grids (on the same level) that border this grid
integer, dimension(maxlv) iregsz
integer, parameter ndihi
global i index of right border of this grid
integer, dimension(maxlv) jregsz
integer, dimension(nsize, maxgr) node
subroutine intfil(val, mi, mj, time, flaguse, nrowst, ncolst, ilo, ihi, jlo, jhi, level, nvar, naux, msrc)
Fill values for a patch at the specified level and location, using values from grids at level level O...
integer, dimension(0:maxlv+1) liststart
integer, dimension(maxlv) numgrids
real(kind=8), dimension(rsize, maxgr) rnode
integer pure function iadnew(ivar, i, j)
integer, dimension(maxgr) listofgrids
integer pure function iadd(ivar, i, j)
integer, parameter gridnbor
integer, parameter bndlistst
pointer (actually it's an index in the bndList array) to the first node of a linked list...
integer, parameter ndilo
global i index of left border of this grid
integer, parameter ndjlo
global j index of lower border of this grid
integer pure function iadold(ivar, i, j)
integer, parameter store1
pointer to the address of memory storing the first copy of solution data on this grid, usually for storing new solution
integer, dimension(maxlv) lstart
integer, parameter store2
pointer to the address of memory storing the second copy of solution data on this grid...
integer, dimension(bndlistsize, 2) bndlist
integer, parameter outunit
real(kind=8), dimension(maxlv) possk
integer, parameter nextfree
integer, parameter ndjhi
global j index of upper border of this grid
subroutine outtre(mlev, outgrd, nvar, naux)
Output a subtree of the grids.
integer, parameter levelptr
node number (index) of next grid on the same level
The module contains the definition of a "node descriptor" as well as other global variables used duri...
real(kind=8), dimension(:), allocatable alloc