-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmisc.f90
287 lines (250 loc) · 9.47 KB
/
misc.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
module misc
use constant
use declaration
use grid
implicit none
contains
subroutine prm2cons
integer(prec) :: i, j
!converting the flow variables to conservative variables
!$omp parallel do &
!$omp default(none) &
!$omp private(i, j) &
!$omp shared(rho, u, v, p, q, nx, ny, g2)
do j = gc+1, ny-gc
do i = gc+1, nx-gc
q(1, i, j) = rho(i, j)
q(2, i, j) = rho(i, j)*u(i, j)
q(3, i, j) = rho(i, j)*v(i, j)
q(4, i, j) = p(i, j)/g2 + half*rho(i, j)*(u(i, j)**2 + v(i, j)**2)
end do
end do
!$omp end parallel do
end subroutine
subroutine cons2prm
integer(prec) :: i, j
!$omp parallel do &
!$omp default(none) &
!$omp private(i, j) &
!$omp shared(rho, u, v, p, c, h, t, q, nx, ny, xc, yc, g2, g1, r, counter)
do j = gc+1, ny-gc
do i = gc+1, nx-gc
rho(i, j) = q(1, i, j)
u(i, j) = q(2, i, j)/rho(i, j)
v(i, j) = q(3, i, j)/rho(i, j)
p(i, j) = (q(4, i, j) - half*rho(i, j)*(u(i, j)**2 + v(i, j)**2))*g2
c(i, j) = dsqrt(g1*p(i, j)/rho(i, j))
h(i, j) = (p(i, j)/g2 + p(i, j) + half*rho(i, j)*(u(i, j)**2 + v(i, j)**2))/rho(i,j)
t(i,j) = p(i,j)/(R*rho(i,j))
if(isnan(rho(i,j))) then
print *, counter, i, j, xc(i, j), yc(i,j)
print *, q(:, i, j)
STOP
endif
! if(isnan(p(i,j)) .or. isnan(rho(i,j)) .or. isnan(c(i,j))) then
! print *, q(:,i,j), i, j, xc(i, j), yc(i, j), p(i,j), rho(i,j)
! stop
! endif
end do
end do
!$omp end parallel do
end subroutine
!function to calculate the roe averages
real(prec) function roeavg(r1, r2, u1, u2)
real(prec), intent(in) :: r1, r2, u1, u2
real(prec) :: d1, d2
! d = sqrt(r2/r1)
! roeavg = (u1 + d*u2)/(1 + d)
d1 = dsqrt(r1)
d2 = dsqrt(r2)
roeavg = (u1*d1 + u2*d2)/(d1 + d2)
end function
!calculation of leigenvec given nx, ny, and roe averages
subroutine lefteigenvector(u_, v_, c_, normalx, normaly, leigenvec)
real(prec), intent(in) :: u_, v_, c_, normalx, normaly
real(prec) :: un, ut, ek ! normal values of velocities will be used in calc
real(prec), dimension(4, 4), intent(out) :: leigenvec
!calculating the normal velocities
un = u_*normalx + v_*normaly
ut = v_*normalx - u_*normaly
ek = half*(u_**2 + v_**2) ! kinetic energy
leigenvec(1, 1) = half*(g2*ek + c_*un)/c_**2
leigenvec(1, 2) = -half*(g2*u_ + c_*normalx)/c_**2
leigenvec(1, 3) = -half*(g2*v_ + c_*normaly)/c_**2
leigenvec(1, 4) = half*g2/c_**2
leigenvec(2, 1) = (c_**2 - g2*ek)/c_**2
leigenvec(2, 2) = (g2*u_)/c_**2
leigenvec(2, 3) = (g2*v_)/c_**2
leigenvec(2, 4) = -g2/c_**2
leigenvec(3, 1) = ut
leigenvec(3, 2) = normaly
leigenvec(3, 3) = -normalx
leigenvec(3, 4) = zero
leigenvec(4, 1) = half*(g2*ek - c_*un)/c_**2
leigenvec(4, 2) = half*(c_*normalx - g2*u_)/c_**2
leigenvec(4, 3) = half*(c_*normaly - g2*v_)/c_**2
leigenvec(4, 4) = half*g2/c_**2
end subroutine
subroutine righteigenvector(u_, v_, c_, h_, normalx, normaly, reigenvec)
real(prec), intent(in) :: u_, v_, c_, h_, normalx, normaly
real(prec) :: un, ut, h0, ek
real(prec), dimension(4, 4) :: reigenvec
!calculation of normal velocities
un = u_*normalx + v_*normaly
ut = v_*normalx - u_*normaly
ek = half*(u_**2 + v_**2)
h0 = c_**2/g2 + ek
reigenvec(1, 1) = 1
reigenvec(2, 1) = u_ - c_*normalx
reigenvec(3, 1) = v_ - c_*normaly
reigenvec(4, 1) = h0 - c_*un
reigenvec(1, 2) = 1
reigenvec(2, 2) = u_
reigenvec(3, 2) = v_
reigenvec(4, 2) = ek
reigenvec(1, 3) = 0
reigenvec(2, 3) = normaly
reigenvec(3, 3) = -normalx
reigenvec(4, 3) = -ut
reigenvec(1, 4) = 1
reigenvec(2, 4) = u_ + c_*normalx
reigenvec(3, 4) = v_ + c_*normaly
reigenvec(4, 4) = h0 + c_*un
end subroutine
subroutine fluxcalc(rhof, uf, vf, pf, normalx, normaly, f)
real(prec), intent(in) :: rhof, uf, vf, pf, normalx, normaly
real(prec) :: un
real(prec), dimension(4), intent(out) :: f
un = uf*normalx + vf*normaly
f(1) = rhof*un
f(2) = rhof*uf*un + pf*normalx
f(3) = rhof*vf*un + pf*normaly
f(4) = (pf/g2 + half*rhof*(uf**2 + vf**2) + pf)*un
end subroutine
subroutine eigenval(ur, vr, cr, normalx, normaly, el1, el2, el3, el4)
real(prec), intent(in) :: ur, vr, cr, normalx, normaly
real(prec) :: un
real(prec), intent(out) :: el1, el2, el3, el4
un = ur*normalx + vr*normaly
el1 = un - cr
el2 = un
el3 = un
el4 = un + cr
end subroutine
subroutine timestep
integer :: i, j
real(prec) :: qplus, asound
if(cfl) then
dtmin = 1.0d8
!$omp parallel do reduction (min:dtmin) &
!$omp default(none) &
!$omp private(i, j, qplus, asound) &
!$omp shared(u, v, p, rho, g1, cfl_no, dl, dt_cell, nx, ny)
do j = gc+1, ny-gc
do i = gc+1, nx-gc
asound = dsqrt(g1*p(i, j)/rho(i, j))
qplus = dsqrt(u(i,j)*u(i,j) + v(i,j)*v(i,j)) + asound
dt_cell(i, j) = cfl_no*dl(i,j)/qplus
dtmin = dmin1(dtmin, dt_cell(i, j))
enddo
enddo
!$omp end parallel do
if(timeaccurate) then
dt_cell = dtmin
endif
else
dt_cell = dt
dtmin = dt
endif
end subroutine timestep
real(prec) function minmod(a, b)
real(prec), intent(in) :: a, b
!real(prec), intent(out) :: minmod
if((a > zero .and. b > zero) .or. (a < zero .and. b < zero)) then
if(abs(a) > abs(b)) then
minmod = b
else
minmod = a
end if
else
minmod = zero
end if
end function
! real(prec) function psi(x)
! real(prec), intent(in) :: x
! !real(prec), intent(out) :: psi
! psi = sqrt(dd + x**2)
! end function
subroutine output()
integer(prec) :: i, j
open (unit = unit_id, status = 'replace', file = filename, form = 'formatted')
! !headers for tecplot file
! write (unit_id, *) "variables = x, y, rho, u ,v, p, m"
! write (unit_id, *) "zone j=", ny-2*gc, "i=",nx-2*gc, "f=point"
! !writing the data
! do j = gc+1, ny-gc
! do i = gc+1, nx-gc
! asound = dsqrt(g1*p(i, j)/rho(i, j))
! M_num = dsqrt(u(i, j)*u(i, j) + v(i, j)*v(i, j))/asound
! write(unit_id, '(6f15.5)') xc(i, j), yc(i, j), rho(i, j), u(i, j), v(i, j), p(i, j), M_num
! end do
! end do
!headers for tecplot files
write (unit_id, *) "VARIABLES = x, y, rho, u, v, p, t"
write (unit_id, *) "ZONE I =", imax, ",J =", jmax, ", DATAPACKING=BLOCK, VARLOCATION=([3,4,5,6,7]=CELLCENTERED)"
write (unit_id, *) "SOLUTIONTIME=", time
do j = 1, jmax
do i = 1, imax
write(unit_id, *) x(gc+i, gc+j)
enddo
enddo
do j = 1, jmax
do i = 1, imax
write(unit_id, *) y(gc+i, gc+j)
enddo
enddo
do j = 1, jmax-1
do i = 1, imax-1
write(unit_id, *) rho(gc+i, gc+j)
enddo
enddo
do j = 1, jmax-1
do i = 1, imax-1
write(unit_id, *) u(gc+i, gc+j)
enddo
enddo
do j = 1, jmax-1
do i = 1, imax-1
write(unit_id, *) v(gc+i, gc+j)
enddo
enddo
do j = 1, jmax-1
do i = 1, imax-1
write(unit_id, *) p(gc+i, gc+j)
enddo
enddo
do j = 1, jmax-1
do i = 1, imax-1
write(unit_id, *) t(gc+i, gc+j)
enddo
enddo
! write(unit_id, *) ((x(gc+i,gc+j), i = 1, imax), j= 1, jmax)
! write(unit_id, *) ((y(gc+i,gc+j), i = 1, imax), j= 1, jmax)
! write(unit_id, *) ((rho(gc+i,gc+j), i = 1, imax-1), j= 1, jmax-1)
! write(unit_id, *) ((u(gc+i,gc+j), i = 1, imax-1), j= 1, jmax-1)
! write(unit_id, *) ((v(gc+i,gc+j), i = 1, imax-1), j= 1, jmax-1)
! write(unit_id, *) ((p(gc+i,gc+j), i = 1, imax-1), j= 1, jmax-1)
! write(unit_id, *) ((t(gc+i,gc+j), i = 1, imax-1), j= 1, jmax-1)
! ! write(unit_id, *) ((t(gc+i,gc+j), i = 1, imax-1), j= 1, jmax-1)
close (unit_id)
end subroutine output
subroutine restart()
integer :: i, j
open(unit = 11, status='replace', file=restartfile, form='formatted')
write(11, *) ((q(1,i,j), i = gc+1, nx-gc), j= gc+1, ny-gc)
write(11, *) ((q(2,i,j), i = gc+1, nx-gc), j= gc+1, ny-gc)
write(11, *) ((q(3,i,j), i = gc+1, nx-gc), j= gc+1, ny-gc)
write(11, *) ((q(4,i,j), i = gc+1, nx-gc), j= gc+1, ny-gc)
close(11)
end subroutine restart
end module misc