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