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_4.f90 in vendors/PPR/example – NEMO

source: vendors/PPR/example/ex_4.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: 3.5 KB
Line 
1
2!   gfortran -cpp -O3 -flto ex_4.f90 -o ex_4
3!   ./ex_4
4
5!   Test for multi-tracer remapping: remap a set of profiles
6!   between randomly perturbed grids.
7!
8
9#   include "../src/ppr_1d.f90"
10
11    program ex
12
13        use ppr_1d
14
15        implicit none
16
17        integer, parameter :: npos = 53 ! no. edge (old grid)
18        integer, parameter :: ntmp = 43 ! no. edge (new grid)
19        integer, parameter :: nvar = 3  ! no. variables to remap
20        integer, parameter :: ndof = 1  ! no. FV DoF per cell
21        integer :: ipos
22
23    !------------------------------ position of cell edges !
24        real*8  :: xpos(npos),xtmp(ntmp)
25        real*8  :: xdel,xmid
26       
27    !-------------------------------- finite-volume arrays !
28
29    !   Arrays represent a "block" of finite-volume tracers
30    !   to remap. The 1st dim. is the no. of DoF per cell,
31    !   NDOF=1 is a standard finite-volume scheme where the
32    !   data is specified as cell means. NDOF>1 is reserved
33    !   for future use with DG-style schemes. NVAR is the
34    !   number of tracers to remap. Processing tracers in a
35    !   batch is typically more efficient than one-by-one.
36    !   The last dim. is the no. cells (layers) in the grid.
37
38        real*8  :: init(ndof,nvar,npos-1)
39        real*8  :: ftmp(ndof,nvar,ntmp-1)
40        real*8  :: fdat(ndof,nvar,npos-1)
41
42    !------------------------------ method data-structures !
43        type(rmap_work) :: work
44        type(rmap_opts) :: opts
45        type(rcon_ends) :: bc_l(nvar)
46        type(rcon_ends) :: bc_r(nvar)
47
48    !------------------------------ define a simple domain !
49
50        call linspace(0.d0,1.d0,npos,xpos)
51        call rndspace(0.d0,1.d0,ntmp,xtmp)
52
53    !------------------------------ setup some simple data !
54
55        do ipos = +1, npos-1
56
57            xmid = xpos(ipos+0) * 0.5d+0 &
58    &            + xpos(ipos+1) * 0.5d+0
59
60            init(1,1,ipos) = &
61    &   .8d+0 * exp( -75.0d+0 * (xmid - 0.275d+0) ** 2 )
62           
63            init(1,2,ipos) = &
64    & + .9d+0 * exp(-100.0d+0 * (xmid - 0.500d+0) ** 2 )
65           
66            init(1,3,ipos) = &
67    & + 1.d+0 * exp(-125.0d+0 * (xmid - 0.725d+0) ** 2 )
68           
69        end do
70
71    !------------------------------ specify method options !
72
73        opts%edge_meth = p5e_method     ! 5th-order edge interp.
74        opts%cell_meth = pqm_method     ! PQM method in cells
75        opts%cell_lims = mono_limit     ! monotone limiter
76       
77    !------------------------------ set BC.'s at endpoints !
78
79        bc_l%bcopt = bcon_loose         ! "loose" = extrapolate
80        bc_r%bcopt = bcon_loose
81
82    !------------------------------ init. method workspace !
83
84        call work%init(npos,nvar,opts)
85
86    !------------------------------ re-map back-and-forth: !
87
88        fdat = init
89
90        do ipos = +1, +1000
91
92    !------------------------------ re-map from dat-to-tmp !
93
94        call rmap1d(npos,ntmp,nvar,ndof, &
95        &           xpos,xtmp,fdat,ftmp, &
96        &           bc_l,bc_r,work,opts)
97
98    !------------------------------ re-map from tmp-to-dat !
99
100        call rmap1d(ntmp,npos,nvar,ndof, &
101        &           xtmp,xpos,ftmp,fdat, &
102        &           bc_l,bc_r,work,opts)
103
104        end do
105
106    !------------------------------ clear method workspace !
107
108        call work%free()
109
110    !------------------------------ dump results to stdout !
111
112        print*,"Cell data: [INIT] [RMAP] "
113
114        do ipos = +1, npos-1
115
116            print *, init(1,:,ipos) &
117        &          , fdat(1,:,ipos)
118
119        end do
120
121        print*,"Conservation defect := " &
122        &     , sum(init) - sum(fdat)
123
124    end program
125
126
127
Note: See TracBrowser for help on using the repository browser.