|
| 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 | + module ppr_1d |
| 30 | + |
| 31 | + ! |
| 32 | + ! PPR-1D.f90: 1-d piecewise polynomial reconstructions. |
| 33 | + ! |
| 34 | + ! Darren Engwirda |
| 35 | + ! 31-Mar-2019 |
| 36 | + ! de2363 [at] columbia [dot] edu |
| 37 | + ! |
| 38 | + ! |
| 39 | + |
| 40 | + implicit none |
| 41 | + |
| 42 | + !------------------------------------ compile-time def ! |
| 43 | + |
| 44 | +! define __PPR_WARNMAT__ |
| 45 | +! define __PPR_DOTIMER__ |
| 46 | + |
| 47 | +# ifdef __PPR_DOTIMER__ |
| 48 | + |
| 49 | +# define __TIC__ \ |
| 50 | + call system_clock (ttic,rate) |
| 51 | + |
| 52 | +# define __TOC__(time, mark) \ |
| 53 | + call system_clock (ttoc,rate) ; \ |
| 54 | + if ( present(time)) \ |
| 55 | + time%mark=time%mark+(ttoc-ttic) |
| 56 | + |
| 57 | +# else |
| 58 | + |
| 59 | +# define __TIC__ |
| 60 | +# define __TOC__(time, mark) |
| 61 | + |
| 62 | +# endif |
| 63 | + |
| 64 | + !------------------------------------ method selection ! |
| 65 | + |
| 66 | + integer, parameter :: p1e_method = +100 |
| 67 | + integer, parameter :: p3e_method = +101 |
| 68 | + integer, parameter :: p5e_method = +102 |
| 69 | + |
| 70 | + integer, parameter :: pcm_method = +200 |
| 71 | + integer, parameter :: plm_method = +201 |
| 72 | + integer, parameter :: ppm_method = +202 |
| 73 | + integer, parameter :: pqm_method = +203 |
| 74 | + |
| 75 | + integer, parameter :: null_limit = +300 |
| 76 | + integer, parameter :: mono_limit = +301 |
| 77 | + integer, parameter :: weno_limit = +302 |
| 78 | + |
| 79 | + integer, parameter :: bcon_loose = +400 |
| 80 | + integer, parameter :: bcon_value = +401 |
| 81 | + integer, parameter :: bcon_slope = +402 |
| 82 | + |
| 83 | + type rmap_tics |
| 84 | + !------------------------------- tCPU timer for RCON1D ! |
| 85 | + integer :: rmap_time |
| 86 | + integer :: edge_time |
| 87 | + integer :: cell_time |
| 88 | + integer :: oscl_time |
| 89 | + end type rmap_tics |
| 90 | + |
| 91 | + type rcon_opts |
| 92 | + !------------------------------- parameters for RCON1D ! |
| 93 | + integer :: edge_meth |
| 94 | + integer :: cell_meth |
| 95 | + integer :: cell_lims |
| 96 | + integer :: wall_lims |
| 97 | + end type rcon_opts |
| 98 | + |
| 99 | + type rcon_ends |
| 100 | + !------------------------------- end-conditions struct ! |
| 101 | + integer :: bcopt |
| 102 | + real*8 :: value |
| 103 | + real*8 :: slope |
| 104 | + end type rcon_ends |
| 105 | + |
| 106 | + type rcon_work |
| 107 | + !------------------------------- work-space for RCON1D ! |
| 108 | + real*8, allocatable :: edge_func(:,:) |
| 109 | + real*8, allocatable :: edge_dfdx(:,:) |
| 110 | + real*8, allocatable :: cell_oscl(:,:,:) |
| 111 | + contains |
| 112 | + procedure :: init => init_rcon_work |
| 113 | + procedure :: free => free_rcon_work |
| 114 | + end type rcon_work |
| 115 | + |
| 116 | + type, extends(rcon_opts) :: rmap_opts |
| 117 | + !------------------------------- parameters for RMAP1D ! |
| 118 | + end type rmap_opts |
| 119 | + |
| 120 | + type, extends(rcon_work) :: rmap_work |
| 121 | + !------------------------------- work-space for RMAP1D ! |
| 122 | + real*8, allocatable :: cell_spac(:) |
| 123 | + real*8, allocatable :: cell_func(:,:,:) |
| 124 | + contains |
| 125 | + procedure :: init => init_rmap_work |
| 126 | + procedure :: free => free_rmap_work |
| 127 | + end type rmap_work |
| 128 | + |
| 129 | + contains |
| 130 | + |
| 131 | + !------------------------------------------------------! |
| 132 | + ! INIT-RCON-WORK: init. work-space for RCON1D. ! |
| 133 | + !------------------------------------------------------! |
| 134 | + |
| 135 | + subroutine init_rcon_work(this,npos,nvar,opts) |
| 136 | + |
| 137 | + ! |
| 138 | + ! THIS work-space structure for RCON1D . |
| 139 | + ! NPOS no. edges over grid. |
| 140 | + ! NVAR no. state variables. |
| 141 | + ! OPTS parameters structure for RCON1D . |
| 142 | + ! |
| 143 | + |
| 144 | + implicit none |
| 145 | + |
| 146 | + !------------------------------------------- arguments ! |
| 147 | + class(rcon_work) , intent(inout) :: this |
| 148 | + integer, intent(in):: npos |
| 149 | + integer, intent(in):: nvar |
| 150 | + class(rcon_opts) , optional :: opts |
| 151 | + |
| 152 | + !------------------------------------------- variables ! |
| 153 | + integer :: okay |
| 154 | + |
| 155 | + allocate(this% & |
| 156 | + & edge_func( nvar,npos), & |
| 157 | + & this% & |
| 158 | + & edge_dfdx( nvar,npos), & |
| 159 | + & this% & |
| 160 | + & cell_oscl(2,nvar,npos), & |
| 161 | + & stat=okay) |
| 162 | + |
| 163 | + end subroutine |
| 164 | + |
| 165 | + !------------------------------------------------------! |
| 166 | + ! INIT-RMAP-WORK: init. work-space for RMAP1D. ! |
| 167 | + !------------------------------------------------------! |
| 168 | + |
| 169 | + subroutine init_rmap_work(this,npos,nvar,opts) |
| 170 | + |
| 171 | + ! |
| 172 | + ! THIS work-space structure for RMAP1D . |
| 173 | + ! NPOS no. edges over grid. |
| 174 | + ! NVAR no. state variables. |
| 175 | + ! OPTS parameters structure for RMAP1D . |
| 176 | + ! |
| 177 | + |
| 178 | + implicit none |
| 179 | + |
| 180 | + !------------------------------------------- arguments ! |
| 181 | + class(rmap_work) , intent(inout) :: this |
| 182 | + integer, intent(in) :: npos |
| 183 | + integer, intent(in) :: nvar |
| 184 | + class(rcon_opts) , optional :: opts |
| 185 | + |
| 186 | + !------------------------------------------- variables ! |
| 187 | + integer :: okay,ndof |
| 188 | + |
| 189 | + ndof = ndof1d(opts%cell_meth) |
| 190 | + |
| 191 | + allocate(this% & |
| 192 | + & edge_func( nvar,npos), & |
| 193 | + & this% & |
| 194 | + & edge_dfdx( nvar,npos), & |
| 195 | + & this% & |
| 196 | + & cell_oscl(2,nvar,npos), & |
| 197 | + & this% & |
| 198 | + & cell_spac( npos), & |
| 199 | + & this% & |
| 200 | + & cell_func(ndof,nvar,npos) , & |
| 201 | + & stat=okay) |
| 202 | + |
| 203 | + end subroutine |
| 204 | + |
| 205 | + !------------------------------------------------------! |
| 206 | + ! FREE-RCON-WORK: free work-space for RCON1D . ! |
| 207 | + !------------------------------------------------------! |
| 208 | + |
| 209 | + subroutine free_rcon_work(this) |
| 210 | + |
| 211 | + implicit none |
| 212 | + |
| 213 | + !------------------------------------------- arguments ! |
| 214 | + class(rcon_work), intent(inout) :: this |
| 215 | + |
| 216 | + deallocate(this%edge_func, & |
| 217 | + & this%edge_dfdx, & |
| 218 | + & this%cell_oscl) |
| 219 | + |
| 220 | + end subroutine |
| 221 | + |
| 222 | + !------------------------------------------------------! |
| 223 | + ! FREE-RMAP-WORK: free work-space for RMAP1D . ! |
| 224 | + !------------------------------------------------------! |
| 225 | + |
| 226 | + subroutine free_rmap_work(this) |
| 227 | + |
| 228 | + implicit none |
| 229 | + |
| 230 | + !------------------------------------------- arguments ! |
| 231 | + class(rmap_work), intent(inout) :: this |
| 232 | + |
| 233 | + |
| 234 | + deallocate(this%edge_func, & |
| 235 | + & this%edge_dfdx, & |
| 236 | + & this%cell_oscl, & |
| 237 | + & this%cell_func, & |
| 238 | + & this%cell_spac) |
| 239 | + |
| 240 | + end subroutine |
| 241 | + |
| 242 | + !------------------------------------------------------! |
| 243 | + ! NDOF1D : no. degrees-of-freedom per polynomial . ! |
| 244 | + !------------------------------------------------------! |
| 245 | + |
| 246 | + pure function ndof1d(meth) result(rdof) |
| 247 | + |
| 248 | + implicit none |
| 249 | + |
| 250 | + !------------------------------------------- arguments ! |
| 251 | + integer, intent( in) :: meth |
| 252 | + |
| 253 | + !------------------------------------------- variables ! |
| 254 | + integer :: rdof |
| 255 | + |
| 256 | + select case(meth) |
| 257 | + !-------------------------------- edge reconstructions ! |
| 258 | + case (p1e_method) |
| 259 | + rdof = +2 |
| 260 | + case (p3e_method) |
| 261 | + rdof = +4 |
| 262 | + case (p5e_method) |
| 263 | + rdof = +6 |
| 264 | + !-------------------------------- cell reconstructions ! |
| 265 | + case (pcm_method) |
| 266 | + rdof = +1 |
| 267 | + case (plm_method) |
| 268 | + rdof = +2 |
| 269 | + case (ppm_method) |
| 270 | + rdof = +3 |
| 271 | + case (pqm_method) |
| 272 | + rdof = +5 |
| 273 | + |
| 274 | + case default |
| 275 | + rdof = +0 |
| 276 | + |
| 277 | + end select |
| 278 | + |
| 279 | + end function ndof1d |
| 280 | + |
| 281 | + !------------------------------------------------------! |
| 282 | + ! BFUN1D : one-dimensional poly. basis-functions . ! |
| 283 | + !------------------------------------------------------! |
| 284 | + |
| 285 | +# include "bfun1d.f90" |
| 286 | + |
| 287 | + !------------------------------------------------------! |
| 288 | + ! UTIL1D : one-dimensional grid manip. utilities . ! |
| 289 | + !------------------------------------------------------! |
| 290 | + |
| 291 | +# include "util1d.f90" |
| 292 | + |
| 293 | + !------------------------------------------------------! |
| 294 | + ! WENO1D : "essentially" non-oscillatory limiter . ! |
| 295 | + !------------------------------------------------------! |
| 296 | + |
| 297 | +# include "weno1d.f90" |
| 298 | + |
| 299 | +# include "oscl1d.f90" |
| 300 | + |
| 301 | + !------------------------------------------------------! |
| 302 | + ! RCON1D : one-dimensional poly. reconstructions . ! |
| 303 | + !------------------------------------------------------! |
| 304 | + |
| 305 | +# include "rcon1d.f90" |
| 306 | + |
| 307 | +# include "inv.f90" |
| 308 | + |
| 309 | +# include "pbc.f90" |
| 310 | +# include "p1e.f90" |
| 311 | +# include "p3e.f90" |
| 312 | +# include "p5e.f90" |
| 313 | + |
| 314 | +# include "root1d.f90" |
| 315 | + |
| 316 | +# include "pcm.f90" |
| 317 | +# include "plm.f90" |
| 318 | +# include "ppm.f90" |
| 319 | +# include "pqm.f90" |
| 320 | + |
| 321 | + !------------------------------------------------------! |
| 322 | + ! RMAP1D : one-dimensional conservative "re-map" . ! |
| 323 | + !------------------------------------------------------! |
| 324 | + |
| 325 | +# include "rmap1d.f90" |
| 326 | + |
| 327 | + !------------------------------------------------------! |
| 328 | + ! FFSL1D : one-dimensional FFSL scalar transport . ! |
| 329 | + !------------------------------------------------------! |
| 330 | + |
| 331 | +# include "ffsl1d.f90" |
| 332 | + |
| 333 | + |
| 334 | + !------------------------------------------ end ppr_1d ! |
| 335 | + |
| 336 | +# undef __TIC__ |
| 337 | +# undef __TOC__ |
| 338 | + |
| 339 | + end module |
| 340 | + |
| 341 | + |
| 342 | + |
0 commit comments