2D AMRCLAW
step2.f90
Go to the documentation of this file.
1 
7 subroutine step2(maxm,meqn,maux,mbc,mx,my,qold,aux,dx,dy,dt,cflgrid,fm,fp,gm,gp,rpn2,rpt2)
8 !
9 ! clawpack routine ... modified for AMRCLAW
10 !
11 ! Take one time step, updating q.
12 ! On entry, qold gives
13 ! initial data for this step
14 ! and is unchanged in this version.
15 !
16 ! fm, fp are fluxes to left and right of single cell edge
17 ! See the flux2 documentation for more information.
18 !
19 ! Converted to f90 2012-1-04 (KTM)
20 !
21 
22  use amr_module
23 
24  implicit none
25 
26  external rpn2, rpt2
27 
28  ! Arguments
29  integer, intent(in) :: maxm,meqn,maux,mbc,mx,my
30  real(kind=8), intent(in) :: dx,dy,dt
31  real(kind=8), intent(inout) :: cflgrid
32  real(kind=8), intent(inout) :: qold(meqn, 1-mbc:mx+mbc, 1-mbc:my+mbc)
33  real(kind=8), intent(inout) :: aux(maux,1-mbc:mx+mbc, 1-mbc:my+mbc)
34  real(kind=8), intent(inout) :: fm(meqn, 1-mbc:mx+mbc, 1-mbc:my+mbc)
35  real(kind=8), intent(inout) :: fp(meqn,1-mbc:mx+mbc, 1-mbc:my+mbc)
36  real(kind=8), intent(inout) :: gm(meqn,1-mbc:mx+mbc, 1-mbc:my+mbc)
37  real(kind=8), intent(inout) :: gp(meqn,1-mbc:mx+mbc, 1-mbc:my+mbc)
38 
39  ! Local storage for flux accumulation
40  real(kind=8) :: faddm(meqn,1-mbc:maxm+mbc)
41  real(kind=8) :: faddp(meqn,1-mbc:maxm+mbc)
42  real(kind=8) :: gaddm(meqn,1-mbc:maxm+mbc,2)
43  real(kind=8) :: gaddp(meqn,1-mbc:maxm+mbc,2)
44 
45  ! Scratch storage for Sweeps and Riemann problems
46  real(kind=8) :: q1d(meqn,1-mbc:maxm+mbc)
47  real(kind=8) :: aux1(maux,1-mbc:maxm+mbc)
48  real(kind=8) :: aux2(maux,1-mbc:maxm+mbc)
49  real(kind=8) :: aux3(maux,1-mbc:maxm+mbc)
50  real(kind=8) :: dtdx1d(1-mbc:maxm+mbc)
51  real(kind=8) :: dtdy1d(1-mbc:maxm+mbc)
52 
53  real(kind=8) :: wave(meqn, mwaves, 1-mbc:maxm+mbc)
54  real(kind=8) :: s(mwaves, 1-mbc:maxm + mbc)
55  real(kind=8) :: amdq(meqn,1-mbc:maxm + mbc)
56  real(kind=8) :: apdq(meqn,1-mbc:maxm + mbc)
57  real(kind=8) :: cqxx(meqn,1-mbc:maxm + mbc)
58  real(kind=8) :: bmadq(meqn,1-mbc:maxm + mbc)
59  real(kind=8) :: bpadq(meqn,1-mbc:maxm + mbc)
60 
61  ! Looping scalar storage
62  integer :: i,j,thread_num
63  real(kind=8) :: dtdx,dtdy,cfl1d
64 
65  ! Common block storage
66  integer :: icom,jcom
67  real(kind=8) :: dtcom,dxcom,dycom,tcom
68  common /comxyt/ dtcom,dxcom,dycom,tcom,icom,jcom
69 
70  ! Store mesh parameters in common block
71  dxcom = dx
72  dycom = dy
73  dtcom = dt
74 
75  cflgrid = 0.d0
76  dtdx = dt/dx
77  dtdy = dt/dy
78 
79  fm = 0.d0
80  fp = 0.d0
81  gm = 0.d0
82  gp = 0.d0
83 
84  ! ============================================================================
85  ! Perform X-Sweeps
86  do j = 0,my+1
87  ! Copy old q into 1d slice
88  q1d(:,1-mbc:mx+mbc) = qold(:,1-mbc:mx+mbc,j)
89 
90  ! Set dtdx slice if a capacity array exists
91  if (mcapa > 0) then
92  dtdx1d(1-mbc:mx+mbc) = dtdx / aux(mcapa,1-mbc:mx+mbc,j)
93  else
94  dtdx1d = dtdx
95  endif
96 
97  ! Copy aux array into slices
98  if (maux > 0) then
99  aux1(:,1-mbc:mx+mbc) = aux(:,1-mbc:mx+mbc,j-1)
100  aux2(:,1-mbc:mx+mbc) = aux(:,1-mbc:mx+mbc,j )
101  aux3(:,1-mbc:mx+mbc) = aux(:,1-mbc:mx+mbc,j+1)
102  endif
103 
104  ! Store value of j along the slice into common block
105  jcom = j
106 
107  ! Compute modifications fadd and gadd to fluxes along this slice:
108  call flux2(1,maxm,meqn,maux,mbc,mx,q1d,dtdx1d,aux1,aux2,aux3, &
109  faddm,faddp,gaddm,gaddp,cfl1d,wave,s, &
110  amdq,apdq,cqxx,bmadq,bpadq,rpn2,rpt2)
111 
112  cflgrid = max(cflgrid,cfl1d)
113 
114  ! Update fluxes
115  ! here gm(:,i,j) and gp(:,i,j) are the same since
116  ! they are both \tilde{G}_{i-1/2,j} in the textbook
117  fm(:,1:mx+1,j) = fm(:,1:mx+1,j) + faddm(:,1:mx+1)
118  fp(:,1:mx+1,j) = fp(:,1:mx+1,j) + faddp(:,1:mx+1)
119  gm(:,1:mx+1,j) = gm(:,1:mx+1,j) + gaddm(:,1:mx+1,1)
120  gp(:,1:mx+1,j) = gp(:,1:mx+1,j) + gaddp(:,1:mx+1,1)
121  gm(:,1:mx+1,j+1) = gm(:,1:mx+1,j+1) + gaddm(:,1:mx+1,2)
122  gp(:,1:mx+1,j+1) = gp(:,1:mx+1,j+1) + gaddp(:,1:mx+1,2)
123 
124  enddo
125 
126  ! ============================================================================
127  ! y-sweeps
128  !
129  do i = 0,mx+1
130 
131  ! Copy data along a slice into 1d arrays:
132  q1d(:,1-mbc:my+mbc) = qold(:,i,1-mbc:my+mbc)
133 
134  ! Set dt/dy ratio in slice
135  if (mcapa > 0) then
136  dtdy1d(1-mbc:my+mbc) = dtdy / aux(mcapa,i,1-mbc:my+mbc)
137  else
138  dtdy1d = dtdy
139  endif
140 
141  ! Copy aux slices
142  if (maux .gt. 0) then
143  aux1(:,1-mbc:my+mbc) = aux(:,i-1,1-mbc:my+mbc)
144  aux2(:,1-mbc:my+mbc) = aux(:,i,1-mbc:my+mbc)
145  aux3(:,1-mbc:my+mbc) = aux(:,i+1,1-mbc:my+mbc)
146  endif
147 
148  ! Store the value of i along this slice in the common block
149  icom = i
150 
151  ! Compute modifications fadd and gadd to fluxes along this slice
152  call flux2(2,maxm,meqn,maux,mbc,my,q1d,dtdy1d,aux1,aux2,aux3, &
153  faddm,faddp,gaddm,gaddp,cfl1d,wave,s,amdq,apdq,cqxx, &
154  bmadq,bpadq,rpn2,rpt2)
155 
156  cflgrid = max(cflgrid,cfl1d)
157 
158  ! Update fluxes
159  gm(:,i,1:my+1) = gm(:,i,1:my+1) + faddm(:,1:my+1)
160  gp(:,i,1:my+1) = gp(:,i,1:my+1) + faddp(:,1:my+1)
161  fm(:,i,1:my+1) = fm(:,i,1:my+1) + gaddm(:,1:my+1,1)
162  fp(:,i,1:my+1) = fp(:,i,1:my+1) + gaddp(:,1:my+1,1)
163  fm(:,i+1,1:my+1) = fm(:,i+1,1:my+1) + gaddm(:,1:my+1,2)
164  fp(:,i+1,1:my+1) = fp(:,i+1,1:my+1) + gaddp(:,1:my+1,2)
165 
166  enddo
167 
168 
169 end subroutine step2
integer mcapa
Definition: amr_module.f90:253
subroutine rpn2(ixy, maxm, meqn, mwaves, maux, mbc, mx, ql, qr, auxl, auxr, wave, s, amdq, apdq)
Definition: rpn2from1d.f:7
subroutine flux2(ixy, maxm, meqn, maux, mbc, mx, q1d, dtdx1d, aux1, aux2, aux3, faddm, faddp, gaddm, gaddp, cfl1d, wave, s, amdq, apdq, cqxx, bmasdq, bpasdq, rpn2, rpt2)
Solve Riemann problems based on 1D slice of the 2D grid.
Definition: flux2.f:29
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
integer mwaves
Definition: amr_module.f90:253