4 subroutine regrid (nvar,lbase,cut,naux,start_time)
7 implicit double precision (a-h,o-z)
8 integer newnumgrids(
maxlv)
9 integer clock_start2, clock_finish, clock_rate
41 20
if (lcheck .lt. lbase)
go to 50
42 call grdfit(lbase,lcheck,nvar,naux,cut,time,start_time)
43 if (
newstl(lcheck+1) .eq. 0)
go to 40
44 lfnew = max0(lcheck + 1,lfnew)
56 do lev = lbase+1,lfnew
59 52
if (mptr .eq. 0)
go to 55
60 ngridcount = ngridcount + 1
64 55 newnumgrids(lev) = ngridcount
65 maxnumnewgrids = max(maxnumnewgrids,ngridcount)
71 call system_clock(clock_start2,clock_rate)
72 call gfixup(lbase,lfnew,nvar,naux,newnumgrids,maxnumnewgrids)
73 call system_clock(clock_finish,clock_rate)
80 do 60 level = lbase,
lfine-1
81 call prepf(level+1,nvar,naux)
82 call prepc(level,nvar)
89 do 72 levnew = lbase+1,
mxnest 93 do while (mptr .gt. 0)
94 ngridcount = ngridcount + 1
105 if (ngridcount .gt. 1)
call arrangegrids(levnew,ngridcount)
108 write(*,100) ngridcount,ncells,levnew
109 write(
outunit,100) ngridcount,ncells,levnew
110 100
format(
"there are ",i6,
" grids with ",i10,
111 &
" cells at level ", i3)
117 do levnew = lbase+1,
lfine 135 implicit double precision (a-h,o-z)
136 integer listgrids(numg), cost(numg), index(numg), prevptr
153 call qsorti(index, numg, cost)
160 lstart(level) = listgrids(index(numg))
161 prevptr = listgrids(index(numg))
164 prevptr = listgrids(index(numg-i))
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
subroutine regrid(nvar, lbase, cut, naux, start_time)
integer, dimension(maxlv) numcells
subroutine prepc(level, nvar)
This routine is called because regridding just changed the fine grids.
subroutine qsorti(ORD, N, A)
subroutine makegridlist(lbase)
Make (modify) array of grid numbers, listOfGrids, (after sorting them in the linked list so they are ...
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
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, 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...
subroutine gfixup(lbase, lfnew, nvar, naux, newnumgrids, maxnumnewgrids)
Interpolate initial values for the newly created grids, whose levels start from level lbase+1...
integer, parameter outunit
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 prepf(level, nvar, naux)
For new fine grids (grids on level level), allocate space for saving fluxes at boundary (even if it's...