2D AMRCLAW
Functions/Subroutines
stepgrid_dimSplit.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine stepgrid_dimsplit (q, fm, fp, gm, gp, mitot, mjtot, mbc, dt, dtnew, dx, dy, nvar, xlow, ylow, time, mptr, maux, aux)
 

Function/Subroutine Documentation

◆ stepgrid_dimsplit()

subroutine stepgrid_dimsplit ( dimension(nvar,mitot,mjtot)  q,
dimension(nvar,mitot,mjtot)  fm,
dimension(nvar,mitot,mjtot)  fp,
dimension(nvar,mitot,mjtot)  gm,
dimension(nvar,mitot,mjtot)  gp,
  mitot,
  mjtot,
  mbc,
  dt,
  dtnew,
  dx,
  dy,
  nvar,
  xlow,
  ylow,
  time,
  mptr,
  maux,
dimension(maux,mitot,mjtot)  aux 
)

Definition at line 7 of file stepgrid_dimSplit.f.

References b4step2(), amr_module::cfl, amr_module::cfl_level, amr_module::cflv1, amr_module::dbugunit, amr_module::max1d, amr_module::maxaux, amr_module::maxvar, amr_module::mcapa, amr_module::method, amr_module::outunit, amr_module::rinfinity, rpn2(), src2(), step2x(), and step2y().

Referenced by par_advanc().

7 c
8 c
9 c ::::::::::::::::::: STEPGRID_DIMSPLIT ::::::::::::::::::::::::::::::::::::
10 c dimensionally split version of stepgrid
11 c take a step in x, then y. for now only godunov splitting
12 c not strang splitting.
13 c
14 c take a time step on a single grid. overwrite solution array q.
15 c A modified version of the clawpack routine step2 is used.
16 c
17 c return fluxes in fm,fp and gm,gp.
18 c patch has room for ghost cells (mbc of them) around the grid.
19 c everything is the enlarged size (mitot by mjtot).
20 c
21 c mbc = number of ghost cells (= lwidth)
22 c mptr = grid number (for debugging)
23 c xlow,ylow = lower left corner of enlarged grid (including ghost cells).
24 c dt = incoming time step
25 c dx,dy = mesh widths for this grid
26 c dtnew = return suggested new time step for this grid's soln.
27 c
28 c :::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
29 
30  use amr_module
31  implicit double precision (a-h,o-z)
32  external rpn2
33 
34  parameter(msize=max1d+4)
35  parameter(mwork=msize*(maxvar*maxvar + 13*maxvar + 3*maxaux +2))
36 
37  dimension q(nvar,mitot,mjtot)
38  dimension fp(nvar,mitot,mjtot),gp(nvar,mitot,mjtot)
39  dimension fm(nvar,mitot,mjtot),gm(nvar,mitot,mjtot)
40  dimension aux(maux,mitot,mjtot)
41 C dimension work(mwork)
42 
43  logical debug, dump
44  data debug/.false./, dump/.false./
45 
46 c
47  if (dump) then
48  write(outunit,*) "dumping grid ",mptr," at time ",time
49  do i = 1, mitot
50  do j = 1, mjtot
51  write(outunit,545) i,j,(q(ivar,i,j),ivar=1,nvar),
52  . (aux(ivar,i,j),ivar=1,maux)
53  545 format(2i4,4e15.7)
54  end do
55  end do
56  endif
57 c
58  meqn = nvar
59  mx = mitot - 2*mbc
60  my = mjtot - 2*mbc
61  maxm = max(mx,my) !# size for 1d scratch array
62  mbig = maxm
63  xlowmbc = xlow + mbc*dx
64  ylowmbc = ylow + mbc*dy
65 
66 c
67 c old work array code now removed.
68 c
69 c
70  call b4step2(mbc,mx,my,nvar,q,
71  & xlowmbc,ylowmbc,dx,dy,time,dt,maux,aux)
72 c
73 c
74 c # take one step on the conservation law: Godunov splitting
75 c
76 c 1. First step in x direction:
77 
78  call step2x(mbig,nvar,maux, mbc,mx,my,q,aux,dx,dt,cflgrid,
79  & fm,fp,rpn2)
80 c
81 !$OMP CRITICAL (cflm)
82 
83  cfl_level = dmax1(cfl_level,cflgrid)
84 
85 !$OMP END CRITICAL (cflm)
86 
87 c # update q with x fluxes first. all rows in y get updated
88  dtdx = dt/dx
89  do 50 j=1,mjtot
90  do 50 i=mbc+1,mitot-mbc
91  do 50 m=1,nvar
92  if (mcapa.eq.0) then
93 c
94 c # no capa array. Standard flux differencing:
95 
96  q(m,i,j) = q(m,i,j) - dtdx * (fm(m,i+1,j) - fp(m,i,j))
97 
98  else
99 c # with capa array.
100  q(m,i,j) = q(m,i,j)
101  & - dtdx * (fm(m,i+1,j) - fp(m,i,j)) / aux(mcapa,i,j)
102  endif
103 
104  50 continue
105 
106 c 2. Second step in x direction:
107 
108  call step2y(mbig,nvar,maux, mbc,mx,my,q,aux,dy,dt,cflgrid,
109  & gm,gp,rpn2)
110 c
111 !$OMP CRITICAL (cflm)
112 
113  cfl_level = dmax1(cfl_level,cflgrid)
114 
115 !$OMP END CRITICAL (cflm)
116 c
117 c # update q
118  dtdy = dt/dy
119  do 51 j=mbc+1,mjtot-mbc
120  do 51 i=mbc+1,mitot-mbc
121  do 51 m=1,nvar
122  if (mcapa.eq.0) then
123 c
124 c # no capa array. Standard flux differencing:
125 
126  q(m,i,j) = q(m,i,j) - dtdy * (gm(m,i,j+1) - gp(m,i,j))
127  else
128 c # with capa array.
129  q(m,i,j) = q(m,i,j)
130  & - dtdy * (gm(m,i,j+1) - gp(m,i,j)) / aux(mcapa,i,j)
131  endif
132 
133  51 continue
134 c
135 !-- write(outunit,1001) mptr, node(nestlevel,mptr),cflgrid
136 !-- 1001 format(' Courant # of grid', i4,
137 !-- & ' on level', i3, ' is ', e10.3)
138 
139 c
140  if (method(5).eq.1) then
141 c # with source term: use Godunov splitting for this too
142  call src2(nvar,mbc,mx,my,xlowmbc,ylowmbc,dx,dy,
143  & q,maux,aux,time,dt)
144  endif
145 c
146 c
147 c
148 c # output fluxes for debugging purposes:
149  if (debug) then
150  write(dbugunit,*)" fluxes for grid ",mptr
151 c do 830 j = mbc+1, mjtot-1
152  do 830 i = mbc+1, mitot-1
153  do 830 j = mbc+1, mjtot-1
154  write(dbugunit,831) i,j,fm(1,i,j),fp(1,i,j),
155  . gm(1,i,j),gp(1,i,j)
156  do 830 m = 2, meqn
157  write(dbugunit,832) fm(m,i,j),fp(m,i,j),
158  . gm(m,i,j),gp(m,i,j)
159  831 format(2i4,4d16.6)
160  832 format(8x,4d16.6)
161  830 continue
162  endif
163 
164 c
165 c
166 c For variable time stepping, use max speed seen on this grid to
167 c choose the allowable new time step dtnew. This will later be
168 c compared to values seen on other grids.
169 c NB: orig code tested cflgrid, but with dimensional splitting
170 c there are 2 - one for each x and y steps, so changed 12/23/14 (mjb)
171 c to use level
172 c
173  if (cfl_level .gt. 0.d0) then
174  dtnew = dt*cfl/cfl_level
175  else
176 c # velocities are all zero on this grid so there's no
177 c # time step restriction coming from this grid.
178  dtnew = rinfinity
179  endif
180 
181 c # give a warning if Courant number too large...
182 c
183  if (cflgrid .gt. cflv1) then
184  write(*,810) cflgrid
185  write(outunit,810) cflgrid, cflv1
186  810 format('*** WARNING *** Courant number =', d12.4,
187  & ' is larger than input cfl_max = ', d12.4)
188  endif
189 c
190  if (dump) then
191  write(outunit,*) "dumping grid ",mptr," after stepgrid"
192  do i = 1, mitot
193  do j = 1, mjtot
194  write(outunit,545) i,j,(q(ivar,i,j),ivar=1,nvar)
195  end do
196  end do
197  endif
198  return
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 step2y(maxm, meqn, maux, mbc, mx, my, qold, aux, dy, dt, cflgrid, gm, gp, rpn2)
Definition: step2y.f90:5
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 step2x(maxm, meqn, maux, mbc, mx, my, qold, aux, dx, dt, cflgrid, fm, fp, rpn2)
Definition: step2x.f90:5
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
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21
Here is the call graph for this function:
Here is the caller graph for this function: