/
rpn2_geoclaw.f
300 lines (258 loc) · 9.8 KB
/
rpn2_geoclaw.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
c======================================================================
subroutine rpn2(ixy,maxm,meqn,mwaves,maux,mbc,mx,
& ql,qr,auxl,auxr,fwave,s,amdq,apdq)
c======================================================================
c
c Solves normal Riemann problems for the 2D SHALLOW WATER equations
c with topography:
c # h_t + (hu)_x + (hv)_y = 0 #
c # (hu)_t + (hu^2 + 0.5gh^2)_x + (huv)_y = -ghb_x #
c # (hv)_t + (huv)_x + (hv^2 + 0.5gh^2)_y = -ghb_y #
c On input, ql contains the state vector at the left edge of each cell
c qr contains the state vector at the right edge of each cell
c
c This data is along a slice in the x-direction if ixy=1
c or the y-direction if ixy=2.
c Note that the i'th Riemann problem has left state qr(i-1,:)
c and right state ql(i,:)
c From the basic clawpack routines, this routine is called with
c ql = qr
c
c
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! !
! # This Riemann solver is for the shallow water equations. !
! !
! It allows the user to easily select a Riemann solver in !
! riemannsolvers_geo.f. this routine initializes all the variables !
! for the shallow water equations, accounting for wet dry boundary !
! dry cells, wave speeds etc. !
! !
! David George, Vancouver WA, Feb. 2009 !
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
use geoclaw_module, only: g => grav, drytol => dry_tolerance, rho
use geoclaw_module, only: earth_radius, deg2rad
use amr_module, only: mcapa
use storm_module, only: pressure_forcing, pressure_index
implicit none
!input
integer maxm,meqn,maux,mwaves,mbc,mx,ixy
double precision fwave(meqn, mwaves, 1-mbc:maxm+mbc)
double precision s(mwaves, 1-mbc:maxm+mbc)
double precision ql(meqn, 1-mbc:maxm+mbc)
double precision qr(meqn, 1-mbc:maxm+mbc)
double precision apdq(meqn,1-mbc:maxm+mbc)
double precision amdq(meqn,1-mbc:maxm+mbc)
double precision auxl(maux,1-mbc:maxm+mbc)
double precision auxr(maux,1-mbc:maxm+mbc)
!local only
integer m,i,mw,maxiter,mu,nv
double precision wall(3)
double precision fw(3,3)
double precision sw(3)
double precision hR,hL,huR,huL,uR,uL,hvR,hvL,vR,vL,phiR,phiL,pL,pR
double precision bR,bL,sL,sR,sRoe1,sRoe2,sE1,sE2,uhat,chat
double precision s1m,s2m
double precision hstar,hstartest,hstarHLL,sLtest,sRtest
double precision tw,dxdc
logical rare1,rare2
! In case there is no pressure forcing
pL = 0.d0
pR = 0.d0
!loop through Riemann problems at each grid cell
do i=2-mbc,mx+mbc
!-----------------------Initializing-----------------------------------
!inform of a bad riemann problem from the start
if((qr(1,i-1).lt.0.d0).or.(ql(1,i) .lt. 0.d0)) then
write(*,*) 'Negative input: hl,hr,i=',qr(1,i-1),ql(1,i),i
endif
!Initialize Riemann problem for grid interface
do mw=1,mwaves
s(mw,i)=0.d0
fwave(1,mw,i)=0.d0
fwave(2,mw,i)=0.d0
fwave(3,mw,i)=0.d0
enddo
c !set normal direction
if (ixy.eq.1) then
mu=2
nv=3
else
mu=3
nv=2
endif
!zero (small) negative values if they exist
if (qr(1,i-1).lt.0.d0) then
qr(1,i-1)=0.d0
qr(2,i-1)=0.d0
qr(3,i-1)=0.d0
endif
if (ql(1,i).lt.0.d0) then
ql(1,i)=0.d0
ql(2,i)=0.d0
ql(3,i)=0.d0
endif
!skip problem if in a completely dry area
if (qr(1,i-1) <= drytol .and. ql(1,i) <= drytol) then
go to 30
endif
!Riemann problem variables
hL = qr(1,i-1)
hR = ql(1,i)
huL = qr(mu,i-1)
huR = ql(mu,i)
bL = auxr(1,i-1)
bR = auxl(1,i)
if (pressure_forcing) then
pL = auxr(pressure_index, i-1)
pR = auxl(pressure_index, i)
end if
hvL=qr(nv,i-1)
hvR=ql(nv,i)
!check for wet/dry boundary
if (hR.gt.drytol) then
uR=huR/hR
vR=hvR/hR
phiR = 0.5d0*g*hR**2 + huR**2/hR
else
hR = 0.d0
huR = 0.d0
hvR = 0.d0
uR = 0.d0
vR = 0.d0
phiR = 0.d0
endif
if (hL.gt.drytol) then
uL=huL/hL
vL=hvL/hL
phiL = 0.5d0*g*hL**2 + huL**2/hL
else
hL=0.d0
huL=0.d0
hvL=0.d0
uL=0.d0
vL=0.d0
phiL = 0.d0
endif
wall(1) = 1.d0
wall(2) = 1.d0
wall(3) = 1.d0
if (hR.le.drytol) then
call riemanntype(hL,hL,uL,-uL,hstar,s1m,s2m,
& rare1,rare2,1,drytol,g)
hstartest=max(hL,hstar)
if (hstartest+bL.lt.bR) then !right state should become ghost values that mirror left for wall problem
c bR=hstartest+bL
wall(2)=0.d0
wall(3)=0.d0
hR=hL
huR=-huL
bR=bL
phiR=phiL
uR=-uL
vR=vL
elseif (hL+bL.lt.bR) then
bR=hL+bL
endif
elseif (hL.le.drytol) then ! right surface is lower than left topo
call riemanntype(hR,hR,-uR,uR,hstar,s1m,s2m,
& rare1,rare2,1,drytol,g)
hstartest=max(hR,hstar)
if (hstartest+bR.lt.bL) then !left state should become ghost values that mirror right
c bL=hstartest+bR
wall(1)=0.d0
wall(2)=0.d0
hL=hR
huL=-huR
bL=bR
phiL=phiR
uL=-uR
vL=vR
elseif (hR+bR.lt.bL) then
bL=hR+bR
endif
endif
!determine wave speeds
sL=uL-sqrt(g*hL) ! 1 wave speed of left state
sR=uR+sqrt(g*hR) ! 2 wave speed of right state
uhat=(sqrt(g*hL)*uL + sqrt(g*hR)*uR)/(sqrt(g*hR)+sqrt(g*hL)) ! Roe average
chat=sqrt(g*0.5d0*(hR+hL)) ! Roe average
sRoe1=uhat-chat ! Roe wave speed 1 wave
sRoe2=uhat+chat ! Roe wave speed 2 wave
sE1 = min(sL,sRoe1) ! Eindfeldt speed 1 wave
sE2 = max(sR,sRoe2) ! Eindfeldt speed 2 wave
!--------------------end initializing...finally----------
!solve Riemann problem.
maxiter = 1
call riemann_aug_JCP(maxiter,3,3,hL,hR,huL,
& huR,hvL,hvR,bL,bR,uL,uR,vL,vR,phiL,phiR,pL,pR,sE1,sE2,
& drytol,g,rho,sw,fw)
C call riemann_ssqfwave(maxiter,meqn,mwaves,hL,hR,huL,huR,
C & hvL,hvR,bL,bR,uL,uR,vL,vR,phiL,phiR,pL,pR,sE1,sE2,drytol,g,
C & rho,sw,fw)
C call riemann_fwave(meqn,mwaves,hL,hR,huL,huR,hvL,hvR,
C & bL,bR,uL,uR,vL,vR,phiL,phiR,pL,pR,sE1,sE2,drytol,g,rho,sw,
C & fw)
c !eliminate ghost fluxes for wall
do mw=1,3
sw(mw)=sw(mw)*wall(mw)
fw(1,mw)=fw(1,mw)*wall(mw)
fw(2,mw)=fw(2,mw)*wall(mw)
fw(3,mw)=fw(3,mw)*wall(mw)
enddo
do mw=1,mwaves
s(mw,i)=sw(mw)
fwave(1,mw,i)=fw(1,mw)
fwave(mu,mw,i)=fw(2,mw)
fwave(nv,mw,i)=fw(3,mw)
! write(51,515) sw(mw),fw(1,mw),fw(2,mw),fw(3,mw)
!515 format("++sw",4e25.16)
enddo
30 continue
enddo
c==========Capacity for mapping from latitude longitude to physical space====
if (mcapa.gt.0) then
do i=2-mbc,mx+mbc
if (ixy.eq.1) then
dxdc=(earth_radius*deg2rad)
else
dxdc=earth_radius*cos(auxl(3,i))*deg2rad
endif
do mw=1,mwaves
c if (s(mw,i) .gt. 316.d0) then
c # shouldn't happen unless h > 10 km!
c write(6,*) 'speed > 316: i,mw,s(mw,i): ',i,mw,s(mw,i)
c endif
s(mw,i)=dxdc*s(mw,i)
fwave(1,mw,i)=dxdc*fwave(1,mw,i)
fwave(2,mw,i)=dxdc*fwave(2,mw,i)
fwave(3,mw,i)=dxdc*fwave(3,mw,i)
enddo
enddo
endif
c===============================================================================
c============= compute fluctuations=============================================
amdq(1:3,:) = 0.d0
apdq(1:3,:) = 0.d0
do i=2-mbc,mx+mbc
do mw=1,mwaves
if (s(mw,i) < 0.d0) then
amdq(1:3,i) = amdq(1:3,i) + fwave(1:3,mw,i)
else if (s(mw,i) > 0.d0) then
apdq(1:3,i) = apdq(1:3,i) + fwave(1:3,mw,i)
else
amdq(1:3,i) = amdq(1:3,i) + 0.5d0 * fwave(1:3,mw,i)
apdq(1:3,i) = apdq(1:3,i) + 0.5d0 * fwave(1:3,mw,i)
endif
enddo
enddo
!-- do i=2-mbc,mx+mbc
!-- do m=1,meqn
!-- write(51,151) m,i,amdq(m,i),apdq(m,i)
!-- write(51,152) fwave(m,1,i),fwave(m,2,i),fwave(m,3,i)
!--151 format("++3 ampdq ",2i4,2e25.15)
!--152 format("++3 fwave ",8x,3e25.15)
!-- enddo
!-- enddo
return
end subroutine