mom_horizontal_regridding module reference¶
Horizontal interpolation.
Functions/Subroutines¶
Write to the terminal some basic statistics about the k-th level of an array. |
|
Use ICE-9 algorithm to populate points (fill=1) with valid data (good=1). |
|
Extrapolate and interpolate from a file record. |
|
Extrapolate and interpolate using a FMS time interpolation handle. |
|
Create a 2d-mesh of grid coordinates from 1-d arrays. |
|
Fill grid edges for integer data. |
|
Fill grid edges for real data. |
|
Solve del2 (zi) = 0 using successive iterations with a 5 point stencil. |
Detailed Description¶
Horizontal interpolation.
Function/Subroutine Documentation¶
-
subroutine
mom_horizontal_regridding/
mystats
(array, missing, is, ie, js, je, k, mesg)¶ Write to the terminal some basic statistics about the k-th level of an array.
- Parameters
array :: [in] input array (ND)
missing :: [in] missing value (ND)
is :: Start index in i
ie :: End index in i
js :: Start index in j
je :: End index in j
k :: Level to calculate statistics for
mesg :: Label to use in message
- Called from
horiz_interp_and_extrap_tracer_fms_id
horiz_interp_and_extrap_tracer_record
mom_tracer_initialization_from_z::mom_initialize_tracer_from_z
-
subroutine
mom_horizontal_regridding/
fill_miss_2d
(aout, good, fill, prev, G, smooth, num_pass, relc, crit, debug, answers_2018)¶ Use ICE-9 algorithm to populate points (fill=1) with valid data (good=1). If no information is available, Then use a previous guess (prev). Optionally (smooth) blend the filled points to achieve a more desirable result.
- Parameters
g :: [inout] The ocean’s grid structure.
aout :: [inout] The array with missing values to fill
good :: [in] Valid data mask for incoming array
fill :: [in] Same shape array of points which need
prev :: [in] First guess where isolated holes exist.
smooth :: [in] If present and true, apply a number of Laplacian iterations to the interpolated data
num_pass :: [in] The maximum number of iterations
relc :: [in] A relaxation coefficient for Laplacian (ND)
crit :: [in] A minimal value for deltas between iterations.
debug :: [in] If true, write verbose debugging messages.
answers_2018 :: [in] If true, use expressions that give the same answers as the code did in late 2018. Otherwise add parentheses for rotational symmetry.
- Called from
horiz_interp_and_extrap_tracer_fms_id
horiz_interp_and_extrap_tracer_record
-
subroutine
mom_horizontal_regridding/
horiz_interp_and_extrap_tracer_record
(filename, varnam, conversion, recnum, G, tr_z, mask_z, z_in, z_edges_in, missing_value, reentrant_x, tripolar_n, homogenize, m_to_Z, answers_2018, ongrid)¶ Extrapolate and interpolate from a file record.
- Parameters
filename :: [in] Path to file containing tracer to be interpolated.
varnam :: [in] Name of tracer in filee.
conversion :: [in] Conversion factor for tracer.
recnum :: [in] Record number of tracer to be read.
g :: [inout] Grid object
tr_z :: pointer to allocatable tracer array on local model grid and input-file vertical levels.
mask_z :: pointer to allocatable tracer mask array on local model grid and input-file vertical levels.
z_in :: Cell grid values for input data.
z_edges_in :: Cell grid edge values for input data.
missing_value :: [out] The missing value in the returned array.
reentrant_x :: [in] If true, this grid is reentrant in the x-direction
tripolar_n :: [in] If true, this is a northern tripolar grid
homogenize :: [in] If present and true, horizontally homogenize data to produce perfectly “flat” initial conditions
m_to_z :: [in] A conversion factor from meters to the units of depth. If missing, GbathyT must be in m.
answers_2018 :: [in] If true, use expressions that give the same answers as the code did in late 2018. Otherwise add parentheses for rotational symmetry.
ongrid :: [in] If true, then data are assumed to have been interpolated to the model horizontal grid. In this case, only extrapolation is performed by this routine
- Call to
-
subroutine
mom_horizontal_regridding/
horiz_interp_and_extrap_tracer_fms_id
(fms_id, Time, conversion, G, tr_z, mask_z, z_in, z_edges_in, missing_value, reentrant_x, tripolar_n, homogenize, spongeOngrid, m_to_Z, answers_2018)¶ Extrapolate and interpolate using a FMS time interpolation handle.
- Parameters
fms_id :: [in] A unique id used by the FMS time interpolator
time :: [in] A FMS time type
conversion :: [in] Conversion factor for tracer.
g :: [inout] Grid object
tr_z :: pointer to allocatable tracer array on local model grid and native vertical levels.
mask_z :: pointer to allocatable tracer mask array on local model grid and native vertical levels.
z_in :: Cell grid values for input data.
z_edges_in :: Cell grid edge values for input data. (Intent out)
missing_value :: [out] The missing value in the returned array.
reentrant_x :: [in] If true, this grid is reentrant in the x-direction
tripolar_n :: [in] If true, this is a northern tripolar grid
homogenize :: [in] If present and true, horizontally homogenize data to produce perfectly “flat” initial conditions
spongeongrid :: [in] If present and true, the sponge data are on the model grid
m_to_z :: [in] A conversion factor from meters to the units of depth. If missing, GbathyT must be in m.
answers_2018 :: [in] If true, use expressions that give the same answers as the code did in late 2018. Otherwise add parentheses for rotational symmetry.
- Call to
-
subroutine
mom_horizontal_regridding/
meshgrid
(x, y, x_T, y_T)¶ Create a 2d-mesh of grid coordinates from 1-d arrays.
- Parameters
x :: [in] input 1-dimensional vector
y :: [in] input 1-dimensional vector
x_t :: [inout] output 2-dimensional array
y_t :: [inout] output 2-dimensional array
- Called from
horiz_interp_and_extrap_tracer_fms_id
horiz_interp_and_extrap_tracer_record
-
function
mom_horizontal_regridding/
fill_boundaries_int
(m, cyclic_x, tripolar_n) [integer]¶ Fill grid edges for integer data.
- Parameters
m :: [in] input array (ND)
cyclic_x :: [in] True if domain is zonally re-entrant
tripolar_n :: [in] True if domain has an Arctic fold
- Call to
-
function
mom_horizontal_regridding/
fill_boundaries_real
(m, cyclic_x, tripolar_n) [real]¶ Fill grid edges for real data.
- Parameters
m :: [in] input array (ND)
cyclic_x :: [in] True if domain is zonally re-entrant
tripolar_n :: [in] True if domain has an Arctic fold
- Called from
-
subroutine
mom_horizontal_regridding/
smooth_heights
(zi, fill, bad, sor, niter, cyclic_x, tripolar_n)¶ Solve del2 (zi) = 0 using successive iterations with a 5 point stencil. Only points fill==1 are modified. Except where bad==1, information propagates isotropically in index space. The resulting solution in each region is an approximation to del2(zi)=0 subject to boundary conditions along the valid points curve bounding this region.
- Parameters
zi :: [inout] input and output array (ND)
fill :: [in] same shape as zi, 1=fill
bad :: [in] same shape as zi, 1=bad data
sor :: [in] relaxation coefficient (ND)
niter :: [in] maximum number of iterations
cyclic_x :: [in] true if domain is zonally reentrant
tripolar_n :: [in] true if domain has an Arctic fold