-
Notifications
You must be signed in to change notification settings - Fork 12
/
Copy pathppm.f90
385 lines (280 loc) · 10.9 KB
/
ppm.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
!
! This program may be freely redistributed under the
! condition that the copyright notices (including this
! entire header) are not removed, and no compensation
! is received through use of the software. Private,
! research, and institutional use is free. You may
! distribute modified versions of this code UNDER THE
! CONDITION THAT THIS CODE AND ANY MODIFICATIONS MADE
! TO IT IN THE SAME FILE REMAIN UNDER COPYRIGHT OF THE
! ORIGINAL AUTHOR, BOTH SOURCE AND OBJECT CODE ARE
! MADE FREELY AVAILABLE WITHOUT CHARGE, AND CLEAR
! NOTICE IS GIVEN OF THE MODIFICATIONS. Distribution
! of this code as part of a commercial system is
! permissible ONLY BY DIRECT ARRANGEMENT WITH THE
! AUTHOR. (If you are not directly supplying this
! code to a customer, and you are instead telling them
! how they can obtain it for free, then you are not
! required to make any arrangement with me.)
!
! Disclaimer: Neither I nor: Columbia University, the
! National Aeronautics and Space Administration, nor
! the Massachusetts Institute of Technology warrant
! or certify this code in any way whatsoever. This
! code is provided "as-is" to be used at your own risk.
!
!
!
! PPM.f90: 1d slope-limited, piecewise parabolic recon.
!
! Darren Engwirda
! 03-Jul-2020
! d [dot] engwirda [at] gmail [dot] com
!
!
! P. Colella and PR. Woodward, The Piecewise Parabolic
! Method (PPM) for gas-dynamical simulations., J. Comp.
! Phys., 54 (1), 1984, 174-201,
! https://doi.org/10.1016/0021-9991(84)90143-8
!
pure subroutine ppm(npos,nvar,ndof,delx, &
& fdat,fhat,edge,oscl, &
& dmin,ilim,wlim,halo)
!
! NPOS no. edges over grid.
! NVAR no. state variables.
! NDOF no. degrees-of-freedom per grid-cell.
! DELX grid-cell spacing array. LENGTH(DELX) == +1 if
! spacing is uniform .
! FDAT grid-cell moments array. FDAT is an array with
! SIZE = NDOF-by-NVAR-by-NPOS-1 .
! FHAT grid-cell re-con. array. FHAT is an array with
! SIZE = MDOF-by-NVAR-by-NPOS-1 .
! EDGE edge-centred interp. for function-value. EDGE
! is an array with SIZE = NVAR-by-NPOS .
! OSCL grid-cell oscil. dof.'s. OSCL is an array with
! SIZE = +2 -by-NVAR-by-NPOS-1 .
! DMIN min. grid-cell spacing thresh .
! ILIM cell slope-limiting selection .
! WLIM wall slope-limiting selection .
! HALO width of re-con. stencil, symmetric about mid. .
!
implicit none
!------------------------------------------- arguments !
integer , intent(in) :: npos,nvar,ndof
real(kind=dp), intent(in) :: dmin
real(kind=dp), intent(out) :: fhat(:,:,:)
real(kind=dp), intent(in) :: oscl(:,:,:)
real(kind=dp), intent(in) :: delx(:)
real(kind=dp), intent(in) :: fdat(:,:,:)
real(kind=dp), intent(in) :: edge(:,:)
integer , intent(in) :: ilim,wlim,halo
!------------------------------------------- variables !
integer :: ipos,ivar,iill,iirr,head,tail
real(kind=dp) :: ff00,ffll,ffrr,hh00,hhll,hhrr
integer :: mono
real(kind=dp) :: fell,ferr
real(kind=dp) :: dfds(-1:+1)
real(kind=dp) :: wval(+1:+2)
real(kind=dp) :: uhat(+1:+3)
real(kind=dp) :: lhat(+1:+3)
head = +1; tail = npos - 1
if (npos.eq.2) then
!----- default to reduced order if insufficient points !
do ivar = +1, nvar
fhat(1,ivar,+1) = &
& fdat(1,ivar,+1)
fhat(2,ivar,+1) = 0.d0
fhat(3,ivar,+1) = 0.d0
end do
end if
if (npos.le.2) return
!------------------- reconstruct function on each cell !
uhat = +0.d+0
lhat = +0.d+0
do ipos = +1 , npos-1
iill = max(head,ipos-1)
iirr = min(tail,ipos+1)
do ivar = +1 , nvar-0
!----------------------------- cell mean + edge values !
ff00 = fdat(1,ivar,ipos)
ffll = fdat(1,ivar,iill)
ffrr = fdat(1,ivar,iirr)
fell = edge(ivar,ipos+0)
ferr = edge(ivar,ipos+1)
!----------------------------- calc. LL/00/RR gradient !
if (size(delx).gt.+1) then
hh00 = delx(ipos)
hhll = delx(iill)
hhrr = delx(iirr)
call plsv (ffll,hhll,ff00, &
& hh00,ffrr,hhrr, &
& dfds)
else
call plsc (ffll,ff00,ffrr, &
& dfds)
end if
!----------------------------- calc. cell-wise profile !
select case(ilim)
case (null_limit)
!----------------------------- calc. unlimited profile !
call ppmfn(ff00,ffll,ffrr, &
& fell,ferr,dfds, &
& uhat,lhat,mono)
!----------------------------- pref. unlimited profile !
wval(1) = +1.d+0
wval(2) = +0.d+0
case (mono_limit)
!----------------------------- calc. monotonic profile !
call ppmfn(ff00,ffll,ffrr, &
& fell,ferr,dfds, &
& uhat,lhat,mono)
!----------------------------- pref. monotonic profile !
wval(1) = +0.d+0
wval(2) = +1.d+0
case (weno_limit)
!----------------------------- calc. unlimited profile !
call ppmfn(ff00,ffll,ffrr, &
& fell,ferr,dfds, &
& uhat,lhat,mono)
if (mono.gt.+0) then
!----------------------------- calc. WENO-type weights !
call wenoi(npos,delx,oscl, &
& ipos,ivar,halo, &
& wlim,wval)
else
!----------------------------- pref. unlimited profile !
wval(1) = +1.d+0
wval(2) = +0.d+0
end if
end select
!----------------------------- blend "null" and "mono" !
fhat(1,ivar,ipos) = &
& wval(1) * uhat(1) + &
& wval(2) * lhat(1)
fhat(2,ivar,ipos) = &
& wval(1) * uhat(2) + &
& wval(2) * lhat(2)
fhat(3,ivar,ipos) = &
& wval(1) * uhat(3) + &
& wval(2) * lhat(3)
end do
end do
return
end subroutine
!--------- assemble piecewise parabolic reconstruction !
pure subroutine ppmfn(ff00,ffll,ffrr,fell,&
& ferr,dfds,uhat,lhat,&
& mono)
!
! FF00 centred grid-cell mean.
! FFLL left -biased grid-cell mean.
! FFRR right-biased grid-cell mean.
! FELL left -biased edge interp.
! FERR right-biased edge interp.
! DFDS piecewise linear gradients in local co-ord.'s.
! DFDS(+0) is a centred, slope-limited estimate,
! DFDS(-1), DFDS(+1) are left- and right-biased
! estimates (unlimited).
! UHAT unlimited PPM reconstruction coefficients .
! LHAT monotonic PPM reconstruction coefficients .
! MONO slope-limiting indicator, MONO > +0 if some
! limiting has occured .
!
implicit none
!------------------------------------------- arguments !
real(kind=dp), intent(in) :: ff00
real(kind=dp), intent(in) :: ffll,ffrr
real(kind=dp), intent(inout) :: fell,ferr
real(kind=dp), intent(in) :: dfds(-1:+1)
real(kind=dp), intent(out) :: uhat(+1:+3)
real(kind=dp), intent(out) :: lhat(+1:+3)
integer , intent(out) :: mono
!------------------------------------------- variables !
real(kind=dp) :: turn,fmid,fmin,fmax
mono = 0
!-------------------------------- "null" slope-limiter !
uhat( 1 ) = &
& + (3.0d+0 / 2.0d+0) * ff00 &
& - (1.0d+0 / 4.0d+0) *(ferr+fell)
uhat( 2 ) = &
& + (1.0d+0 / 2.0d+0) *(ferr-fell)
uhat( 3 ) = &
& - (3.0d+0 / 2.0d+0) * ff00 &
& + (3.0d+0 / 4.0d+0) *(ferr+fell)
!-------------------------------- "mono" slope-limiter !
if((ffrr - ff00) * &
& (ff00 - ffll) .lt. 0.d+0) then
!----------------------------------- "flatten" extrema !
mono = +1
lhat(1) = ff00
lhat(2) = 0.d0
lhat(3) = 0.d0
return
end if
!----------------------------------- limit edge values !
if((ffll - fell) * &
& (fell - ff00) .le. 0.d+0) then
mono = +1
fell = ff00 - dfds(0)
end if
if((ffrr - ferr) * &
& (ferr - ff00) .le. 0.d+0) then
mono = +1
ferr = ff00 + dfds(0)
end if
!----------------------------------- update ppm coeff. !
lhat( 1 ) = &
& + (3.0d+0 / 2.0d+0) * ff00 &
& - (1.0d+0 / 4.0d+0) *(ferr+fell)
lhat( 2 ) = &
& + (1.0d+0 / 2.0d+0) *(ferr-fell)
lhat( 3 ) = &
& - (3.0d+0 / 2.0d+0) * ff00 &
& + (3.0d+0 / 4.0d+0) *(ferr+fell)
!----------------------------------- limit cell values !
if (abs(lhat(3)) .gt. &
& abs(lhat(2))*.5d+0) then
turn = -0.5d+0 * lhat(2) &
& / lhat(3)
fmid = lhat(1) + lhat(2) * turn**1 &
& + lhat(3) * turn**2
fmin = min(ffll,ffrr)
fmax = max(ffll,ffrr)
if ((fmid .le. fmax)&
& .and.(fmid .ge. fmin)) return
if ((turn .ge. -1.d+0)&
& .and.(turn .le. +0.d+0)) then
mono = +2
!--------------------------- push TURN onto lower edge !
fell = ffll
ferr = +3.0d+0 * ff00 &
& -2.0d+0 * fell
lhat( 1 ) = &
& + (3.0d+0 / 2.0d+0) * ff00 &
& - (1.0d+0 / 4.0d+0) *(ferr+fell)
lhat( 2 ) = &
& + (1.0d+0 / 2.0d+0) *(ferr-fell)
lhat( 3 ) = &
& - (3.0d+0 / 2.0d+0) * ff00 &
& + (3.0d+0 / 4.0d+0) *(ferr+fell)
else &
& if ((turn .gt. +0.d+0)&
& .and.(turn .le. +1.d+0)) then
mono = +2
!--------------------------- push TURN onto upper edge !
ferr = ffrr
fell = +3.0d+0 * ff00 &
& -2.0d+0 * ferr
lhat( 1 ) = &
& + (3.0d+0 / 2.0d+0) * ff00 &
& - (1.0d+0 / 4.0d+0) *(ferr+fell)
lhat( 2 ) = &
& + (1.0d+0 / 2.0d+0) *(ferr-fell)
lhat( 3 ) = &
& - (3.0d+0 / 2.0d+0) * ff00 &
& + (3.0d+0 / 4.0d+0) *(ferr+fell)
end if
end if
return
end subroutine