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

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