forked from clawpack/riemann
-
Notifications
You must be signed in to change notification settings - Fork 0
/
rp1_layered_shallow_water.f90
766 lines (632 loc) · 28.9 KB
/
rp1_layered_shallow_water.f90
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
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
! ============================================================================
! Reimann Solver for the two-layer shallow water equations
! ============================================================================
subroutine rp1(maxmx,meqn,mwaves,maux,mbc,mx,ql,qr,auxl,auxr,fwave,s,amdq,apdq)
implicit none
! Input arguments
integer, intent(in) :: maxmx,meqn,mwaves,mbc,mx,maux
double precision, intent(in), dimension(meqn, 1-mbc:maxmx+mbc) :: ql,qr
double precision, intent(in), dimension(maux, 1-mbc:maxmx+mbc) :: auxl,auxr
! Output arguments
double precision, intent(out) :: s(mwaves, 1-mbc:maxmx+mbc)
double precision, intent(out) :: fwave(meqn, mwaves, 1-mbc:maxmx+mbc)
double precision, intent(out), dimension(meqn, 1-mbc:maxmx+mbc) :: amdq,apdq
! Local storage
integer :: i,j,m,mw
double precision :: eig_vec(4,4),delta(4),alpha(4),beta(4)
double precision, dimension(4) :: flux_l,flux_r
double precision, dimension(2) :: h_l,u_l,hu_l,h_r,u_r,hu_r,h_ave,u_ave
double precision :: b_l,b_r,gamma_l,gamma_r,tau,w_l,w_r
double precision :: wind_speed,lambda(4),eta_l(2),eta_r(2),h_hat_l(2),h_hat_r(2)
logical :: dry_state_l(2), dry_state_r(2), inundation
integer :: layer_index
double precision :: momentum_transfer(2),flux_transfer_r,flux_transfer_l
double precision :: inundation_height(2),lambda_l,lambda_r,lambda_hat
integer :: trans_wave(2-mbc:mx+mbc)
double precision :: wave_correction(2-mbc:mx+mbc)
logical :: stop_computation
! Function interfaces
interface
function eval_lapack_solve(eig_vec,delta) result(beta)
implicit none
double precision, intent(in) :: eig_vec(4,4), delta(4)
double precision :: beta(4)
end function eval_lapack_solve
end interface
interface
double precision function wind_drag(wind_speed)
implicit none
double precision, intent(in) :: wind_speed
end function wind_drag
end interface
! Common blocks
double precision :: rho(2),rho_air,g,r,one_minus_r
double precision :: dry_tolerance
integer :: eigen_method,inundation_method
logical :: entropy_fix
common /cparam/ rho,rho_air,g,r,one_minus_r, &
dry_tolerance,eigen_method,inundation_method,entropy_fix
! Initialize return variables
amdq = 0.d0
apdq = 0.d0
! | |
! | |
! | |
!----------|----------|----------
! i-1 i i+1
! qr(i-1) ql(i)
! Riemann problem at i-1/2 is between qr(i-1) (left state) and ql(i)
! (right state)
do i=2-mbc,mx+mbc
inundation = .false.
dry_state_l = .false.
dry_state_r = .false.
! These state variables are the actual depth and momenta
do j=1,2
layer_index = 2*(j-1)
h_l(j) = qr(layer_index+1,i-1) / rho(j)
h_r(j) = ql(layer_index+1,i) / rho(j)
hu_l(j) = qr(layer_index+2,i-1) / rho(j)
hu_r(j) = ql(layer_index+2,i) / rho(j)
h_hat_l(j) = auxr(j+2,i-1)
h_hat_r(j) = auxl(j+2,i)
! Check for dry states in this layer
if (h_l(j) < dry_tolerance) then
dry_state_l(j) = .true.
h_l(j) = 0.d0
u_l(j) = 0.d0
else
u_l(j) = qr(layer_index+2,i-1) / qr(layer_index+1,i-1)
endif
if (h_r(j) < dry_tolerance) then
dry_state_r(j) = .true.
h_r(j) = 0.d0
u_r(j) = 0.d0
else
u_r(j) = ql(layer_index+2,i) / ql(layer_index+1,i)
endif
enddo
h_ave = 0.5d0 * (h_l(:) + h_r(:))
b_l = auxr(1,i-1)
b_r = auxl(1,i)
w_l = auxr(2,i-1)
w_r = auxl(2,i)
! ====================================================================
! Solve Single layer problem seperately
if (dry_state_r(2).and.dry_state_l(2)) then
call single_layer_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
delta(1) = rho(1) * (hu_r(1) - hu_l(1))
flux_r(2) = rho(1) * (h_r(1) * u_r(1)**2 + 0.5d0 * g * h_r(1)**2)
flux_l(2) = rho(1) * (h_l(1) * u_l(1)**2 + 0.5d0 * g * h_l(1)**2)
delta(2) = flux_r(2) - flux_l(2) + g * rho(1) * h_ave(1) * (b_r - b_l)
beta(1) = (delta(1)*s(4,i)-delta(2)) / (s(4,i)-s(1,i))
beta(2) = 0.d0
beta(3) = 0.d0
beta(4) = (delta(2) - s(1,i)*beta(1)) / s(4,i)
! Calculate waves
forall(mw=1:4)
fwave(:,mw,i) = eig_vec(:,mw) * beta(mw)
end forall
cycle
endif
! ====================================================================
! Calculate eigen-space values
! ====================================================================
! ====================================================================
! Inundation cases
if (dry_state_r(2).and.(.not.dry_state_l(2)).and.(h_l(2) + b_l > b_r)) then
! print *,"Right inundation problem at i=",i
inundation = .true.
if (inundation_method == 0) then
stop "Inundation not allowed."
else if (inundation_method == 1) then
! Linear eigensystem
inundation_height = [h_r(1),0.d0]
call linear_eigen(h_l,inundation_height,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! Corrective wave
s(3,i) = u_l(2) + 2.d0 * sqrt(g*(1.d0-r)*h_l(2))
s(4,i) = u_r(1) + sqrt(g*h_r(1))
alpha(3) = r * g * h_l(2) / ((s(3,i) - u_l(2))**2 - g * h_l(2))
alpha(4) = 0.d0
eig_vec(1,3:4) = 1.d0
eig_vec(2,3:4) = s(3:4,i)
eig_vec(3,3:4) = alpha(3:4)
eig_vec(4,3:4) = s(3:4,i)*alpha(3:4)
else if (inundation_method == 2) then
inundation_height = [h_r(1),dry_tolerance]
h_r(2) = dry_tolerance
call linear_eigen(h_l,inundation_height,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! s(4,i) = u_r(1) + sqrt(g*h_r(1))
! eig_vec(:,4) = [1.d0,s(4,i),0.d0,0.d0]
else if (inundation_method == 3) then
inundation_height = [h_r(1),dry_tolerance]
call velocity_eigen(h_l,inundation_height,u_l,u_r,b_l,b_r,lambda,eig_vec)
s(:,i) = lambda
! Correction for fast wave
s(4,i) = u_r(1) + sqrt(g*h_r(1))
eig_vec(:,4) = [1.d0,s(4,i),0.d0,0.d0]
else if (inundation_method == 4) then
inundation_height = [h_r(1),dry_tolerance]
call lapack_eigen(h_l,inundation_height,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! Correction for the fast waves
s(2,i) = u_r(1) + sqrt(g*h_r(1))
eig_vec(:,2) = [1.d0,s(2,i),0.d0,0.d0]
else if (inundation_method == 5) then
inundation_height = [h_r(1),0.d0]
call lapack_eigen(h_l,inundation_height,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
endif
else if (dry_state_l(2).and.(.not.dry_state_r(2)).and.(h_r(2) + b_r > b_l)) then
! print *,"Left inundation problem at i=",i
inundation = .true.
! Inundation problem eigen
if (inundation_method == 0) then
stop "Inundation not allowed."
else if (inundation_method == 1) then
! Linear eigensystem
inundation_height = [h_l(1),dry_tolerance]
call linear_eigen(inundation_height,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! Corrections to internal wave
s(2,i) = u_r(2) - 2.d0 * sqrt(g*(1.d0-r)*h_r(2))
alpha(2) = r * g * h_r(2) / ((s(2,i) - u_r(2))**2 - g * h_r(2))
eig_vec(:,2) = [1.d0,s(2,i),alpha(2),alpha(2)*s(2,i)]
! Correction for the fast waves
s(1,i) = u_l(1) - sqrt(g*h_l(1))
eig_vec(:,1) = [1.d0,s(1,i),0.d0,0.d0]
else if (inundation_method == 2) then
! Use linearized eigensystem
inundation_height = [h_l(1),dry_tolerance]
call linear_eigen(inundation_height,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! Correction for the fast waves
s(1,i) = u_l(1) - sqrt(g*h_l(1))
eig_vec(:,1) = [1.d0,s(1,i),0.d0,0.d0]
else if (inundation_method == 3) then
! Use velocity difference expansion eigensystems
inundation_height = [h_l(1),dry_tolerance]
call velocity_eigen(inundation_height,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! Correction for the fast waves
s(1,i) = u_r(1) - sqrt(g*h_r(1))
eig_vec(:,1) = [1.d0,s(1,i),0.d0,0.d0]
else if (inundation_method == 4) then
! LAPACK solver with corrective wave and small wet layer
inundation_height = [h_r(1),dry_tolerance]
call lapack_eigen(h_l,inundation_height,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
! Correction for the fast waves
s(1,i) = u_l(1) - sqrt(g*h_l(1))
eig_vec(:,1) = [1.d0,s(1,i),0.d0,0.d0]
else if (inundation_method == 5) then
! Use the LAPACK solver with no correction
inundation_height = [h_l(1),0.d0]
call lapack_eigen(inundation_height,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
endif
! ====================================================================
! Wall or wet case
else
! Wall dry state or completely wet case
if (eigen_method == 1) then
call linear_eigen(h_hat_l,h_hat_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
else if (eigen_method == 2) then
call linear_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
else if (eigen_method == 3) then
call velocity_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
else if (eigen_method == 4) then
if (dry_state_r(2).and.(.not.dry_state_l(2)).or. &
dry_state_l(2).and.(.not.dry_state_r(2))) then
call linear_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
else
call lapack_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& trans_wave(i),wave_correction(i),lambda,eig_vec)
s(:,i) = lambda
endif
else
stop "Invalid eigensystem method requested, method = (1,4)."
endif
endif
! end of eigenspace calculation
! ====================================================================
! ====================================================================
! Calculate flux vector to be projected onto e-space
! Calculate jumps in fluxes
if (dry_state_r(2).and.(.not.dry_state_l(2)).and.(.not.inundation)) then
! Wall boundary conditions
h_r(2) = h_l(2)
hu_r(2) = -hu_l(2)
u_r(2) = -u_l(2)
! Top layer eta(2) = b_r - h_l(2) - b_l
momentum_transfer(1) = g * rho(1) * h_ave(1) * (b_r - h_l(2) - b_l)
momentum_transfer(2) = 0.d0
flux_transfer_r = 0.d0
flux_transfer_l = 0.d0
! Left state dry, right wet
else if (dry_state_l(2).and.(.not.dry_state_r(2)).and.(.not.inundation)) then
! Wall boundary conditions
h_l(2) = h_r(2)
hu_l(2) = -hu_r(2)
u_l(2) = -u_r(2)
! Top layer eta(2) = h_r(2) + b_r - b_l
momentum_transfer(1) = g * rho(1) * h_ave(1) * (b_r + h_r(2) - b_l)
momentum_transfer(2) = 0.d0
flux_transfer_r = 0.d0
flux_transfer_l = 0.d0
! Fully wet bottom layer or inundation
else
momentum_transfer(1) = g * rho(1) * h_ave(1) * (h_r(2) - h_l(2) + b_r - b_l)
momentum_transfer(2) = - g * rho(1) * h_ave(1) * (h_r(2) - h_l(2)) + g * rho(2) * h_ave(2) * (b_r - b_l)
! Bottom layer momentum transfer flux
flux_transfer_r = g * rho(1) * product(h_r)
flux_transfer_l = g * rho(1) * product(h_l)
endif
! Flux jumps
do j=1,2
layer_index = 2*(j-1)
flux_r(layer_index+1) = rho(j) * hu_r(j)
flux_r(layer_index+2) = rho(j) * (h_r(j) * u_r(j)**2 + 0.5d0 * g * h_r(j)**2)
flux_l(layer_index+1) = rho(j) * hu_l(j)
flux_l(layer_index+2) = rho(j) * (h_l(j) * u_l(j)**2 + 0.5d0 * g * h_l(j)**2)
enddo
delta = flux_r - flux_l
! Bottom layer additional flux
delta(4) = delta(4) + flux_transfer_r - flux_transfer_l
! Momentum source term from layers
delta(2) = delta(2) + momentum_transfer(1)
delta(4) = delta(4) + momentum_transfer(2)
! Wind forcing
wind_speed = 0.5d0 * (w_l + w_r)
tau = wind_drag(wind_speed) * rho_air * wind_speed
delta(2) = delta(2) - tau * wind_speed
! ====================================================================
! Solve system, solution is stored in delta
beta = eval_lapack_solve(eig_vec,delta)
! Calculate waves
forall(mw=1:4)
fwave(:,mw,i) = eig_vec(:,mw) * beta(mw)
end forall
enddo
! Calculate amdq and apdq
! No entropy fix requested
if (.not.entropy_fix) then
stop_computation = .false.
do i=2-mbc,mx+mbc
do mw=1,mwaves
if (s(mw,i) > 0.d0) then
apdq(:,i) = apdq(:,i) + fwave(:,mw,i)
else
amdq(:,i) = amdq(:,i) + fwave(:,mw,i)
endif
enddo
enddo
else
do i=2-mbc,mx+mbc
! Check to see if an entropy fix is necessary
if (trans_wave(i) /= 0) then
print *,"Entropy fix needed: i=",i
print *," s(:,i)=",(s(m,i),m=1,4)
print *," trans_wave=",trans_wave(i)
print *," beta =",wave_correction(i)
do mw=1,trans_wave(i)-1
amdq(:,i) = amdq(:,i) + fwave(:,mw,i)
enddo
amdq(:,i) = amdq(:,i) + wave_correction(i) * fwave(:,trans_wave(i),i)
apdq(:,i) = apdq(:,i) + (1.d0 - wave_correction(i)) * fwave(:,trans_wave(i),i)
do mw=trans_wave(i)+1,mwaves
apdq(:,i) = apdq(:,i) + fwave(:,mw,i)
enddo
else
! No entropy fix needed
do mw=1,mwaves
if (s(mw,i) > 0.d0) then
apdq(:,i) = apdq(:,i) + fwave(:,mw,i)
else
amdq(:,i) = amdq(:,i) + fwave(:,mw,i)
endif
enddo
endif
enddo
endif
end subroutine rp1
! ============================================================================
! ============================================================================
! Eigenvalue routines
! ============================================================================
! ============================================================================
! Eigenspace reconstruction using linear approximation
subroutine linear_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& transonic_wave,wave_correction,s,eig_vec)
implicit none
! I/O
double precision, intent(in) :: h_l(2),h_r(2),u_l(2),u_r(2),b_l,b_r
integer, intent(inout) :: transonic_wave
double precision, intent(inout) :: wave_correction
double precision, intent(inout) :: s(4),eig_vec(4,4)
! Locals
integer :: mw
double precision :: alpha(4),speeds(4,2),gamma_l,gamma_r
double precision :: h_ave(2),u_ave(2)
double precision :: s_l(4),s_r(4),s_ave(4),work_vec(4,4)
! Common blocks
double precision :: rho(2),rho_air,g,r,one_minus_r
double precision :: dry_tolerance
integer :: eigen_method,inundation_method
logical :: entropy_fix
common /cparam/ rho,rho_air,g,r,one_minus_r, &
dry_tolerance,eigen_method,inundation_method,entropy_fix
gamma_l = h_l(2) / h_l(1)
gamma_r = h_r(2) / h_r(1)
! Left state alphas
alpha(1) = 0.5d0*(gamma_l-1.d0+sqrt((gamma_l-1.d0)**2+4.d0*r*gamma_l))
alpha(2) = 0.5d0*(gamma_l-1.d0-sqrt((gamma_l-1.d0)**2+4.d0*r*gamma_l))
! Right state alphas
alpha(3) = 0.5d0*(gamma_r-1.d0-sqrt((gamma_r-1.d0)**2+4.d0*r*gamma_r))
alpha(4) = 0.5d0*(gamma_r-1.d0+sqrt((gamma_r-1.d0)**2+4.d0*r*gamma_r))
! These are calculated seperately due to entropy corrections
! Left state speeds
speeds(1,1) = u_l(1) - sqrt(g*h_l(1)*(1+alpha(1)))
speeds(2,1) = u_l(2) - sqrt(g*h_l(1)*(1+alpha(2)))
speeds(3,1) = u_l(2) + sqrt(g*h_l(1)*(1+alpha(2)))
speeds(4,1) = u_l(1) + sqrt(g*h_l(1)*(1+alpha(1)))
! Right state speeds
speeds(1,2) = u_r(1) - sqrt(g*h_r(1)*(1+alpha(4)))
speeds(2,2) = u_r(2) - sqrt(g*h_r(1)*(1+alpha(3)))
speeds(3,2) = u_r(2) + sqrt(g*h_r(1)*(1+alpha(3)))
speeds(4,2) = u_r(1) + sqrt(g*h_r(1)*(1+alpha(4)))
! Determine wave speeds
transonic_wave = 0
wave_correction = 0.d0
if (entropy_fix) then
transonic_wave = 0
wave_correction = 0.d0
do mw=1,4
! Both speeds right going
if (speeds(mw,1) > 0.d0) then
s(mw) = speeds(mw,2)
! Transonic rarefaction
! s_l < 0.0 < s_r
else if (speeds(mw,2) > 0.d0) then
! Assign base speed for rarefaction, should approximate true speed
s(mw) = 0.5d0 * sum(speeds(mw,:))
transonic_wave = mw
wave_correction = abs(speeds(mw,1)) / (abs(speeds(mw,1)) + abs(speeds(mw,2)))
else
s(mw) = speeds(mw,1)
endif
enddo
else
s(:) = [speeds(1,1),speeds(2,1),speeds(3,2),speeds(4,2)]
endif
eig_vec(1,:) = 1.d0
eig_vec(2,:) = s(:)
eig_vec(3,:) = alpha
eig_vec(4,:) = s(:)*alpha(:)
end subroutine linear_eigen
! ============================================================================
! ============================================================================
! Eigenspace reconstruction using the velocity difference expansion
subroutine velocity_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& transonic_wave,wave_correction,s,eig_vec)
implicit none
! I/O
double precision, intent(in) :: h_l(2),h_r(2),u_l(2),u_r(2),b_l,b_r
integer, intent(inout) :: transonic_wave
double precision, intent(inout) :: wave_correction
double precision, intent(inout) :: s(4),eig_vec(4,4)
! Locals
double precision :: total_depth_l,total_depth_r,mult_depth_l,mult_depth_r
double precision :: alpha(4)
! Common blocks
double precision :: rho(2),rho_air,g,r,one_minus_r
double precision :: dry_tolerance
integer :: eigen_method,inundation_method
logical :: entropy_fix
common /cparam/ rho,rho_air,g,r,one_minus_r, &
dry_tolerance,eigen_method,inundation_method,entropy_fix
total_depth_l = sum(h_l)
total_depth_r = sum(h_r)
mult_depth_l = product(h_l)
mult_depth_r = product(h_r)
s(1) = (h_l(1)*u_l(1)+h_l(2)*u_l(2)) / total_depth_l &
- sqrt(g*total_depth_l)
s(2) = (h_l(2)*u_l(1)+h_l(1)*u_l(2)) / total_depth_l &
- sqrt(g*one_minus_r*mult_depth_l/total_depth_l &
* (1-(u_l(1)-u_l(2))**2/(g*one_minus_r*total_depth_l)))
s(3) = (h_r(2)*u_r(1)+h_r(1)*u_r(2)) / total_depth_r &
+ sqrt(g*one_minus_r*mult_depth_r/total_depth_r &
* (1-(u_r(1)-u_r(2))**2/(g*one_minus_r*total_depth_r)))
s(4) = (h_r(1)*u_r(1)+h_r(2)*u_r(2)) / total_depth_r &
+ sqrt(g*total_depth_r)
! Determine wave speeds
transonic_wave = 0
wave_correction = 0.d0
alpha(1:2) = ((s(1:2) - u_l(1))**2 - g * h_l(1)) / (g*h_l(1))
alpha(3:4) = ((s(3:4) - u_r(1))**2 - g * h_r(1)) / (g*h_r(1))
eig_vec(1,:) = 1.d0
eig_vec(2,:) = s(:)
eig_vec(3,:) = alpha
eig_vec(4,:) = s(:)*alpha(:)
end subroutine velocity_eigen
! ============================================================================
! ============================================================================
! Eigenspace reconstruction using LAPACK
subroutine lapack_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& transonic_wave,wave_correction,s,eig_vec)
implicit none
! I/O
double precision, intent(in) :: h_l(2),h_r(2),u_l(2),u_r(2),b_l,b_r
integer, intent(inout) :: transonic_wave
double precision, intent(inout) :: wave_correction
double precision, intent(inout) :: s(4),eig_vec(4,4)
! Local
integer :: j
double precision :: h_ave(2),u_ave(2)
double precision :: s_l(4),s_r(4),vec_work(4,4)
! Common blocks
double precision :: rho(2),rho_air,g,r,one_minus_r
double precision :: dry_tolerance
integer :: eigen_method,inundation_method
logical :: entropy_fix
common /cparam/ rho,rho_air,g,r,one_minus_r, &
dry_tolerance,eigen_method,inundation_method,entropy_fix
! Solve eigenvalue problem
h_ave(:) = 0.5d0 * (h_l(:) + h_r(:))
u_ave(:) = 0.5d0 * (u_l(:) + u_r(:))
call eval_lapack_eigen(h_ave,u_ave,s,eig_vec)
transonic_wave = 0
wave_correction = 0.d0
if (entropy_fix) then
! Check to see if we may be at a transonic rarefaction
call eval_lapack_eigen(h_l,u_l,s_l,vec_work)
call eval_lapack_eigen(h_r,u_r,s_r,vec_work)
! Check each wave for a transonic problem
do j=1,4
if (s_l(j) < 0.d0 .and. 0.d0 < s_r(j)) then
print *,"Transonic wave detected in wave family ",j,"."
transonic_wave = j
wave_correction = (s_l(j) - s(j)) / (s_l(j) - s_r(j))
endif
enddo
endif
end subroutine lapack_eigen
! ============================================================================
! ============================================================================
! Single layer eigenspace reconstruction
subroutine single_layer_eigen(h_l,h_r,u_l,u_r,b_l,b_r, &
& transonic_wave,wave_correction,s,eig_vec)
implicit none
! I/O
double precision, intent(in) :: h_l(2),h_r(2),u_l(2),u_r(2),b_l,b_r
integer, intent(inout) :: transonic_wave
double precision, intent(inout) :: wave_correction
double precision, intent(inout) :: s(4),eig_vec(4,4)
! Common blocks
double precision :: rho(2),rho_air,g,r,one_minus_r
double precision :: dry_tolerance
integer :: eigen_method,inundation_method
logical :: entropy_fix
common /cparam/ rho,rho_air,g,r,one_minus_r, &
dry_tolerance,eigen_method,inundation_method,entropy_fix
transonic_wave = 0
wave_correction = 0.d0
s(1) = u_l(1) - sqrt(g*h_l(1))
s(2) = 0.d0
s(3) = 0.d0
s(4) = u_r(1) + sqrt(g*h_r(1))
eig_vec(1,:) = 1.d0
eig_vec(2,:) = s(:)
eig_vec(3,:) = 0.d0
eig_vec(4,:) = 0.d0
end subroutine single_layer_eigen
! ============================================================================
! ============================================================================
! Helper routines
! ============================================================================
! ============================================================================
! Evaluate eigenvalues using LAPACK's DGEEV function
subroutine eval_lapack_eigen(h,u,lambda,vec)
implicit none
double precision, intent(in) :: h(2),u(2)
double precision, intent(inout) :: lambda(4),vec(4,4)
integer, parameter :: lwork = 16
integer :: info
double precision :: A(4,4),real_lambda(4),imaginary_lambda(4)
double precision :: empty,work(1,lwork)
! Common blocks
double precision :: rho(2),rho_air,g,r,one_minus_r
double precision :: dry_tolerance
integer :: eigen_method,inundation_method
logical :: entropy_fix
common /cparam/ rho,rho_air,g,r,one_minus_r, &
dry_tolerance,eigen_method,inundation_method,entropy_fix
! Quasi-linear matrix
A(1,:) = [0.d0,1.d0,0.d0,0.d0]
A(2,:) = [-u(1)**2 + g*h(1),2.d0*u(1),g*h(1),0.d0]
A(3,:) = [0.d0,0.d0,0.d0,1.d0]
A(4,:) = [g*r*h(2),0.d0,-u(2)**2 + g*h(2),2.d0*u(2)]
! Call LAPACK eigen solver
call dgeev('N','V',4,A,4,lambda,imaginary_lambda,empty,1,vec,4,work,lwork,info)
if (info < 0) then
info = -info
print "(a,i1,a)","The ",info,"th argument had an illegal value."
stop
else if (info > 0) then
print "(a)","The QR algorithm failed to compute all the"
print "(a)","eigenvalues, and no eigenvectors have been"
print "(a,i1,a)","computed; elements",info,"+1:4 of WR and WI"
print "(a)","contain eigenvalues which have converged."
stop
endif
! Only need to check if there is a positive imaginary eigenvalue as they
! should come in pairs
if (any(imaginary_lambda > 0.d0)) then
print "(a,4d16.8)","Imaginary eigenvalues computed: ",imaginary_lambda
stop
endif
end subroutine eval_lapack_eigen
! ============================================================================
! ============================================================================
! Solve R beta = delta using LAPACK's DGESV function
function eval_lapack_solve(eig_vec,delta) result(beta)
implicit none
! Matrix and RHS
double precision, intent(in) :: eig_vec(4,4), delta(4)
! Output
double precision :: beta(4)
! Local storage
double precision :: R(4,4)
integer :: ipiv(4),info,i,m
! Move input data into work arrays
R = eig_vec
beta = delta
! Call LAPACK dense linear solver
call dgesv(4,1,R,4,ipiv,beta,4,info)
if (.not.(info == 0)) then
print *, "Error solving R beta = delta, ",info
print *, " R=",(eig_vec(1,m),m=1,4)
do i=2,4
print *, " ",(eig_vec(i,m),m=1,4)
enddo
print *, " delta=",(delta(m),m=1,4)
stop
endif
end function eval_lapack_solve
! ============================================================================
! ============================================================================
! Calculate wind drag coefficient
double precision function wind_drag(wind_speed)
implicit none
! Input
double precision, intent(in) :: wind_speed
if (wind_speed <= 11.d0) then
wind_drag = 1.2d0
else if ((wind_speed > 11.d0).and.(wind_speed <= 25.d0)) then
wind_drag = 0.49d0 + 0.065d0 * wind_speed
else
wind_drag = 0.49 + 0.065d0 * 25.d0
endif
wind_drag = wind_drag * 10.d-3
end function wind_drag