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