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