2D AMRCLAW
step2y.f90
Go to the documentation of this file.
1 !
2 ! -------------------------------------------------------------
3 !
4 subroutine step2y(maxm,meqn,maux,mbc,mx,my,qold,aux,dy,dt,cflgrid,gm,gp,rpn2)
5 !
6 ! clawpack routine ... modified for AMRCLAW
7 !
8 ! Take one time step, updating q.
9 ! On entry, qold gives
10 ! initial data for this step
11 ! and is unchanged in this version.
12 !
13 ! fm, fp are fluxes to left and right of single cell edge
14 ! See the flux2 documentation for more information.
15 !
16 ! Converted to f90 2012-1-04 (KTM)
17 !
18 
19  use amr_module
20 
21  implicit none
22 
23  external rpn2
24 
25  ! Arguments
26  integer, intent(in) :: maxm,meqn,maux,mbc,mx,my
27  real(kind=8), intent(in) :: dy,dt
28  real(kind=8), intent(inout) :: cflgrid
29  real(kind=8), intent(inout) :: qold(meqn, 1-mbc:mx+mbc, 1-mbc:my+mbc)
30  real(kind=8), intent(inout) :: aux(maux,1-mbc:mx+mbc, 1-mbc:my+mbc)
31  real(kind=8), intent(inout) :: gm(meqn,1-mbc:mx+mbc, 1-mbc:my+mbc)
32  real(kind=8), intent(inout) :: gp(meqn,1-mbc:mx+mbc, 1-mbc:my+mbc)
33 
34  ! Local storage for flux accumulation (these help compute normal g flux, based on step2 usage)
35  real(kind=8) :: faddm(meqn,1-mbc:maxm+mbc)
36  real(kind=8) :: faddp(meqn,1-mbc:maxm+mbc)
37 
38 
39  ! Scratch storage for Sweeps and Riemann problems
40  real(kind=8) :: q1d(meqn,1-mbc:maxm+mbc)
41  real(kind=8) :: aux2(maux,1-mbc:maxm+mbc)
42  real(kind=8) :: dtdy1d(1-mbc:maxm+mbc)
43 
44  real(kind=8) :: wave(meqn, mwaves, 1-mbc:maxm+mbc)
45  real(kind=8) :: s(mwaves, 1-mbc:maxm + mbc)
46  real(kind=8) :: cqxx(meqn,1-mbc:maxm + mbc)
47 
48  ! Looping scalar storage
49  integer :: i
50  real(kind=8) :: dtdy,cfl1d
51 
52 
53  cflgrid = 0.d0
54  dtdy = dt/dy
55 
56  gm = 0.d0
57  gp = 0.d0
58 
59 
60  ! ============================================================================
61  ! y-sweeps (for Godunov split method, called from stepgrid_dimSplit)
62  ! loop indices assume that x sweep has gone first, this is second sweep.
63  do i = 1,mx
64 
65  ! Copy data along a slice into 1d arrays:
66  q1d(:,1-mbc:my+mbc) = qold(:,i,1-mbc:my+mbc)
67 
68  ! Set dt/dy ratio in slice
69  if (mcapa > 0) then
70  dtdy1d(1-mbc:my+mbc) = dtdy / aux(mcapa,i,1-mbc:my+mbc)
71  else
72  dtdy1d = dtdy
73  endif
74 
75  ! Copy aux slices
76  if (maux .gt. 0) then
77 
78  aux2(:,1-mbc:my+mbc) = aux(:,i,1-mbc:my+mbc)
79 
80  endif
81 
82  ! Compute modifications fadd and gadd to fluxes along this slice
83 !!! ::: dont need aux1 and aux3 for dimensionally split method
84 !!! ::: but flux2 needs 3 aux arrays, so passing in aux2 3 times
85  call flux2_dimsplit(2,maxm,meqn,maux,mbc,my,q1d,dtdy1d,aux2, &
86  faddm,faddp,cfl1d,wave,s,cqxx,rpn2)
87 
88  cflgrid = max(cflgrid,cfl1d)
89 
90  ! Update fluxes
91  gm(:,i,1:my+1) = gm(:,i,1:my+1) + faddm(:,1:my+1)
92  gp(:,i,1:my+1) = gp(:,i,1:my+1) + faddp(:,1:my+1)
93 
94  enddo
95 
96 
97 end subroutine step2y
subroutine flux2_dimsplit(ixy, maxm, meqn, maux, mbc, mx, q1d, dtdx1d, aux2, faddm, faddp, cfl1d, wave, s, cqxx, rpn2)
Definition: flux2_dimSplit.f:8
subroutine step2y(maxm, meqn, maux, mbc, mx, my, qold, aux, dy, dt, cflgrid, gm, gp, rpn2)
Definition: step2y.f90:5
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
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