4 subroutine advanc (level,nvar,dtlevnew,vtime,naux)
7 implicit double precision (a-h,o-z)
11 integer omp_get_thread_num, omp_get_max_threads
12 integer mythread/0/, maxthreads/1/
13 integer listgrids(numgrids(level))
14 integer clock_start, clock_finish, clock_rate
15 integer clock_startstepgrid,clock_startbound,clock_finishbound
16 real(kind=8) cpu_start, cpu_finish
17 real(kind=8) cpu_startbound, cpu_finishbound
18 real(kind=8) cpu_startstepgrid, cpu_finishstepgrid
35 call system_clock(clock_start,clock_rate)
36 call cpu_time(cpu_start)
46 call system_clock(clock_startbound,clock_rate)
47 call cpu_time(cpu_startbound)
61 do j = 1, numgrids(level)
63 levst = liststart(level)
64 mptr = listofgrids(levst+j-1)
66 nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
67 ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
70 locnew = node(store1,mptr)
71 locaux = node(storeaux,mptr)
72 time = rnode(timemult,mptr)
74 call
bound(time,nvar,nghost,alloc(locnew),mitot,mjtot,mptr,
79 call system_clock(clock_finishbound,clock_rate)
80 call cpu_time(cpu_finishbound)
81 timebound = timebound + clock_finishbound - clock_startbound
82 timeboundcpu=timeboundcpu+cpu_finishbound-cpu_startbound
86 if (level+1 .le. mxnest)
then
87 if (lstart(level+1) .ne. null)
then
88 call
saveqc(level+1,nvar,naux)
96 call system_clock(clock_startstepgrid,clock_rate)
97 call cpu_time(cpu_startstepgrid)
108 do j = 1, numgrids(level)
110 levst = liststart(level)
111 mptr = listofgrids(levst+j-1)
112 nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
113 ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
114 mitot = nx + 2*nghost
115 mjtot = ny + 2*nghost
117 call
par_advanc(mptr,mitot,mjtot,nvar,naux,dtnew)
119 dtlevnew = dmin1(dtlevnew,dtnew)
125 call system_clock(clock_finish,clock_rate)
126 call cpu_time(cpu_finish)
127 tvoll(level) = tvoll(level) + clock_finish - clock_start
128 tvollcpu(level) = tvollcpu(level) + cpu_finish - cpu_start
129 timestepgrid = timestepgrid +clock_finish-clock_startstepgrid
130 timestepgridcpu=timestepgridcpu+cpu_finish-cpu_startstepgrid
132 cflmax = dmax1(cflmax, cfl_level)
143 implicit double precision (a-h,o-z)
144 integer listgrids(num)
149 mptr = node(levelptr, mptr)
152 if (mptr .ne. 0)
then
153 write(*,*)
" Error in routine setting up grid array "
167 implicit double precision (a-h,o-z)
170 integer omp_get_thread_num, omp_get_max_threads
171 integer mythread/0/, maxthreads/1/
173 double precision fp(nvar,mitot,mjtot),fm(nvar,mitot,mjtot)
174 double precision gp(nvar,mitot,mjtot),gm(nvar,mitot,mjtot)
186 level = node(nestlevel,mptr)
190 nx = node(ndihi,mptr) - node(ndilo,mptr) + 1
191 ny = node(ndjhi,mptr) - node(ndjlo,mptr) + 1
192 time = rnode(timemult,mptr)
196 locold = node(store2, mptr)
197 locnew = node(store1, mptr)
205 if (level .lt. mxnest)
then
206 ntot = mitot * mjtot * nvar
209 10 alloc(locold + i - 1) = alloc(locnew + i - 1)
212 xlow = rnode(cornxlo,mptr) - nghost*hx
213 ylow = rnode(cornylo,mptr) - nghost*hy
216 rvol = rvol + nx * ny
217 rvoll(level) = rvoll(level) + nx * ny
221 locaux = node(storeaux,mptr)
223 if (node(ffluxptr,mptr) .ne. 0)
then
224 lenbc = 2*(nx/intratx(level-1)+ny/intraty(level-1))
225 locsvf = node(ffluxptr,mptr)
226 locsvq = locsvf + nvar*lenbc
227 locx1d = locsvq + nvar*lenbc
228 call
qad(alloc(locnew),mitot,mjtot,nvar,
229 1 alloc(locsvf),alloc(locsvq),lenbc,
230 2 intratx(level-1),intraty(level-1),hx,hy,
231 3 naux,alloc(locaux),alloc(locx1d),delt,mptr)
242 if (num_gauges > 0)
then
244 . alloc(locaux:locaux+nvar*mitot*mjtot),
245 . xlow,ylow,nvar,mitot,mjtot,naux,mptr)
249 if (dimensional_split .eq. 0)
then
251 call
stepgrid(alloc(locnew),fm,fp,gm,gp,
252 2 mitot,mjtot,nghost,
253 3 delt,dtnew,hx,hy,nvar,
254 4 xlow,ylow,time,mptr,naux,alloc(locaux))
255 else if (dimensional_split .eq. 1)
then
258 2 mitot,mjtot,nghost,
259 3 delt,dtnew,hx,hy,nvar,
260 4 xlow,ylow,time,mptr,naux,alloc(locaux))
263 write(6,*)
'*** Strang splitting not supported'
267 if (node(cfluxptr,mptr) .ne. 0)
268 2 call
fluxsv(mptr,fm,fp,gm,gp,
269 3 alloc(node(cfluxptr,mptr)),mitot,mjtot,
270 4 nvar,listsp(level),delt,hx,hy)
271 if (node(ffluxptr,mptr) .ne. 0)
then
272 lenbc = 2*(nx/intratx(level-1)+ny/intraty(level-1))
273 locsvf = node(ffluxptr,mptr)
275 2 alloc(locsvf),mptr,mitot,mjtot,nvar,
276 4 lenbc,intratx(level-1),intraty(level-1),
282 rnode(timemult,mptr) = rnode(timemult,mptr)+delt
subroutine saveqc(level, nvar, naux)
subroutine fluxad(xfluxm, xfluxp, yfluxm, yfluxp, svdflx, mptr, mitot, mjtot, nvar, lenbc, lratiox, lratioy, ng, dtf, dx, dy)
subroutine qad(valbig, mitot, mjtot, nvar, svdflx, qc1d, lenbc, lratiox, lratioy, hx, hy, maux, aux, auxc1d, delt, mptr)
subroutine update_gauges(q, aux, xlow, ylow, num_eqn, mitot, mjtot, num_aux, mptr)
subroutine advanc(level, nvar, dtlevnew, vtime, naux)
subroutine prepgrids(listgrids, num, level)
subroutine stepgrid(q, fm, fp, gm, gp, mitot, mjtot, mbc, dt, dtnew, dx, dy, nvar, xlow, ylow, time, mptr, maux, aux)
subroutine bound(time, nvar, ng, valbig, mitot, mjtot, mptr, aux, naux)
This routine sets the boundary values for a given grid at level level.
subroutine par_advanc(mptr, mitot, mjtot, nvar, naux, dtnew)
subroutine fluxsv(mptr, xfluxm, xfluxp, yfluxm, yfluxp, listbc, ndimx, ndimy, nvar, maxsp, dtc, hx, hy)
subroutine stepgrid_dimsplit(q, fm, fp, gm, gp, mitot, mjtot, mbc, dt, dtnew, dx, dy, nvar, xlow, ylow, time, mptr, maux, aux)