[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 | ! |
---|
| 30 | ! WENO1D.f90: WENO-style slope-limiting for 1d reconst. |
---|
| 31 | ! |
---|
| 32 | ! Darren Engwirda |
---|
| 33 | ! 08-Sep-2016 |
---|
| 34 | ! de2363 [at] columbia [dot] edu |
---|
| 35 | ! |
---|
| 36 | ! |
---|
| 37 | |
---|
| 38 | pure subroutine wenoi (npos,delx,oscl,ipos, & |
---|
| 39 | & ivar,halo,& |
---|
| 40 | & wlim,wval ) |
---|
| 41 | |
---|
| 42 | ! |
---|
| 43 | ! NPOS no. edges over grid. |
---|
| 44 | ! DELX grid-cell spacing array. SIZE(DELX) == +1 if |
---|
| 45 | ! the grid is uniformly spaced . |
---|
| 46 | ! OSCL cell-centred oscillation-detectors, where OSCL |
---|
| 47 | ! has SIZE = +2-by-NVAR-by-NPOS-1. OSCL is given |
---|
| 48 | ! by calls to OSCLI(). |
---|
| 49 | ! IPOS grid-cell index for which to calc. weights . |
---|
| 50 | ! IVAR state-var index for which to calc/ weights . |
---|
| 51 | ! HALO width of recon. stencil, symmetric about IPOS . |
---|
| 52 | ! WLIM limiter treatment at endpoints, monotonic or |
---|
| 53 | ! otherwise . |
---|
| 54 | ! WVAL WENO weights vector, such that FHAT = WVAL(1) * |
---|
| 55 | ! UHAT + WVAL(2) * LHAT, where UHAT and LHAT are |
---|
| 56 | ! the unlimited and monotonic grid-cell profiles |
---|
| 57 | ! respectively . |
---|
| 58 | ! |
---|
| 59 | |
---|
| 60 | implicit none |
---|
| 61 | |
---|
| 62 | !------------------------------------------- arguments ! |
---|
| 63 | integer, intent(in) :: npos,halo |
---|
| 64 | integer, intent(in) :: ipos,ivar |
---|
| 65 | integer, intent(in) :: wlim |
---|
| 66 | real*8 , intent(in) :: delx(:) |
---|
| 67 | real*8 , intent(in) :: oscl(:,:,:) |
---|
| 68 | real*8 , intent(out) :: wval(2) |
---|
| 69 | |
---|
| 70 | !------------------------------------------- variables ! |
---|
| 71 | real*8 :: omin,omax,wsum |
---|
| 72 | |
---|
| 73 | real*8 , parameter :: ZERO = +1.d-16 |
---|
| 74 | |
---|
| 75 | if (size(delx).gt.+1) then |
---|
| 76 | |
---|
| 77 | !------------------- use variable grid spacing variant ! |
---|
| 78 | |
---|
| 79 | call wenov(npos,delx,oscl, & |
---|
| 80 | & ipos,ivar,halo, & |
---|
| 81 | & wlim,omin,omax) |
---|
| 82 | |
---|
| 83 | else |
---|
| 84 | |
---|
| 85 | !------------------- use constant grid spacing variant ! |
---|
| 86 | |
---|
| 87 | call wenoc(npos,delx,oscl, & |
---|
| 88 | & ipos,ivar,halo, & |
---|
| 89 | & wlim,omin,omax) |
---|
| 90 | |
---|
| 91 | end if |
---|
| 92 | |
---|
| 93 | !------------------ compute WENO-style profile weights ! |
---|
| 94 | |
---|
| 95 | omax = omax + ZERO |
---|
| 96 | omin = omin + ZERO |
---|
| 97 | |
---|
| 98 | if (halo .ge. +3) then |
---|
| 99 | |
---|
| 100 | wval(1) = +1.0d+7 / omax ** 3 |
---|
| 101 | wval(2) = +1.0d+0 / omin ** 3 |
---|
| 102 | |
---|
| 103 | else & |
---|
| 104 | & if (halo .le. +2) then |
---|
| 105 | |
---|
| 106 | wval(1) = +1.0d+5 / omax ** 3 |
---|
| 107 | wval(2) = +1.0d+0 / omin ** 3 |
---|
| 108 | |
---|
| 109 | end if |
---|
| 110 | |
---|
| 111 | wsum = wval(1) + wval(2) + ZERO |
---|
| 112 | wval(1) = wval(1) / wsum |
---|
| 113 | ! wval(2) = wval(2) / wsum |
---|
| 114 | wval(2) =-wval(1) + 1.d0 ! wval(2)/wsum but robust ! |
---|
| 115 | |
---|
| 116 | return |
---|
| 117 | |
---|
| 118 | end subroutine |
---|
| 119 | |
---|
| 120 | pure subroutine wenov (npos,delx,oscl,ipos, & |
---|
| 121 | & ivar,halo,& |
---|
| 122 | & wlim,omin,omax) |
---|
| 123 | |
---|
| 124 | ! |
---|
| 125 | ! *this is the variable grid-spacing variant . |
---|
| 126 | ! |
---|
| 127 | ! NPOS no. edges over grid. |
---|
| 128 | ! DELX grid-cell spacing array. SIZE(DELX) == +1 if |
---|
| 129 | ! the grid is uniformly spaced . |
---|
| 130 | ! OSCL cell-centred oscillation-detectors, where OSCL |
---|
| 131 | ! has SIZE = +2-by-NVAR-by-NPOS-1. OSCL is given |
---|
| 132 | ! by calls to OSCLI(). |
---|
| 133 | ! IPOS grid-cell index for which to calc. weights . |
---|
| 134 | ! IVAR state-var index for which to calc/ weights . |
---|
| 135 | ! HALO width of recon. stencil, symmetric about IPOS . |
---|
| 136 | ! WLIM limiter treatment at endpoints, monotonic or |
---|
| 137 | ! otherwise . |
---|
| 138 | ! OMIN min. and max. oscillation indicators over the |
---|
| 139 | ! OMAX local re-con. stencil . |
---|
| 140 | ! |
---|
| 141 | |
---|
| 142 | implicit none |
---|
| 143 | |
---|
| 144 | !------------------------------------------- arguments ! |
---|
| 145 | integer, intent(in) :: npos,halo |
---|
| 146 | integer, intent(in) :: ipos,ivar |
---|
| 147 | integer, intent(in) :: wlim |
---|
| 148 | real*8 , intent(in) :: delx(:) |
---|
| 149 | real*8 , intent(in) :: oscl(:,:,:) |
---|
| 150 | real*8 , intent(out) :: omin,omax |
---|
| 151 | |
---|
| 152 | !------------------------------------------- variables ! |
---|
| 153 | integer :: hpos |
---|
| 154 | integer :: head,tail |
---|
| 155 | integer :: imin,imax |
---|
| 156 | real*8 :: deli,delh |
---|
| 157 | real*8 :: hh00,hsqr |
---|
| 158 | real*8 :: dfx1,dfx2 |
---|
| 159 | real*8 :: oval |
---|
| 160 | |
---|
| 161 | !------------------- calc. lower//upper stencil bounds ! |
---|
| 162 | |
---|
| 163 | head = 1; tail = npos - 1 |
---|
| 164 | |
---|
| 165 | if(wlim.eq.mono_limit) then |
---|
| 166 | |
---|
| 167 | !---------------------- deactivate WENO at boundaries ! |
---|
| 168 | |
---|
| 169 | if (ipos-halo.lt.head) then |
---|
| 170 | |
---|
| 171 | omax = 1.d0 |
---|
| 172 | omin = 0.d0 ; return |
---|
| 173 | |
---|
| 174 | end if |
---|
| 175 | |
---|
| 176 | if (ipos+halo.gt.tail) then |
---|
| 177 | |
---|
| 178 | omax = 1.d0 |
---|
| 179 | omin = 0.d0 ; return |
---|
| 180 | |
---|
| 181 | end if |
---|
| 182 | |
---|
| 183 | end if |
---|
| 184 | |
---|
| 185 | !---------------------- truncate stencil at boundaries ! |
---|
| 186 | |
---|
| 187 | imin = max(ipos-halo,head) |
---|
| 188 | imax = min(ipos+halo,tail) |
---|
| 189 | |
---|
| 190 | !------------------ find min/max indicators on stencil ! |
---|
| 191 | |
---|
| 192 | dfx1 = oscl(1,ivar,ipos) |
---|
| 193 | dfx2 = oscl(2,ivar,ipos) |
---|
| 194 | |
---|
| 195 | hh00 = delx(ipos+0)**1 |
---|
| 196 | hsqr = delx(ipos+0)**2 |
---|
| 197 | |
---|
| 198 | oval =(hh00 * dfx1)**2 & |
---|
| 199 | & +(hsqr * dfx2)**2 |
---|
| 200 | |
---|
| 201 | omin = oval |
---|
| 202 | omax = oval |
---|
| 203 | |
---|
| 204 | !---------------------------------------- "lower" part ! |
---|
| 205 | |
---|
| 206 | delh = 0.d0 |
---|
| 207 | |
---|
| 208 | do hpos = ipos-1, imin, -1 |
---|
| 209 | |
---|
| 210 | !------------------ calc. derivatives centred on IPOS. ! |
---|
| 211 | |
---|
| 212 | deli = delx(hpos+0) & |
---|
| 213 | & + delx(hpos+1) |
---|
| 214 | |
---|
| 215 | delh = delh + deli*.5d0 |
---|
| 216 | |
---|
| 217 | dfx1 = oscl(1,ivar,hpos) |
---|
| 218 | dfx2 = oscl(2,ivar,hpos) |
---|
| 219 | |
---|
| 220 | dfx1 = dfx1 + dfx2*delh |
---|
| 221 | |
---|
| 222 | !------------------ indicator: NORM(H^N * D^N/DX^N(F)) ! |
---|
| 223 | |
---|
| 224 | oval = (hh00 * dfx1)**2 & |
---|
| 225 | & + (hsqr * dfx2)**2 |
---|
| 226 | |
---|
| 227 | if (oval .lt. omin) then |
---|
| 228 | omin = oval |
---|
| 229 | else & |
---|
| 230 | & if (oval .gt. omax) then |
---|
| 231 | omax = oval |
---|
| 232 | end if |
---|
| 233 | |
---|
| 234 | end do |
---|
| 235 | |
---|
| 236 | !---------------------------------------- "upper" part ! |
---|
| 237 | |
---|
| 238 | delh = 0.d0 |
---|
| 239 | |
---|
| 240 | do hpos = ipos+1, imax, +1 |
---|
| 241 | |
---|
| 242 | !------------------ calc. derivatives centred on IPOS. ! |
---|
| 243 | |
---|
| 244 | deli = delx(hpos+0) & |
---|
| 245 | & + delx(hpos-1) |
---|
| 246 | |
---|
| 247 | delh = delh - deli*.5d0 |
---|
| 248 | |
---|
| 249 | dfx1 = oscl(1,ivar,hpos) |
---|
| 250 | dfx2 = oscl(2,ivar,hpos) |
---|
| 251 | |
---|
| 252 | dfx1 = dfx1 + dfx2*delh |
---|
| 253 | |
---|
| 254 | !------------------ indicator: NORM(H^N * D^N/DX^N(F)) ! |
---|
| 255 | |
---|
| 256 | oval = (hh00 * dfx1)**2 & |
---|
| 257 | & + (hsqr * dfx2)**2 |
---|
| 258 | |
---|
| 259 | if (oval .lt. omin) then |
---|
| 260 | omin = oval |
---|
| 261 | else & |
---|
| 262 | & if (oval .gt. omax) then |
---|
| 263 | omax = oval |
---|
| 264 | end if |
---|
| 265 | |
---|
| 266 | end do |
---|
| 267 | |
---|
| 268 | return |
---|
| 269 | |
---|
| 270 | end subroutine |
---|
| 271 | |
---|
| 272 | pure subroutine wenoc (npos,delx,oscl,ipos, & |
---|
| 273 | & ivar,halo,& |
---|
| 274 | & wlim,omin,omax) |
---|
| 275 | |
---|
| 276 | ! |
---|
| 277 | ! *this is the constant grid-spacing variant . |
---|
| 278 | ! |
---|
| 279 | ! NPOS no. edges over grid. |
---|
| 280 | ! DELX grid-cell spacing array. SIZE(DELX) == +1 if |
---|
| 281 | ! the grid is uniformly spaced . |
---|
| 282 | ! OSCL cell-centred oscillation-detectors, where OSCL |
---|
| 283 | ! has SIZE = +2-by-NVAR-by-NPOS-1. OSCL is given |
---|
| 284 | ! by calls to OSCLI(). |
---|
| 285 | ! IPOS grid-cell index for which to calc. weights . |
---|
| 286 | ! IVAR state-var index for which to calc/ weights . |
---|
| 287 | ! HALO width of recon. stencil, symmetric about IPOS . |
---|
| 288 | ! WLIM limiter treatment at endpoints, monotonic or |
---|
| 289 | ! otherwise . |
---|
| 290 | ! OMIN min. and max. oscillation indicators over the |
---|
| 291 | ! OMAX local re-con. stencil . |
---|
| 292 | ! |
---|
| 293 | |
---|
| 294 | implicit none |
---|
| 295 | |
---|
| 296 | !------------------------------------------- arguments ! |
---|
| 297 | integer, intent(in) :: npos,halo |
---|
| 298 | integer, intent(in) :: ipos,ivar |
---|
| 299 | integer, intent(in) :: wlim |
---|
| 300 | real*8 , intent(in) :: delx(1) |
---|
| 301 | real*8 , intent(in) :: oscl(:,:,:) |
---|
| 302 | real*8 , intent(out) :: omin,omax |
---|
| 303 | |
---|
| 304 | !------------------------------------------- variables ! |
---|
| 305 | integer :: hpos |
---|
| 306 | integer :: head,tail |
---|
| 307 | integer :: imin,imax |
---|
| 308 | real*8 :: delh |
---|
| 309 | real*8 :: dfx1,dfx2 |
---|
| 310 | real*8 :: oval |
---|
| 311 | |
---|
| 312 | !------------------- calc. lower//upper stencil bounds ! |
---|
| 313 | |
---|
| 314 | head = 1; tail = npos - 1 |
---|
| 315 | |
---|
| 316 | if(wlim.eq.mono_limit) then |
---|
| 317 | |
---|
| 318 | !---------------------- deactivate WENO at boundaries ! |
---|
| 319 | |
---|
| 320 | if (ipos-halo.lt.head) then |
---|
| 321 | |
---|
| 322 | omax = 1.d0 |
---|
| 323 | omin = 0.d0 ; return |
---|
| 324 | |
---|
| 325 | end if |
---|
| 326 | |
---|
| 327 | if (ipos+halo.gt.tail) then |
---|
| 328 | |
---|
| 329 | omax = 1.d0 |
---|
| 330 | omin = 0.d0 ; return |
---|
| 331 | |
---|
| 332 | end if |
---|
| 333 | |
---|
| 334 | end if |
---|
| 335 | |
---|
| 336 | !---------------------- truncate stencil at boundaries ! |
---|
| 337 | |
---|
| 338 | imin = max(ipos-halo,head) |
---|
| 339 | imax = min(ipos+halo,tail) |
---|
| 340 | |
---|
| 341 | !------------------ find min/max indicators on stencil ! |
---|
| 342 | |
---|
| 343 | dfx1 = oscl(1,ivar,ipos) |
---|
| 344 | dfx2 = oscl(2,ivar,ipos) |
---|
| 345 | |
---|
| 346 | oval = (2.d0**1*dfx1)**2 & |
---|
| 347 | & + (2.d0**2*dfx2)**2 |
---|
| 348 | |
---|
| 349 | omin = oval |
---|
| 350 | omax = oval |
---|
| 351 | |
---|
| 352 | !---------------------------------------- "lower" part ! |
---|
| 353 | |
---|
| 354 | delh = 0.d0 |
---|
| 355 | |
---|
| 356 | do hpos = ipos-1, imin, -1 |
---|
| 357 | |
---|
| 358 | !------------------ calc. derivatives centred on IPOS. ! |
---|
| 359 | |
---|
| 360 | delh = delh + 2.d0 |
---|
| 361 | |
---|
| 362 | dfx1 = oscl(1,ivar,hpos) |
---|
| 363 | dfx2 = oscl(2,ivar,hpos) |
---|
| 364 | |
---|
| 365 | dfx1 = dfx1 + dfx2*delh |
---|
| 366 | |
---|
| 367 | !------------------ indicator: NORM(H^N * D^N/DX^N(F)) ! |
---|
| 368 | |
---|
| 369 | oval = (2.d0**1*dfx1)**2 & |
---|
| 370 | & + (2.d0**2*dfx2)**2 |
---|
| 371 | |
---|
| 372 | if (oval .lt. omin) then |
---|
| 373 | omin = oval |
---|
| 374 | else & |
---|
| 375 | & if (oval .gt. omax) then |
---|
| 376 | omax = oval |
---|
| 377 | end if |
---|
| 378 | |
---|
| 379 | end do |
---|
| 380 | |
---|
| 381 | !---------------------------------------- "upper" part ! |
---|
| 382 | |
---|
| 383 | delh = 0.d0 |
---|
| 384 | |
---|
| 385 | do hpos = ipos+1, imax, +1 |
---|
| 386 | |
---|
| 387 | !------------------ calc. derivatives centred on IPOS. ! |
---|
| 388 | |
---|
| 389 | delh = delh - 2.d0 |
---|
| 390 | |
---|
| 391 | dfx1 = oscl(1,ivar,hpos) |
---|
| 392 | dfx2 = oscl(2,ivar,hpos) |
---|
| 393 | |
---|
| 394 | dfx1 = dfx1 + dfx2*delh |
---|
| 395 | |
---|
| 396 | !------------------ indicator: NORM(H^N * D^N/DX^N(F)) ! |
---|
| 397 | |
---|
| 398 | oval = (2.d0**1*dfx1)**2 & |
---|
| 399 | & + (2.d0**2*dfx2)**2 |
---|
| 400 | |
---|
| 401 | if (oval .lt. omin) then |
---|
| 402 | omin = oval |
---|
| 403 | else & |
---|
| 404 | & if (oval .gt. omax) then |
---|
| 405 | omax = oval |
---|
| 406 | end if |
---|
| 407 | |
---|
| 408 | end do |
---|
| 409 | |
---|
| 410 | return |
---|
| 411 | |
---|
| 412 | end subroutine |
---|
| 413 | |
---|
| 414 | |
---|
| 415 | |
---|