4 subroutine setgrd (nvar,cut,naux,dtinit,start_time)
7 implicit double precision (a-h,o-z)
24 if (
mxnest .eq. 1)
go to 99
30 10
if (levnew .gt.
mxnest)
go to 30
32 if (
lstart(levold) .eq. 0)
go to 30
42 call advanc(levold,nvar,dtlev,vtime,naux)
49 dtinit = min(dtinit, dtlev*kfac)
53 20
rvoll(level) = 0.d0
60 call grdfit(lbase,levold,nvar,naux,cut,time,start_time)
61 if (
newstl(levnew) .ne. 0) lfnew = levnew
67 call ginit(
newstl(levnew),.true., nvar, naux,start_time)
70 call ginit(
lstart(levold),.false., nvar, naux,start_time)
79 do while (mptr .gt. 0)
91 write(*,100) ngrids,ncells,levnew
92 100
format(
"there are ",i6,
" grids with ",i10,
93 &
" cells at level ", i3)
111 40
if ((lev .eq.
mxnest) .or. (lev .gt.
lfine))
go to 60
118 if (mptr .ne. 0)
go to 50
126 do 65 level = 1,
lfine-1
132 do 70 level = 1,
lfine-1
133 call prepf(level+1,nvar,naux)
134 call prepc(level,nvar)
subroutine makebndrylist(level)
Preprocess each grid on level level to have a linked list of other grids at the same level that suppl...
integer, dimension(7) method
integer, parameter timemult
current simulation time on this grid
integer, dimension(maxlv) kratio
integer, dimension(maxlv) numcells
subroutine prepc(level, nvar)
This routine is called because regridding just changed the fine grids.
subroutine makegridlist(lbase)
Make (modify) array of grid numbers, listOfGrids, (after sorting them in the linked list so they are ...
subroutine ginit(msave, first, nvar, naux, start_time)
Initializes solution on all grids at level of grid msave, by calling qinit().
integer, dimension(maxlv) newstl
integer, parameter ndihi
global i index of right border of this grid
integer, dimension(nsize, maxgr) node
integer, dimension(maxlv) numgrids
real(kind=8), dimension(rsize, maxgr) rnode
integer, dimension(maxlv) iregridcount
real(kind=8), dimension(maxlv) avenumgrids
subroutine update(level, nvar, naux)
Synchronize between all grids on level level and grids on level level+1.
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, parameter store1
pointer to the address of memory storing the first copy of solution data on this grid, usually for storing new solution
subroutine advanc(level, nvar, dtlevnew, vtime, naux)
Integrate all grids at the input level by one step of its delta(t)
integer, dimension(maxlv) lstart
subroutine grdfit(lbase, lcheck, nvar, naux, cut, time, start_time)
This subroutine is called by setgrd and regrid to actually fit the new grids on each level...
integer, parameter store2
pointer to the address of memory storing the second copy of solution data on this grid...
real(kind=8), dimension(maxlv) rvoll
real(kind=8), dimension(maxlv) possk
integer, parameter ndjhi
global j index of upper border of this grid
subroutine arrangegrids(level, numg)
Sort all grids at level level.
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...
subroutine setgrd(nvar, cut, naux, dtinit, start_time)
subroutine prepf(level, nvar, naux)
For new fine grids (grids on level level), allocate space for saving fluxes at boundary (even if it's...