[13926] | 1 | |
---|
| 2 | ! |
---|
| 3 | ! This program may be freely redistributed under the |
---|
| 4 | ! condition that the copyright notices (including this |
---|
| 5 | ! entire header) are not removed, and no compensation |
---|
| 6 | ! is received through use of the software. Private, |
---|
| 7 | ! research, and institutional use is free. You may |
---|
| 8 | ! distribute modified versions of this code UNDER THE |
---|
| 9 | ! CONDITION THAT THIS CODE AND ANY MODIFICATIONS MADE |
---|
| 10 | ! TO IT IN THE SAME FILE REMAIN UNDER COPYRIGHT OF THE |
---|
| 11 | ! ORIGINAL AUTHOR, BOTH SOURCE AND OBJECT CODE ARE |
---|
| 12 | ! MADE FREELY AVAILABLE WITHOUT CHARGE, AND CLEAR |
---|
| 13 | ! NOTICE IS GIVEN OF THE MODIFICATIONS. Distribution |
---|
| 14 | ! of this code as part of a commercial system is |
---|
| 15 | ! permissible ONLY BY DIRECT ARRANGEMENT WITH THE |
---|
| 16 | ! AUTHOR. (If you are not directly supplying this |
---|
| 17 | ! code to a customer, and you are instead telling them |
---|
| 18 | ! how they can obtain it for free, then you are not |
---|
| 19 | ! required to make any arrangement with me.) |
---|
| 20 | ! |
---|
| 21 | ! Disclaimer: Neither I nor: Columbia University, the |
---|
| 22 | ! National Aeronautics and Space Administration, nor |
---|
| 23 | ! the Massachusetts Institute of Technology warrant |
---|
| 24 | ! or certify this code in any way whatsoever. This |
---|
| 25 | ! code is provided "as-is" to be used at your own risk. |
---|
| 26 | ! |
---|
| 27 | ! |
---|
| 28 | |
---|
| 29 | ! |
---|
[14212] | 30 | ! PPM.h90: 1d slope-limited, piecewise parabolic recon. |
---|
[13926] | 31 | ! |
---|
| 32 | ! Darren Engwirda |
---|
| 33 | ! 08-Sep-2016 |
---|
| 34 | ! de2363 [at] columbia [dot] edu |
---|
| 35 | ! |
---|
| 36 | ! |
---|
| 37 | |
---|
| 38 | ! P. Colella and PR. Woodward, The Piecewise Parabolic |
---|
| 39 | ! Method (PPM) for gas-dynamical simulations., J. Comp. |
---|
| 40 | ! Phys., 54 (1), 1984, 174-201, |
---|
| 41 | ! https://doi.org/10.1016/0021-9991(84)90143-8 |
---|
| 42 | ! |
---|
| 43 | |
---|
| 44 | pure subroutine ppm(npos,nvar,ndof,delx, & |
---|
| 45 | & fdat,fhat,edge,oscl, & |
---|
| 46 | & dmin,ilim,wlim,halo) |
---|
| 47 | |
---|
| 48 | ! |
---|
| 49 | ! NPOS no. edges over grid. |
---|
| 50 | ! NVAR no. state variables. |
---|
| 51 | ! NDOF no. degrees-of-freedom per grid-cell. |
---|
| 52 | ! DELX grid-cell spacing array. LENGTH(DELX) == +1 if |
---|
| 53 | ! spacing is uniform . |
---|
| 54 | ! FDAT grid-cell moments array. FDAT is an array with |
---|
| 55 | ! SIZE = NDOF-by-NVAR-by-NPOS-1 . |
---|
| 56 | ! FHAT grid-cell re-con. array. FHAT is an array with |
---|
| 57 | ! SIZE = MDOF-by-NVAR-by-NPOS-1 . |
---|
| 58 | ! EDGE edge-centred interp. for function-value. EDGE |
---|
| 59 | ! is an array with SIZE = NVAR-by-NPOS . |
---|
| 60 | ! OSCL grid-cell oscil. dof.'s. OSCL is an array with |
---|
| 61 | ! SIZE = +2 -by-NVAR-by-NPOS-1 . |
---|
| 62 | ! DMIN min. grid-cell spacing thresh . |
---|
| 63 | ! ILIM cell slope-limiting selection . |
---|
| 64 | ! WLIM wall slope-limiting selection . |
---|
| 65 | ! HALO width of re-con. stencil, symmetric about mid. . |
---|
| 66 | ! |
---|
| 67 | |
---|
| 68 | implicit none |
---|
| 69 | |
---|
| 70 | !------------------------------------------- arguments ! |
---|
| 71 | integer, intent(in) :: npos,nvar,ndof |
---|
| 72 | real*8 , intent(in) :: dmin |
---|
| 73 | real*8 , intent(out) :: fhat(:,:,:) |
---|
| 74 | real*8 , intent(in) :: oscl(:,:,:) |
---|
| 75 | real*8 , intent(in) :: delx(:) |
---|
| 76 | real*8 , intent(in) :: fdat(:,:,:) |
---|
| 77 | real*8 , intent(in) :: edge(:,:) |
---|
| 78 | integer, intent(in) :: ilim,wlim,halo |
---|
| 79 | |
---|
| 80 | !------------------------------------------- variables ! |
---|
| 81 | integer :: ipos,ivar,iill,iirr,head,tail |
---|
| 82 | real*8 :: ff00,ffll,ffrr,hh00,hhll,hhrr |
---|
| 83 | integer :: mono |
---|
| 84 | real*8 :: fell,ferr |
---|
| 85 | real*8 :: dfds(-1:+1) |
---|
| 86 | real*8 :: wval(+1:+2) |
---|
| 87 | real*8 :: uhat(+1:+3) |
---|
| 88 | real*8 :: lhat(+1:+3) |
---|
| 89 | |
---|
| 90 | head = +1; tail = npos - 1 |
---|
| 91 | |
---|
| 92 | if (npos.eq.2) then |
---|
| 93 | !----- default to reduced order if insufficient points ! |
---|
| 94 | do ivar = +1, nvar |
---|
| 95 | fhat(1,ivar,+1) = & |
---|
| 96 | & fdat(1,ivar,+1) |
---|
[14387] | 97 | fhat(2,ivar,+1) = 0.e0 |
---|
| 98 | fhat(3,ivar,+1) = 0.e0 |
---|
[13926] | 99 | end do |
---|
| 100 | end if |
---|
| 101 | |
---|
| 102 | if (npos.le.2) return |
---|
| 103 | |
---|
| 104 | !------------------- reconstruct function on each cell ! |
---|
| 105 | |
---|
[14387] | 106 | uhat = +0.e+0 |
---|
| 107 | lhat = +0.e+0 |
---|
[13926] | 108 | |
---|
| 109 | do ipos = +1 , npos-1 |
---|
| 110 | |
---|
| 111 | iill = max(head,ipos-1) |
---|
| 112 | iirr = min(tail,ipos+1) |
---|
| 113 | |
---|
| 114 | do ivar = +1 , nvar-0 |
---|
| 115 | |
---|
| 116 | !----------------------------- cell mean + edge values ! |
---|
| 117 | |
---|
| 118 | ff00 = fdat(1,ivar,ipos) |
---|
| 119 | ffll = fdat(1,ivar,iill) |
---|
| 120 | ffrr = fdat(1,ivar,iirr) |
---|
| 121 | |
---|
| 122 | fell = edge(ivar,ipos+0) |
---|
| 123 | ferr = edge(ivar,ipos+1) |
---|
| 124 | |
---|
| 125 | !----------------------------- calc. LL/00/RR gradient ! |
---|
| 126 | |
---|
| 127 | if (size(delx).gt.+1) then |
---|
| 128 | |
---|
| 129 | hh00 = delx(ipos) |
---|
| 130 | hhll = delx(iill) |
---|
| 131 | hhrr = delx(iirr) |
---|
| 132 | |
---|
| 133 | call plsv (ffll,hhll,ff00, & |
---|
| 134 | & hh00,ffrr,hhrr, & |
---|
| 135 | & dfds) |
---|
| 136 | else |
---|
| 137 | |
---|
| 138 | call plsc (ffll,ff00,ffrr, & |
---|
| 139 | & dfds) |
---|
| 140 | |
---|
| 141 | end if |
---|
| 142 | |
---|
| 143 | !----------------------------- calc. cell-wise profile ! |
---|
| 144 | |
---|
| 145 | select case(ilim) |
---|
| 146 | case (null_limit) |
---|
| 147 | |
---|
| 148 | !----------------------------- calc. unlimited profile ! |
---|
| 149 | |
---|
| 150 | call ppmfn(ff00,ffll,ffrr, & |
---|
| 151 | & fell,ferr,dfds, & |
---|
| 152 | & uhat,lhat,mono) |
---|
| 153 | |
---|
| 154 | !----------------------------- pref. unlimited profile ! |
---|
| 155 | |
---|
[14387] | 156 | wval(1) = +1.e+0 |
---|
| 157 | wval(2) = +0.e+0 |
---|
[13926] | 158 | |
---|
| 159 | case (mono_limit) |
---|
| 160 | |
---|
| 161 | !----------------------------- calc. monotonic profile ! |
---|
| 162 | |
---|
| 163 | call ppmfn(ff00,ffll,ffrr, & |
---|
| 164 | & fell,ferr,dfds, & |
---|
| 165 | & uhat,lhat,mono) |
---|
| 166 | |
---|
| 167 | !----------------------------- pref. monotonic profile ! |
---|
| 168 | |
---|
[14387] | 169 | wval(1) = +0.e+0 |
---|
| 170 | wval(2) = +1.e+0 |
---|
[13926] | 171 | |
---|
| 172 | case (weno_limit) |
---|
| 173 | |
---|
| 174 | !----------------------------- calc. unlimited profile ! |
---|
| 175 | |
---|
| 176 | call ppmfn(ff00,ffll,ffrr, & |
---|
| 177 | & fell,ferr,dfds, & |
---|
| 178 | & uhat,lhat,mono) |
---|
| 179 | |
---|
| 180 | if (mono.gt.+0) then |
---|
| 181 | |
---|
| 182 | !----------------------------- calc. WENO-type weights ! |
---|
| 183 | |
---|
| 184 | call wenoi(npos,delx,oscl, & |
---|
| 185 | & ipos,ivar,halo, & |
---|
| 186 | & wlim,wval) |
---|
| 187 | |
---|
| 188 | else |
---|
| 189 | |
---|
| 190 | !----------------------------- pref. unlimited profile ! |
---|
| 191 | |
---|
[14387] | 192 | wval(1) = +1.e+0 |
---|
| 193 | wval(2) = +0.e+0 |
---|
[13926] | 194 | |
---|
| 195 | end if |
---|
| 196 | |
---|
| 197 | end select |
---|
| 198 | |
---|
| 199 | !----------------------------- blend "null" and "mono" ! |
---|
| 200 | |
---|
| 201 | fhat(1,ivar,ipos) = & |
---|
| 202 | & wval(1) * uhat(1) + & |
---|
| 203 | & wval(2) * lhat(1) |
---|
| 204 | fhat(2,ivar,ipos) = & |
---|
| 205 | & wval(1) * uhat(2) + & |
---|
| 206 | & wval(2) * lhat(2) |
---|
| 207 | fhat(3,ivar,ipos) = & |
---|
| 208 | & wval(1) * uhat(3) + & |
---|
| 209 | & wval(2) * lhat(3) |
---|
| 210 | |
---|
| 211 | end do |
---|
| 212 | |
---|
| 213 | end do |
---|
| 214 | |
---|
| 215 | return |
---|
| 216 | |
---|
| 217 | end subroutine |
---|
| 218 | |
---|
| 219 | !--------- assemble piecewise parabolic reconstruction ! |
---|
| 220 | |
---|
| 221 | pure subroutine ppmfn(ff00,ffll,ffrr,fell,& |
---|
| 222 | & ferr,dfds,uhat,lhat,& |
---|
| 223 | & mono) |
---|
| 224 | |
---|
| 225 | ! |
---|
| 226 | ! FF00 centred grid-cell mean. |
---|
| 227 | ! FFLL left -biased grid-cell mean. |
---|
| 228 | ! FFRR right-biased grid-cell mean. |
---|
| 229 | ! FELL left -biased edge interp. |
---|
| 230 | ! FERR right-biased edge interp. |
---|
| 231 | ! DFDS piecewise linear gradients in local co-ord.'s. |
---|
| 232 | ! DFDS(+0) is a centred, slope-limited estimate, |
---|
| 233 | ! DFDS(-1), DFDS(+1) are left- and right-biased |
---|
| 234 | ! estimates (unlimited). |
---|
| 235 | ! UHAT unlimited PPM reconstruction coefficients . |
---|
| 236 | ! LHAT monotonic PPM reconstruction coefficients . |
---|
| 237 | ! MONO slope-limiting indicator, MONO > +0 if some |
---|
| 238 | ! limiting has occured . |
---|
| 239 | ! |
---|
| 240 | |
---|
| 241 | implicit none |
---|
| 242 | |
---|
| 243 | !------------------------------------------- arguments ! |
---|
| 244 | real*8 , intent(in) :: ff00 |
---|
| 245 | real*8 , intent(in) :: ffll,ffrr |
---|
| 246 | real*8 , intent(inout) :: fell,ferr |
---|
| 247 | real*8 , intent(in) :: dfds(-1:+1) |
---|
| 248 | real*8 , intent(out) :: uhat(+1:+3) |
---|
| 249 | real*8 , intent(out) :: lhat(+1:+3) |
---|
| 250 | integer, intent(out) :: mono |
---|
| 251 | |
---|
| 252 | !------------------------------------------- variables ! |
---|
| 253 | real*8 :: turn |
---|
| 254 | |
---|
| 255 | mono = 0 |
---|
| 256 | |
---|
| 257 | !-------------------------------- "null" slope-limiter ! |
---|
| 258 | |
---|
| 259 | uhat( 1 ) = & |
---|
| 260 | & + (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 261 | & - (1.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 262 | uhat( 2 ) = & |
---|
| 263 | & + (1.0d+0 / 2.0d+0) *(ferr-fell) |
---|
| 264 | uhat( 3 ) = & |
---|
| 265 | & - (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 266 | & + (3.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 267 | |
---|
| 268 | !-------------------------------- "mono" slope-limiter ! |
---|
| 269 | |
---|
| 270 | if((ffrr - ff00) * & |
---|
[14387] | 271 | & (ff00 - ffll) .lt. 0.e+0) then |
---|
[13926] | 272 | |
---|
| 273 | !----------------------------------- "flatten" extrema ! |
---|
| 274 | |
---|
| 275 | mono = +1 |
---|
| 276 | |
---|
| 277 | lhat(1) = ff00 |
---|
[14387] | 278 | lhat(2) = 0.e0 |
---|
| 279 | lhat(3) = 0.e0 |
---|
[13926] | 280 | |
---|
| 281 | return |
---|
| 282 | |
---|
| 283 | end if |
---|
| 284 | |
---|
| 285 | !----------------------------------- limit edge values ! |
---|
| 286 | |
---|
| 287 | if((ffll - fell) * & |
---|
[14387] | 288 | & (fell - ff00) .le. 0.e+0) then |
---|
[13926] | 289 | |
---|
| 290 | mono = +1 |
---|
| 291 | |
---|
| 292 | fell = ff00 - dfds(0) |
---|
| 293 | |
---|
| 294 | end if |
---|
| 295 | |
---|
| 296 | if((ffrr - ferr) * & |
---|
[14387] | 297 | & (ferr - ff00) .le. 0.e+0) then |
---|
[13926] | 298 | |
---|
| 299 | mono = +1 |
---|
| 300 | |
---|
| 301 | ferr = ff00 + dfds(0) |
---|
| 302 | |
---|
| 303 | end if |
---|
| 304 | |
---|
| 305 | !----------------------------------- update ppm coeff. ! |
---|
| 306 | |
---|
| 307 | lhat( 1 ) = & |
---|
| 308 | & + (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 309 | & - (1.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 310 | lhat( 2 ) = & |
---|
| 311 | & + (1.0d+0 / 2.0d+0) *(ferr-fell) |
---|
| 312 | lhat( 3 ) = & |
---|
| 313 | & - (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 314 | & + (3.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 315 | |
---|
| 316 | !----------------------------------- limit cell values ! |
---|
| 317 | |
---|
| 318 | if (abs(lhat(3)) .gt. & |
---|
| 319 | & abs(lhat(2))*.5d+0) then |
---|
| 320 | |
---|
| 321 | turn = -0.5d+0 * lhat(2) & |
---|
| 322 | & / lhat(3) |
---|
| 323 | |
---|
[14387] | 324 | if ((turn .ge. -1.e+0)& |
---|
| 325 | & .and.(turn .le. +0.e+0)) then |
---|
[13926] | 326 | |
---|
| 327 | mono = +2 |
---|
| 328 | |
---|
| 329 | !--------------------------- push TURN onto lower edge ! |
---|
| 330 | |
---|
| 331 | ferr = +3.0d+0 * ff00 & |
---|
| 332 | & -2.0d+0 * fell |
---|
| 333 | |
---|
| 334 | lhat( 1 ) = & |
---|
| 335 | & + (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 336 | & - (1.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 337 | lhat( 2 ) = & |
---|
| 338 | & + (1.0d+0 / 2.0d+0) *(ferr-fell) |
---|
| 339 | lhat( 3 ) = & |
---|
| 340 | & - (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 341 | & + (3.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 342 | |
---|
| 343 | else & |
---|
[14387] | 344 | & if ((turn .gt. +0.e+0)& |
---|
| 345 | & .and.(turn .le. +1.e+0)) then |
---|
[13926] | 346 | |
---|
| 347 | mono = +2 |
---|
| 348 | |
---|
| 349 | !--------------------------- push TURN onto upper edge ! |
---|
| 350 | |
---|
| 351 | fell = +3.0d+0 * ff00 & |
---|
| 352 | & -2.0d+0 * ferr |
---|
| 353 | |
---|
| 354 | lhat( 1 ) = & |
---|
| 355 | & + (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 356 | & - (1.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 357 | lhat( 2 ) = & |
---|
| 358 | & + (1.0d+0 / 2.0d+0) *(ferr-fell) |
---|
| 359 | lhat( 3 ) = & |
---|
| 360 | & - (3.0d+0 / 2.0d+0) * ff00 & |
---|
| 361 | & + (3.0d+0 / 4.0d+0) *(ferr+fell) |
---|
| 362 | |
---|
| 363 | end if |
---|
| 364 | |
---|
| 365 | end if |
---|
| 366 | |
---|
| 367 | return |
---|
| 368 | |
---|
| 369 | end subroutine |
---|
| 370 | |
---|
| 371 | |
---|
| 372 | |
---|