New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
ex_6.f90 in vendors/PPR/example – NEMO

source: vendors/PPR/example/ex_6.f90 @ 15787

Last change on this file since 15787 was 13926, checked in by jchanut, 4 years ago

#2222, add Piecewise Polynomial Reconstruction library

File size: 4.6 KB
Line 
1
2!   gfortran -cpp -O3 -flto ex_6.f90 -o ex_6
3!   ./ex_6
4
5!   1d scalar transport in a periodic domain. Fluxes are
6!   computed in a flux-form semi-lagrangian sense --
7!   integrating over the upwind regions "swept" by edges
8!   in each time-step. This implementation requires CFL<1,
9!   with the upwind regions covering adjacent cells only.
10!   A CFL>=1 variant could be constructed using RMAP1D().
11!
12
13#   include "../src/ppr_1d.f90"
14
15    program ex
16
17        use ppr_1d
18
19        implicit none
20
21        integer, parameter :: halo = 4  ! halo cells at boundary
22        integer, parameter :: npos = 51 ! no. edge
23        integer, parameter :: nvar = 1  ! no. tracers to trnsprt
24        integer, parameter :: ndof = 1  ! no. FV DoF per cell
25        integer :: ipos,il,ir,step
26
27    !------------------------------- domain discretisation !
28        real*8  :: xpos(1-halo:npos+halo)
29        real*8  :: xmid,xdel(1),tDEL
30       
31    !-------------------------------- finite-volume arrays !
32
33    !   INIT: initial cell-wise finite-volume profile
34    !   QBAR: dynamic cell-wise finite-volume profile
35    !   MASK: cell-wise "land" mask (all TRUE here)
36    !   UVEL: edge-wise velocity distribution
37    !   FLUX: edge-wise distribution of upwind fluxes
38    !   QDIV: cell-wise distribution of divergence
39
40        real*8  :: init(ndof,nvar,1-halo:npos+halo-1)
41        real*8  :: qbar(ndof,nvar,1-halo:npos+halo-1)
42        logical :: mask(          1-halo:npos+halo-1)
43        real*8  :: uvel(          1-halo:npos+halo)       
44        real*8  :: flux(     nvar,1-halo:npos+halo)       
45        real*8  :: qdiv(     nvar,1-halo:npos+halo-1)
46       
47
48    !------------------------------ method data-structures !
49        type(rmap_work) :: work
50        type(rmap_opts) :: opts
51        type(rcon_ends) :: bc_l(nvar)
52        type(rcon_ends) :: bc_r(nvar)
53
54    !------------------------------ define a simple domain !
55
56        il =      + 1                   ! 1st real interior cell
57        ir = npos - 1                   ! Nth real interior cell
58
59        xpos(il+0) = 0.0d+00
60        xpos(ir+1) = 1.0d+00
61       
62        xdel(1) = (xpos(ir+1)-xpos(il+0))/(npos-1)
63
64        do ipos = il+1, ir-0
65
66            xpos(ipos) = (ipos-1) * xdel(1)         
67
68        end do
69
70    !------------------------------ setup some simple data !
71
72        uvel       = +1.0d+0
73        mask       = .true.
74
75        tDEL       = +1.0d-2
76
77        do ipos = +1, npos-1
78
79            xmid = xpos(ipos+0)* 0.5d+00 &
80    &            + xpos(ipos+1)* 0.5d+00
81
82            init(1,1,ipos) = &
83    &   .8d+0 * exp( -75.0d+0 * (xmid - 0.275d+0) ** 2 ) &
84    & + .9d+0 * exp(-100.0d+0 * (xmid - 0.500d+0) ** 2 ) &
85    & + 1.d+0 * exp(-125.0d+0 * (xmid - 0.725d+0) ** 2 )
86           
87        end do
88
89    !------------------------------ specify method options !
90
91        opts%edge_meth = p3e_method     ! 3rd-order edge interp.
92        opts%cell_meth = ppm_method     ! PPM method in cells
93        opts%cell_lims = weno_limit     ! "non-oscillatory" lim.
94       
95    !------------------------------ set BC.'s at endpoints !
96
97        bc_l%bcopt = bcon_loose         ! "loose" = extrapolate
98        bc_r%bcopt = bcon_loose
99
100    !------------------------------ init. method workspace !
101
102        call work%init(npos+2*halo,nvar,opts)
103
104    !------------------------------ calc. scalar transport !
105
106        qbar = init
107
108        do step = +1, +100        ! 100 steps => full loop
109
110    !------------------------------ periodicity via halo's !
111
112            qbar(1,:,il-halo:il-1) = & 
113    &       qbar(1,:,ir-halo+1:ir)
114            qbar(1,:,ir+1:ir+halo) = &
115    &       qbar(1,:,il:il+halo-1)
116
117    !------------------------------ form lagrangian fluxes !
118
119            call ffsl1d (npos+2*halo,nvar, &
120    &            ndof,xdel,tDEL,mask,uvel, &
121    &            qbar,flux,bc_l,bc_r,work, &
122    &            opts)
123
124 
125    !------------------------------ flux divergence eval's !
126            qdiv(  1,il:ir) = &
127    &              flux(1,il+1:ir+1) &
128    &            - flux(1,il+0:ir+0)
129
130    !------------------------------ take a single timestep !
131
132            qbar(1,1,il:ir) = &
133    &           qbar(1,1,il:ir) - &
134    &               qdiv(1,il:ir) / xdel(1)
135
136        end do
137
138    !------------------------------ clear method workspace !
139
140        call work%free()
141
142    !------------------------------ dump results to stdout !
143
144        print*,"End timestep profile : "
145
146        do ipos = il+0, ir-0
147
148            print *, init(1,:,ipos)  &
149    &              , qbar(1,:,ipos)
150
151        end do
152
153        print*,"Conservation defect := " &
154    &         , sum(init(1,:,il:ir)) &
155    &         - sum(qbar(1,:,il:ir))
156       
157    end program
158
159
160
Note: See TracBrowser for help on using the repository browser.