Recon1d_PLM_WLS.F90
1!> Piecewise Linear Method using Weighted Conservative Least Squares 1D reconstruction
2module recon1d_plm_wls
3
4! This file is part of MOM6. See LICENSE.md for the license.
5
6use recon1d_type, only : recon1d, testing
7
8implicit none ; private
9
10public plm_wls, testing
11
12!> PLM reconstruction using Weighted Least Squares constrained to conserve for central cell
13!!
14!! The source for the methods ultimately used by this class are:
15!! - init() *locally defined
16!! - reconstruct() *locally defined
17!! - average() *locally defined
18!! - f() *locally defined
19!! - dfdx() *locally defined
20!! - x() -> recon1d_type.x()
21!! - check_reconstruction() *locally defined
22!! - unit_tests() *locally defined
23!! - destroy() *locally defined
24!! - remap_to_sub_grid() -> recon1d_type.remap_to_sub_grid()
25!! - init_parent() -> init()
26!! - reconstruct_parent() -> reconstruct()
28
29 real, allocatable :: ul(:) !< Left edge value [A]
30 real, allocatable :: ur(:) !< Right edge value [A]
31 real, allocatable, private :: slp(:) !< Difference across cell, ur - ul [A].
32 !! This is redundant with ul and ur and not used
33 !! in any evaluations, but is needed for testing.
34
35contains
36 !> Implementation of the PLM_WLS initialization
38 !> Implementation of the PLM_WLS reconstruction
39 procedure :: reconstruct => reconstruct
40 !> Implementation of the PLM_WLS average over an interval [A]
42 !> Implementation of evaluating the PLM_WLS reconstruction at a point [A]
44 !> Implementation of the derivative of the PLM_WLS reconstruction at a point [A]
46 !> Implementation of deallocation for PLM_WLS
48 !> Implementation of check reconstruction for the PLM_WLS reconstruction
49 procedure :: check_reconstruction => check_reconstruction
50 !> Implementation of unit tests for the PLM_WLS reconstruction
51 procedure :: unit_tests => unit_tests
52
53 !> Duplicate interface to init()
54 procedure :: init_parent => init
55 !> Duplicate interface to reconstruct()
56 procedure :: reconstruct_parent => reconstruct
57
58end type plm_wls
59
60contains
61
62!> Initialize a 1D PLM reconstruction for n cells
63subroutine init(this, n, h_neglect, check)
64 class(plm_wls), intent(out) :: this !< This reconstruction
65 integer, intent(in) :: n !< Number of cells in this column
66 real, optional, intent(in) :: h_neglect !< A negligibly small width used in cell reconstructions [H]
67 logical, optional, intent(in) :: check !< If true, enable some consistency checking
68
69 this%n = n
70
71 allocate( this%u_mean(n) )
72 allocate( this%ul(n) )
73 allocate( this%ur(n) )
74 allocate( this%slp(n) )
75
76 this%h_neglect = tiny( this%u_mean(1) )
77 if (present(h_neglect)) this%h_neglect = h_neglect
78 this%check = .false.
79 if (present(check)) this%check = check
80
81end subroutine init
82
83!> Calculate a 1D PLM_WLS reconstruction based on h(:) and u(:)
84subroutine reconstruct(this, h, u)
85 class(plm_wls), intent(inout) :: this !< This reconstruction
86 real, intent(in) :: h(*) !< Grid spacing (thickness) [typically H]
87 real, intent(in) :: u(*) !< Cell mean values [A]
88 ! Local variables
89 real :: slp ! The PLM slopes (difference across cell) [A]
90 real :: u_l, u_r, u_c ! Left, right, and center values [A]
91 real :: h_l, h_c, h_r ! Thickness of left, center and right cells [H]
92 real :: h_l0, h_r0 ! Thickness of left and right cells with h_neglect added [H]
93 real :: hx2l, hx2r ! Contributions to denominator, <h x^2> [H3]
94 real :: hxyl, hxyr ! Contributions to numerator, <h x y> [H2 A]
95 integer :: n, km1, k, kp1
96
97 n = this%n
98
99 ! Loop over all cells
100 do k = 1, n
101 km1 = max(1, k-1)
102 kp1 = min(n, k+1)
103 u_l = u(km1)
104 u_c = u(k)
105 u_r = u(kp1)
106
107 h_l = h(km1) * real( k - km1 ) ! This zeroes h_l at k==1
108 h_c = h(k)
109 h_r = h(kp1) * real( kp1 - k ) ! This zeroes h_r at k==n
110
111 ! This is the slope that minimizes the error
112 ! sum_l={-1,1} h(k+l) * [ u(k+l) - u(k) + slp * ( z(k+l) - z(k) ) ]
113 ! i.e. volume weighted least squares
114 h_l0 = h_l + this%h_neglect
115 h_r0 = h_r + this%h_neglect
116 hxyl = ( h_l * ( h_c + h_l ) ) * ( u_c - u_l )
117 hxyr = ( h_r * ( h_c + h_r ) ) * ( u_r - u_c )
118 hx2l = h_l0 * ( h_c + h_l0 )**2
119 hx2r = h_r0 * ( h_c + h_r0 )**2
120 slp = 2. * h_c * ( hxyr + hxyl ) / ( hx2l + hx2r )
121
122 ! Mean value
123 this%u_mean(k) = u_c
124
125 ! Left edge
126 this%ul(k) = u_c - 0.5 * slp
127
128 ! Right edge
129 this%ur(k) = u_c + 0.5 * slp
130
131 ! Store slope
132 this%slp(k) = slp
133 enddo
134
135end subroutine reconstruct
136
137!> Value of PLM_WLS reconstruction at a point in cell k [A]
138real function f(this, k, x)
139 class(plm_wls), intent(in) :: this !< This reconstruction
140 integer, intent(in) :: k !< Cell number
141 real, intent(in) :: x !< Non-dimensional position within element [nondim]
142 real :: du ! Difference across cell [A]
143
144 du = this%ur(k) - this%ul(k)
145
146 ! This expression might be used beyond the element to evaluate
147 ! LS errors. In other PLM implementations x is bounded to the
148 ! element and the expressions are constructed to not exceed
149 ! bounds. There are no such constraints for PLM_WLS.
150 f = this%u_mean(k) + du * ( x - 0.5)
151 !f = this%u_mean(k) + this%slp(k) * ( x - 0.5)
152
153end function f
154
155!> Derivative of PLM_WLS reconstruction at a point in cell k [A]
156real function dfdx(this, k, x)
157 class(plm_wls), intent(in) :: this !< This reconstruction
158 integer, intent(in) :: k !< Cell number
159 real, intent(in) :: x !< Non-dimensional position within element [nondim]
160
161 dfdx = this%ur(k) - this%ul(k)
162
163end function dfdx
164
165!> Average between xa and xb for cell k of a 1D PLM reconstruction [A]
166real function average(this, k, xa, xb)
167 class(plm_wls), intent(in) :: this !< This reconstruction
168 integer, intent(in) :: k !< Cell number
169 real, intent(in) :: xa !< Start of averaging interval on element (0 to 1)
170 real, intent(in) :: xb !< End of averaging interval on element (0 to 1)
171 real :: xmab ! Mid-point between xa and xb (0 to 1)
172 real :: u_a, u_b ! Values at xa and xb [A]
173
174 ! Mid-point between xa and xb
175 xmab = 0.5 * ( xa + xb )
176
177 ! This expression for u_a can overshoot u_r but is good for xmab<<1
178 u_a = this%ul(k) + ( this%ur(k) - this%ul(k) ) * xmab
179 ! This expression for u_b can overshoot u_l but is good for 1-xmab<<1
180 u_b = this%ur(k) + ( this%ul(k) - this%ur(k) ) * ( 1. - xmab )
181
182 ! Since u_a and u_b are both bounded, this will perserve uniformity but will the
183 ! sum be bounded? Emperically it seems to work...
184 average = 0.5 * ( u_a + u_b )
185
186end function average
187
188!> Deallocate the PLM reconstruction
189subroutine destroy(this)
190 class(plm_wls), intent(inout) :: this !< This reconstruction
191
192 deallocate( this%u_mean, this%ul, this%ur )
193
194end subroutine destroy
195
196!> Checks the PLM_WLS reconstruction for consistency
197logical function check_reconstruction(this, h, u)
198 class(plm_wls), intent(in) :: this !< This reconstruction
199 real, intent(in) :: h(*) !< Grid spacing (thickness) [typically H]
200 real, intent(in) :: u(*) !< Cell mean values [A]
201 ! Local variables
202 integer :: k
203 real :: slp ! Cell slope [A]
204 type(plm_wls) :: perturbed !< A perturbed reconstruction
205 real :: u_l, u_r, u_c ! Left, right, and center values [A]
206 real :: h_l, h_c, h_r ! Thickness of left, center and right cells [H]
207 real :: h_l0, h_r0, h_c0 ! Thickness of left, right, center cells with h_neglect added [H]
208 real :: x_l, x_r ! Positions of left and right cells [H]
209 real :: hx2l, hx2r ! Contributions to denominator, <h x^2> [H3]
210 real :: hxyl, hxyr ! Contributions to numerator, <h x y> [H2 A]
211 real :: hy2l, hy2r ! Contributions to error, <h y^2> [H3]
212 real :: y_l, y_r ! Left, right, value differencess [A]
213 real :: b_h, bp_h ! slp / h_c [A H-1]
214 integer :: km1, kp1
215
216 check_reconstruction = .false.
217
218 do k = 1, this%n
219 if ( abs( this%u_mean(k) - u(k) ) > 0. ) check_reconstruction = .true.
220 enddo
221
222 ! Check the cell reconstruction is monotonic within each cell (it should be as a straight line)
223 do k = 1, this%n
224 if ( ( this%u_mean(k) - this%ul(k) ) * ( this%ur(k) - this%u_mean(k) ) < 0. ) check_reconstruction = .true.
225 enddo
226
227 ! Check the cell is a straight line (to within machine precision)
228 do k = 1, this%n
229 if ( abs(2. * this%u_mean(k) - ( this%ul(k) + this%ur(k) )) > epsilon(this%u_mean(1)) * &
230 max(abs(2. * this%u_mean(k)), abs(this%ul(k)), abs(this%ur(k))) ) check_reconstruction = .true.
231 enddo
232
233 ! Create a perturbable reconstruction
234 call perturbed%init( this%n, h_neglect=this%h_neglect )
235 call perturbed%reconstruct( h, u ) ! Should reproduce "this"
236 ! Check the copy is identical
237 do k = 1, this%n
238 if ( abs( perturbed%u_mean(k) - this%u_mean(k) ) > 0. ) check_reconstruction = .true.
239 if ( abs( perturbed%ul(k) - this%ul(k) ) > 0. ) check_reconstruction = .true.
240 if ( abs( perturbed%ur(k) - this%ur(k) ) > 0. ) check_reconstruction = .true.
241 if ( abs( perturbed%slp(k) - this%slp(k) ) > 0. ) check_reconstruction = .true.
242 enddo
243 ! Now perturb the slope. The local error should not decrease.
244 do k = 1, this%n
245 slp = this%slp(k) * ( 1.0 + 1. * epsilon(slp) )
246 perturbed%slp(k) = slp
247 perturbed%ul(k) = u(k) - 0.5 * slp
248 perturbed%ur(k) = u(k) + 0.5 * slp
249 if ( ls_error(perturbed, k, h, u) < ls_error(this, k, h, u) ) check_reconstruction = .true.
250
251 slp = this%slp(k) * ( 1.0 - 1. * epsilon(slp) )
252 perturbed%slp(k) = slp
253 perturbed%ul(k) = u(k) - 0.5 * slp
254 perturbed%ur(k) = u(k) + 0.5 * slp
255 if ( ls_error(perturbed, k, h, u) < ls_error(this, k, h, u) ) check_reconstruction = .true.
256 enddo
257
258end function check_reconstruction
259
260!> Returns local least squares error for a particular cell
261!!
262!! Note that this is the error relative to the minimum of the loss function so that at the
263!! true solution this function returns zero. See module documentation.
264real function ls_error(this, k, h, u)
265 type(plm_wls), intent(in) :: this !< This reconstruction
266 integer, intent(in) :: k !< Cell number
267 real, intent(in) :: h(*) !< Grid spacing (thickness) [typically H]
268 real, intent(in) :: u(*) !< Cell mean values [A]
269 ! Local variables
270 real :: u_l, u_r, u_c ! Left, right, and center values [A]
271 real :: h_l, h_c, h_r ! Thickness of left, center and right cells [H]
272 real :: h_l0, h_r0, hc0 ! Thickness of left, right, center cells with h_neglect added [H]
273 real :: hx2l, hx2r ! Contributions to denominator, <h x^2> [H3]
274 real :: hxyl, hxyr ! Contributions to numerator, <h x y> [H2 A]
275 integer :: km1, kp1
276
277 km1 = max(1, k-1)
278 kp1 = min(this%n, k+1)
279 u_l = u(km1)
280 u_c = u(k)
281 u_r = u(kp1)
282
283 h_l = h(km1) * real( k - km1 ) ! This zeroes h_l at k==1
284 h_r = h(kp1) * real( kp1 - k ) ! This zeroes h_r at k==n
285 h_c = h(k)
286 hc0 = h_c + this%h_neglect
287
288 h_l0 = h_l + this%h_neglect
289 h_r0 = h_r + this%h_neglect
290 hxyl = ( h_l * 0.5 * ( h_c + h_l ) ) * ( u_c - u_l )
291 hxyr = ( h_r * 0.5 * ( h_c + h_r ) ) * ( u_r - u_c )
292 hx2l = h_l0 * 0.25 * ( h_c + h_l0 )**2
293 hx2r = h_r0 * 0.25 * ( h_c + h_r0 )**2
294 ls_error = h_c * ( ( hx2l + hx2r ) * this%slp(k) - h(k) * ( hxyl + hxyr ) )**2
296end function ls_error
297
298!> Runs PLM_WLS reconstruction unit tests and returns True for any fails, False otherwise
299logical function unit_tests(this, verbose, stdout, stderr)
300 class(plm_wls), intent(inout) :: this !< This reconstruction
301 logical, intent(in) :: verbose !< True, if verbose
302 integer, intent(in) :: stdout !< I/O channel for stdout
303 integer, intent(in) :: stderr !< I/O channel for stderr
304 ! Local variables
305 real, allocatable :: ul(:), ur(:), um(:) ! test values [A]
306 real, allocatable :: ull(:), urr(:) ! test values [A]
307 type(testing) :: test ! convenience functions
308 integer :: k
309
310 call test%set( stdout=stdout ) ! Sets the stdout channel in test
311 call test%set( stderr=stderr ) ! Sets the stderr channel in test
312 call test%set( verbose=verbose ) ! Sets the verbosity flag in test
313
314 call this%init(3, h_neglect=1.e-20)
315 call test%test( this%n /= 3, "Setting number of levels")
316 allocate( um(3), ul(3), ur(3), ull(3), urr(3) )
317
318 call this%reconstruct( (/1.,1.,1./), (/-1.,0.,2./) )
319 call test%real_arr(3, this%slp, (/1.,1.5,2./), "(1,1,1)(-1,0,2) slope")
320
321 do k = 1, 3
322 um(k) = ls_error(this, k, (/1.,1.,1./), (/-1.,0.,2./) )
323 enddo
324 call test%real_arr(3, um, (/0.,0.,0./), "(1,1,1)(-1,0,2) LS' rel error")
325
326 call this%reconstruct( (/0.,1.,1./), (/-1.,0.,2./) )
327 call test%real_arr(3, this%slp, (/0.,2.,2./), "(0,1,1)(-1,0,2) slope")
328
329 do k = 1, 3
330 um(k) = ls_error(this, k, (/0.,1.,1./), (/-1.,0.,2./) )
331 enddo
332 call test%real_arr(3, um, (/0.,0.,0./), "(0,1,1)(-1,0,2) LS' rel error")
333
334 call this%reconstruct( (/1.,1.,1./), (/-2.,0.,1./) )
335 call test%real_arr(3, this%slp, (/2.,1.5,1./), "(1,1,1)(-2,0,1) slope")
336
337 call this%reconstruct( (/1.,1.,0./), (/-2.,0.,1./) )
338 call test%real_arr(3, this%slp, (/2.,2.,0./), "(1,1,0)(-2,0,1) slope")
339
340 call this%destroy()
341 call this%init(3) ! Reset to defaults
342
343 ! Straight line data on uniform grid
344 call this%reconstruct( (/2.,2.,2./), (/1.,3.,5./) )
345 call test%real_arr(3, this%u_mean, (/1.,3.,5./), "Straight line data")
346
347 do k = 1, 3
348 ul(k) = this%f(k, 0.)
349 um(k) = this%f(k, 0.5)
350 ur(k) = this%f(k, 1.)
351 enddo
352 call test%real_arr(3, ul, (/0.,2.,4./), "Evaluation on left edge")
353 call test%real_arr(3, um, (/1.,3.,5./), "Evaluation in center")
354 call test%real_arr(3, ur, (/2.,4.,6./), "Evaluation on right edge")
355
356 do k = 1, 3
357 ul(k) = this%dfdx(k, 0.)
358 um(k) = this%dfdx(k, 0.5)
359 ur(k) = this%dfdx(k, 1.)
360 enddo
361 call test%real_arr(3, ul, (/2.,2.,2./), "dfdx on left edge")
362 call test%real_arr(3, um, (/2.,2.,2./), "dfdx in center")
363 call test%real_arr(3, ur, (/2.,2.,2./), "dfdx on right edge")
364
365 do k = 1, 3
366 um(k) = ls_error(this, k, (/2.,2.,2./), (/1.,3.,5./) )
367 enddo
368 call test%real_arr(3, um, (/0.,0.,0./), "Rel error is 0")
369
370 do k = 1, 3
371 um(k) = this%average(k, 0.5, 0.75) ! Average from x=0.5 to 0.75 in each cell
372 enddo
373 call test%real_arr(3, um, (/1.25,3.25,5.25/), "Return interval average")
374
375 call test%real_scalar( this%x(1,0.), 0., 'f-1(1,0)=0')
376 call test%real_scalar( this%x(1,1.), 0.5, 'f-1(1,1)=0.5')
377 call test%real_scalar( this%x(1,3.), 1., 'f-1(1,3)=1')
378 call test%real_scalar( this%x(2,1.), 0., 'f-1(2,1)=0')
379 call test%real_scalar( this%x(2,3.), 0.5, 'f-1(2,3)=0.5')
380 call test%real_scalar( this%x(2,5.), 1., 'f-1(2,5)=1')
381 call test%real_scalar( this%x(3,3.), 0., 'f-1(3,3)=0')
382 call test%real_scalar( this%x(3,5.), 0.5, 'f-1(3,5)=0.5')
383 call test%real_scalar( this%x(3,7.), 1., 'f-1(3,7)=1')
384
385 call this%destroy()
386 deallocate( um, ul, ur, ull, urr )
387
388 allocate( um(4), ul(4), ur(4) )
389 call this%init(4)
390
391 deallocate( um, ul, ur )
392
393 unit_tests = test%summarize("PLM_WLS:unit_tests")
394
395end function unit_tests
396
397!> \namespace recon1d_plm_wls
398!!
399!! This implementation of PLM fits the slope using least squares, but retains conservation
400!! for the central cell by passing through the central value.
401!! Cell-wise reconstructions are NOT limited by neighbours.
402!! Thus, this reconstruction does not yield monotonic profiles needed for the general remapping problem.
403!!
404!! The algorithm solves the least squares problem of fitting a straight line through
405!! the neighboring data. The line is constained to pass through the center cell,
406!! \f$ (x_{k}, y_{k}) \f$, so that the construction is conservative. The more general
407!! function \f$ f(x) = a_{k} + b_{k} x \f$ would not conserve for arbitrary data.
408!!
409!! The unknown parameter \f$ b_{k} \f$ in the line
410!! \f[
411!! f(x) = y_{k} + b_{k} ( x - x_{k} ) / h_{k}
412!! \f]
413!! is fit to neighbors \f$ x_{k-1}, y_{k-1} \f$ and \f$ x_{k+1}, y_{k+1} \f$.
414!!
415!! Denoting \f$ y'_{k+j} = y_{k+j} - y_{k} \f$ and \f$ x'_{k+j} = x_{k+j} - x_{k} \f$
416!! the local error is
417!! \f{align}{
418!! e_{k+j} &= b_k \frac{ x_{k+j} - x_{k} }{ h_{k} } + y_{k} - y_{k+j} \\\\
419!! &= b_k \frac{ x'_{k+j} }{ h_{k} } - y'_{k+j}
420!! \;\; . \f}
421!!
422!! We use volume weighting in the loss
423!! \f[
424!! G(b) = h_{k-1} e_{k-1}^2 + h_{k+1} e_{k+1}^2
425!! \;\; . \f]
426!!
427!! When solving for \f$ b_k \f$, we solve \f$ dG/db = 0 \f$ where
428!! \f{align}{
429!! dG/db &= 2 h_{k-1} e_{k-1} \frac{ de_{k-1} }{db} + 2 h_{k+1} e_{k+1} \frac{ de_{k+1} }{db} \\\\
430!! &= 2 h_{k-1} ( b_k \frac{ x'_{k-1} }{ h_{k} } - \frac{ y'_{k-1} ) x'_{k-1} }{ h_{k} } +
431!! 2 h_{k+1} ( b_k \frac{ x'_{k+1} }{ h_{k} } - \frac{ y'_{k+1} ) x'_{k+1} }{ h_{k} } \\\\
432!! &= 4 b_k \frac{ < h x'^2 > }{ h_{k}^2 } - 4 \frac{ < h x' y' > }{ h_{k} }
433!! \f}
434!! and where \f$ < a > = \frac{1}{2} ( a_{k-1} + a_{k+1} ) \f$.
435!! Thus
436!! \f[
437!! b_k = \frac{ h_{k} < h x' y' > }{ < h x'^2 > } \;\; .
438!! \f]
439!!
440!! When evaluating the loss, \f$ G \f$, some rearrangement is necessary to reduce truncation
441!! errors. Since
442!! \f{align}{
443!! e_{k+j}^2 &= \left( b \frac{ x'_{k+j} }{ h_{k} } - y'_{k+j} \right)^2 \\\\
444!! &= b^2 \frac{ {x'}_{k+j}^2 }{ h_{k}^2 } - 2 b \frac{ x'_{k+j} y'_{k+j} }{ h_{k} } + {y'}_{k+j}^2
445!! \f}
446!! then
447!! \f{align}{
448!! G(b) &= 2 < h e^2 > \\\\
449!! &= 2 b^2 \frac{ < h {x'}^2 > }{ h_{k}^2 } - 4 b \frac{ < h x' y' > }{ h_{k} } + 2 < h' {y'}^2 >
450!! \;\; .
451!! \f}
452!!
453!! If we denote the value of b that yields the minimum value as \f$ b^* \f$ then
454!! \f[
455!! G(b^*) = < h {y'}^2 > - \frac{ < h x' y' >^2 }{ < h {x'}^2 > }
456!! \;\; .
457!! \f]
458!!
459!! Let
460!! \f{align}{
461!! G''(b) &= G(b) - G(b^*) \\\\
462!! &= b^2 \frac{ < h {x'}^2 > }{ h_{k}^2 } - 2 b \frac{ < h x' y' > }{ h_{k} }
463!! + \frac{ < h x' y' > }{ < h {x'}^2 > } \\\\
464!! &= \frac{ \left( b < h {x'}^2 > - h_{k} < h x' y' > \right)^2 }{ h_{k} < h {x'}^2 > }
465!! \;\; .
466!! \f}
467!! Minimizing \f$ G''(b) \f$ is equivalent to minimizing \f$ G(b) \f$ for the same data.
468!! \f$ G''(b^*)=0 \f$ so evaluation with the last form, in the vicinity of \f$ b^* \f$, avoids
469!! large cancelling terms.
470
471end module recon1d_plm_wls