mom_tracer_z_init module reference¶
Used to initialize tracers from a depth- (or z*-) space file.
Functions/Subroutines¶
This function initializes a tracer by reading a Z-space file, returning .true. |
|
Layer model routine for remapping tracers from pseudo-z coordinates into layers defined by target interface positions. |
|
This subroutine reads the vertical coordinate data for a field from a NetCDF file. |
|
Determines the layers bounded by interfaces e that overlap with the depth range between Z_top and Z_bot, and the fractional weights of each layer. |
|
This subroutine determines a limited slope for val to be advected with a piecewise limited scheme. |
|
This subroutine determines the potential temperature and salinity that is consistent with the target density using provided initial guess. |
Detailed Description¶
Used to initialize tracers from a depth- (or z*-) space file.
Function/Subroutine Documentation¶
-
function
mom_tracer_z_init/
tracer_z_init
(tr, h, filename, tr_name, G, GV, US, missing_val, land_val, scale) [logical]¶ This function initializes a tracer by reading a Z-space file, returning .true. if this appears to have been successful, and false otherwise.
- Return:
undefined :: A return code indicating if the initialization has been successful
- Parameters:
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure.
us :: [in] A dimensional unit scaling type
tr :: [out] The tracer to initialize [CU ~> conc]
h :: [in] Layer thicknesses [H ~> m or kg m-2] or other
filename :: [in] The name of the file to read from
tr_name :: [in] The name of the tracer in the file
missing_val :: [in] The missing value for the tracer [CU ~> conc]
land_val :: [in] A value to use to fill in land points [CU ~> conc]
scale :: [in] A factor by which to scale the output tracers from the their units in the file [CU conc-1 ~> 1]
- Call to:
find_limited_slope
find_overlap
mom_error_handler::mom_error
read_z_edges
- Called from:
mom_cfc_cap::init_tracer_cfc
mom_ocmip2_cfc::init_tracer_cfc
ideal_age_example::initialize_ideal_age_tracer
mom_generic_tracer::initialize_mom_generic_tracer
oil_tracer::initialize_oil_tracer
-
subroutine
mom_tracer_z_init/
tracer_z_init_array
(tr_in, z_edges, nk_data, e, land_fill, G, nlay, nlevs, eps_z, tr, scale)¶ Layer model routine for remapping tracers from pseudo-z coordinates into layers defined by target interface positions.
- Parameters:
g :: [in] The ocean’s grid structure
nk_data :: [in] The number of levels in the input data
tr_in :: [in] The z-space array of tracer concentrations
z_edges :: [in] The depths of the cell edges in the input z* data [Z ~> m] or [m]
nlay :: [in] The number of vertical layers in the target grid
e :: [in] The depths of the target layer interfaces [Z ~> m] or [m]
land_fill :: [in] fill in data over land [B]
nlevs :: [in] The number of input levels with valid data
eps_z :: [in] A negligibly thin layer thickness [Z ~> m].
tr :: [out] tracers in model space [B]
scale :: [in] A factor by which to scale the output tracers from the input tracers [B A-1 ~> 1]
- Call to:
-
subroutine
mom_tracer_z_init/
read_z_edges
(filename, tr_name, z_edges, nz_out, has_edges, use_missing, missing, scale, missing_scale)¶ This subroutine reads the vertical coordinate data for a field from a NetCDF file. It also might read the missing value attribute for that same field.
- Parameters:
filename :: [in] The name of the file to read from.
tr_name :: [in] The name of the tracer in the file.
z_edges :: [out] The depths of the vertical edges of the tracer array [Z ~> m]
nz_out :: [out] The number of vertical layers in the tracer array
has_edges :: [out] If true the values in z_edges are the edges of the tracer cells, otherwise they are the cell centers
use_missing :: [inout] If false on input, see whether the tracer has a missing value, and if so return true
missing :: [inout] The missing value, if one has been found [CU ~> conc]
scale :: [in] A scaling factor for z_edges into new units [Z m-1 ~> 1]
missing_scale :: [in] A scaling factor to use to convert the tracers and their missing value from the units in the file into their internal units [CU conc-1 ~> 1]
- Call to:
mom_io::close_file_to_read
mom_error_handler::mom_error
mom_io::open_file_to_read
- Called from:
-
subroutine
mom_tracer_z_init/
find_overlap
(e, Z_top, Z_bot, k_max, k_start, k_top, k_bot, wt, z1, z2)¶ Determines the layers bounded by interfaces e that overlap with the depth range between Z_top and Z_bot, and the fractional weights of each layer. It also calculates the normalized relative depths of the range of each layer that overlaps that depth range.
- Parameters:
e :: [in] Column interface heights, [Z ~> m] or other units.
z_top :: [in] Top of range being mapped to, in the units of e [Z ~> m].
z_bot :: [in] Bottom of range being mapped to, in the units of e [Z ~> m].
k_max :: [in] Number of valid layers.
k_start :: [in] Layer at which to start searching.
k_top :: [out] Indices of top layers that overlap with the depth range.
k_bot :: [out] Indices of bottom layers that overlap with the depth range.
wt :: [out] Relative weights of each layer from k_top to k_bot [nondim].
z1 :: [out] Depth of the top limits of the part of a layer that contributes to a depth level, relative to the cell center and normalized by the cell thickness [nondim]. Note that -1/2 <= z1 < z2 <= 1/2.
z2 :: [out] Depths of the bottom limit of the part of a layer that contributes to a depth level, relative to the cell center and normalized by the cell thickness [nondim]. Note that -1/2 <= z1 < z2 <= 1/2.
- Called from:
-
function
mom_tracer_z_init/
find_limited_slope
(val, e, k) [real]¶ This subroutine determines a limited slope for val to be advected with a piecewise limited scheme.
- Parameters:
val :: [in] A column of the values that are being interpolated, in arbitrary units [A]
e :: [in] A column’s interface heights [Z ~> m] or other units.
k :: [in] The layer whose slope is being determined.
- Return:
undefined :: The normalized slope in the intracell distribution of val [A]
- Called from:
-
subroutine
mom_tracer_z_init/
determine_temperature
(temp, salt, R_tgt, EOS, p_ref, niter, k_start, G, GV, US, PF, just_read)¶ This subroutine determines the potential temperature and salinity that is consistent with the target density using provided initial guess.
- Parameters:
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure.
temp :: [inout] potential temperature [C ~> degC]
salt :: [inout] salinity [S ~> ppt]
r_tgt :: [in] desired potential density [R ~> kg m-3].
eos :: [in] seawater equation of state control structure
p_ref :: [in] reference pressure [R L2 T-2 ~> Pa].
niter :: [in] maximum number of iterations
k_start :: [in] starting index (i.e. below the buffer layer)
us :: [in] A dimensional unit scaling type
pf :: [in] A structure indicating the open file to parse for model parameter values.
just_read :: [in] If true, this call will only read parameters without changing T or S.
- Call to: