2D AMRCLAW
stepgrid.f
Go to the documentation of this file.
1 c
2 c -------------------------------------------------------------
3 c
24  subroutine stepgrid(q,fm,fp,gm,gp,mitot,mjtot,mbc,dt,dtnew,dx,dy,
25  & nvar,xlow,ylow,time,mptr,maux,aux)
26 c
27 c
28 c ::::::::::::::::::: STEPGRID ::::::::::::::::::::::::::::::::::::
29 c take a time step on a single grid. overwrite solution array q.
30 c A modified version of the clawpack routine step2 is used.
31 c
32 c return fluxes in fm,fp and gm,gp.
33 c patch has room for ghost cells (mbc of them) around the grid.
34 c everything is the enlarged size (mitot by mjtot).
35 c
36 c mbc = number of ghost cells (= lwidth)
37 c mptr = grid number (for debugging)
38 c xlow,ylow = lower left corner of enlarged grid (including ghost cells).
39 c dt = incoming time step
40 c dx,dy = mesh widths for this grid
41 c dtnew = return suggested new time step for this grid's soln.
42 c
43 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
44 
45  use amr_module
46  implicit double precision (a-h,o-z)
47  external rpn2,rpt2
48 
49  common /comxyt/ dtcom,dxcom,dycom,tcom,icom,jcom
50 
51  parameter(msize=max1d+4)
52  parameter(mwork=msize*(maxvar*maxvar + 13*maxvar + 3*maxaux +2))
53 
54  dimension q(nvar,mitot,mjtot)
55  dimension fp(nvar,mitot,mjtot),gp(nvar,mitot,mjtot)
56  dimension fm(nvar,mitot,mjtot),gm(nvar,mitot,mjtot)
57  dimension aux(maux,mitot,mjtot)
58 c dimension work(mwork)
59 
60  logical debug, dump
61  data debug/.false./, dump/.false./
62 
63 c
64 c # set tcom = time. This is in the common block comxyt that could
65 c # be included in the Riemann solver, for example, if t is explicitly
66 c # needed there.
67 
68  tcom = time
69 
70  if (dump) then
71  write(outunit,*) "dumping grid ",mptr," at time ",time
72  do i = 1, mitot
73  do j = 1, mjtot
74  write(outunit,545) i,j,(q(ivar,i,j),ivar=1,nvar)
75 c . ,(aux(ivar,i,j),ivar=1,maux)
76  545 format(2i4,5e15.7)
77  end do
78  end do
79  endif
80 c
81  meqn = nvar
82  mx = mitot - 2*mbc
83  my = mjtot - 2*mbc
84  maxm = max(mx,my) !# size for 1d scratch array
85  mbig = maxm
86  xlowmbc = xlow + mbc*dx
87  ylowmbc = ylow + mbc*dy
88 
89 c # method(2:7) and mthlim
90 c # are set in the amr2ez file (read by amr)
91 c
92  method(1) = 0
93 c
94 c
95 c # fluxes initialized in step2
96 c
97 C mwork0 = (maxm+2*mbc)*(12*meqn + mwaves + meqn*mwaves + 2)
98 c
99 C if (mwork .lt. mwork0) then
100 C write(outunit,*) 'CLAW2 ERROR... mwork must be increased to ',
101 C & mwork0
102 C write(* ,*) 'CLAW2 ERROR... mwork must be increased to ',
103 C & mwork0
104 C stop
105 C endif
106 c
107 c
108 c # partition work array into pieces needed for local storage in
109 c # step2 routine. Find starting index of each piece:
110 c
111 C i0faddm = 1
112 C i0faddp = i0faddm + (maxm+2*mbc)*meqn
113 C i0gaddm = i0faddp + (maxm+2*mbc)*meqn
114 C i0gaddp = i0gaddm + 2*(maxm+2*mbc)*meqn
115 C i0q1d = i0gaddp + 2*(maxm+2*mbc)*meqn
116 C i0dtdx1 = i0q1d + (maxm+2*mbc)*meqn
117 C i0dtdy1 = i0dtdx1 + (maxm+2*mbc)
118 C i0aux1 = i0dtdy1 + (maxm+2*mbc)
119 C i0aux2 = i0aux1 + (maxm+2*mbc)*maux
120 C i0aux3 = i0aux2 + (maxm+2*mbc)*maux
121 c
122 c
123 C i0next = i0aux3 + (maxm+2*mbc)*maux !# next free space
124 C mused = i0next - 1 !# space already used
125 C mwork1 = mwork - mused !# remaining space (passed to step2)
126 
127 c
128 c
129  call b4step2(mbc,mx,my,nvar,q,
130  & xlowmbc,ylowmbc,dx,dy,time,dt,maux,aux)
131 c
132 c
133 c # take one step on the conservation law:
134 c
135  call step2(mbig,nvar,maux,
136  & mbc,mx,my,
137  & q,aux,dx,dy,dt,cflgrid,
138  & fm,fp,gm,gp,rpn2,rpt2)
139 c
140 c
141 c write(outunit,1001) mptr, node(nestlevel,mptr),cflgrid
142 c1001 format(' Courant # of grid', i4,
143 c & ' on level', i3, ' is ', e10.3)
144 c
145 
146 !$OMP CRITICAL (cflm)
147 
148  cfl_level = dmax1(cfl_level,cflgrid)
149 
150 !$OMP END CRITICAL (cflm)
151 
152 c
153 c # update q
154  dtdx = dt/dx
155  dtdy = dt/dy
156  do 50 j=mbc+1,mjtot-mbc
157  do 50 i=mbc+1,mitot-mbc
158  do 50 m=1,nvar
159  if (mcapa.eq.0) then
160 c
161 c # no capa array. Standard flux differencing:
162 
163  q(m,i,j) = q(m,i,j)
164  & - dtdx * (fm(m,i+1,j) - fp(m,i,j))
165  & - dtdy * (gm(m,i,j+1) - gp(m,i,j))
166  else
167 c # with capa array.
168  q(m,i,j) = q(m,i,j)
169  & - (dtdx * (fm(m,i+1,j) - fp(m,i,j))
170  & + dtdy * (gm(m,i,j+1) - gp(m,i,j))) / aux(mcapa,i,j)
171  endif
172 
173  50 continue
174 c
175 c
176  if (method(5).eq.1) then
177 c # with source term: use Godunov splitting
178  call src2(nvar,mbc,mx,my,xlowmbc,ylowmbc,dx,dy,
179  & q,maux,aux,time,dt)
180  endif
181 c
182 c
183 c
184 c # output fluxes for debugging purposes:
185  if (debug) then
186  write(dbugunit,*)" fluxes for grid ",mptr
187 c do 830 j = mbc+1, mjtot-1
188  do 830 i = mbc+1, mitot-1
189  do 830 j = mbc+1, mjtot-1
190  write(dbugunit,831) i,j,fm(1,i,j),fp(1,i,j),
191  . gm(1,i,j),gp(1,i,j)
192  do 830 m = 2, meqn
193  write(dbugunit,832) fm(m,i,j),fp(m,i,j),
194  . gm(m,i,j),gp(m,i,j)
195  831 format(2i4,4d16.6)
196  832 format(8x,4d16.6)
197  830 continue
198  endif
199 
200 c
201 c
202 c For variable time stepping, use max speed seen on this grid to
203 c choose the allowable new time step dtnew. This will later be
204 c compared to values seen on other grids.
205 c
206  if (cflgrid .gt. 0.d0) then
207  dtnew = dt*cfl/cflgrid
208  else
209 c # velocities are all zero on this grid so there's no
210 c # time step restriction coming from this grid.
211  dtnew = rinfinity
212  endif
213 
214 c # give a warning if Courant number too large...
215 c
216  if (cflgrid .gt. cflv1) then
217  write(*,810) cflgrid
218  write(outunit,810) cflgrid, cflv1
219  810 format('*** WARNING *** Courant number =', d12.4,
220  & ' is larger than input cfl_max = ', d12.4)
221  endif
222 c
223  if (dump) then
224  write(outunit,*) "dumping grid ",mptr," after stepgrid"
225  do i = 1, mitot
226  do j = 1, mjtot
227  write(outunit,545) i,j,(q(ivar,i,j),ivar=1,nvar)
228  end do
229  end do
230  endif
231  return
232  end
233 
234 
integer, parameter dbugunit
Definition: amr_module.f90:293
integer, dimension(7) method
Definition: amr_module.f90:253
integer, parameter maxvar
Definition: amr_module.f90:183
real(kind=8), parameter rinfinity
Definition: amr_module.f90:169
real(kind=8) cfl
Definition: amr_module.f90:255
integer, parameter max1d
Definition: amr_module.f90:181
real(kind=8) cflv1
Definition: amr_module.f90:255
subroutine b4step2(mbc, mx, my, meqn, q, xlower, ylower, dx, dy, t, dt, maux, aux)
Called before each call to step2.
Definition: b4step2.f90:6
real(kind=8) cfl_level
Definition: amr_module.f90:255
integer, parameter maxaux
Definition: amr_module.f90:184
subroutine src2(meqn, mbc, mx, my, xlower, ylower, dx, dy, q, maux, aux, t, dt)
Definition: src2.f90:2
integer mcapa
Definition: amr_module.f90:253
subroutine stepgrid(q, fm, fp, gm, gp, mitot, mjtot, mbc, dt, dtnew, dx, dy, nvar, xlow, ylow, time, mptr, maux, aux)
Take a time step on a single grid mptr and overwrite solution array q.
Definition: stepgrid.f:26
integer, parameter outunit
Definition: amr_module.f90:290
subroutine rpn2(ixy, maxm, meqn, mwaves, maux, mbc, mx, ql, qr, auxl, auxr, wave, s, amdq, apdq)
Definition: rpn2from1d.f:7
subroutine step2(maxm, meqn, maux, mbc, mx, my, qold, aux, dx, dy, dt, cflgrid, fm, fp, gm, gp, rpn2, rpt2)
Compute all fluxes at cell edges.
Definition: step2.f90:8
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21