1 | |
---|
2 | ! gfortran -cpp -O3 -flto ex_3.f90 -o ex_3 |
---|
3 | ! ./ex_3 |
---|
4 | |
---|
5 | ! "Convergence" testing: remap a smooth (Gaussian) profile |
---|
6 | ! between perturbed grids and compute error in the L2-nrm. |
---|
7 | ! Re-run with N = 50, 100, 200, ... etc (and N-tmp at 80%) |
---|
8 | ! to assess the order of spatial accuracy. With NUll/WENO- |
---|
9 | ! LIMIT, all methods should achieve their best-case O(h^p) |
---|
10 | ! scaling. |
---|
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 :: npos = 50 ! no. edge (old grid) |
---|
22 | integer, parameter :: ntmp = 40 ! no. edge (new grid) |
---|
23 | integer, parameter :: nvar = 1 ! no. variables to remap |
---|
24 | integer, parameter :: ndof = 1 ! no. FV DoF per cell |
---|
25 | integer :: ipos |
---|
26 | |
---|
27 | !------------------------------ position of cell edges ! |
---|
28 | real*8 :: xpos(npos),xtmp(ntmp) |
---|
29 | real*8 :: xdel,xmid |
---|
30 | |
---|
31 | !-------------------------------- finite-volume arrays ! |
---|
32 | |
---|
33 | ! Arrays represent a "block" of finite-volume tracers |
---|
34 | ! to remap. The 1st dim. is the no. of DoF per cell, |
---|
35 | ! NDOF=1 is a standard finite-volume scheme where the |
---|
36 | ! data is specified as cell means. NDOF>1 is reserved |
---|
37 | ! for future use with DG-style schemes. NVAR is the |
---|
38 | ! number of tracers to remap. Processing tracers in a |
---|
39 | ! batch is typically more efficient than one-by-one. |
---|
40 | ! The last dim. is the no. cells (layers) in the grid. |
---|
41 | |
---|
42 | real*8 :: init(ndof,nvar,npos-1) |
---|
43 | real*8 :: ftmp(ndof,nvar,ntmp-1) |
---|
44 | real*8 :: fdat(ndof,nvar,npos-1) |
---|
45 | |
---|
46 | !------------------------------ method data-structures ! |
---|
47 | type(rmap_work) :: work |
---|
48 | type(rmap_opts) :: opts |
---|
49 | type(rcon_ends) :: bc_l(nvar) |
---|
50 | type(rcon_ends) :: bc_r(nvar) |
---|
51 | |
---|
52 | !------------------------------ define a simple domain ! |
---|
53 | |
---|
54 | call linspace(0.d0,1.d0,npos,xpos) |
---|
55 | call linspace(0.d0,1.d0,ntmp,xtmp) |
---|
56 | |
---|
57 | !------------------------------ setup some simple data ! |
---|
58 | |
---|
59 | do ipos = +1, npos-1 |
---|
60 | |
---|
61 | xmid = xpos(ipos+0) * 0.5d+0 & |
---|
62 | & + xpos(ipos+1) * 0.5d+0 |
---|
63 | |
---|
64 | init(1,1,ipos) = & |
---|
65 | & .8d+0 * exp( -75.0d+0 * (xmid - 0.275d+0) ** 2 ) & |
---|
66 | & + .9d+0 * exp(-100.0d+0 * (xmid - 0.500d+0) ** 2 ) & |
---|
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 ! 3rd-order edge interp. |
---|
74 | opts%cell_meth = pqm_method ! PPM method in cells |
---|
75 | opts%cell_lims = null_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,1,ipos) & |
---|
117 | & , fdat(1,1,ipos) |
---|
118 | |
---|
119 | end do |
---|
120 | |
---|
121 | print*,"Conservation defect := " & |
---|
122 | & , sum(init) - sum(fdat) |
---|
123 | |
---|
124 | !------------------------------ calc. errors in L2 nrm ! |
---|
125 | |
---|
126 | xdel = (xpos(npos)-xpos(1))/(npos-1) |
---|
127 | |
---|
128 | print*,"Discrete (L2) Error := ",& |
---|
129 | sqrt(sum(xdel * (init(1,1,:) - fdat(1,1,:)) ** 2)) |
---|
130 | |
---|
131 | end program |
---|
132 | |
---|
133 | |
---|
134 | |
---|