[13926] | 1 | |
---|
| 2 | ! gfortran -cpp -O3 -flto ex_5.f90 -o ex_5 |
---|
| 3 | ! ./ex_5 |
---|
| 4 | |
---|
| 5 | ! Assemble high-order interpolants over a uniform domain. |
---|
| 6 | ! |
---|
| 7 | |
---|
| 8 | # include "../src/ppr_1d.f90" |
---|
| 9 | |
---|
| 10 | program ex |
---|
| 11 | |
---|
| 12 | use ppr_1d |
---|
| 13 | |
---|
| 14 | implicit none |
---|
| 15 | |
---|
| 16 | integer, parameter :: npos = 43 ! no. edge |
---|
| 17 | integer, parameter :: nvar = 1 ! no. variables to build |
---|
| 18 | integer, parameter :: ndof = 1 ! no. FV DoF per cell |
---|
| 19 | integer :: ipos,jpos,mdof |
---|
| 20 | |
---|
| 21 | !------------------------------- domain discretisation ! |
---|
| 22 | real*8 :: xpos(npos),xdel(1) |
---|
| 23 | real*8 :: xmid,xhat,xloc,floc |
---|
| 24 | |
---|
| 25 | !-------------------------------- finite-volume arrays ! |
---|
| 26 | |
---|
| 27 | ! Arrays represent a "block" of finite-volume tracers |
---|
| 28 | ! to remap. The 1st dim. is the no. of DoF per cell, |
---|
| 29 | ! NDOF=1 is a standard finite-volume scheme where the |
---|
| 30 | ! data is specified as cell means. NDOF>1 is reserved |
---|
| 31 | ! for future use with DG-style schemes. NVAR is the |
---|
| 32 | ! number of tracers to remap. Processing tracers in a |
---|
| 33 | ! batch is typically more efficient than one-by-one. |
---|
| 34 | ! The last dim. is the no. cells (layers) in the grid. |
---|
| 35 | |
---|
| 36 | real*8 :: fdat(ndof,nvar,npos-1) |
---|
| 37 | |
---|
| 38 | !-------------------------------- reconstruction coeff ! |
---|
| 39 | |
---|
| 40 | ! Coeff. for the piecewise polynomial reconstruction. |
---|
| 41 | ! A polynomial is assembled for each cell w.r.t. a |
---|
| 42 | ! "local" cell coordinate system: each cell is mapped |
---|
| 43 | ! onto [-1,+1]. The interpolants can be evaluated by |
---|
| 44 | ! taking the product FHAT*BVEC, where BVEC is a basis |
---|
| 45 | ! vector assembled at the interpolation points. Basis |
---|
| 46 | ! vectors can eb assembled via calls to BFUN1D(). |
---|
| 47 | |
---|
| 48 | real*8 :: fhat( 5,nvar,npos-1) |
---|
| 49 | real*8 :: bvec( 5) |
---|
| 50 | real*8 :: spos( 5) |
---|
| 51 | |
---|
| 52 | !------------------------------ method data-structures ! |
---|
| 53 | type(rcon_work) :: work |
---|
| 54 | type(rcon_opts) :: opts |
---|
| 55 | type(rcon_ends) :: bc_l(nvar) |
---|
| 56 | type(rcon_ends) :: bc_r(nvar) |
---|
| 57 | |
---|
| 58 | !------------------------------ define a simple domain ! |
---|
| 59 | |
---|
| 60 | call linspace(0.d0,1.d0,npos,xpos) |
---|
| 61 | |
---|
| 62 | xdel(1) = (xpos(npos)& |
---|
| 63 | & - xpos( 1)) / (npos- 1) |
---|
| 64 | |
---|
| 65 | !------------------------------ setup some simple data ! |
---|
| 66 | |
---|
| 67 | do ipos = +1, npos-1 |
---|
| 68 | |
---|
| 69 | xmid = xpos(ipos+0) * 0.5d+0 & |
---|
| 70 | & + xpos(ipos+1) * 0.5d+0 |
---|
| 71 | |
---|
| 72 | fdat(1,1,ipos) = & |
---|
| 73 | & .8d+0 * exp( -75.0d+0 * (xmid - 0.275d+0) ** 2 ) & |
---|
| 74 | & + .9d+0 * exp(-100.0d+0 * (xmid - 0.500d+0) ** 2 ) & |
---|
| 75 | & + 1.d+0 * exp(-125.0d+0 * (xmid - 0.725d+0) ** 2 ) |
---|
| 76 | |
---|
| 77 | end do |
---|
| 78 | |
---|
| 79 | !------------------------------ specify method options ! |
---|
| 80 | |
---|
| 81 | opts%edge_meth = p5e_method ! 5th-order edge interp. |
---|
| 82 | opts%cell_meth = pqm_method ! PPM method in cells |
---|
| 83 | opts%cell_lims = mono_limit ! monotone limiter |
---|
| 84 | |
---|
| 85 | !------------------------------ set BC.'s at endpoints ! |
---|
| 86 | |
---|
| 87 | bc_l%bcopt = bcon_loose ! "loose" = extrapolate |
---|
| 88 | bc_r%bcopt = bcon_loose |
---|
| 89 | |
---|
| 90 | !------------------------------ init. method workspace ! |
---|
| 91 | |
---|
| 92 | call work%init(npos,nvar,opts) |
---|
| 93 | |
---|
| 94 | !------------------------------ build cell polynomials ! |
---|
| 95 | |
---|
| 96 | fhat = 0.d+0 |
---|
| 97 | |
---|
| 98 | mdof = ndof1d (opts%cell_meth) |
---|
| 99 | |
---|
| 100 | call rcon1d(npos,nvar,ndof,xdel, & |
---|
| 101 | & fdat,bc_l,bc_r,fhat, & |
---|
| 102 | & work,opts) |
---|
| 103 | |
---|
| 104 | !------------------------------ clear method workspace ! |
---|
| 105 | |
---|
| 106 | call work%free() |
---|
| 107 | |
---|
| 108 | !------------------------------ dump results to stdout ! |
---|
| 109 | |
---|
| 110 | print*,"Eval. PPR interpolant: " |
---|
| 111 | |
---|
| 112 | spos(1) = -1.0d+0 ! eval. at local points |
---|
| 113 | spos(2) = -0.5d+0 |
---|
| 114 | spos(3) = +0.0d+0 |
---|
| 115 | spos(4) = +0.5d+0 |
---|
| 116 | spos(5) = +1.0d+0 |
---|
| 117 | |
---|
| 118 | do ipos = +1, npos-1 |
---|
| 119 | do jpos = +1, +5 |
---|
| 120 | |
---|
| 121 | xmid = xpos(ipos+1)* 0.5d+0 & |
---|
| 122 | & + xpos(ipos+0)* 0.5d+0 |
---|
| 123 | |
---|
| 124 | xhat = xpos(ipos+1)* 0.5d+0 & |
---|
| 125 | & - xpos(ipos+0)* 0.5d+0 |
---|
| 126 | |
---|
| 127 | xloc = xmid + spos(jpos)*xhat |
---|
| 128 | |
---|
| 129 | call bfun1d(0,mdof,spos(jpos),bvec) |
---|
| 130 | |
---|
| 131 | floc = dot_product( & |
---|
| 132 | & fhat(+1:mdof,1,ipos),bvec(+1:mdof)) |
---|
| 133 | |
---|
| 134 | print *, xloc, floc |
---|
| 135 | |
---|
| 136 | end do |
---|
| 137 | end do |
---|
| 138 | |
---|
| 139 | end program |
---|
| 140 | |
---|
| 141 | |
---|
| 142 | |
---|