-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathrpn2_riemann_solver_adjoint.f90
executable file
·327 lines (286 loc) · 12.7 KB
/
rpn2_riemann_solver_adjoint.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
! -------------------------------- NUMERICAL FLUXES --------------------------------
subroutine rpn2(ixy,maxm,meqn,mwaves,maux,mbc,mx,ql,qr,auxl,auxr,wave,s,amdq,apdq)
! =====================================================
! INPUT
! ixy: dimension, ixy == 1 -> x-dimension
! ixy == 2 -> y-dimension
! maxm: TOdo: ??
! meqn: the number of equations in the system
! mwaves: the number of waves in each Riemann solution
! maux: the number of auxiliary variables
! mbc: number of boundary(ghost) cells
! mx: the number of grid cells
! ql: contains the state vector at the left edge of each cell
! qr: contains the state vector at the right edge of each cell
! NOTE: the i'th Riemann problem has left state qr(:, i-1) and right state ql(:,i)
! NOTE: as we do not use reconstruction the left and right state of each cell should be the same
! auxl: contains the aux variables at the left edge of each cell
! auxr: contains the aux variables at the right edge of each cell
! NOTE: as we do not use reconstruction the left and right aux variables of each cell should be the same
! OUTPUT
! wave: contains the waves
! s: contains the speeds
! amdq: flux difference leftgoing part (A- Delta q)
! apdq: flux difference rightgoing part (A+ Delta q)
! NOTE: flux difference: f(qr(i-1)) - f(ql(i))
! ====================================================
implicit none
double precision :: flux_precision, condition
common /cparam/ flux_precision
integer, intent(in) :: ixy
integer, intent(in) :: maxm
integer, intent(in) :: meqn
integer, intent(in) :: mwaves
integer, intent(in) :: maux
integer, intent(in) :: mbc
integer, intent(in) :: mx
double precision, dimension(meqn, mwaves, 1-mbc:maxm+mbc), intent(out) :: wave
double precision, dimension(mwaves, 1-mbc:maxm+mbc), intent(out) :: s
double precision, dimension(meqn, 1-mbc:maxm+mbc), intent(in) :: ql
double precision, dimension(meqn, 1-mbc:maxm+mbc), intent(in) :: qr
double precision, dimension(meqn, 1-mbc:maxm+mbc), intent(out) :: apdq
double precision, dimension(meqn, 1-mbc:maxm+mbc), intent(out) :: amdq
double precision, dimension(maux, 1-mbc:maxm+mbc), intent(in) :: auxl
double precision, dimension(maux, 1-mbc:maxm+mbc), intent(in) :: auxr
double precision, dimension(3,3) :: jacobian, eigvecs, inveigvecs
double precision, dimension(3) :: eigvals, q_mean, dlambda, alpha
logical :: cplx
integer :: i, m
do i = 1, mx+2
q_mean = 0.5d0*(auxl(4:6, i) + auxr(4:6, i-1))
dlambda = ql(1:3, i) - qr(1:3, i-1)
if ((abs(q_mean(1)) > flux_precision)) THEN
if (ixy == 1) THEN
call calcJacobianX(q_mean, jacobian)
ELSE
call calcJacobianY(q_mean, jacobian)
end if
jacobian = -transpose(jacobian)
call calcEigen(jacobian, eigvecs, eigvals, cplx)
if (cplx) then
read *, q_mean
end if
! jacobian is overwritten!
call calcInverse(eigvecs, inveigvecs)
call decomposeJump(inveigvecs, dlambda, alpha)
call calcWaves(alpha, eigvecs, wave(:, :, i))
s(:, i) = eigvals
call calcAMDQ(inveigvecs, eigvals, eigvecs, dlambda, amdq(:, i))
call calcAPDQ(inveigvecs, eigvals, eigvecs, dlambda, apdq(:, i))
ELSE
s(:,i) = 1.d0
wave(:, :, i) = 0.d0
amdq(:,i) = 0.d0
apdq(:,i) = 0.d0
end if
end do
end subroutine rpn2
! -------------------------------- XI APPROXIMATION --------------------------------
subroutine calcChi(abs_a, chi)
implicit none
double precision, intent(in) :: abs_a
double precision, intent(out) :: chi
double precision :: a_0, a_2, a_4, a_6
double precision :: b_0, b_2
a_0 = 0.621529d+0
a_2 = 0.348509d+0
a_4 = -0.139318d+0
a_6 = 0.720371d+0
b_0 = 1.87095d+0
b_2 = -1.32002d+0
chi = (abs_a**6.d0*a_6 + abs_a**4.d0*a_4 + abs_a**2.d0*a_2 + a_0)/(abs_a**4.d0 + abs_a**2.d0*b_2 + b_0)
chi = max(min(chi, 1.d0), 1.d0/3.d0)
end subroutine calcChi
subroutine calcDChi(abs_a, dchi)
implicit none
double precision, intent(in) :: abs_a
double precision, intent(out) :: dchi
double precision :: a_0, a_2, a_4, a_6
double precision :: b_0, b_2
a_0 = 0.621529d+0
a_2 = 0.348509d+0
a_4 = -0.139318d+0
a_6 = 0.720371d+0
b_0 = 1.87095d+0
b_2 = -1.32002d+0
dchi = (-4.d0*abs_a**3.d0 - 2.d0*abs_a*b_2)*(a_0 + a_2*abs_a**2.d0 + a_4*abs_a**4.d0 + a_6*abs_a**6.d0)/(abs_a**4.d0 + abs_a**2.d0*b_2 + b_0)**2.d0 + (2.d0*a_2*abs_a + 4.d0*a_4*abs_a**3.d0 + 6.d0*a_6*abs_a**5.d0)/(abs_a**4.d0 + abs_a**2.d0*b_2 + b_0)
dchi = max(min(dchi, 2.d0), 0.d0)
end subroutine calcDChi
subroutine calcChiTotal(q, chi)
implicit none
double precision, dimension(3), intent(in) :: q
double precision, intent(out) :: chi
double precision :: abs_a
abs_a = sqrt(q(2)**2.d0 + q(3)**2.d0)/q(1)
call calcChi(abs_a, chi)
end subroutine calcChiTotal
subroutine calcDChiTotal(q, dchi)
implicit none
double precision, dimension(3), intent(in) :: q
double precision, dimension(3), intent(out) :: dchi
double precision :: dchi_, abs_a
integer :: m
abs_a = sqrt(q(2)**2.d0 + q(3)**2.d0)/q(1)
call calcDChi(abs_a, dchi_)
dchi(1) = dchi_ * (-sqrt(q(2)**2.d0 + q(3)**2.d0)/q(1)**2.d0)
dchi(2) = dchi_ * (q(2)/(q(1)*sqrt(q(2)**2.d0 + q(3)**2.d0)))
dchi(3) = dchi_ * (q(3)/(q(1)*sqrt(q(2)**2.d0 + q(3)**2.d0)))
do m = 1,3
if (isnan(dchi(m))) then
! print *, "dchi(", m, ") = nan"
dchi(m) = 0.d0
end if
end do
end subroutine calcDChiTotal
! -------------------------------- JACOBIAN --------------------------------
subroutine calcJacobianX(q, jacF)
implicit none
double precision, dimension(3), intent(in) :: q
double precision, dimension(3,3), intent(out) :: jacF
double precision, dimension(3) :: dchi
double precision :: chi
integer :: i,m
call calcChiTotal(q, chi)
call calcDChiTotal(q, dchi)
jacF(1,1) = 0.d0
jacF(1,2) = 1.d0
jacF(1,3) = 0.d0
jacF(2,1) = q(1)*(3.d0*q(2)**2.d0*dchi(1)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - dchi(1)/2.d0) + q(2)**2.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - chi/2.d0 + 1.d0/2
jacF(2,2) = q(1)*(-4.d0*q(2)**3.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(2)**2.d0*dchi(2)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + 2.d0*q(2)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - dchi(2)/2.d0)
jacF(2,3) = q(1)*(-4.d0*q(2)**2.d0*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(2)**2.d0*dchi(3)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - dchi(3)/2.d0)
jacF(3,1) = 3.d0*q(1)*q(2)*q(3)*dchi(1)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + q(2)*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)
jacF(3,2) = -4.d0*q(1)*q(2)**2.d0*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(1)*q(2)*q(3)*dchi(2)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + q(1)*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)
jacF(3,3) = -4.d0*q(1)*q(2)*q(3)**2.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(1)*q(2)*q(3)*dchi(3)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + q(1)*q(2)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)
do i = 1,3
do m = 1,3
if (isnan(jacF(i,m))) then
! print *, "nan", i, m
! print *, q
jacF(i,m) = 0.d0
end if
end do
end do
end subroutine calcJacobianX
subroutine calcJacobianY(q, jacF)
implicit none
double precision, dimension(3), intent(in) :: q
double precision, dimension(3,3), intent(out) :: jacF
double precision, dimension(3) :: dchi
double precision :: chi
integer :: i,m
call calcChiTotal(q, chi)
call calcDChiTotal(q, dchi)
jacF(1,1) = 0.d0
jacF(1,2) = 0.d0
jacF(1,3) = 1.d0
jacF(2,1) = 3.d0*q(1)*q(2)*q(3)*dchi(1)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + q(2)*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)
jacF(2,2) = -4.d0*q(1)*q(2)**2.d0*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(1)*q(2)*q(3)*dchi(2)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + q(1)*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)
jacF(2,3) = -4.d0*q(1)*q(2)*q(3)**2.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(1)*q(2)*q(3)*dchi(3)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + q(1)*q(2)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)
jacF(3,1) = q(1)*(3.d0*q(3)**2.d0*dchi(1)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - dchi(1)/2.d0) + q(3)**2.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - chi/2.d0 + 1.d0/2
jacF(3,2) = q(1)*(-4.d0*q(2)*q(3)**2.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(3)**2.d0*dchi(2)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - dchi(2)/2.d0)
jacF(3,3) = q(1)*(-4.d0*q(3)**3.d0*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0)**2.d0 + 3.d0*q(3)**2.d0*dchi(3)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) + 2.d0*q(3)*(3.d0*chi - 1.d0)/(2.d0*q(2)**2.d0 + 2.d0*q(3)**2.d0) - dchi(3)/2.d0)
do i = 1,3
do m = 1,3
if (isnan(jacF(i,m))) then
! print *, "nan"
! print *, q
jacF(i,m) = 0.d0
end if
end do
end do
end subroutine calcJacobianY
subroutine calcEigen(mat, evcR, evlR, cplx)
double precision, dimension(3,3), intent(in) :: mat
double precision, dimension(3,3), intent(out) :: evcR
logical, intent(out) :: cplx
double precision, dimension(3,3) :: dummy
double precision, dimension(3), intent(out) :: evlR
double precision, dimension(3) :: evlI
integer :: LWORK, INFO, N
double precision, dimension(102) :: WORK
N = 3
LWORK = 102
cplx = .false.
call DGEEV( 'N', 'V', N, mat, N, evlR, evlI, dummy, N, evcR, N, WORK, LWORK, INFO )
if (ANY(abs(evlI) > 1.d-30)) THEN
! print *, "complex eigenvalue"
cplx = .true.
!read *, s
end if
! if (INFO /= 0) THEN
! print *, "DGEEV", INFO
! end if
end subroutine calcEigen
subroutine calcInverse(mat, inv)
implicit none
double precision, dimension(3,3), intent(in) :: mat
double precision, dimension(3,3), intent(out) :: inv
integer, dimension(3) :: IPIV
integer :: LWORK, INFO, N
double precision, dimension(102) :: WORK
N = 3
LWORK = 102
inv = mat
call DGETRF(N, N, inv, N, IPIV, INFO)
! if (INFO /= 0) THEN
! print *, "DGETRF", INFO
! end if
call DGETRI(N, inv, N, IPIV, WORK, LWORK, INFO)
! if (INFO /= 0) THEN
! print *, "DGETRI", INFO
! end if
end subroutine calcInverse
subroutine calcAP(eigvals, Aplus)
implicit none
double precision, dimension(3, 3), intent(out) :: Aplus
double precision, dimension(3), intent(in) :: eigvals
integer :: i
Aplus = 0.d0
do i = 1,3
Aplus(i,i) = max(0.d0, eigvals(i))
end do
end subroutine calcAP
subroutine calcAM(eigvals, Aminus)
implicit none
double precision, dimension(3, 3), intent(out) :: Aminus
double precision, dimension(3), intent(in) :: eigvals
integer :: i
Aminus = 0.d0
do i = 1,3
Aminus(i,i) = min(0.d0, eigvals(i))
end do
end subroutine calcAM
subroutine decomposeJump(inveigvecs, dq, alpha)
implicit none
double precision, dimension(3,3), intent(in) :: inveigvecs
double precision, dimension(3), intent(in) :: dq
double precision, dimension(3), intent(out) :: alpha
alpha = matmul(inveigvecs, dq)
end subroutine decomposeJump
subroutine calcWaves(alpha, eigvecs, wave)
implicit none
double precision, dimension(3), intent(in) :: alpha
double precision, dimension(3, 3), intent(in) :: eigvecs
double precision, dimension(3, 3), intent(out) :: wave
integer :: i
do i = 1, 3
wave(:, i) = eigvecs(:, i)*alpha(i)
!print *, i, "th wave", eigvecs(:, i)*alpha(i)
end do
end subroutine calcWaves
subroutine calcAPDQ(inveigvecs, eigvals, eigvecs, dq, apdq)
implicit none
double precision, dimension(3, 3), intent(in) :: inveigvecs, eigvecs
double precision, dimension(3), intent(in) :: dq, eigvals
double precision, dimension(3), intent(out) :: apdq
!amdq= matmul(matmul(inveigvecs, matmul(APM, eigvecs)), dq)
apdq = matmul(eigvecs, max(eigvals, 0.d0)*matmul(inveigvecs, dq))
end subroutine calcAPDQ
subroutine calcAMDQ(inveigvecs, eigvals, eigvecs, dq, amdq)
implicit none
double precision, dimension(3, 3), intent(in) :: inveigvecs, eigvecs
double precision, dimension(3), intent(in) :: dq, eigvals
double precision, dimension(3), intent(out) :: amdq
!amdq= matmul(matmul(inveigvecs, matmul(APM, eigvecs)), dq)
amdq = matmul(eigvecs, min(eigvals, 0.d0)*matmul(inveigvecs, dq))
end subroutine calcAMDQ