MOM6
MOM_remapping.F90
1 !> Provides column-wise vertical remapping functions
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 ! Original module written by Laurent White, 2008.06.09
6 
7 use mom_error_handler, only : mom_error, fatal
8 use mom_string_functions, only : uppercase
9 use regrid_edge_values, only : edge_values_explicit_h4, edge_values_implicit_h4
10 use regrid_edge_values, only : edge_values_implicit_h4, edge_values_implicit_h6
11 use regrid_edge_values, only : edge_slopes_implicit_h3, edge_slopes_implicit_h5
12 use pcm_functions, only : pcm_reconstruction
13 use plm_functions, only : plm_reconstruction, plm_boundary_extrapolation
14 use ppm_functions, only : ppm_reconstruction, ppm_boundary_extrapolation
15 use pqm_functions, only : pqm_reconstruction, pqm_boundary_extrapolation_v1
16 
17 use iso_fortran_env, only : stdout=>output_unit, stderr=>error_unit
18 
19 implicit none ; private
20 
21 #include <MOM_memory.h>
22 
23 !> Container for remapping parameters
24 type, public :: remapping_cs
25  private
26  !> Determines which reconstruction to use
27  integer :: remapping_scheme = -911
28  !> Degree of polynomial reconstruction
29  integer :: degree = 0
30  !> If true, extrapolate boundaries
31  logical :: boundary_extrapolation = .true.
32  !> If true, reconstructions are checked for consistency.
33  logical :: check_reconstruction = .false.
34  !> If true, the result of remapping are checked for conservation and bounds.
35  logical :: check_remapping = .false.
36  !> If true, the intermediate values used in remapping are forced to be bounded.
37  logical :: force_bounds_in_subcell = .false.
38  !> If true use older, less acccurate expressions.
39  logical :: answers_2018 = .true.
40 end type
41 
42 ! The following routines are visible to the outside world
43 public remapping_core_h, remapping_core_w
44 public initialize_remapping, end_remapping, remapping_set_param, extract_member_remapping_cs
45 public remapping_unit_tests, build_reconstructions_1d, average_value_ppoly
46 public dzfromh1h2
47 
48 ! The following are private parameter constants
49 integer, parameter :: remapping_pcm = 0 !< O(h^1) remapping scheme
50 integer, parameter :: remapping_plm = 1 !< O(h^2) remapping scheme
51 integer, parameter :: remapping_ppm_h4 = 2 !< O(h^3) remapping scheme
52 integer, parameter :: remapping_ppm_ih4 = 3 !< O(h^3) remapping scheme
53 integer, parameter :: remapping_pqm_ih4ih3 = 4 !< O(h^4) remapping scheme
54 integer, parameter :: remapping_pqm_ih6ih5 = 5 !< O(h^5) remapping scheme
55 
56 integer, parameter :: integration_pcm = 0 !< Piecewise Constant Method
57 integer, parameter :: integration_plm = 1 !< Piecewise Linear Method
58 integer, parameter :: integration_ppm = 3 !< Piecewise Parabolic Method
59 integer, parameter :: integration_pqm = 5 !< Piecewise Quartic Method
60 
61 character(len=40) :: mdl = "MOM_remapping" !< This module's name.
62 
63 !> Documentation for external callers
64 character(len=256), public :: remappingschemesdoc = &
65  "PCM (1st-order accurate)\n"//&
66  "PLM (2nd-order accurate)\n"//&
67  "PPM_H4 (3rd-order accurate)\n"//&
68  "PPM_IH4 (3rd-order accurate)\n"//&
69  "PQM_IH4IH3 (4th-order accurate)\n"//&
70  "PQM_IH6IH5 (5th-order accurate)\n"
71 character(len=3), public :: remappingdefaultscheme = "PLM" !< Default remapping method
72 
73 ! This CPP macro turns on/off bounding of integrations limits so that they are
74 ! always within the cell. Roundoff can lead to the non-dimensional bounds being
75 ! outside of the range 0 to 1.
76 #define __USE_ROUNDOFF_SAFE_ADJUSTMENTS__
77 
78 real, parameter :: hneglect_dflt = 1.e-30 !< A thickness [H ~> m or kg m-2] that can be
79  !! added to thicknesses in a denominator without
80  !! changing the numerical result, except where
81  !! a division by zero would otherwise occur.
82 
83 logical, parameter :: old_algorithm = .false. !< Use the old "broken" algorithm.
84  !! This is a temporary measure to assist
85  !! debugging until we delete the old algorithm.
86 
87 contains
88 
89 !> Set parameters within remapping object
90 subroutine remapping_set_param(CS, remapping_scheme, boundary_extrapolation, &
91  check_reconstruction, check_remapping, force_bounds_in_subcell, answers_2018)
92  type(remapping_cs), intent(inout) :: cs !< Remapping control structure
93  character(len=*), optional, intent(in) :: remapping_scheme !< Remapping scheme to use
94  logical, optional, intent(in) :: boundary_extrapolation !< Indicate to extrapolate in boundary cells
95  logical, optional, intent(in) :: check_reconstruction !< Indicate to check reconstructions
96  logical, optional, intent(in) :: check_remapping !< Indicate to check results of remapping
97  logical, optional, intent(in) :: force_bounds_in_subcell !< Force subcells values to be bounded
98  logical, optional, intent(in) :: answers_2018 !< If true use older, less acccurate expressions.
99 
100  if (present(remapping_scheme)) then
101  call setreconstructiontype( remapping_scheme, cs )
102  endif
103  if (present(boundary_extrapolation)) then
104  cs%boundary_extrapolation = boundary_extrapolation
105  endif
106  if (present(check_reconstruction)) then
107  cs%check_reconstruction = check_reconstruction
108  endif
109  if (present(check_remapping)) then
110  cs%check_remapping = check_remapping
111  endif
112  if (present(force_bounds_in_subcell)) then
113  cs%force_bounds_in_subcell = force_bounds_in_subcell
114  endif
115  if (present(answers_2018)) then
116  cs%answers_2018 = answers_2018
117  endif
118 end subroutine remapping_set_param
119 
120 subroutine extract_member_remapping_cs(CS, remapping_scheme, degree, boundary_extrapolation, check_reconstruction, &
121  check_remapping, force_bounds_in_subcell)
122  type(remapping_cs), intent(in) :: cs !< Control structure for remapping module
123  integer, optional, intent(out) :: remapping_scheme !< Determines which reconstruction scheme to use
124  integer, optional, intent(out) :: degree !< Degree of polynomial reconstruction
125  logical, optional, intent(out) :: boundary_extrapolation !< If true, extrapolate boundaries
126  logical, optional, intent(out) :: check_reconstruction !< If true, reconstructions are checked for consistency.
127  logical, optional, intent(out) :: check_remapping !< If true, the result of remapping are checked
128  !! for conservation and bounds.
129  logical, optional, intent(out) :: force_bounds_in_subcell !< If true, the intermediate values used in
130  !! remapping are forced to be bounded.
131 
132  if (present(remapping_scheme)) remapping_scheme = cs%remapping_scheme
133  if (present(degree)) degree = cs%degree
134  if (present(boundary_extrapolation)) boundary_extrapolation = cs%boundary_extrapolation
135  if (present(check_reconstruction)) check_reconstruction = cs%check_reconstruction
136  if (present(check_remapping)) check_remapping = cs%check_remapping
137  if (present(force_bounds_in_subcell)) force_bounds_in_subcell = cs%force_bounds_in_subcell
138 
139 end subroutine extract_member_remapping_cs
140 !> Calculate edge coordinate x from cell width h
141 subroutine buildgridfromh(nz, h, x)
142  integer, intent(in) :: nz !< Number of cells
143  real, dimension(nz), intent(in) :: h !< Cell widths
144  real, dimension(nz+1), intent(inout) :: x !< Edge coordiantes starting at x(1)=0
145  ! Local variables
146  integer :: k
147 
148  x(1) = 0.0
149  do k = 1,nz
150  x(k+1) = x(k) + h(k)
151  enddo
152 
153 end subroutine buildgridfromh
154 
155 !> Compare two summation estimates of positive data and judge if due to more
156 !! than round-off.
157 !! When two sums are calculated from different vectors that should add up to
158 !! the same value, the results can differ by round off. The round off error
159 !! can be bounded to be proportional to the number of operations.
160 !! This function returns true if the difference between sum1 and sum2 is
161 !! larger than than the estimated round off bound.
162 !! \note This estimate/function is only valid for summation of positive data.
163 function ispossumerrsignificant(n1, sum1, n2, sum2)
164  integer, intent(in) :: n1 !< Number of values in sum1
165  integer, intent(in) :: n2 !< Number of values in sum2
166  real, intent(in) :: sum1 !< Sum of n1 values
167  real, intent(in) :: sum2 !< Sum of n2 values
168  logical :: ispossumerrsignificant !< True if difference in sums is large
169  ! Local variables
170  real :: sumerr, allowederr, eps
171 
172  if (sum1<0.) call mom_error(fatal,'isPosSumErrSignificant: sum1<0 is not allowed!')
173  if (sum2<0.) call mom_error(fatal,'isPosSumErrSignificant: sum2<0 is not allowed!')
174  sumerr = abs(sum1-sum2)
175  eps = epsilon(sum1)
176  allowederr = eps*0.5*(real(n1-1)*sum1+real(n2-1)*sum2)
177  if (sumerr>allowederr) then
178  write(0,*) 'isPosSumErrSignificant: sum1,sum2=',sum1,sum2
179  write(0,*) 'isPosSumErrSignificant: eps=',eps
180  write(0,*) 'isPosSumErrSignificant: err,n*eps=',sumerr,allowederr
181  write(0,*) 'isPosSumErrSignificant: err/eps,n1,n2,n1+n2=',sumerr/eps,n1,n2,n1+n2
182  ispossumerrsignificant = .true.
183  else
184  ispossumerrsignificant = .false.
185  endif
186 end function ispossumerrsignificant
187 
188 !> Remaps column of values u0 on grid h0 to grid h1 assuming the top edge is aligned.
189 subroutine remapping_core_h(CS, n0, h0, u0, n1, h1, u1, h_neglect, h_neglect_edge)
190  type(remapping_cs), intent(in) :: cs !< Remapping control structure
191  integer, intent(in) :: n0 !< Number of cells on source grid
192  real, dimension(n0), intent(in) :: h0 !< Cell widths on source grid
193  real, dimension(n0), intent(in) :: u0 !< Cell averages on source grid
194  integer, intent(in) :: n1 !< Number of cells on target grid
195  real, dimension(n1), intent(in) :: h1 !< Cell widths on target grid
196  real, dimension(n1), intent(out) :: u1 !< Cell averages on target grid
197  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
198  !! purpose of cell reconstructions
199  !! in the same units as h0.
200  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
201  !! for the purpose of edge value
202  !! calculations in the same units as h0.
203  ! Local variables
204  integer :: imethod
205  real, dimension(n0,2) :: ppoly_r_e !Edge value of polynomial
206  real, dimension(n0,2) :: ppoly_r_s !Edge slope of polynomial
207  real, dimension(n0,CS%degree+1) :: ppoly_r_coefs !Coefficients of polynomial
208  integer :: k
209  real :: eps, h0tot, h0err, h1tot, h1err, u0tot, u0err, u0min, u0max, u1tot, u1err, u1min, u1max, uh_err
210  real :: hneglect, hneglect_edge
211 
212  hneglect = 1.0e-30 ; if (present(h_neglect)) hneglect = h_neglect
213  hneglect_edge = 1.0e-10 ; if (present(h_neglect_edge)) hneglect_edge = h_neglect_edge
214 
215  call build_reconstructions_1d( cs, n0, h0, u0, ppoly_r_coefs, ppoly_r_e, ppoly_r_s, imethod, &
216  hneglect, hneglect_edge )
217 
218  if (cs%check_reconstruction) call check_reconstructions_1d(n0, h0, u0, cs%degree, &
219  cs%boundary_extrapolation, ppoly_r_coefs, ppoly_r_e, ppoly_r_s)
220 
221 
222  call remap_via_sub_cells( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, n1, h1, imethod, &
223  cs%force_bounds_in_subcell, u1, uh_err )
224 
225  if (cs%check_remapping) then
226  ! Check errors and bounds
227  call measure_input_bounds( n0, h0, u0, ppoly_r_e, h0tot, h0err, u0tot, u0err, u0min, u0max )
228  call measure_output_bounds( n1, h1, u1, h1tot, h1err, u1tot, u1err, u1min, u1max )
229  if (imethod<5) then ! We except PQM until we've debugged it
230  if ( (abs(u1tot-u0tot)>(u0err+u1err)+uh_err .and. abs(h1tot-h0tot)<h0err+h1err) &
231  .or. (u1min<u0min .or. u1max>u0max) ) then
232  write(0,*) 'iMethod = ',imethod
233  write(0,*) 'H: h0tot=',h0tot,'h1tot=',h1tot,'dh=',h1tot-h0tot,'h0err=',h0err,'h1err=',h1err
234  if (abs(h1tot-h0tot)>h0err+h1err) &
235  write(0,*) 'H non-conservation difference=',h1tot-h0tot,'allowed err=',h0err+h1err,' <-----!'
236  write(0,*) 'UH: u0tot=',u0tot,'u1tot=',u1tot,'duh=',u1tot-u0tot,'u0err=',u0err,'u1err=',u1err,'uh_err=',uh_err
237  if (abs(u1tot-u0tot)>(u0err+u1err)+uh_err) &
238  write(0,*) 'U non-conservation difference=',u1tot-u0tot,'allowed err=',u0err+u1err+uh_err,' <-----!'
239  write(0,*) 'U: u0min=',u0min,'u1min=',u1min
240  if (u1min<u0min) write(0,*) 'U minimum overshoot=',u1min-u0min,' <-----!'
241  write(0,*) 'U: u0max=',u0max,'u1max=',u1max
242  if (u1max>u0max) write(0,*) 'U maximum overshoot=',u1max-u0max,' <-----!'
243  write(0,'(a3,6a24)') 'k','h0','left edge','u0','right edge','h1','u1'
244  do k = 1, max(n0,n1)
245  if (k<=min(n0,n1)) then
246  write(0,'(i3,1p6e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2),h1(k),u1(k)
247  elseif (k>n0) then
248  write(0,'(i3,96x,1p2e24.16)') k,h1(k),u1(k)
249  else
250  write(0,'(i3,1p4e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2)
251  endif
252  enddo
253  write(0,'(a3,2a24)') 'k','u0','Polynomial coefficients'
254  do k = 1, n0
255  write(0,'(i3,1p6e24.16)') k,u0(k),ppoly_r_coefs(k,:)
256  enddo
257  call mom_error( fatal, 'MOM_remapping, remapping_core_h: '//&
258  'Remapping result is inconsistent!' )
259  endif
260  endif ! method<5
261  endif
262 
263 end subroutine remapping_core_h
264 
265 !> Remaps column of values u0 on grid h0 to implied grid h1
266 !! where the interfaces of h1 differ from those of h0 by dx.
267 subroutine remapping_core_w( CS, n0, h0, u0, n1, dx, u1, h_neglect, h_neglect_edge )
268  type(remapping_cs), intent(in) :: cs !< Remapping control structure
269  integer, intent(in) :: n0 !< Number of cells on source grid
270  real, dimension(n0), intent(in) :: h0 !< Cell widths on source grid
271  real, dimension(n0), intent(in) :: u0 !< Cell averages on source grid
272  integer, intent(in) :: n1 !< Number of cells on target grid
273  real, dimension(n1+1), intent(in) :: dx !< Cell widths on target grid
274  real, dimension(n1), intent(out) :: u1 !< Cell averages on target grid
275  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
276  !! purpose of cell reconstructions
277  !! in the same units as h0.
278  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
279  !! for the purpose of edge value
280  !! calculations in the same units as h0.
281  ! Local variables
282  integer :: imethod
283  real, dimension(n0,2) :: ppoly_r_e !Edge value of polynomial
284  real, dimension(n0,2) :: ppoly_r_s !Edge slope of polynomial
285  real, dimension(n0,CS%degree+1) :: ppoly_r_coefs !Coefficients of polynomial
286  integer :: k
287  real :: eps, h0tot, h0err, h1tot, h1err
288  real :: u0tot, u0err, u0min, u0max, u1tot, u1err, u1min, u1max, uh_err
289  real, dimension(n1) :: h1 !< Cell widths on target grid
290  real :: hneglect, hneglect_edge
291 
292  hneglect = 1.0e-30 ; if (present(h_neglect)) hneglect = h_neglect
293  hneglect_edge = 1.0e-10 ; if (present(h_neglect_edge)) hneglect_edge = h_neglect_edge
294 
295  call build_reconstructions_1d( cs, n0, h0, u0, ppoly_r_coefs, ppoly_r_e, ppoly_r_s, imethod,&
296  hneglect, hneglect_edge )
297 
298  if (cs%check_reconstruction) call check_reconstructions_1d(n0, h0, u0, cs%degree, &
299  cs%boundary_extrapolation, ppoly_r_coefs, ppoly_r_e, ppoly_r_s)
300 
301  ! This is a temporary step prior to switching to remapping_core_h()
302  do k = 1, n1
303  if (k<=n0) then
304  h1(k) = max( 0., h0(k) + ( dx(k+1) - dx(k) ) )
305  else
306  h1(k) = max( 0., dx(k+1) - dx(k) )
307  endif
308  enddo
309  call remap_via_sub_cells( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, n1, h1, imethod, &
310  cs%force_bounds_in_subcell,u1, uh_err )
311 ! call remapByDeltaZ( n0, h0, u0, ppoly_r_E, ppoly_r_coefs, n1, dx, iMethod, u1, hNeglect )
312 ! call remapByProjection( n0, h0, u0, CS%ppoly_r, n1, h1, iMethod, u1, hNeglect )
313 
314  if (cs%check_remapping) then
315  ! Check errors and bounds
316  call measure_input_bounds( n0, h0, u0, ppoly_r_e, h0tot, h0err, u0tot, u0err, u0min, u0max )
317  call measure_output_bounds( n1, h1, u1, h1tot, h1err, u1tot, u1err, u1min, u1max )
318  if (imethod<5) then ! We except PQM until we've debugged it
319  if ( (abs(u1tot-u0tot)>(u0err+u1err)+uh_err .and. abs(h1tot-h0tot)<h0err+h1err) &
320  .or. (u1min<u0min .or. u1max>u0max) ) then
321  write(0,*) 'iMethod = ',imethod
322  write(0,*) 'H: h0tot=',h0tot,'h1tot=',h1tot,'dh=',h1tot-h0tot,'h0err=',h0err,'h1err=',h1err
323  if (abs(h1tot-h0tot)>h0err+h1err) &
324  write(0,*) 'H non-conservation difference=',h1tot-h0tot,'allowed err=',h0err+h1err,' <-----!'
325  write(0,*) 'UH: u0tot=',u0tot,'u1tot=',u1tot,'duh=',u1tot-u0tot,'u0err=',u0err,'u1err=',u1err,'uh_err=',uh_err
326  if (abs(u1tot-u0tot)>(u0err+u1err)+uh_err) &
327  write(0,*) 'U non-conservation difference=',u1tot-u0tot,'allowed err=',u0err+u1err+uh_err,' <-----!'
328  write(0,*) 'U: u0min=',u0min,'u1min=',u1min
329  if (u1min<u0min) write(0,*) 'U minimum overshoot=',u1min-u0min,' <-----!'
330  write(0,*) 'U: u0max=',u0max,'u1max=',u1max
331  if (u1max>u0max) write(0,*) 'U maximum overshoot=',u1max-u0max,' <-----!'
332  write(0,'(a3,6a24)') 'k','h0','left edge','u0','right edge','h1','u1'
333  do k = 1, max(n0,n1)
334  if (k<=min(n0,n1)) then
335  write(0,'(i3,1p6e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2),h1(k),u1(k)
336  elseif (k>n0) then
337  write(0,'(i3,96x,1p2e24.16)') k,h1(k),u1(k)
338  else
339  write(0,'(i3,1p4e24.16)') k,h0(k),ppoly_r_e(k,1),u0(k),ppoly_r_e(k,2)
340  endif
341  enddo
342  write(0,'(a3,2a24)') 'k','u0','Polynomial coefficients'
343  do k = 1, n0
344  write(0,'(i3,1p6e24.16)') k,u0(k),ppoly_r_coefs(k,:)
345  enddo
346  call mom_error( fatal, 'MOM_remapping, remapping_core_w: '//&
347  'Remapping result is inconsistent!' )
348  endif
349  endif ! method<5
350  endif
351 
352 end subroutine remapping_core_w
353 
354 !> Creates polynomial reconstructions of u0 on the source grid h0.
355 subroutine build_reconstructions_1d( CS, n0, h0, u0, ppoly_r_coefs, &
356  ppoly_r_E, ppoly_r_S, iMethod, h_neglect, &
357  h_neglect_edge )
358  type(remapping_cs), intent(in) :: cs !< Remapping control structure
359  integer, intent(in) :: n0 !< Number of cells on source grid
360  real, dimension(n0), intent(in) :: h0 !< Cell widths on source grid
361  real, dimension(n0), intent(in) :: u0 !< Cell averages on source grid
362  real, dimension(n0,CS%degree+1), &
363  intent(out) :: ppoly_r_coefs !< Coefficients of polynomial
364  real, dimension(n0,2), intent(out) :: ppoly_r_e !< Edge value of polynomial
365  real, dimension(n0,2), intent(out) :: ppoly_r_s !< Edge slope of polynomial
366  integer, intent(out) :: imethod !< Integration method
367  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
368  !! purpose of cell reconstructions
369  !! in the same units as h0.
370  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
371  !! for the purpose of edge value
372  !! calculations in the same units as h0.
373  ! Local variables
374  integer :: local_remapping_scheme
375  integer :: remapping_scheme !< Remapping scheme
376  logical :: boundary_extrapolation !< Extrapolate at boundaries if true
377 
378  ! Reset polynomial
379  ppoly_r_e(:,:) = 0.0
380  ppoly_r_s(:,:) = 0.0
381  ppoly_r_coefs(:,:) = 0.0
382  imethod = -999
383 
384  local_remapping_scheme = cs%remapping_scheme
385  if (n0<=1) then
386  local_remapping_scheme = remapping_pcm
387  elseif (n0<=3) then
388  local_remapping_scheme = min( local_remapping_scheme, remapping_plm )
389  elseif (n0<=4) then
390  local_remapping_scheme = min( local_remapping_scheme, remapping_ppm_h4 )
391  endif
392  select case ( local_remapping_scheme )
393  case ( remapping_pcm )
394  call pcm_reconstruction( n0, u0, ppoly_r_e, ppoly_r_coefs)
395  imethod = integration_pcm
396  case ( remapping_plm )
397  call plm_reconstruction( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, h_neglect )
398  if ( cs%boundary_extrapolation ) then
399  call plm_boundary_extrapolation( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, h_neglect)
400  endif
401  imethod = integration_plm
402  case ( remapping_ppm_h4 )
403  call edge_values_explicit_h4( n0, h0, u0, ppoly_r_e, h_neglect_edge, answers_2018=cs%answers_2018 )
404  call ppm_reconstruction( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, h_neglect, answers_2018=cs%answers_2018 )
405  if ( cs%boundary_extrapolation ) then
406  call ppm_boundary_extrapolation( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, h_neglect )
407  endif
408  imethod = integration_ppm
409  case ( remapping_ppm_ih4 )
410  call edge_values_implicit_h4( n0, h0, u0, ppoly_r_e, h_neglect_edge, answers_2018=cs%answers_2018 )
411  call ppm_reconstruction( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, h_neglect, answers_2018=cs%answers_2018 )
412  if ( cs%boundary_extrapolation ) then
413  call ppm_boundary_extrapolation( n0, h0, u0, ppoly_r_e, ppoly_r_coefs, h_neglect )
414  endif
415  imethod = integration_ppm
416  case ( remapping_pqm_ih4ih3 )
417  call edge_values_implicit_h4( n0, h0, u0, ppoly_r_e, h_neglect_edge, answers_2018=cs%answers_2018 )
418  call edge_slopes_implicit_h3( n0, h0, u0, ppoly_r_s, h_neglect, answers_2018=cs%answers_2018 )
419  call pqm_reconstruction( n0, h0, u0, ppoly_r_e, ppoly_r_s, ppoly_r_coefs, h_neglect, &
420  answers_2018=cs%answers_2018 )
421  if ( cs%boundary_extrapolation ) then
422  call pqm_boundary_extrapolation_v1( n0, h0, u0, ppoly_r_e, ppoly_r_s, &
423  ppoly_r_coefs, h_neglect )
424  endif
425  imethod = integration_pqm
426  case ( remapping_pqm_ih6ih5 )
427  call edge_values_implicit_h6( n0, h0, u0, ppoly_r_e, h_neglect_edge, answers_2018=cs%answers_2018 )
428  call edge_slopes_implicit_h5( n0, h0, u0, ppoly_r_s, h_neglect, answers_2018=cs%answers_2018 )
429  call pqm_reconstruction( n0, h0, u0, ppoly_r_e, ppoly_r_s, ppoly_r_coefs, h_neglect, &
430  answers_2018=cs%answers_2018 )
431  if ( cs%boundary_extrapolation ) then
432  call pqm_boundary_extrapolation_v1( n0, h0, u0, ppoly_r_e, ppoly_r_s, &
433  ppoly_r_coefs, h_neglect )
434  endif
435  imethod = integration_pqm
436  case default
437  call mom_error( fatal, 'MOM_remapping, build_reconstructions_1d: '//&
438  'The selected remapping method is invalid' )
439  end select
440 
441 end subroutine build_reconstructions_1d
442 
443 !> Checks that edge values and reconstructions satisfy bounds
444 subroutine check_reconstructions_1d(n0, h0, u0, deg, boundary_extrapolation, &
445  ppoly_r_coefs, ppoly_r_E, ppoly_r_S)
446  integer, intent(in) :: n0 !< Number of cells on source grid
447  real, dimension(n0), intent(in) :: h0 !< Cell widths on source grid
448  real, dimension(n0), intent(in) :: u0 !< Cell averages on source grid
449  integer, intent(in) :: deg !< Degree of polynomial reconstruction
450  logical, intent(in) :: boundary_extrapolation !< Extrapolate at boundaries if true
451  real, dimension(n0,deg+1),intent(out) :: ppoly_r_coefs !< Coefficients of polynomial
452  real, dimension(n0,2), intent(out) :: ppoly_r_E !< Edge value of polynomial
453  real, dimension(n0,2), intent(out) :: ppoly_r_S !< Edge slope of polynomial
454  ! Local variables
455  integer :: i0, n
456  real :: u_l, u_c, u_r ! Cell averages
457  real :: u_min, u_max
458  logical :: problem_detected
459 
460  problem_detected = .false.
461  do i0 = 1, n0
462  u_l = u0(max(1,i0-1))
463  u_c = u0(i0)
464  u_r = u0(min(n0,i0+1))
465  if (i0 > 1 .or. .not. boundary_extrapolation) then
466  u_min = min(u_l, u_c)
467  u_max = max(u_l, u_c)
468  if (ppoly_r_e(i0,1) < u_min) then
469  write(0,'(a,i4,5(x,a,1pe24.16))') 'Left edge undershoot at',i0,'u(i0-1)=',u_l,'u(i0)=',u_c, &
470  'edge=',ppoly_r_e(i0,1),'err=',ppoly_r_e(i0,1)-u_min
471  problem_detected = .true.
472  endif
473  if (ppoly_r_e(i0,1) > u_max) then
474  write(0,'(a,i4,5(x,a,1pe24.16))') 'Left edge overshoot at',i0,'u(i0-1)=',u_l,'u(i0)=',u_c, &
475  'edge=',ppoly_r_e(i0,1),'err=',ppoly_r_e(i0,1)-u_max
476  problem_detected = .true.
477  endif
478  endif
479  if (i0 < n0 .or. .not. boundary_extrapolation) then
480  u_min = min(u_c, u_r)
481  u_max = max(u_c, u_r)
482  if (ppoly_r_e(i0,2) < u_min) then
483  write(0,'(a,i4,5(x,a,1pe24.16))') 'Right edge undershoot at',i0,'u(i0)=',u_c,'u(i0+1)=',u_r, &
484  'edge=',ppoly_r_e(i0,2),'err=',ppoly_r_e(i0,2)-u_min
485  problem_detected = .true.
486  endif
487  if (ppoly_r_e(i0,2) > u_max) then
488  write(0,'(a,i4,5(x,a,1pe24.16))') 'Right edge overshoot at',i0,'u(i0)=',u_c,'u(i0+1)=',u_r, &
489  'edge=',ppoly_r_e(i0,2),'err=',ppoly_r_e(i0,2)-u_max
490  problem_detected = .true.
491  endif
492  endif
493  if (i0 > 1) then
494  if ( (u_c-u_l)*(ppoly_r_e(i0,1)-ppoly_r_e(i0-1,2)) < 0.) then
495  write(0,'(a,i4,5(x,a,1pe24.16))') 'Non-monotonic edges at',i0,'u(i0-1)=',u_l,'u(i0)=',u_c, &
496  'right edge=',ppoly_r_e(i0-1,2),'left edge=',ppoly_r_e(i0,1)
497  write(0,'(5(a,1pe24.16,x))') 'u(i0)-u(i0-1)',u_c-u_l,'edge diff=',ppoly_r_e(i0,1)-ppoly_r_e(i0-1,2)
498  problem_detected = .true.
499  endif
500  endif
501  if (problem_detected) then
502  write(0,'(a,1p9e24.16)') 'Polynomial coeffs:',ppoly_r_coefs(i0,:)
503  write(0,'(3(a,1pe24.16,x))') 'u_l=',u_l,'u_c=',u_c,'u_r=',u_r
504  write(0,'(a4,10a24)') 'i0','h0(i0)','u0(i0)','left edge','right edge','Polynomial coefficients'
505  do n = 1, n0
506  write(0,'(i4,1p10e24.16)') n,h0(n),u0(n),ppoly_r_e(n,1),ppoly_r_e(n,2),ppoly_r_coefs(n,:)
507  enddo
508  call mom_error(fatal, 'MOM_remapping, check_reconstructions_1d: '// &
509  'Edge values or polynomial coefficients were inconsistent!')
510  endif
511  enddo
512 
513 end subroutine check_reconstructions_1d
514 
515 !> Remaps column of n0 values u0 on grid h0 to grid h1 with n1 cells by calculating
516 !! the n0+n1+1 sub-integrals of the intersection of h0 and h1, and the summing the
517 !! appropriate integrals into the h1*u1 values. h0 and h1 must have the same units.
518 subroutine remap_via_sub_cells( n0, h0, u0, ppoly0_E, ppoly0_coefs, n1, h1, method, &
519  force_bounds_in_subcell, u1, uh_err, ah_sub, aisub_src, aiss, aise )
520  integer, intent(in) :: n0 !< Number of cells in source grid
521  real, intent(in) :: h0(n0) !< Source grid widths (size n0)
522  real, intent(in) :: u0(n0) !< Source cell averages (size n0)
523  real, intent(in) :: ppoly0_E(n0,2) !< Edge value of polynomial
524  real, intent(in) :: ppoly0_coefs(:,:) !< Coefficients of polynomial
525  integer, intent(in) :: n1 !< Number of cells in target grid
526  real, intent(in) :: h1(n1) !< Target grid widths (size n1)
527  integer, intent(in) :: method !< Remapping scheme to use
528  logical, intent(in) :: force_bounds_in_subcell !< Force sub-cell values to be bounded
529  real, intent(out) :: u1(n1) !< Target cell averages (size n1)
530  real, intent(out) :: uh_err !< Estimate of bound on error in sum of u*h
531  real, optional, intent(out) :: ah_sub(n0+n1+1) !< h_sub
532  integer, optional, intent(out) :: aisub_src(n0+n1+1) !< i_sub_src
533  integer, optional, intent(out) :: aiss(n0) !< isrc_start
534  integer, optional, intent(out) :: aise(n0) !< isrc_ens
535  ! Local variables
536  integer :: i_sub ! Index of sub-cell
537  integer :: i0 ! Index into h0(1:n0), source column
538  integer :: i1 ! Index into h1(1:n1), target column
539  integer :: i_start0 ! Used to record which sub-cells map to source cells
540  integer :: i_start1 ! Used to record which sub-cells map to target cells
541  integer :: i_max ! Used to record which sub-cell is the largest contribution of a source cell
542  real :: dh_max ! Used to record which sub-cell is the largest contribution of a source cell
543  real, dimension(n0+n1+1) :: h_sub ! Width of each each sub-cell
544  real, dimension(n0+n1+1) :: uh_sub ! Integral of u*h over each sub-cell
545  real, dimension(n0+n1+1) :: u_sub ! Average of u over each sub-cell
546  integer, dimension(n0+n1+1) :: isub_src ! Index of source cell for each sub-cell
547  integer, dimension(n0) :: isrc_start ! Index of first sub-cell within each source cell
548  integer, dimension(n0) :: isrc_end ! Index of last sub-cell within each source cell
549  integer, dimension(n0) :: isrc_max ! Index of thickest sub-cell within each source cell
550  real, dimension(n0) :: h0_eff ! Effective thickness of source cells
551  real, dimension(n0) :: u0_min ! Minimum value of reconstructions in source cell
552  real, dimension(n0) :: u0_max ! Minimum value of reconstructions in source cell
553  integer, dimension(n1) :: itgt_start ! Index of first sub-cell within each target cell
554  integer, dimension(n1) :: itgt_end ! Index of last sub-cell within each target cell
555  real :: xa, xb ! Non-dimensional position within a source cell (0..1)
556  real :: h0_supply, h1_supply ! The amount of width available for constructing sub-cells
557  real :: dh ! The width of the sub-cell
558  real :: duh ! The total amount of accumulated stuff (u*h)
559  real :: dh0_eff ! Running sum of source cell thickness
560  ! For error checking/debugging
561  logical, parameter :: force_bounds_in_target = .true. ! To fix round-off issues
562  logical, parameter :: adjust_thickest_subcell = .true. ! To fix round-off conservation issues
563  logical, parameter :: debug_bounds = .false. ! For debugging overshoots etc.
564  integer :: k, i0_last_thick_cell
565  real :: h0tot, h0err, h1tot, h1err, h2tot, h2err, u02_err
566  real :: u0tot, u0err, u0min, u0max, u1tot, u1err, u1min, u1max, u2tot, u2err, u2min, u2max, u_orig
567  logical :: src_has_volume !< True if h0 has not been consumed
568  logical :: tgt_has_volume !< True if h1 has not been consumed
569 
570  if (old_algorithm) isrc_max(:)=1
571 
572  i0_last_thick_cell = 0
573  do i0 = 1, n0
574  u0_min(i0) = min(ppoly0_e(i0,1), ppoly0_e(i0,2))
575  u0_max(i0) = max(ppoly0_e(i0,1), ppoly0_e(i0,2))
576  if (h0(i0)>0.) i0_last_thick_cell = i0
577  enddo
578 
579  ! Initialize algorithm
580  h0_supply = h0(1)
581  h1_supply = h1(1)
582  src_has_volume = .true.
583  tgt_has_volume = .true.
584  i0 = 1 ; i1 = 1
585  i_start0 = 1 ; i_start1 = 1
586  i_max = 1
587  dh_max = 0.
588  dh0_eff = 0.
589 
590  ! First sub-cell is always vanished
591  h_sub(1) = 0.
592  isrc_start(1) = 1
593  isrc_end(1) = 1
594  isrc_max(1) = 1
595  isub_src(1) = 1
596 
597  ! Loop over each sub-cell to calculate intersections with source and target grids
598  do i_sub = 2, n0+n1+1
599 
600  ! This is the width of the sub-cell, determined by which ever column has the least
601  ! supply available to consume.
602  dh = min(h0_supply, h1_supply)
603 
604  ! This is the running sum of the source cell thickness. After summing over each
605  ! sub-cell, the sum of sub-cell thickness might differ from the original source
606  ! cell thickness due to round off.
607  dh0_eff = dh0_eff + min(dh, h0_supply)
608 
609  ! Record the source index (i0) that this sub-cell integral belongs to. This
610  ! is needed to index the reconstruction coefficients for the source cell
611  ! used in the integrals of the sub-cell width.
612  isub_src(i_sub) = i0
613  h_sub(i_sub) = dh
614 
615  ! For recording the largest sub-cell within a source cell.
616  if (dh >= dh_max) then
617  i_max = i_sub
618  dh_max = dh
619  endif
620 
621  ! Which ever column (source or target) has the least width left to consume determined
622  ! the width, dh, of sub-cell i_sub in the expression for dh above.
623  if (h0_supply <= h1_supply .and. src_has_volume) then
624  ! h0_supply is smaller than h1_supply) so we consume h0_supply and increment the
625  ! source cell index.
626  h1_supply = h1_supply - dh ! Although this is a difference the result will
627  ! be non-negative because of the conditional.
628  ! Record the sub-cell start/end index that span the source cell i0.
629  isrc_start(i0) = i_start0
630  isrc_end(i0) = i_sub
631  i_start0 = i_sub + 1
632  ! Record the sub-cell that is the largest fraction of the source cell.
633  isrc_max(i0) = i_max
634  i_max = i_sub + 1
635  dh_max = 0.
636  ! Record the source cell thickness found by summing the sub-cell thicknesses.
637  h0_eff(i0) = dh0_eff
638  ! Move the source index.
639  if (old_algorithm) then
640  if (i0 < i0_last_thick_cell) then
641  i0 = i0 + 1
642  h0_supply = h0(i0)
643  dh0_eff = 0.
644  do while (h0_supply==0. .and. i0<i0_last_thick_cell)
645  ! This loop skips over vanished source cells
646  i0 = i0 + 1
647  h0_supply = h0(i0)
648  enddo
649  else
650  h0_supply = 1.e30
651  endif
652  else
653  if (i0 < n0) then
654  i0 = i0 + 1
655  h0_supply = h0(i0)
656  dh0_eff = 0.
657  else
658  h0_supply = 0.
659  src_has_volume = .false.
660  endif
661  endif
662  elseif (h0_supply >= h1_supply .and. tgt_has_volume) then
663  ! h1_supply is smaller than h0_supply) so we consume h1_supply and increment the
664  ! target cell index.
665  h0_supply = h0_supply - dh ! Although this is a difference the result will
666  ! be non-negative because of the conditional.
667  ! Record the sub-cell start/end index that span the target cell i1.
668  itgt_start(i1) = i_start1
669  itgt_end(i1) = i_sub
670  i_start1 = i_sub + 1
671  ! Move the target index.
672  if (i1 < n1) then
673  i1 = i1 + 1
674  h1_supply = h1(i1)
675  else
676  if (old_algorithm) then
677  h1_supply = 1.e30
678  else
679  h1_supply = 0.
680  tgt_has_volume = .false.
681  endif
682  endif
683  elseif (src_has_volume) then
684  ! We ran out of target volume but still have source cells to consume
685  h_sub(i_sub) = h0_supply
686  ! Record the sub-cell start/end index that span the source cell i0.
687  isrc_start(i0) = i_start0
688  isrc_end(i0) = i_sub
689  i_start0 = i_sub + 1
690  ! Record the sub-cell that is the largest fraction of the source cell.
691  isrc_max(i0) = i_max
692  i_max = i_sub + 1
693  dh_max = 0.
694  ! Record the source cell thickness found by summing the sub-cell thicknesses.
695  h0_eff(i0) = dh0_eff
696  if (i0 < n0) then
697  i0 = i0 + 1
698  h0_supply = h0(i0)
699  dh0_eff = 0.
700  else
701  h0_supply = 0.
702  src_has_volume = .false.
703  endif
704  elseif (tgt_has_volume) then
705  ! We ran out of source volume but still have target cells to consume
706  h_sub(i_sub) = h1_supply
707  ! Record the sub-cell start/end index that span the target cell i1.
708  itgt_start(i1) = i_start1
709  itgt_end(i1) = i_sub
710  i_start1 = i_sub + 1
711  ! Move the target index.
712  if (i1 < n1) then
713  i1 = i1 + 1
714  h1_supply = h1(i1)
715  else
716  h1_supply = 0.
717  tgt_has_volume = .false.
718  endif
719  else
720  stop 'remap_via_sub_cells: THIS SHOULD NEVER HAPPEN!'
721  endif
722 
723  enddo
724 
725  ! Loop over each sub-cell to calculate average/integral values within each sub-cell.
726  xa = 0.
727  dh0_eff = 0.
728  uh_sub(1) = 0.
729  u_sub(1) = ppoly0_e(1,1)
730  u02_err = 0.
731  do i_sub = 2, n0+n1
732 
733  ! Sub-cell thickness from loop above
734  dh = h_sub(i_sub)
735 
736  ! Source cell
737  i0 = isub_src(i_sub)
738 
739  ! Evaluate average and integral for sub-cell i_sub.
740  ! Integral is over distance dh but expressed in terms of non-dimensional
741  ! positions with source cell from xa to xb (0 <= xa <= xb <= 1).
742  dh0_eff = dh0_eff + dh ! Cumulative thickness within the source cell
743  if (h0_eff(i0)>0.) then
744  xb = dh0_eff / h0_eff(i0) ! This expression yields xa <= xb <= 1.0
745  xb = min(1., xb) ! This is only needed when the total target column is wider than the source column
746  u_sub(i_sub) = average_value_ppoly( n0, u0, ppoly0_e, ppoly0_coefs, method, i0, xa, xb)
747  else ! Vanished cell
748  xb = 1.
749  u_sub(i_sub) = u0(i0)
750  endif
751  if (debug_bounds) then
752  if (method<5 .and.(u_sub(i_sub)<u0_min(i0) .or. u_sub(i_sub)>u0_max(i0))) then
753  write(0,*) 'Sub cell average is out of bounds',i_sub,'method=',method
754  write(0,*) 'xa,xb: ',xa,xb
755  write(0,*) 'Edge values: ',ppoly0_e(i0,:),'mean',u0(i0)
756  write(0,*) 'a_c: ',(u0(i0)-ppoly0_e(i0,1))+(u0(i0)-ppoly0_e(i0,2))
757  write(0,*) 'Polynomial coeffs: ',ppoly0_coefs(i0,:)
758  write(0,*) 'Bounds min=',u0_min(i0),'max=',u0_max(i0)
759  write(0,*) 'Average: ',u_sub(i_sub),'rel to min=',u_sub(i_sub)-u0_min(i0),'rel to max=',u_sub(i_sub)-u0_max(i0)
760  call mom_error( fatal, 'MOM_remapping, remap_via_sub_cells: '//&
761  'Sub-cell average is out of bounds!' )
762  endif
763  endif
764  if (force_bounds_in_subcell) then
765  ! These next two lines should not be needed but when using PQM we found roundoff
766  ! can lead to overshoots. These lines sweep issues under the rug which need to be
767  ! properly .. later. -AJA
768  u_orig = u_sub(i_sub)
769  u_sub(i_sub) = max( u_sub(i_sub), u0_min(i0) )
770  u_sub(i_sub) = min( u_sub(i_sub), u0_max(i0) )
771  u02_err = u02_err + dh*abs( u_sub(i_sub) - u_orig )
772  endif
773  uh_sub(i_sub) = dh * u_sub(i_sub)
774 
775  if (isub_src(i_sub+1) /= i0) then
776  ! If the next sub-cell is in a different source cell, reset the position counters
777  dh0_eff = 0.
778  xa = 0.
779  else
780  xa = xb ! Next integral will start at end of last
781  endif
782 
783  enddo
784  u_sub(n0+n1+1) = ppoly0_e(n0,2) ! This value is only needed when total target column
785  uh_sub(n0+n1+1) = ppoly0_e(n0,2) * h_sub(n0+n1+1) ! is wider than the source column
786 
787  if (adjust_thickest_subcell) then
788  ! Loop over each source cell substituting the integral/average for the thickest sub-cell (within
789  ! the source cell) with the residual of the source cell integral minus the other sub-cell integrals
790  ! aka a genius algorithm for accurate conservation when remapping from Robert Hallberg (@Hallberg-NOAA).
791  do i0 = 1, i0_last_thick_cell
792  i_max = isrc_max(i0)
793  dh_max = h_sub(i_max)
794  if (dh_max > 0.) then
795  ! duh will be the sum of sub-cell integrals within the source cell except for the thickest sub-cell.
796  duh = 0.
797  do i_sub = isrc_start(i0), isrc_end(i0)
798  if (i_sub /= i_max) duh = duh + uh_sub(i_sub)
799  enddo
800  uh_sub(i_max) = u0(i0)*h0(i0) - duh
801  u02_err = u02_err + max( abs(uh_sub(i_max)), abs(u0(i0)*h0(i0)), abs(duh) )
802  endif
803  enddo
804  endif
805 
806  ! Loop over each target cell summing the integrals from sub-cells within the target cell.
807  uh_err = 0.
808  do i1 = 1, n1
809  if (h1(i1) > 0.) then
810  duh = 0. ; dh = 0.
811  i_sub = itgt_start(i1)
812  if (force_bounds_in_target) then
813  u1min = u_sub(i_sub)
814  u1max = u_sub(i_sub)
815  endif
816  do i_sub = itgt_start(i1), itgt_end(i1)
817  if (force_bounds_in_target) then
818  u1min = min(u1min, u_sub(i_sub))
819  u1max = max(u1max, u_sub(i_sub))
820  endif
821  dh = dh + h_sub(i_sub)
822  duh = duh + uh_sub(i_sub)
823  ! This accumulates the contribution to the error bound for the sum of u*h
824  uh_err = uh_err + max(abs(duh),abs(uh_sub(i_sub)))*epsilon(duh)
825  enddo
826  u1(i1) = duh / dh
827  ! This is the contribution from the division to the error bound for the sum of u*h
828  uh_err = uh_err + abs(duh)*epsilon(duh)
829  if (force_bounds_in_target) then
830  u_orig = u1(i1)
831  u1(i1) = max(u1min, min(u1max, u1(i1)))
832  ! Adjusting to be bounded contributes to the error for the sum of u*h
833  uh_err = uh_err + dh*abs( u1(i1)-u_orig )
834  endif
835  else
836  u1(i1) = u_sub(itgt_start(i1))
837  endif
838  enddo
839 
840  ! Check errors and bounds
841  if (debug_bounds) then
842  call measure_input_bounds( n0, h0, u0, ppoly0_e, h0tot, h0err, u0tot, u0err, u0min, u0max )
843  call measure_output_bounds( n1, h1, u1, h1tot, h1err, u1tot, u1err, u1min, u1max )
844  call measure_output_bounds( n0+n1+1, h_sub, u_sub, h2tot, h2err, u2tot, u2err, u2min, u2max )
845  if (method<5) then ! We except PQM until we've debugged it
846  if ( (abs(u1tot-u0tot)>(u0err+u1err)+uh_err+u02_err .and. abs(h1tot-h0tot)<h0err+h1err) &
847  .or. (abs(u2tot-u0tot)>u0err+u2err+u02_err .and. abs(h2tot-h0tot)<h0err+h2err) &
848  .or. (u1min<u0min .or. u1max>u0max) ) then
849  write(0,*) 'method = ',method
850  write(0,*) 'Source to sub-cells:'
851  write(0,*) 'H: h0tot=',h0tot,'h2tot=',h2tot,'dh=',h2tot-h0tot,'h0err=',h0err,'h2err=',h2err
852  if (abs(h2tot-h0tot)>h0err+h2err) &
853  write(0,*) 'H non-conservation difference=',h2tot-h0tot,'allowed err=',h0err+h2err,' <-----!'
854  write(0,*) 'UH: u0tot=',u0tot,'u2tot=',u2tot,'duh=',u2tot-u0tot,'u0err=',u0err,'u2err=',u2err,&
855  'adjustment err=',u02_err
856  if (abs(u2tot-u0tot)>u0err+u2err) &
857  write(0,*) 'U non-conservation difference=',u2tot-u0tot,'allowed err=',u0err+u2err,' <-----!'
858  write(0,*) 'Sub-cells to target:'
859  write(0,*) 'H: h2tot=',h2tot,'h1tot=',h1tot,'dh=',h1tot-h2tot,'h2err=',h2err,'h1err=',h1err
860  if (abs(h1tot-h2tot)>h2err+h1err) &
861  write(0,*) 'H non-conservation difference=',h1tot-h2tot,'allowed err=',h2err+h1err,' <-----!'
862  write(0,*) 'UH: u2tot=',u2tot,'u1tot=',u1tot,'duh=',u1tot-u2tot,'u2err=',u2err,'u1err=',u1err,'uh_err=',uh_err
863  if (abs(u1tot-u2tot)>u2err+u1err) &
864  write(0,*) 'U non-conservation difference=',u1tot-u2tot,'allowed err=',u2err+u1err,' <-----!'
865  write(0,*) 'Source to target:'
866  write(0,*) 'H: h0tot=',h0tot,'h1tot=',h1tot,'dh=',h1tot-h0tot,'h0err=',h0err,'h1err=',h1err
867  if (abs(h1tot-h0tot)>h0err+h1err) &
868  write(0,*) 'H non-conservation difference=',h1tot-h0tot,'allowed err=',h0err+h1err,' <-----!'
869  write(0,*) 'UH: u0tot=',u0tot,'u1tot=',u1tot,'duh=',u1tot-u0tot,'u0err=',u0err,'u1err=',u1err,'uh_err=',uh_err
870  if (abs(u1tot-u0tot)>(u0err+u1err)+uh_err) &
871  write(0,*) 'U non-conservation difference=',u1tot-u0tot,'allowed err=',u0err+u1err+uh_err,' <-----!'
872  write(0,*) 'U: u0min=',u0min,'u1min=',u1min,'u2min=',u2min
873  if (u1min<u0min) write(0,*) 'U minimum overshoot=',u1min-u0min,' <-----!'
874  if (u2min<u0min) write(0,*) 'U2 minimum overshoot=',u2min-u0min,' <-----!'
875  write(0,*) 'U: u0max=',u0max,'u1max=',u1max,'u2max=',u2max
876  if (u1max>u0max) write(0,*) 'U maximum overshoot=',u1max-u0max,' <-----!'
877  if (u2max>u0max) write(0,*) 'U2 maximum overshoot=',u2max-u0max,' <-----!'
878  write(0,'(a3,6a24,2a3)') 'k','h0','left edge','u0','right edge','h1','u1','is','ie'
879  do k = 1, max(n0,n1)
880  if (k<=min(n0,n1)) then
881  write(0,'(i3,1p6e24.16,2i3)') k,h0(k),ppoly0_e(k,1),u0(k),ppoly0_e(k,2),h1(k),u1(k),itgt_start(k),itgt_end(k)
882  elseif (k>n0) then
883  write(0,'(i3,96x,1p2e24.16,2i3)') k,h1(k),u1(k),itgt_start(k),itgt_end(k)
884  else
885  write(0,'(i3,1p4e24.16)') k,h0(k),ppoly0_e(k,1),u0(k),ppoly0_e(k,2)
886  endif
887  enddo
888  write(0,'(a3,2a24)') 'k','u0','Polynomial coefficients'
889  do k = 1, n0
890  write(0,'(i3,1p6e24.16)') k,u0(k),ppoly0_coefs(k,:)
891  enddo
892  write(0,'(a3,3a24,a3,2a24)') 'k','Sub-cell h','Sub-cell u','Sub-cell hu','i0','xa','xb'
893  xa = 0.
894  dh0_eff = 0.
895  do k = 1, n0+n1+1
896  dh = h_sub(k)
897  i0 = isub_src(k)
898  dh0_eff = dh0_eff + dh ! Cumulative thickness within the source cell
899  xb = dh0_eff / h0_eff(i0) ! This expression yields xa <= xb <= 1.0
900  xb = min(1., xb) ! This is only needed when the total target column is wider than the source column
901  write(0,'(i3,1p3e24.16,i3,1p2e24.16)') k,h_sub(k),u_sub(k),uh_sub(k),i0,xa,xb
902  if (k<=n0+n1) then
903  if (isub_src(k+1) /= i0) then
904  dh0_eff = 0.; xa = 0.
905  else
906  xa = xb
907  endif
908  endif
909  enddo
910  call mom_error( fatal, 'MOM_remapping, remap_via_sub_cells: '//&
911  'Remapping result is inconsistent!' )
912  endif
913  endif ! method<5
914  endif ! debug_bounds
915 
916  ! Include the error remapping from source to sub-cells in the estimate of total remapping error
917  uh_err = uh_err + u02_err
918 
919  if (present(ah_sub)) ah_sub(1:n0+n1+1) = h_sub(1:n0+n1+1)
920  if (present(aisub_src)) aisub_src(1:n0+n1+1) = isub_src(1:n0+n1+1)
921  if (present(aiss)) aiss(1:n0) = isrc_start(1:n0)
922  if (present(aise)) aise(1:n0) = isrc_end(1:n0)
923 
924 end subroutine remap_via_sub_cells
925 
926 !> Returns the average value of a reconstruction within a single source cell, i0,
927 !! between the non-dimensional positions xa and xb (xa<=xb) with dimensional
928 !! separation dh.
929 real function average_value_ppoly( n0, u0, ppoly0_E, ppoly0_coefs, method, i0, xa, xb)
930  integer, intent(in) :: n0 !< Number of cells in source grid
931  real, intent(in) :: u0(:) !< Cell means
932  real, intent(in) :: ppoly0_e(:,:) !< Edge value of polynomial
933  real, intent(in) :: ppoly0_coefs(:,:) !< Coefficients of polynomial
934  integer, intent(in) :: method !< Remapping scheme to use
935  integer, intent(in) :: i0 !< Source cell index
936  real, intent(in) :: xa !< Non-dimensional start position within source cell
937  real, intent(in) :: xb !< Non-dimensional end position within source cell
938  ! Local variables
939  real :: u_ave, xa_2, xb_2, xa2pxb2, xapxb
940  real, parameter :: r_3 = 1.0/3.0 ! Used in evaluation of integrated polynomials
941 
942  real :: mx, a_l, a_r, u_c, ya, yb, my, xa2b2ab, ya2b2ab, a_c
943 
944  if (xb > xa) then
945  select case ( method )
946  case ( integration_pcm )
947  u_ave = u0(i0)
948  case ( integration_plm )
949  u_ave = ( &
950  ppoly0_coefs(i0,1) &
951  + ppoly0_coefs(i0,2) * 0.5 * ( xb + xa ) )
952  case ( integration_ppm )
953  mx = 0.5 * ( xa + xb )
954  a_l = ppoly0_e(i0, 1)
955  a_r = ppoly0_e(i0, 2)
956  u_c = u0(i0)
957  a_c = 0.5 * ( ( u_c - a_l ) + ( u_c - a_r ) ) ! a_6 / 6
958  if (mx<0.5) then
959  ! This integration of the PPM reconstruction is expressed in distances from the left edge
960  xa2b2ab = (xa*xa+xb*xb)+xa*xb
961  u_ave = a_l + ( ( a_r - a_l ) * mx &
962  + a_c * ( 3. * ( xb + xa ) - 2.*xa2b2ab ) )
963  else
964  ! This integration of the PPM reconstruction is expressed in distances from the right edge
965  ya = 1. - xa
966  yb = 1. - xb
967  my = 0.5 * ( ya + yb )
968  ya2b2ab = (ya*ya+yb*yb)+ya*yb
969  u_ave = a_r + ( ( a_l - a_r ) * my &
970  + a_c * ( 3. * ( yb + ya ) - 2.*ya2b2ab ) )
971  endif
972  case ( integration_pqm )
973  xa_2 = xa*xa
974  xb_2 = xb*xb
975  xa2pxb2 = xa_2 + xb_2
976  xapxb = xa + xb
977  u_ave = ( &
978  ppoly0_coefs(i0,1) &
979  + ( ppoly0_coefs(i0,2) * 0.5 * ( xapxb ) &
980  + ( ppoly0_coefs(i0,3) * r_3 * ( xa2pxb2 + xa*xb ) &
981  + ( ppoly0_coefs(i0,4) * 0.25* ( xa2pxb2 * xapxb ) &
982  + ppoly0_coefs(i0,5) * 0.2 * ( ( xb*xb_2 + xa*xa_2 ) * xapxb + xa_2*xb_2 ) ) ) ) )
983  case default
984  call mom_error( fatal,'The selected integration method is invalid' )
985  end select
986  else ! dh == 0.
987  select case ( method )
988  case ( integration_pcm )
989  u_ave = ppoly0_coefs(i0,1)
990  case ( integration_plm )
991  !u_ave = ppoly0_coefs(i0,1) &
992  ! + xa * ppoly0_coefs(i0,2)
993  a_l = ppoly0_e(i0, 1)
994  a_r = ppoly0_e(i0, 2)
995  ya = 1. - xa
996  if (xa < 0.5) then
997  u_ave = a_l + xa * ( a_r - a_l )
998  else
999  u_ave = a_r + ya * ( a_l - a_r )
1000  endif
1001  case ( integration_ppm )
1002  !u_ave = ppoly0_coefs(i0,1) &
1003  ! + xa * ( ppoly0_coefs(i0,2) &
1004  ! + xa * ppoly0_coefs(i0,3) )
1005  a_l = ppoly0_e(i0, 1)
1006  a_r = ppoly0_e(i0, 2)
1007  u_c = u0(i0)
1008  a_c = 3. * ( ( u_c - a_l ) + ( u_c - a_r ) ) ! a_6
1009  ya = 1. - xa
1010  if (xa < 0.5) then
1011  u_ave = a_l + xa * ( ( a_r - a_l ) + a_c * ya )
1012  else
1013  u_ave = a_r + ya * ( ( a_l - a_r ) + a_c * xa )
1014  endif
1015  case ( integration_pqm )
1016  u_ave = ppoly0_coefs(i0,1) &
1017  + xa * ( ppoly0_coefs(i0,2) &
1018  + xa * ( ppoly0_coefs(i0,3) &
1019  + xa * ( ppoly0_coefs(i0,4) &
1020  + xa * ppoly0_coefs(i0,5) ) ) )
1021  case default
1022  call mom_error( fatal,'The selected integration method is invalid' )
1023  end select
1024  endif
1025  average_value_ppoly = u_ave
1026 
1027 end function average_value_ppoly
1028 
1029 !> Measure totals and bounds on source grid
1030 subroutine measure_input_bounds( n0, h0, u0, edge_values, h0tot, h0err, u0tot, u0err, u0min, u0max )
1031  integer, intent(in) :: n0 !< Number of cells on source grid
1032  real, dimension(n0), intent(in) :: h0 !< Cell widths on source grid
1033  real, dimension(n0), intent(in) :: u0 !< Cell averages on source grid
1034  real, dimension(n0,2), intent(in) :: edge_values !< Cell edge values on source grid
1035  real, intent(out) :: h0tot !< Sum of cell widths
1036  real, intent(out) :: h0err !< Magnitude of round-off error in h0tot
1037  real, intent(out) :: u0tot !< Sum of cell widths times values
1038  real, intent(out) :: u0err !< Magnitude of round-off error in u0tot
1039  real, intent(out) :: u0min !< Minimum value in reconstructions of u0
1040  real, intent(out) :: u0max !< Maximum value in reconstructions of u0
1041  ! Local variables
1042  integer :: k
1043  real :: eps
1044 
1045  eps = epsilon(h0(1))
1046  h0tot = h0(1)
1047  h0err = 0.
1048  u0tot = h0(1) * u0(1)
1049  u0err = 0.
1050  u0min = min( edge_values(1,1), edge_values(1,2) )
1051  u0max = max( edge_values(1,1), edge_values(1,2) )
1052  do k = 2, n0
1053  h0tot = h0tot + h0(k)
1054  h0err = h0err + eps * max(h0tot, h0(k))
1055  u0tot = u0tot + h0(k) * u0(k)
1056  u0err = u0err + eps * max(abs(u0tot), abs(h0(k) * u0(k)))
1057  u0min = min( u0min, edge_values(k,1), edge_values(k,2) )
1058  u0max = max( u0max, edge_values(k,1), edge_values(k,2) )
1059  enddo
1060 
1061 end subroutine measure_input_bounds
1062 
1063 !> Measure totals and bounds on destination grid
1064 subroutine measure_output_bounds( n1, h1, u1, h1tot, h1err, u1tot, u1err, u1min, u1max )
1065  integer, intent(in) :: n1 !< Number of cells on destination grid
1066  real, dimension(n1), intent(in) :: h1 !< Cell widths on destination grid
1067  real, dimension(n1), intent(in) :: u1 !< Cell averages on destination grid
1068  real, intent(out) :: h1tot !< Sum of cell widths
1069  real, intent(out) :: h1err !< Magnitude of round-off error in h1tot
1070  real, intent(out) :: u1tot !< Sum of cell widths times values
1071  real, intent(out) :: u1err !< Magnitude of round-off error in u1tot
1072  real, intent(out) :: u1min !< Minimum value in reconstructions of u1
1073  real, intent(out) :: u1max !< Maximum value in reconstructions of u1
1074  ! Local variables
1075  integer :: k
1076  real :: eps
1077 
1078  eps = epsilon(h1(1))
1079  h1tot = h1(1)
1080  h1err = 0.
1081  u1tot = h1(1) * u1(1)
1082  u1err = 0.
1083  u1min = u1(1)
1084  u1max = u1(1)
1085  do k = 2, n1
1086  h1tot = h1tot + h1(k)
1087  h1err = h1err + eps * max(h1tot, h1(k))
1088  u1tot = u1tot + h1(k) * u1(k)
1089  u1err = u1err + eps * max(abs(u1tot), abs(h1(k) * u1(k)))
1090  u1min = min(u1min, u1(k))
1091  u1max = max(u1max, u1(k))
1092  enddo
1093 
1094 end subroutine measure_output_bounds
1095 
1096 !> Remaps column of values u0 on grid h0 to grid h1 by integrating
1097 !! over the projection of each h1 cell onto the h0 grid.
1098 subroutine remapbyprojection( n0, h0, u0, ppoly0_E, ppoly0_coefs, &
1099  n1, h1, method, u1, h_neglect )
1100  integer, intent(in) :: n0 !< Number of cells in source grid
1101  real, intent(in) :: h0(:) !< Source grid widths (size n0)
1102  real, intent(in) :: u0(:) !< Source cell averages (size n0)
1103  real, intent(in) :: ppoly0_E(:,:) !< Edge value of polynomial
1104  real, intent(in) :: ppoly0_coefs(:,:) !< Coefficients of polynomial
1105  integer, intent(in) :: n1 !< Number of cells in target grid
1106  real, intent(in) :: h1(:) !< Target grid widths (size n1)
1107  integer, intent(in) :: method !< Remapping scheme to use
1108  real, intent(out) :: u1(:) !< Target cell averages (size n1)
1109  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
1110  !! purpose of cell reconstructions
1111  !! in the same units as h.
1112  ! Local variables
1113  integer :: iTarget
1114  real :: xL, xR ! coordinates of target cell edges
1115  integer :: jStart ! Used by integrateReconOnInterval()
1116  real :: xStart ! Used by integrateReconOnInterval()
1117 
1118  ! Loop on cells in target grid (grid1). For each target cell, we need to find
1119  ! in which source cells the target cell edges lie. The associated indexes are
1120  ! noted j0 and j1.
1121  xr = 0. ! Left boundary is at x=0
1122  jstart = 1
1123  xstart = 0.
1124  do itarget = 1,n1
1125  ! Determine the coordinates of the target cell edges
1126  xl = xr
1127  xr = xl + h1(itarget)
1128 
1129  call integraterecononinterval( n0, h0, u0, ppoly0_e, ppoly0_coefs, method, &
1130  xl, xr, h1(itarget), u1(itarget), jstart, xstart, h_neglect )
1131 
1132  enddo ! end iTarget loop on target grid cells
1133 
1134 end subroutine remapbyprojection
1135 
1136 
1137 !> Remaps column of values u0 on grid h0 to implied grid h1
1138 !! where the interfaces of h1 differ from those of h0 by dx.
1139 !! The new grid is defined relative to the original grid by change
1140 !! dx1(:) = xNew(:) - xOld(:)
1141 !! and the remapping calculated so that
1142 !! hNew(k) qNew(k) = hOld(k) qOld(k) + F(k+1) - F(k)
1143 !! where
1144 !! F(k) = dx1(k) qAverage
1145 !! and where qAverage is the average qOld in the region zOld(k) to zNew(k).
1146 subroutine remapbydeltaz( n0, h0, u0, ppoly0_E, ppoly0_coefs, n1, dx1, &
1147  method, u1, h1, h_neglect )
1148  integer, intent(in) :: n0 !< Number of cells in source grid
1149  real, dimension(:), intent(in) :: h0 !< Source grid sizes (size n0)
1150  real, dimension(:), intent(in) :: u0 !< Source cell averages (size n0)
1151  real, dimension(:,:), intent(in) :: ppoly0_E !< Edge value of polynomial
1152  real, dimension(:,:), intent(in) :: ppoly0_coefs !< Coefficients of polynomial
1153  integer, intent(in) :: n1 !< Number of cells in target grid
1154  real, dimension(:), intent(in) :: dx1 !< Target grid edge positions (size n1+1)
1155  integer, intent(in) :: method !< Remapping scheme to use
1156  real, dimension(:), intent(out) :: u1 !< Target cell averages (size n1)
1157  real, dimension(:), &
1158  optional, intent(out) :: h1 !< Target grid widths (size n1)
1159  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
1160  !! purpose of cell reconstructions
1161  !! in the same units as h.
1162  ! Local variables
1163  integer :: iTarget
1164  real :: xL, xR ! coordinates of target cell edges
1165  real :: xOld, hOld, uOld
1166  real :: xNew, hNew, h_err
1167  real :: uhNew, hFlux, uAve, fluxL, fluxR
1168  integer :: jStart ! Used by integrateReconOnInterval()
1169  real :: xStart ! Used by integrateReconOnInterval()
1170 
1171  ! Loop on cells in target grid. For each cell, iTarget, the left flux is
1172  ! the right flux of the cell to the left, iTarget-1.
1173  ! The left flux is initialized by started at iTarget=0 to calculate the
1174  ! right flux which can take into account the target left boundary being
1175  ! in the interior of the source domain.
1176  fluxr = 0.
1177  h_err = 0. ! For measuring round-off error
1178  jstart = 1
1179  xstart = 0.
1180  do itarget = 0,n1
1181  fluxl = fluxr ! This does nothing for iTarget=0
1182 
1183  if (itarget == 0) then
1184  xold = 0. ! Left boundary is at x=0
1185  hold = -1.e30 ! Should not be used for iTarget = 0
1186  uold = -1.e30 ! Should not be used for iTarget = 0
1187  elseif (itarget <= n0) then
1188  xold = xold + h0(itarget) ! Position of right edge of cell
1189  hold = h0(itarget)
1190  uold = u0(itarget)
1191  h_err = h_err + epsilon(hold) * max(hold, xold)
1192  else
1193  hold = 0. ! as if for layers>n0, they were vanished
1194  uold = 1.e30 ! and the initial value should not matter
1195  endif
1196  xnew = xold + dx1(itarget+1)
1197  xl = min( xold, xnew )
1198  xr = max( xold, xnew )
1199 
1200  ! hFlux is the positive width of the remapped volume
1201  hflux = abs(dx1(itarget+1))
1202  call integraterecononinterval( n0, h0, u0, ppoly0_e, ppoly0_coefs, method, &
1203  xl, xr, hflux, uave, jstart, xstart )
1204  ! uAve is the average value of u, independent of sign of dx1
1205  fluxr = dx1(itarget+1)*uave ! Includes sign of dx1
1206 
1207  if (itarget>0) then
1208  hnew = hold + ( dx1(itarget+1) - dx1(itarget) )
1209  hnew = max( 0., hnew )
1210  uhnew = ( uold * hold ) + ( fluxr - fluxl )
1211  if (hnew>0.) then
1212  u1(itarget) = uhnew / hnew
1213  else
1214  u1(itarget) = uave
1215  endif
1216  if (present(h1)) h1(itarget) = hnew
1217  endif
1218 
1219  enddo ! end iTarget loop on target grid cells
1220 
1221 end subroutine remapbydeltaz
1222 
1223 
1224 !> Integrate the reconstructed column profile over a single cell
1225 subroutine integraterecononinterval( n0, h0, u0, ppoly0_E, ppoly0_coefs, method, &
1226  xL, xR, hC, uAve, jStart, xStart, h_neglect )
1227  integer, intent(in) :: n0 !< Number of cells in source grid
1228  real, dimension(:), intent(in) :: h0 !< Source grid sizes (size n0)
1229  real, dimension(:), intent(in) :: u0 !< Source cell averages
1230  real, dimension(:,:), intent(in) :: ppoly0_E !< Edge value of polynomial
1231  real, dimension(:,:), intent(in) :: ppoly0_coefs !< Coefficients of polynomial
1232  integer, intent(in) :: method !< Remapping scheme to use
1233  real, intent(in) :: xL !< Left edges of target cell
1234  real, intent(in) :: xR !< Right edges of target cell
1235  real, intent(in) :: hC !< Cell width hC = xR - xL
1236  real, intent(out) :: uAve !< Average value on target cell
1237  integer, intent(inout) :: jStart !< The index of the cell to start searching from
1238  !< On exit, contains index of last cell used
1239  real, intent(inout) :: xStart !< The left edge position of cell jStart
1240  !< On first entry should be 0.
1241  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
1242  !! purpose of cell reconstructions
1243  !! in the same units as h.
1244  ! Local variables
1245  integer :: j, k
1246  integer :: jL, jR ! indexes of source cells containing target
1247  ! cell edges
1248  real :: q ! complete integration
1249  real :: xi0, xi1 ! interval of integration (local -- normalized
1250  ! -- coordinates)
1251  real :: x0jLl, x0jLr ! Left/right position of cell jL
1252  real :: x0jRl, x0jRr ! Left/right position of cell jR
1253  real :: hAct ! The distance actually used in the integration
1254  ! (notionally xR - xL) which differs due to roundoff.
1255  real :: x0_2, x1_2, x02px12, x0px1 ! Used in evaluation of integrated polynomials
1256  real :: hNeglect ! A negligible thicness in the same units as h.
1257  real, parameter :: r_3 = 1.0/3.0 ! Used in evaluation of integrated polynomials
1258 
1259  hneglect = hneglect_dflt ; if (present(h_neglect)) hneglect = h_neglect
1260 
1261  q = -1.e30
1262  x0jll = -1.e30
1263  x0jrl = -1.e30
1264 
1265  ! Find the left most cell in source grid spanned by the target cell
1266  jl = -1
1267  x0jlr = xstart
1268  do j = jstart, n0
1269  x0jll = x0jlr
1270  x0jlr = x0jll + h0(j)
1271  ! Left edge is found in cell j
1272  if ( ( xl >= x0jll ) .AND. ( xl <= x0jlr ) ) then
1273  jl = j
1274  exit ! once target grid cell is found, exit loop
1275  endif
1276  enddo
1277  jstart = jl
1278  xstart = x0jll
1279 
1280 ! ! HACK to handle round-off problems. Need only at j=n0.
1281 ! ! This moves the effective cell boundary outwards a smidgen.
1282 ! if (xL>x0jLr) x0jLr = xL
1283 
1284  ! If, at this point, jL is equal to -1, it means the vanished
1285  ! cell lies outside the source grid. In other words, it means that
1286  ! the source and target grids do not cover the same physical domain
1287  ! and there is something very wrong !
1288  if ( jl == -1 ) call mom_error(fatal, &
1289  'MOM_remapping, integrateReconOnInterval: '//&
1290  'The location of the left-most cell could not be found')
1291 
1292 
1293  ! ============================================================
1294  ! Check whether target cell is vanished. If it is, the cell
1295  ! average is simply the interpolated value at the location
1296  ! of the vanished cell. If it isn't, we need to integrate the
1297  ! quantity within the cell and divide by the cell width to
1298  ! determine the cell average.
1299  ! ============================================================
1300  ! 1. Cell is vanished
1301  !if ( abs(xR - xL) <= epsilon(xR)*max(abs(xR),abs(xL)) ) then
1302  if ( abs(xr - xl) == 0.0 ) then
1303 
1304  ! We check whether the source cell (i.e. the cell in which the
1305  ! vanished target cell lies) is vanished. If it is, the interpolated
1306  ! value is set to be mean of the edge values (which should be the same).
1307  ! If it isn't, we simply interpolate.
1308  if ( h0(jl) == 0.0 ) then
1309  uave = 0.5 * ( ppoly0_e(jl,1) + ppoly0_e(jl,2) )
1310  else
1311  ! WHY IS THIS NOT WRITTEN AS xi0 = ( xL - x0jLl ) / h0(jL) ---AJA
1312  xi0 = xl / ( h0(jl) + hneglect ) - x0jll / ( h0(jl) + hneglect )
1313 
1314  select case ( method )
1315  case ( integration_pcm )
1316  uave = ppoly0_coefs(jl,1)
1317  case ( integration_plm )
1318  uave = ppoly0_coefs(jl,1) &
1319  + xi0 * ppoly0_coefs(jl,2)
1320  case ( integration_ppm )
1321  uave = ppoly0_coefs(jl,1) &
1322  + xi0 * ( ppoly0_coefs(jl,2) &
1323  + xi0 * ppoly0_coefs(jl,3) )
1324  case ( integration_pqm )
1325  uave = ppoly0_coefs(jl,1) &
1326  + xi0 * ( ppoly0_coefs(jl,2) &
1327  + xi0 * ( ppoly0_coefs(jl,3) &
1328  + xi0 * ( ppoly0_coefs(jl,4) &
1329  + xi0 * ppoly0_coefs(jl,5) ) ) )
1330  case default
1331  call mom_error( fatal,'The selected integration method is invalid' )
1332  end select
1333 
1334  endif ! end checking whether source cell is vanished
1335 
1336  ! 2. Cell is not vanished
1337  else
1338 
1339  ! Find the right most cell in source grid spanned by the target cell
1340  jr = -1
1341  x0jrr = xstart
1342  do j = jstart,n0
1343  x0jrl = x0jrr
1344  x0jrr = x0jrl + h0(j)
1345  ! Right edge is found in cell j
1346  if ( ( xr >= x0jrl ) .AND. ( xr <= x0jrr ) ) then
1347  jr = j
1348  exit ! once target grid cell is found, exit loop
1349  endif
1350  enddo ! end loop on source grid cells
1351 
1352  ! If xR>x0jRr then the previous loop reached j=n0 and the target
1353  ! position, xR, was beyond the right edge of the source grid (h0).
1354  ! This can happen due to roundoff, in which case we set jR=n0.
1355  if (xr>x0jrr) jr = n0
1356 
1357  ! To integrate, two cases must be considered: (1) the target cell is
1358  ! entirely contained within a cell of the source grid and (2) the target
1359  ! cell spans at least two cells of the source grid.
1360 
1361  if ( jl == jr ) then
1362  ! The target cell is entirely contained within a cell of the source
1363  ! grid. This situation is represented by the following schematic, where
1364  ! the cell in which xL and xR are located has index jL=jR :
1365  !
1366  ! ----|-----o--------o----------|-------------
1367  ! xL xR
1368  !
1369  ! Determine normalized coordinates
1370 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__
1371  xi0 = max( 0., min( 1., ( xl - x0jll ) / ( h0(jl) + hneglect ) ) )
1372  xi1 = max( 0., min( 1., ( xr - x0jll ) / ( h0(jl) + hneglect ) ) )
1373 #else
1374  xi0 = xl / h0(jl) - x0jll / ( h0(jl) + hneglect )
1375  xi1 = xr / h0(jl) - x0jll / ( h0(jl) + hneglect )
1376 #endif
1377 
1378  hact = h0(jl) * ( xi1 - xi0 )
1379 
1380  ! Depending on which polynomial is used, integrate quantity
1381  ! between xi0 and xi1. Integration is carried out in normalized
1382  ! coordinates, hence: \int_xL^xR p(x) dx = h \int_xi0^xi1 p(xi) dxi
1383  select case ( method )
1384  case ( integration_pcm )
1385  q = ( xr - xl ) * ppoly0_coefs(jl,1)
1386  case ( integration_plm )
1387  q = ( xr - xl ) * ( &
1388  ppoly0_coefs(jl,1) &
1389  + ppoly0_coefs(jl,2) * 0.5 * ( xi1 + xi0 ) )
1390  case ( integration_ppm )
1391  q = ( xr - xl ) * ( &
1392  ppoly0_coefs(jl,1) &
1393  + ( ppoly0_coefs(jl,2) * 0.5 * ( xi1 + xi0 ) &
1394  + ppoly0_coefs(jl,3) * r_3 * ( ( xi1*xi1 + xi0*xi0 ) + xi0*xi1 ) ) )
1395  case ( integration_pqm )
1396  x0_2 = xi0*xi0
1397  x1_2 = xi1*xi1
1398  x02px12 = x0_2 + x1_2
1399  x0px1 = xi1 + xi0
1400  q = ( xr - xl ) * ( &
1401  ppoly0_coefs(jl,1) &
1402  + ( ppoly0_coefs(jl,2) * 0.5 * ( xi1 + xi0 ) &
1403  + ( ppoly0_coefs(jl,3) * r_3 * ( x02px12 + xi0*xi1 ) &
1404  + ppoly0_coefs(jl,4) * 0.25* ( x02px12 * x0px1 ) &
1405  + ppoly0_coefs(jl,5) * 0.2 * ( ( xi1*x1_2 + xi0*x0_2 ) * x0px1 + x0_2*x1_2 ) ) ) )
1406  case default
1407  call mom_error( fatal,'The selected integration method is invalid' )
1408  end select
1409 
1410  else
1411  ! The target cell spans at least two cells of the source grid.
1412  ! This situation is represented by the following schematic, where
1413  ! the cells in which xL and xR are located have indexes jL and jR,
1414  ! respectively :
1415  !
1416  ! ----|-----o---|--- ... --|---o----------|-------------
1417  ! xL xR
1418  !
1419  ! We first integrate from xL up to the right boundary of cell jL, then
1420  ! add the integrated amounts of cells located between jL and jR and then
1421  ! integrate from the left boundary of cell jR up to xR
1422 
1423  q = 0.0
1424 
1425  ! Integrate from xL up to right boundary of cell jL
1426 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__
1427  xi0 = max( 0., min( 1., ( xl - x0jll ) / ( h0(jl) + hneglect ) ) )
1428 #else
1429  xi0 = (xl - x0jll) / ( h0(jl) + hneglect )
1430 #endif
1431  xi1 = 1.0
1432 
1433  hact = h0(jl) * ( xi1 - xi0 )
1434 
1435  select case ( method )
1436  case ( integration_pcm )
1437  q = q + ( x0jlr - xl ) * ppoly0_coefs(jl,1)
1438  case ( integration_plm )
1439  q = q + ( x0jlr - xl ) * ( &
1440  ppoly0_coefs(jl,1) &
1441  + ppoly0_coefs(jl,2) * 0.5 * ( xi1 + xi0 ) )
1442  case ( integration_ppm )
1443  q = q + ( x0jlr - xl ) * ( &
1444  ppoly0_coefs(jl,1) &
1445  + ( ppoly0_coefs(jl,2) * 0.5 * ( xi1 + xi0 ) &
1446  + ppoly0_coefs(jl,3) * r_3 * ( ( xi1*xi1 + xi0*xi0 ) + xi0*xi1 ) ) )
1447  case ( integration_pqm )
1448  x0_2 = xi0*xi0
1449  x1_2 = xi1*xi1
1450  x02px12 = x0_2 + x1_2
1451  x0px1 = xi1 + xi0
1452  q = q + ( x0jlr - xl ) * ( &
1453  ppoly0_coefs(jl,1) &
1454  + ( ppoly0_coefs(jl,2) * 0.5 * ( xi1 + xi0 ) &
1455  + ( ppoly0_coefs(jl,3) * r_3 * ( x02px12 + xi0*xi1 ) &
1456  + ppoly0_coefs(jl,4) * 0.25* ( x02px12 * x0px1 ) &
1457  + ppoly0_coefs(jl,5) * 0.2 * ( ( xi1*x1_2 + xi0*x0_2 ) * x0px1 + x0_2*x1_2 ) ) ) )
1458  case default
1459  call mom_error( fatal, 'The selected integration method is invalid' )
1460  end select
1461 
1462  ! Integrate contents within cells strictly comprised between jL and jR
1463  if ( jr > (jl+1) ) then
1464  do k = jl+1,jr-1
1465  q = q + h0(k) * u0(k)
1466  hact = hact + h0(k)
1467  enddo
1468  endif
1469 
1470  ! Integrate from left boundary of cell jR up to xR
1471  xi0 = 0.0
1472 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__
1473  xi1 = max( 0., min( 1., ( xr - x0jrl ) / ( h0(jr) + hneglect ) ) )
1474 #else
1475  xi1 = (xr - x0jrl) / ( h0(jr) + hneglect )
1476 #endif
1477 
1478  hact = hact + h0(jr) * ( xi1 - xi0 )
1479 
1480  select case ( method )
1481  case ( integration_pcm )
1482  q = q + ( xr - x0jrl ) * ppoly0_coefs(jr,1)
1483  case ( integration_plm )
1484  q = q + ( xr - x0jrl ) * ( &
1485  ppoly0_coefs(jr,1) &
1486  + ppoly0_coefs(jr,2) * 0.5 * ( xi1 + xi0 ) )
1487  case ( integration_ppm )
1488  q = q + ( xr - x0jrl ) * ( &
1489  ppoly0_coefs(jr,1) &
1490  + ( ppoly0_coefs(jr,2) * 0.5 * ( xi1 + xi0 ) &
1491  + ppoly0_coefs(jr,3) * r_3 * ( ( xi1*xi1 + xi0*xi0 ) + xi0*xi1 ) ) )
1492  case ( integration_pqm )
1493  x0_2 = xi0*xi0
1494  x1_2 = xi1*xi1
1495  x02px12 = x0_2 + x1_2
1496  x0px1 = xi1 + xi0
1497  q = q + ( xr - x0jrl ) * ( &
1498  ppoly0_coefs(jr,1) &
1499  + ( ppoly0_coefs(jr,2) * 0.5 * ( xi1 + xi0 ) &
1500  + ( ppoly0_coefs(jr,3) * r_3 * ( x02px12 + xi0*xi1 ) &
1501  + ppoly0_coefs(jr,4) * 0.25* ( x02px12 * x0px1 ) &
1502  + ppoly0_coefs(jr,5) * 0.2 * ( ( xi1*x1_2 + xi0*x0_2 ) * x0px1 + x0_2*x1_2 ) ) ) )
1503  case default
1504  call mom_error( fatal,'The selected integration method is invalid' )
1505  end select
1506 
1507  endif ! end integration for non-vanished cells
1508 
1509  ! The cell average is the integrated value divided by the cell width
1510 #ifdef __USE_ROUNDOFF_SAFE_ADJUSTMENTS__
1511 if (hact==0.) then
1512  uave = ppoly0_coefs(jl,1)
1513 else
1514  uave = q / hact
1515 endif
1516 #else
1517  uave = q / hc
1518 #endif
1519 
1520  endif ! endif clause to check if cell is vanished
1521 
1522 end subroutine integraterecononinterval
1523 
1524 !> Calculates the change in interface positions based on h1 and h2
1525 subroutine dzfromh1h2( n1, h1, n2, h2, dx )
1526  integer, intent(in) :: n1 !< Number of cells on source grid
1527  real, dimension(:), intent(in) :: h1 !< Cell widths of source grid (size n1)
1528  integer, intent(in) :: n2 !< Number of cells on target grid
1529  real, dimension(:), intent(in) :: h2 !< Cell widths of target grid (size n2)
1530  real, dimension(:), intent(out) :: dx !< Change in interface position (size n2+1)
1531  ! Local variables
1532  integer :: k
1533  real :: x1, x2
1534 
1535  x1 = 0.
1536  x2 = 0.
1537  dx(1) = 0.
1538  do k = 1, max(n1,n2)
1539  if (k <= n1) x1 = x1 + h1(k) ! Interface k+1, right of source cell k
1540  if (k <= n2) then
1541  x2 = x2 + h2(k) ! Interface k+1, right of target cell k
1542  dx(k+1) = x2 - x1 ! Change of interface k+1, target - source
1543  endif
1544  enddo
1545 
1546 end subroutine dzfromh1h2
1547 
1548 !> Constructor for remapping control structure
1549 subroutine initialize_remapping( CS, remapping_scheme, boundary_extrapolation, &
1550  check_reconstruction, check_remapping, force_bounds_in_subcell, answers_2018)
1551  ! Arguments
1552  type(remapping_cs), intent(inout) :: cs !< Remapping control structure
1553  character(len=*), intent(in) :: remapping_scheme !< Remapping scheme to use
1554  logical, optional, intent(in) :: boundary_extrapolation !< Indicate to extrapolate in boundary cells
1555  logical, optional, intent(in) :: check_reconstruction !< Indicate to check reconstructions
1556  logical, optional, intent(in) :: check_remapping !< Indicate to check results of remapping
1557  logical, optional, intent(in) :: force_bounds_in_subcell !< Force subcells values to be bounded
1558  logical, optional, intent(in) :: answers_2018 !< If true use older, less acccurate expressions.
1559 
1560  ! Note that remapping_scheme is mandatory for initialize_remapping()
1561  call remapping_set_param(cs, remapping_scheme=remapping_scheme, boundary_extrapolation=boundary_extrapolation, &
1562  check_reconstruction=check_reconstruction, check_remapping=check_remapping, &
1563  force_bounds_in_subcell=force_bounds_in_subcell, answers_2018=answers_2018)
1564 
1565 end subroutine initialize_remapping
1566 
1567 !> Changes the method of reconstruction
1568 !! Use this routine to parse a string parameter specifying the reconstruction
1569 !! and re-allocates work arrays appropriately. It is called from
1570 !! initialize_remapping but can be called from an external module too.
1571 subroutine setreconstructiontype(string,CS)
1572  character(len=*), intent(in) :: string !< String to parse for method
1573  type(remapping_cs), intent(inout) :: CS !< Remapping control structure
1574  ! Local variables
1575  integer :: degree
1576  degree = -99
1577  select case ( uppercase(trim(string)) )
1578  case ("PCM")
1579  cs%remapping_scheme = remapping_pcm
1580  degree = 0
1581  case ("PLM")
1582  cs%remapping_scheme = remapping_plm
1583  degree = 1
1584  case ("PPM_H4")
1585  cs%remapping_scheme = remapping_ppm_h4
1586  degree = 2
1587  case ("PPM_IH4")
1588  cs%remapping_scheme = remapping_ppm_ih4
1589  degree = 2
1590  case ("PQM_IH4IH3")
1591  cs%remapping_scheme = remapping_pqm_ih4ih3
1592  degree = 4
1593  case ("PQM_IH6IH5")
1594  cs%remapping_scheme = remapping_pqm_ih6ih5
1595  degree = 4
1596  case default
1597  call mom_error(fatal, "setReconstructionType: "//&
1598  "Unrecognized choice for REMAPPING_SCHEME ("//trim(string)//").")
1599  end select
1600 
1601  cs%degree = degree
1602 
1603 end subroutine setreconstructiontype
1604 
1605 !> Destrcutor for remapping control structure
1606 subroutine end_remapping(CS)
1607  type(remapping_cs), intent(inout) :: cs !< Remapping control structure
1608 
1609  cs%degree = 0
1610 
1611 end subroutine end_remapping
1612 
1613 !> Runs unit tests on remapping functions.
1614 !! Should only be called from a single/root thread
1615 !! Returns True if a test fails, otherwise False
1616 logical function remapping_unit_tests(verbose)
1617  logical, intent(in) :: verbose !< If true, write results to stdout
1618  ! Local variables
1619  integer, parameter :: n0 = 4, n1 = 3, n2 = 6
1620  real :: h0(n0), x0(n0+1), u0(n0)
1621  real :: h1(n1), x1(n1+1), u1(n1), hn1(n1), dx1(n1+1)
1622  real :: h2(n2), x2(n2+1), u2(n2), hn2(n2), dx2(n2+1)
1623  data u0 /9., 3., -3., -9./ ! Linear profile, 4 at surface to -4 at bottom
1624  data h0 /4*0.75/ ! 4 uniform layers with total depth of 3
1625  data h1 /3*1./ ! 3 uniform layers with total depth of 3
1626  data h2 /6*0.5/ ! 6 uniform layers with total depth of 3
1627  type(remapping_cs) :: cs !< Remapping control structure
1628  real, allocatable, dimension(:,:) :: ppoly0_e, ppoly0_s, ppoly0_coefs
1629  logical :: answers_2018 ! If true use older, less acccurate expressions.
1630  integer :: i
1631  real :: err, h_neglect, h_neglect_edge
1632  logical :: thistest, v
1633 
1634  v = verbose
1635  answers_2018 = .false. ! .true.
1636  h_neglect = hneglect_dflt
1637  h_neglect_edge = hneglect_dflt ; if (answers_2018) h_neglect_edge = 1.0e-10
1638 
1639  write(*,*) '==== MOM_remapping: remapping_unit_tests ================='
1640  remapping_unit_tests = .false. ! Normally return false
1641 
1642  thistest = .false.
1643  call buildgridfromh(n0, h0, x0)
1644  do i=1,n0+1
1645  err=x0(i)-0.75*real(i-1)
1646  if (abs(err)>real(i-1)*epsilon(err)) thistest = .true.
1647  enddo
1648  if (thistest) write(*,*) 'remapping_unit_tests: Failed buildGridFromH() 1'
1649  remapping_unit_tests = remapping_unit_tests .or. thistest
1650  call buildgridfromh(n1, h1, x1)
1651  do i=1,n1+1
1652  err=x1(i)-real(i-1)
1653  if (abs(err)>real(i-1)*epsilon(err)) thistest = .true.
1654  enddo
1655  if (thistest) write(*,*) 'remapping_unit_tests: Failed buildGridFromH() 2'
1656  remapping_unit_tests = remapping_unit_tests .or. thistest
1657 
1658  thistest = .false.
1659  call initialize_remapping(cs, 'PPM_H4', answers_2018=answers_2018)
1660  if (verbose) write(*,*) 'h0 (test data)'
1661  if (verbose) call dumpgrid(n0,h0,x0,u0)
1662 
1663  call dzfromh1h2( n0, h0, n1, h1, dx1 )
1664  call remapping_core_w( cs, n0, h0, u0, n1, dx1, u1, h_neglect, h_neglect_edge)
1665  do i=1,n1
1666  err=u1(i)-8.*(0.5*real(1+n1)-real(i))
1667  if (abs(err)>real(n1-1)*epsilon(err)) thistest = .true.
1668  enddo
1669  if (verbose) write(*,*) 'h1 (by projection)'
1670  if (verbose) call dumpgrid(n1,h1,x1,u1)
1671  if (thistest) write(*,*) 'remapping_unit_tests: Failed remapping_core_w()'
1672  remapping_unit_tests = remapping_unit_tests .or. thistest
1673 
1674  thistest = .false.
1675  allocate(ppoly0_e(n0,2))
1676  allocate(ppoly0_s(n0,2))
1677  allocate(ppoly0_coefs(n0,cs%degree+1))
1678 
1679  ppoly0_e(:,:) = 0.0
1680  ppoly0_s(:,:) = 0.0
1681  ppoly0_coefs(:,:) = 0.0
1682 
1683  call edge_values_explicit_h4( n0, h0, u0, ppoly0_e, h_neglect=1e-10, answers_2018=answers_2018 )
1684  call ppm_reconstruction( n0, h0, u0, ppoly0_e, ppoly0_coefs, h_neglect, answers_2018=answers_2018 )
1685  call ppm_boundary_extrapolation( n0, h0, u0, ppoly0_e, ppoly0_coefs, h_neglect )
1686  u1(:) = 0.
1687  call remapbyprojection( n0, h0, u0, ppoly0_e, ppoly0_coefs, &
1688  n1, h1, integration_ppm, u1, h_neglect )
1689  do i=1,n1
1690  err=u1(i)-8.*(0.5*real(1+n1)-real(i))
1691  if (abs(err)>2.*epsilon(err)) thistest = .true.
1692  enddo
1693  if (thistest) write(*,*) 'remapping_unit_tests: Failed remapByProjection()'
1694  remapping_unit_tests = remapping_unit_tests .or. thistest
1695 
1696  thistest = .false.
1697  u1(:) = 0.
1698  call remapbydeltaz( n0, h0, u0, ppoly0_e, ppoly0_coefs, &
1699  n1, x1-x0(1:n1+1), &
1700  integration_ppm, u1, hn1, h_neglect )
1701  if (verbose) write(*,*) 'h1 (by delta)'
1702  if (verbose) call dumpgrid(n1,h1,x1,u1)
1703  hn1=hn1-h1
1704  do i=1,n1
1705  err=u1(i)-8.*(0.5*real(1+n1)-real(i))
1706  if (abs(err)>2.*epsilon(err)) thistest = .true.
1707  enddo
1708  if (thistest) write(*,*) 'remapping_unit_tests: Failed remapByDeltaZ() 1'
1709  remapping_unit_tests = remapping_unit_tests .or. thistest
1710 
1711  thistest = .false.
1712  call buildgridfromh(n2, h2, x2)
1713  dx2(1:n0+1) = x2(1:n0+1) - x0
1714  dx2(n0+2:n2+1) = x2(n0+2:n2+1) - x0(n0+1)
1715  call remapbydeltaz( n0, h0, u0, ppoly0_e, ppoly0_coefs, &
1716  n2, dx2, &
1717  integration_ppm, u2, hn2, h_neglect )
1718  if (verbose) write(*,*) 'h2'
1719  if (verbose) call dumpgrid(n2,h2,x2,u2)
1720  if (verbose) write(*,*) 'hn2'
1721  if (verbose) call dumpgrid(n2,hn2,x2,u2)
1722 
1723  do i=1,n2
1724  err=u2(i)-8./2.*(0.5*real(1+n2)-real(i))
1725  if (abs(err)>2.*epsilon(err)) thistest = .true.
1726  enddo
1727  if (thistest) write(*,*) 'remapping_unit_tests: Failed remapByDeltaZ() 2'
1728  remapping_unit_tests = remapping_unit_tests .or. thistest
1729 
1730  if (verbose) write(*,*) 'Via sub-cells'
1731  thistest = .false.
1732  call remap_via_sub_cells( n0, h0, u0, ppoly0_e, ppoly0_coefs, &
1733  n2, h2, integration_ppm, .false., u2, err )
1734  if (verbose) call dumpgrid(n2,h2,x2,u2)
1735 
1736  do i=1,n2
1737  err=u2(i)-8./2.*(0.5*real(1+n2)-real(i))
1738  if (abs(err)>2.*epsilon(err)) thistest = .true.
1739  enddo
1740  if (thistest) write(*,*) 'remapping_unit_tests: Failed remap_via_sub_cells() 2'
1741  remapping_unit_tests = remapping_unit_tests .or. thistest
1742 
1743  call remap_via_sub_cells( n0, h0, u0, ppoly0_e, ppoly0_coefs, &
1744  6, (/.125,.125,.125,.125,.125,.125/), integration_ppm, .false., u2, err )
1745  if (verbose) call dumpgrid(6,h2,x2,u2)
1746 
1747  call remap_via_sub_cells( n0, h0, u0, ppoly0_e, ppoly0_coefs, &
1748  3, (/2.25,1.5,1./), integration_ppm, .false., u2, err )
1749  if (verbose) call dumpgrid(3,h2,x2,u2)
1750 
1751  if (.not. remapping_unit_tests) write(*,*) 'Pass'
1752 
1753  write(*,*) '===== MOM_remapping: new remapping_unit_tests =================='
1754 
1755  deallocate(ppoly0_e, ppoly0_s, ppoly0_coefs)
1756  allocate(ppoly0_coefs(5,6))
1757  allocate(ppoly0_e(5,2))
1758  allocate(ppoly0_s(5,2))
1759 
1760  call pcm_reconstruction(3, (/1.,2.,4./), ppoly0_e(1:3,:), &
1761  ppoly0_coefs(1:3,:) )
1762  remapping_unit_tests = remapping_unit_tests .or. &
1763  test_answer(v, 3, ppoly0_e(:,1), (/1.,2.,4./), 'PCM: left edges')
1764  remapping_unit_tests = remapping_unit_tests .or. &
1765  test_answer(v, 3, ppoly0_e(:,2), (/1.,2.,4./), 'PCM: right edges')
1766  remapping_unit_tests = remapping_unit_tests .or. &
1767  test_answer(v, 3, ppoly0_coefs(:,1), (/1.,2.,4./), 'PCM: P0')
1768 
1769  call plm_reconstruction(3, (/1.,1.,1./), (/1.,3.,5./), ppoly0_e(1:3,:), &
1770  ppoly0_coefs(1:3,:), h_neglect )
1771  remapping_unit_tests = remapping_unit_tests .or. &
1772  test_answer(v, 3, ppoly0_e(:,1), (/1.,2.,5./), 'Unlim PLM: left edges')
1773  remapping_unit_tests = remapping_unit_tests .or. &
1774  test_answer(v, 3, ppoly0_e(:,2), (/1.,4.,5./), 'Unlim PLM: right edges')
1775  remapping_unit_tests = remapping_unit_tests .or. &
1776  test_answer(v, 3, ppoly0_coefs(:,1), (/1.,2.,5./), 'Unlim PLM: P0')
1777  remapping_unit_tests = remapping_unit_tests .or. &
1778  test_answer(v, 3, ppoly0_coefs(:,2), (/0.,2.,0./), 'Unlim PLM: P1')
1779 
1780  call plm_reconstruction(3, (/1.,1.,1./), (/1.,2.,7./), ppoly0_e(1:3,:), &
1781  ppoly0_coefs(1:3,:), h_neglect )
1782  remapping_unit_tests = remapping_unit_tests .or. &
1783  test_answer(v, 3, ppoly0_e(:,1), (/1.,1.,7./), 'Left lim PLM: left edges')
1784  remapping_unit_tests = remapping_unit_tests .or. &
1785  test_answer(v, 3, ppoly0_e(:,2), (/1.,3.,7./), 'Left lim PLM: right edges')
1786  remapping_unit_tests = remapping_unit_tests .or. &
1787  test_answer(v, 3, ppoly0_coefs(:,1), (/1.,1.,7./), 'Left lim PLM: P0')
1788  remapping_unit_tests = remapping_unit_tests .or. &
1789  test_answer(v, 3, ppoly0_coefs(:,2), (/0.,2.,0./), 'Left lim PLM: P1')
1790 
1791  call plm_reconstruction(3, (/1.,1.,1./), (/1.,6.,7./), ppoly0_e(1:3,:), &
1792  ppoly0_coefs(1:3,:), h_neglect )
1793  remapping_unit_tests = remapping_unit_tests .or. &
1794  test_answer(v, 3, ppoly0_e(:,1), (/1.,5.,7./), 'Right lim PLM: left edges')
1795  remapping_unit_tests = remapping_unit_tests .or. &
1796  test_answer(v, 3, ppoly0_e(:,2), (/1.,7.,7./), 'Right lim PLM: right edges')
1797  remapping_unit_tests = remapping_unit_tests .or. &
1798  test_answer(v, 3, ppoly0_coefs(:,1), (/1.,5.,7./), 'Right lim PLM: P0')
1799  remapping_unit_tests = remapping_unit_tests .or. &
1800  test_answer(v, 3, ppoly0_coefs(:,2), (/0.,2.,0./), 'Right lim PLM: P1')
1801 
1802  call plm_reconstruction(3, (/1.,2.,3./), (/1.,4.,9./), ppoly0_e(1:3,:), &
1803  ppoly0_coefs(1:3,:), h_neglect )
1804  remapping_unit_tests = remapping_unit_tests .or. &
1805  test_answer(v, 3, ppoly0_e(:,1), (/1.,2.,9./), 'Non-uniform line PLM: left edges')
1806  remapping_unit_tests = remapping_unit_tests .or. &
1807  test_answer(v, 3, ppoly0_e(:,2), (/1.,6.,9./), 'Non-uniform line PLM: right edges')
1808  remapping_unit_tests = remapping_unit_tests .or. &
1809  test_answer(v, 3, ppoly0_coefs(:,1), (/1.,2.,9./), 'Non-uniform line PLM: P0')
1810  remapping_unit_tests = remapping_unit_tests .or. &
1811  test_answer(v, 3, ppoly0_coefs(:,2), (/0.,4.,0./), 'Non-uniform line PLM: P1')
1812 
1813  call edge_values_explicit_h4( 5, (/1.,1.,1.,1.,1./), (/1.,3.,5.,7.,9./), ppoly0_e, &
1814  h_neglect=1e-10, answers_2018=answers_2018 )
1815  ! The next two tests currently fail due to roundoff, but pass when given a reasonable tolerance.
1816  thistest = test_answer(v, 5, ppoly0_e(:,1), (/0.,2.,4.,6.,8./), 'Line H4: left edges', tol=8.0e-15)
1817  remapping_unit_tests = remapping_unit_tests .or. thistest
1818  thistest = test_answer(v, 5, ppoly0_e(:,2), (/2.,4.,6.,8.,10./), 'Line H4: right edges', tol=1.0e-14)
1819  remapping_unit_tests = remapping_unit_tests .or. thistest
1820  ppoly0_e(:,1) = (/0.,2.,4.,6.,8./)
1821  ppoly0_e(:,2) = (/2.,4.,6.,8.,10./)
1822  call ppm_reconstruction(5, (/1.,1.,1.,1.,1./), (/1.,3.,5.,7.,9./), ppoly0_e(1:5,:), &
1823  ppoly0_coefs(1:5,:), h_neglect, answers_2018=answers_2018 )
1824  remapping_unit_tests = remapping_unit_tests .or. &
1825  test_answer(v, 5, ppoly0_coefs(:,1), (/1.,2.,4.,6.,9./), 'Line PPM: P0')
1826  remapping_unit_tests = remapping_unit_tests .or. &
1827  test_answer(v, 5, ppoly0_coefs(:,2), (/0.,2.,2.,2.,0./), 'Line PPM: P1')
1828  remapping_unit_tests = remapping_unit_tests .or. &
1829  test_answer(v, 5, ppoly0_coefs(:,3), (/0.,0.,0.,0.,0./), 'Line PPM: P2')
1830 
1831  call edge_values_explicit_h4( 5, (/1.,1.,1.,1.,1./), (/1.,1.,7.,19.,37./), ppoly0_e, &
1832  h_neglect=1e-10, answers_2018=answers_2018 )
1833  ! The next two tests are now passing when answers_2018 = .false., but otherwise only work to roundoff.
1834  thistest = test_answer(v, 5, ppoly0_e(:,1), (/3.,0.,3.,12.,27./), 'Parabola H4: left edges', tol=2.7e-14)
1835  remapping_unit_tests = remapping_unit_tests .or. thistest
1836  thistest = test_answer(v, 5, ppoly0_e(:,2), (/0.,3.,12.,27.,48./), 'Parabola H4: right edges', tol=4.8e-14)
1837  remapping_unit_tests = remapping_unit_tests .or. thistest
1838  ppoly0_e(:,1) = (/0.,0.,3.,12.,27./)
1839  ppoly0_e(:,2) = (/0.,3.,12.,27.,48./)
1840  call ppm_reconstruction(5, (/1.,1.,1.,1.,1./), (/0.,1.,7.,19.,37./), ppoly0_e(1:5,:), &
1841  ppoly0_coefs(1:5,:), h_neglect, answers_2018=answers_2018 )
1842  remapping_unit_tests = remapping_unit_tests .or. &
1843  test_answer(v, 5, ppoly0_e(:,1), (/0.,0.,3.,12.,37./), 'Parabola PPM: left edges')
1844  remapping_unit_tests = remapping_unit_tests .or. &
1845  test_answer(v, 5, ppoly0_e(:,2), (/0.,3.,12.,27.,37./), 'Parabola PPM: right edges')
1846  remapping_unit_tests = remapping_unit_tests .or. &
1847  test_answer(v, 5, ppoly0_coefs(:,1), (/0.,0.,3.,12.,37./), 'Parabola PPM: P0')
1848  remapping_unit_tests = remapping_unit_tests .or. &
1849  test_answer(v, 5, ppoly0_coefs(:,2), (/0.,0.,6.,12.,0./), 'Parabola PPM: P1')
1850  remapping_unit_tests = remapping_unit_tests .or. &
1851  test_answer(v, 5, ppoly0_coefs(:,3), (/0.,3.,3.,3.,0./), 'Parabola PPM: P2')
1852 
1853  ppoly0_e(:,1) = (/0.,0.,6.,10.,15./)
1854  ppoly0_e(:,2) = (/0.,6.,12.,17.,15./)
1855  call ppm_reconstruction(5, (/1.,1.,1.,1.,1./), (/0.,5.,7.,16.,15./), ppoly0_e(1:5,:), &
1856  ppoly0_coefs(1:5,:), h_neglect, answers_2018=answers_2018 )
1857  remapping_unit_tests = remapping_unit_tests .or. &
1858  test_answer(v, 5, ppoly0_e(:,1), (/0.,3.,6.,16.,15./), 'Limits PPM: left edges')
1859  remapping_unit_tests = remapping_unit_tests .or. &
1860  test_answer(v, 5, ppoly0_e(:,2), (/0.,6.,9.,16.,15./), 'Limits PPM: right edges')
1861  remapping_unit_tests = remapping_unit_tests .or. &
1862  test_answer(v, 5, ppoly0_coefs(:,1), (/0.,3.,6.,16.,15./), 'Limits PPM: P0')
1863  remapping_unit_tests = remapping_unit_tests .or. &
1864  test_answer(v, 5, ppoly0_coefs(:,2), (/0.,6.,0.,0.,0./), 'Limits PPM: P1')
1865  remapping_unit_tests = remapping_unit_tests .or. &
1866  test_answer(v, 5, ppoly0_coefs(:,3), (/0.,-3.,3.,0.,0./), 'Limits PPM: P2')
1867 
1868  call plm_reconstruction(4, (/0.,1.,1.,0./), (/5.,4.,2.,1./), ppoly0_e(1:4,:), &
1869  ppoly0_coefs(1:4,:), h_neglect )
1870  remapping_unit_tests = remapping_unit_tests .or. &
1871  test_answer(v, 4, ppoly0_e(1:4,1), (/5.,5.,3.,1./), 'PPM: left edges h=0110')
1872  remapping_unit_tests = remapping_unit_tests .or. &
1873  test_answer(v, 4, ppoly0_e(1:4,2), (/5.,3.,1.,1./), 'PPM: right edges h=0110')
1874  call remap_via_sub_cells( 4, (/0.,1.,1.,0./), (/5.,4.,2.,1./), ppoly0_e(1:4,:), &
1875  ppoly0_coefs(1:4,:), &
1876  2, (/1.,1./), integration_plm, .false., u2, err )
1877  remapping_unit_tests = remapping_unit_tests .or. &
1878  test_answer(v, 2, u2, (/4.,2./), 'PLM: remapped h=0110->h=11')
1879 
1880  deallocate(ppoly0_e, ppoly0_s, ppoly0_coefs)
1881 
1882  if (.not. remapping_unit_tests) write(*,*) 'Pass'
1883 
1884 end function remapping_unit_tests
1885 
1886 !> Returns true if any cell of u and u_true are not identical. Returns false otherwise.
1887 logical function test_answer(verbose, n, u, u_true, label, tol)
1888  logical, intent(in) :: verbose !< If true, write results to stdout
1889  integer, intent(in) :: n !< Number of cells in u
1890  real, dimension(n), intent(in) :: u !< Values to test
1891  real, dimension(n), intent(in) :: u_true !< Values to test against (correct answer)
1892  character(len=*), intent(in) :: label !< Message
1893  real, optional, intent(in) :: tol !< The tolerance for differences between u and u_true
1894  ! Local variables
1895  real :: tolerance ! The tolerance for differences between u and u_true
1896  integer :: k
1897 
1898  tolerance = 0.0 ; if (present(tol)) tolerance = tol
1899  test_answer = .false.
1900  do k = 1, n
1901  if (abs(u(k) - u_true(k)) > tolerance) test_answer = .true.
1902  enddo
1903  if (test_answer .or. verbose) then
1904  write(stdout,'(a4,2a24,x,a)') 'k','Calculated value','Correct value',label
1905  do k = 1, n
1906  if (abs(u(k) - u_true(k)) > tolerance) then
1907  write(stdout,'(i4,1p2e24.16,a,1pe24.16,a)') k,u(k),u_true(k),' err=',u(k)-u_true(k),' < wrong'
1908  write(stderr,'(i4,1p2e24.16,a,1pe24.16,a)') k,u(k),u_true(k),' err=',u(k)-u_true(k),' < wrong'
1909  else
1910  write(stdout,'(i4,1p2e24.16)') k,u(k),u_true(k)
1911  endif
1912  enddo
1913  endif
1914 
1915 end function test_answer
1916 
1917 !> Convenience function for printing grid to screen
1918 subroutine dumpgrid(n,h,x,u)
1919  integer, intent(in) :: n !< Number of cells
1920  real, dimension(:), intent(in) :: h !< Cell thickness
1921  real, dimension(:), intent(in) :: x !< Interface delta
1922  real, dimension(:), intent(in) :: u !< Cell average values
1923  integer :: i
1924  write(stdout,'("i=",20i10)') (i,i=1,n+1)
1925  write(stdout,'("x=",20es10.2)') (x(i),i=1,n+1)
1926  write(stdout,'("i=",5x,20i10)') (i,i=1,n)
1927  write(stdout,'("h=",5x,20es10.2)') (h(i),i=1,n)
1928  write(stdout,'("u=",5x,20es10.2)') (u(i),i=1,n)
1929 end subroutine dumpgrid
1930 
1931 end module mom_remapping
pcm_functions
Piecewise constant reconstruction functions.
Definition: PCM_functions.F90:2
ppm_functions
Provides functions used with the Piecewise-Parabolic-Method in the vertical ALE algorithm.
Definition: PPM_functions.F90:2
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_remapping::remapping_cs
Container for remapping parameters.
Definition: MOM_remapping.F90:24
regrid_edge_values
Edge value estimation for high-order resconstruction.
Definition: regrid_edge_values.F90:2
mom_remapping
Provides column-wise vertical remapping functions.
Definition: MOM_remapping.F90:2
pqm_functions
Piecewise quartic reconstruction functions.
Definition: PQM_functions.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
plm_functions
Piecewise linear reconstruction functions.
Definition: PLM_functions.F90:2