MOM6
regrid_interp.F90
1 !> Vertical interpolation for regridding
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, fatal
7 use mom_string_functions, only : uppercase
8 
9 use regrid_edge_values, only : edge_values_explicit_h2, edge_values_explicit_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 
13 use plm_functions, only : plm_reconstruction, plm_boundary_extrapolation
14 use ppm_functions, only : ppm_reconstruction, ppm_boundary_extrapolation
15 use pqm_functions, only : pqm_reconstruction, pqm_boundary_extrapolation_v1
16 
17 use p1m_functions, only : p1m_interpolation, p1m_boundary_extrapolation
18 use p3m_functions, only : p3m_interpolation, p3m_boundary_extrapolation
19 
20 implicit none ; private
21 
22 !> Control structure for regrid_interp module
23 type, public :: interp_cs_type ; private
24 
25  !> The following parameter is only relevant when used with the target
26  !! interface densities regridding scheme. It indicates which interpolation
27  !! to use to determine the grid.
28  integer :: interpolation_scheme = -1
29 
30  !> Indicate whether high-order boundary extrapolation should be used within
31  !! boundary cells
32  logical :: boundary_extrapolation
33 end type interp_cs_type
34 
35 public regridding_set_ppolys, interpolate_grid
36 public build_and_interpolate_grid
37 public set_interp_scheme, set_interp_extrap
38 
39 ! List of interpolation schemes
40 integer, parameter :: interpolation_p1m_h2 = 0 !< O(h^2)
41 integer, parameter :: interpolation_p1m_h4 = 1 !< O(h^2)
42 integer, parameter :: interpolation_p1m_ih4 = 2 !< O(h^2)
43 integer, parameter :: interpolation_plm = 3 !< O(h^2)
44 integer, parameter :: interpolation_ppm_h4 = 4 !< O(h^3)
45 integer, parameter :: interpolation_ppm_ih4 = 5 !< O(h^3)
46 integer, parameter :: interpolation_p3m_ih4ih3 = 6 !< O(h^4)
47 integer, parameter :: interpolation_p3m_ih6ih5 = 7 !< O(h^4)
48 integer, parameter :: interpolation_pqm_ih4ih3 = 8 !< O(h^4)
49 integer, parameter :: interpolation_pqm_ih6ih5 = 9 !< O(h^5)
50 
51 !>@{ Interpolant degrees
52 integer, parameter :: degree_1 = 1, degree_2 = 2, degree_3 = 3, degree_4 = 4
53 integer, public, parameter :: degree_max = 5
54 !!@}
55 
56 !> When the N-R algorithm produces an estimate that lies outside [0,1], the
57 !! estimate is set to be equal to the boundary location, 0 or 1, plus or minus
58 !! an offset, respectively, when the derivative is zero at the boundary.
59 real, public, parameter :: nr_offset = 1e-6
60 !> Maximum number of Newton-Raphson iterations. Newton-Raphson iterations are
61 !! used to build the new grid by finding the coordinates associated with
62 !! target densities and interpolations of degree larger than 1.
63 integer, public, parameter :: nr_iterations = 8
64 !> Tolerance for Newton-Raphson iterations (stop when increment falls below this)
65 real, public, parameter :: nr_tolerance = 1e-12
66 
67 contains
68 
69 !> Builds an interpolated profile for the densities within each grid cell.
70 !!
71 !! It may happen that, given a high-order interpolator, the number of
72 !! available layers is insufficient (e.g., there are two available layers for
73 !! a third-order PPM ih4 scheme). In these cases, we resort to the simplest
74 !! continuous linear scheme (P1M h2).
75 subroutine regridding_set_ppolys(CS, densities, n0, h0, ppoly0_E, ppoly0_S, &
76  ppoly0_coefs, degree, h_neglect, h_neglect_edge)
77  type(interp_cs_type),intent(in) :: cs !< Interpolation control structure
78  real, dimension(:), intent(in) :: densities !< Actual cell densities
79  integer, intent(in) :: n0 !< Number of cells on source grid
80  real, dimension(:), intent(in) :: h0 !< cell widths on source grid
81  real, dimension(:,:),intent(inout) :: ppoly0_e !< Edge value of polynomial
82  real, dimension(:,:),intent(inout) :: ppoly0_s !< Edge slope of polynomial
83  real, dimension(:,:),intent(inout) :: ppoly0_coefs !< Coefficients of polynomial
84  integer, intent(inout) :: degree !< The degree of the polynomials
85  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
86  !! purpose of cell reconstructions
87  !! in the same units as h0.
88  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
89  !! for the purpose of edge value calculations
90  !! in the same units as h0.
91  ! Local variables
92  logical :: extrapolate
93 
94  ! Reset piecewise polynomials
95  ppoly0_e(:,:) = 0.0
96  ppoly0_s(:,:) = 0.0
97  ppoly0_coefs(:,:) = 0.0
98 
99  extrapolate = cs%boundary_extrapolation
100 
101  ! Compute the interpolated profile of the density field and build grid
102  select case (cs%interpolation_scheme)
103 
104  case ( interpolation_p1m_h2 )
105  degree = degree_1
106  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
107  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
108  if (extrapolate) then
109  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
110  endif
111 
112  case ( interpolation_p1m_h4 )
113  degree = degree_1
114  if ( n0 >= 4 ) then
115  call edge_values_explicit_h4( n0, h0, densities, ppoly0_e, h_neglect_edge )
116  else
117  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
118  endif
119  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
120  if (extrapolate) then
121  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
122  endif
123 
124  case ( interpolation_p1m_ih4 )
125  degree = degree_1
126  if ( n0 >= 4 ) then
127  call edge_values_implicit_h4( n0, h0, densities, ppoly0_e, h_neglect_edge )
128  else
129  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
130  endif
131  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
132  if (extrapolate) then
133  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
134  endif
135 
136  case ( interpolation_plm )
137  degree = degree_1
138  call plm_reconstruction( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
139  if (extrapolate) then
140  call plm_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
141  endif
142 
143  case ( interpolation_ppm_h4 )
144  if ( n0 >= 4 ) then
145  degree = degree_2
146  call edge_values_explicit_h4( n0, h0, densities, ppoly0_e, h_neglect_edge )
147  call ppm_reconstruction( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
148  if (extrapolate) then
149  call ppm_boundary_extrapolation( n0, h0, densities, ppoly0_e, &
150  ppoly0_coefs, h_neglect )
151  endif
152  else
153  degree = degree_1
154  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
155  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
156  if (extrapolate) then
157  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
158  endif
159  endif
160 
161  case ( interpolation_ppm_ih4 )
162  if ( n0 >= 4 ) then
163  degree = degree_2
164  call edge_values_implicit_h4( n0, h0, densities, ppoly0_e, h_neglect_edge )
165  call ppm_reconstruction( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
166  if (extrapolate) then
167  call ppm_boundary_extrapolation( n0, h0, densities, ppoly0_e, &
168  ppoly0_coefs, h_neglect )
169  endif
170  else
171  degree = degree_1
172  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
173  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
174  if (extrapolate) then
175  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
176  endif
177  endif
178 
179  case ( interpolation_p3m_ih4ih3 )
180  if ( n0 >= 4 ) then
181  degree = degree_3
182  call edge_values_implicit_h4( n0, h0, densities, ppoly0_e, h_neglect_edge )
183  call edge_slopes_implicit_h3( n0, h0, densities, ppoly0_s, h_neglect )
184  call p3m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_s, &
185  ppoly0_coefs, h_neglect )
186  if (extrapolate) then
187  call p3m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_s, &
188  ppoly0_coefs, h_neglect, h_neglect_edge )
189  endif
190  else
191  degree = degree_1
192  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
193  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
194  if (extrapolate) then
195  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
196  endif
197  endif
198 
199  case ( interpolation_p3m_ih6ih5 )
200  if ( n0 >= 6 ) then
201  degree = degree_3
202  call edge_values_implicit_h6( n0, h0, densities, ppoly0_e, h_neglect_edge )
203  call edge_slopes_implicit_h5( n0, h0, densities, ppoly0_s, h_neglect )
204  call p3m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_s, &
205  ppoly0_coefs, h_neglect )
206  if (extrapolate) then
207  call p3m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_s, &
208  ppoly0_coefs, h_neglect, h_neglect_edge )
209  endif
210  else
211  degree = degree_1
212  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
213  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
214  if (extrapolate) then
215  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
216  endif
217  endif
218 
219  case ( interpolation_pqm_ih4ih3 )
220  if ( n0 >= 4 ) then
221  degree = degree_4
222  call edge_values_implicit_h4( n0, h0, densities, ppoly0_e, h_neglect_edge )
223  call edge_slopes_implicit_h3( n0, h0, densities, ppoly0_s, h_neglect )
224  call pqm_reconstruction( n0, h0, densities, ppoly0_e, ppoly0_s, &
225  ppoly0_coefs, h_neglect )
226  if (extrapolate) then
227  call pqm_boundary_extrapolation_v1( n0, h0, densities, ppoly0_e, ppoly0_s, &
228  ppoly0_coefs, h_neglect )
229  endif
230  else
231  degree = degree_1
232  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
233  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
234  if (extrapolate) then
235  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
236  endif
237  endif
238 
239  case ( interpolation_pqm_ih6ih5 )
240  if ( n0 >= 6 ) then
241  degree = degree_4
242  call edge_values_implicit_h6( n0, h0, densities, ppoly0_e, h_neglect_edge )
243  call edge_slopes_implicit_h5( n0, h0, densities, ppoly0_s, h_neglect )
244  call pqm_reconstruction( n0, h0, densities, ppoly0_e, ppoly0_s, &
245  ppoly0_coefs, h_neglect )
246  if (extrapolate) then
247  call pqm_boundary_extrapolation_v1( n0, h0, densities, ppoly0_e, ppoly0_s, &
248  ppoly0_coefs, h_neglect )
249  endif
250  else
251  degree = degree_1
252  call edge_values_explicit_h2( n0, h0, densities, ppoly0_e, h_neglect_edge )
253  call p1m_interpolation( n0, h0, densities, ppoly0_e, ppoly0_coefs, h_neglect )
254  if (extrapolate) then
255  call p1m_boundary_extrapolation( n0, h0, densities, ppoly0_e, ppoly0_coefs )
256  endif
257  endif
258  end select
259 end subroutine regridding_set_ppolys
260 
261 !> Given target values (e.g., density), build new grid based on polynomial
262 !!
263 !! Given the grid 'grid0' and the piecewise polynomial interpolant
264 !! 'ppoly0' (possibly discontinuous), the coordinates of the new grid 'grid1'
265 !! are determined by finding the corresponding target interface densities.
266 subroutine interpolate_grid( n0, h0, x0, ppoly0_E, ppoly0_coefs, &
267  target_values, degree, n1, h1, x1 )
268  integer, intent(in) :: n0 !< Number of points on source grid
269  real, dimension(:), intent(in) :: h0 !< Thicknesses of source grid cells
270  real, dimension(:), intent(in) :: x0 !< Source interface positions
271  real, dimension(:,:), intent(in) :: ppoly0_e !< Edge values of interpolating polynomials
272  real, dimension(:,:), intent(in) :: ppoly0_coefs !< Coefficients of interpolating polynomials
273  real, dimension(:), intent(in) :: target_values !< Target values of interfaces
274  integer, intent(in) :: degree !< Degree of interpolating polynomials
275  integer, intent(in) :: n1 !< Number of points on target grid
276  real, dimension(:), intent(inout) :: h1 !< Thicknesses of target grid cells
277  real, dimension(:), intent(inout) :: x1 !< Target interface positions
278  ! Local variables
279  integer :: k ! loop index
280  real :: t ! current interface target density
281 
282  ! Make sure boundary coordinates of new grid coincide with boundary
283  ! coordinates of previous grid
284  x1(1) = x0(1)
285  x1(n1+1) = x0(n0+1)
286 
287  ! Find coordinates for interior target values
288  do k = 2,n1
289  t = target_values(k)
290  x1(k) = get_polynomial_coordinate( n0, h0, x0, ppoly0_e, ppoly0_coefs, t, degree )
291  h1(k-1) = x1(k) - x1(k-1)
292  enddo
293  h1(n1) = x1(n1+1) - x1(n1)
294 
295 end subroutine interpolate_grid
296 
297 !> Build a grid by interpolating for target values
298 subroutine build_and_interpolate_grid(CS, densities, n0, h0, x0, target_values, &
299  n1, h1, x1, h_neglect, h_neglect_edge)
300  type(interp_cs_type), intent(in) :: cs !< A control structure for regrid_interp
301  real, dimension(:), intent(in) :: densities !< Input cell densities [kg m-3]
302  real, dimension(:), intent(in) :: target_values !< Target values of interfaces
303  integer, intent(in) :: n0 !< The number of points on the input grid
304  real, dimension(:), intent(in) :: h0 !< Initial cell widths
305  real, dimension(:), intent(in) :: x0 !< Source interface positions
306  integer, intent(in) :: n1 !< The number of points on the output grid
307  real, dimension(:), intent(inout) :: h1 !< Output cell widths
308  real, dimension(:), intent(inout) :: x1 !< Target interface positions
309  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
310  !! purpose of cell reconstructions
311  !! in the same units as h0.
312  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
313  !! for the purpose of edge value calculations
314  !! in the same units as h0.
315 
316  real, dimension(n0,2) :: ppoly0_e, ppoly0_s
317  real, dimension(n0,DEGREE_MAX+1) :: ppoly0_c
318  integer :: degree
319 
320  call regridding_set_ppolys(cs, densities, n0, h0, ppoly0_e, ppoly0_s, ppoly0_c, &
321  degree, h_neglect, h_neglect_edge)
322  call interpolate_grid(n0, h0, x0, ppoly0_e, ppoly0_c, target_values, degree, &
323  n1, h1, x1)
324 end subroutine build_and_interpolate_grid
325 
326 !> Given a target value, find corresponding coordinate for given polynomial
327 !!
328 !! Here, 'ppoly' is assumed to be a piecewise discontinuous polynomial of degree
329 !! 'degree' throughout the domain defined by 'grid'. A target value is given
330 !! and we need to determine the corresponding grid coordinate to define the
331 !! new grid.
332 !!
333 !! If the target value is out of range, the grid coordinate is simply set to
334 !! be equal to one of the boundary coordinates, which results in vanished layers
335 !! near the boundaries.
336 !!
337 !! IT IS ASSUMED THAT THE PIECEWISE POLYNOMIAL IS MONOTONICALLY INCREASING.
338 !! IF THIS IS NOT THE CASE, THE NEW GRID MAY BE ILL-DEFINED.
339 !!
340 !! It is assumed that the number of cells defining 'grid' and 'ppoly' are the
341 !! same.
342 function get_polynomial_coordinate( N, h, x_g, ppoly_E, ppoly_coefs, &
343  target_value, degree ) result ( x_tgt )
344  ! Arguments
345  integer, intent(in) :: n !< Number of grid cells
346  real, dimension(:), intent(in) :: h !< Grid cell thicknesses
347  real, dimension(:), intent(in) :: x_g !< Grid interface locations
348  real, dimension(:,:), intent(in) :: ppoly_e !< Edge values of interpolating polynomials
349  real, dimension(:,:), intent(in) :: ppoly_coefs !< Coefficients of interpolating polynomials
350  real, intent(in) :: target_value !< Target value to find position for
351  integer, intent(in) :: degree !< Degree of the interpolating polynomials
352  real :: x_tgt !< The position of x_g at which target_value is found.
353  ! Local variables
354  integer :: i, k ! loop indices
355  integer :: k_found ! index of target cell
356  integer :: iter
357  real :: xi0 ! normalized target coordinate
358  real, dimension(DEGREE_MAX) :: a ! polynomial coefficients
359  real :: numerator
360  real :: denominator
361  real :: delta ! Newton-Raphson increment
362  real :: x ! global target coordinate
363  real :: eps ! offset used to get away from
364  ! boundaries
365  real :: grad ! gradient during N-R iterations
366 
367  eps = nr_offset
368  k_found = -1
369 
370  ! If the target value is outside the range of all values, we
371  ! force the target coordinate to be equal to the lowest or
372  ! largest value, depending on which bound is overtaken
373  if ( target_value <= ppoly_e(1,1) ) then
374  x_tgt = x_g(1)
375  return ! return because there is no need to look further
376  endif
377 
378  ! Since discontinuous edge values are allowed, we check whether the target
379  ! value lies between two discontinuous edge values at interior interfaces
380  do k = 2,n
381  if ( ( target_value >= ppoly_e(k-1,2) ) .AND. &
382  ( target_value <= ppoly_e(k,1) ) ) then
383  x_tgt = x_g(k)
384  return ! return because there is no need to look further
385  exit
386  endif
387  enddo
388 
389  ! If the target value is outside the range of all values, we
390  ! force the target coordinate to be equal to the lowest or
391  ! largest value, depending on which bound is overtaken
392  if ( target_value >= ppoly_e(n,2) ) then
393  x_tgt = x_g(n+1)
394  return ! return because there is no need to look further
395  endif
396 
397  ! At this point, we know that the target value is bounded and does not
398  ! lie between discontinuous, monotonic edge values. Therefore,
399  ! there is a unique solution. We loop on all cells and find which one
400  ! contains the target value. The variable k_found holds the index value
401  ! of the cell where the taregt value lies.
402  do k = 1,n
403  if ( ( target_value > ppoly_e(k,1) ) .AND. &
404  ( target_value < ppoly_e(k,2) ) ) then
405  k_found = k
406  exit
407  endif
408  enddo
409 
410  ! At this point, 'k_found' should be strictly positive. If not, this is
411  ! a major failure because it means we could not find any target cell
412  ! despite the fact that the target value lies between the extremes. It
413  ! means there is a major problem with the interpolant. This needs to be
414  ! reported.
415  if ( k_found == -1 ) then
416  write(*,*) target_value, ppoly_e(1,1), ppoly_e(n,2)
417  write(*,*) 'Could not find target coordinate in ' //&
418  '"get_polynomial_coordinate". This is caused by an '//&
419  'inconsistent interpolant (perhaps not monotonically '//&
420  'increasing)'
421  call mom_error( fatal, 'Aborting execution' )
422  endif
423 
424  ! Reset all polynomial coefficients to 0 and copy those pertaining to
425  ! the found cell
426  a(:) = 0.0
427  do i = 1,degree+1
428  a(i) = ppoly_coefs(k_found,i)
429  enddo
430 
431  ! Guess value to start Newton-Raphson iterations (middle of cell)
432  xi0 = 0.5
433  iter = 1
434  delta = 1e10
435 
436  ! Newton-Raphson iterations
437  do
438  ! break if converged or too many iterations taken
439  if ( ( iter > nr_iterations ) .OR. &
440  ( abs(delta) < nr_tolerance ) ) then
441  exit
442  endif
443 
444  numerator = a(1) + a(2)*xi0 + a(3)*xi0*xi0 + a(4)*xi0*xi0*xi0 + &
445  a(5)*xi0*xi0*xi0*xi0 - target_value
446 
447  denominator = a(2) + 2*a(3)*xi0 + 3*a(4)*xi0*xi0 + 4*a(5)*xi0*xi0*xi0
448 
449  delta = -numerator / denominator
450 
451  xi0 = xi0 + delta
452 
453  ! Check whether new estimate is out of bounds. If the new estimate is
454  ! indeed out of bounds, we manually set it to be equal to the overtaken
455  ! bound with a small offset towards the interior when the gradient of
456  ! the function at the boundary is zero (in which case, the Newton-Raphson
457  ! algorithm does not converge).
458  if ( xi0 < 0.0 ) then
459  xi0 = 0.0
460  grad = a(2)
461  if ( grad == 0.0 ) xi0 = xi0 + eps
462  endif
463 
464  if ( xi0 > 1.0 ) then
465  xi0 = 1.0
466  grad = a(2) + 2*a(3) + 3*a(4) + 4*a(5)
467  if ( grad == 0.0 ) xi0 = xi0 - eps
468  endif
469 
470  iter = iter + 1
471  enddo ! end Newton-Raphson iterations
472 
473  x_tgt = x_g(k_found) + xi0 * h(k_found)
474 end function get_polynomial_coordinate
475 
476 !> Numeric value of interpolation_scheme corresponding to scheme name
477 integer function interpolation_scheme(interp_scheme)
478  character(len=*), intent(in) :: interp_scheme !< Name of the interpolation scheme
479  !! Valid values include "P1M_H2", "P1M_H4", "P1M_IH2", "PLM", "PPM_H4",
480  !! "PPM_IH4", "P3M_IH4IH3", "P3M_IH6IH5", "PQM_IH4IH3", and "PQM_IH6IH5"
481 
482  select case ( uppercase(trim(interp_scheme)) )
483  case ("P1M_H2"); interpolation_scheme = interpolation_p1m_h2
484  case ("P1M_H4"); interpolation_scheme = interpolation_p1m_h4
485  case ("P1M_IH2"); interpolation_scheme = interpolation_p1m_ih4
486  case ("PLM"); interpolation_scheme = interpolation_plm
487  case ("PPM_H4"); interpolation_scheme = interpolation_ppm_h4
488  case ("PPM_IH4"); interpolation_scheme = interpolation_ppm_ih4
489  case ("P3M_IH4IH3"); interpolation_scheme = interpolation_p3m_ih4ih3
490  case ("P3M_IH6IH5"); interpolation_scheme = interpolation_p3m_ih6ih5
491  case ("PQM_IH4IH3"); interpolation_scheme = interpolation_pqm_ih4ih3
492  case ("PQM_IH6IH5"); interpolation_scheme = interpolation_pqm_ih6ih5
493  case default ; call mom_error(fatal, "regrid_interp: "//&
494  "Unrecognized choice for INTERPOLATION_SCHEME ("//trim(interp_scheme)//").")
495  end select
496 end function interpolation_scheme
497 
498 !> Store the interpolation_scheme value in the interp_CS based on the input string.
499 subroutine set_interp_scheme(CS, interp_scheme)
500  type(interp_cs_type), intent(inout) :: cs !< A control structure for regrid_interp
501  character(len=*), intent(in) :: interp_scheme !< Name of the interpolation scheme
502  !! Valid values include "P1M_H2", "P1M_H4", "P1M_IH2", "PLM", "PPM_H4",
503  !! "PPM_IH4", "P3M_IH4IH3", "P3M_IH6IH5", "PQM_IH4IH3", and "PQM_IH6IH5"
504 
505  cs%interpolation_scheme = interpolation_scheme(interp_scheme)
506 end subroutine set_interp_scheme
507 
508 !> Store the boundary_extrapolation value in the interp_CS
509 subroutine set_interp_extrap(CS, extrap)
510  type(interp_cs_type), intent(inout) :: cs !< A control structure for regrid_interp
511  logical, intent(in) :: extrap !< Indicate whether high-order boundary
512  !! extrapolation should be used in boundary cells
513 
514  cs%boundary_extrapolation = extrap
515 end subroutine set_interp_extrap
516 
517 end module regrid_interp
ppm_functions
Provides functions used with the Piecewise-Parabolic-Method in the vertical ALE algorithm.
Definition: PPM_functions.F90:2
p3m_functions
Cubic interpolation functions.
Definition: P3M_functions.F90:2
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
regrid_edge_values
Edge value estimation for high-order resconstruction.
Definition: regrid_edge_values.F90:2
p1m_functions
Linear interpolation functions.
Definition: P1M_functions.F90:2
regrid_interp
Vertical interpolation for regridding.
Definition: regrid_interp.F90:2
pqm_functions
Piecewise quartic reconstruction functions.
Definition: PQM_functions.F90:2
regrid_interp::interp_cs_type
Control structure for regrid_interp module.
Definition: regrid_interp.F90:23
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