10 subroutine filval(val, mitot, mjtot, dx, dy, level, time, mic, &
11 mjc, xleft, xright, ybot, ytop, nvar, mptr, ilo, ihi, &
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
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
27 real(kind=8),
intent(in out) :: val(nvar,mitot,mjtot), aux(naux,mitot,mjtot)
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)
39 logical :: sticksoutxfine, sticksoutyfine,sticksoutxcrse,sticksoutycrse
42 integer :: clock_start, clock_finish, clock_rate
43 real(kind=8) :: cpu_start, cpu_finish
47 real(kind=8) :: get_max_speed
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
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
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)))
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)
80 call
intcopy(valc,mic,mjc,nvar,iclo,ichi,jclo,jchi,level-1,1,1)
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, &
89 call
icall(valc,auxc,mic,mjc,nvar,naux,iclo,ichi,jclo,jchi,level-1,1,1)
102 call
setaux(ng,mic,mjc,xl,yb,dx_coarse,dy_coarse,naux,auxc)
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)
116 nx = mitot - 2*nghost
117 ny = mjtot - 2*nghost
119 if (naux .gt. 0)
then
125 aux(1,:,:) = needs_to_be_set
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)
131 call
icall(val,aux,mitot,mjtot,nvar,naux,ilo-nghost,ihi+nghost, &
132 jlo-nghost,jhi+nghost,level,1,1)
137 call
setaux(nghost,nx,ny,xleft,ybot,dx,dy,naux,aux)
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)
148 call
intcopy(val,mitot,mjtot,nvar,ilo-nghost,ihi+nghost, &
149 jlo-nghost,jhi+nghost,level,1,1)
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))
166 if ( s1m*s1p <= 0.d0) slopex = 0.d0
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
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
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
183 if (
setflags(ifine,jfine) .eq. needs_to_be_set)
then
184 val(ivar,ifine,jfine) = valc(ivar,i,j) + xoff*slopex + yoff*slopey
196 if (mcapa .ne. 0)
then
197 call
fixcapaq(val,aux,mitot,mjtot,valc,auxc,mic,mjc,nvar,naux,level-1,
setflags)
227 real(kind=8) :: aux(naux,mitot,mjtot)
228 integer :: naux,mitot,mjtot,i,j,iaux
232 write(*,444) i,j,(aux(iaux,i,j),iaux=1,naux)
233 444
format(2i4,5e12.5)
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...
subroutine filval(val, mitot, mjtot, dx, dy, level, time, mic, mjc, xleft, xright, ybot, ytop, nvar, mptr, ilo, ihi, jlo, jhi, aux, naux)
subroutine setflags(iflags, isize, jsize, rctold, idim3, mitot, mjtot, mptr)
subroutine icall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, iputst, jputst)
subroutine intcopy(val, mitot, mjtot, nvar, ilo, ihi, jlo, jhi, level, iputst, jputst)
subroutine fixcapaq(val, aux, mitot, mjtot, valc, auxc, mic, mjc, nvar, naux, levc, setflags)
subroutine preintcopy(val, mitot, mjtot, nvar, ilo, ihi, jlo, jhi, level, fliparray)
subroutine dumpaux(aux, naux, mitot, mjtot)
subroutine preicall(val, aux, nrow, ncol, nvar, naux, ilo, ihi, jlo, jhi, level, fliparray)
subroutine setaux(mbc, mx, my, xlower, ylower, dx, dy, maux, aux)