MOM6
|
Data Types | |
type | diag_to_z_cs |
Functions/Subroutines | |
real function, dimension(cs%nk_zspace) | global_z_mean (var, G, CS, tracer) |
subroutine, public | calculate_z_diag_fields (u, v, h, ssh_in, frac_shelf_h, G, GV, CS) |
This subroutine maps tracers and velocities into depth space for diagnostics. More... | |
subroutine, public | calculate_z_transport (uh_int, vh_int, h, dt, G, GV, CS) |
This subroutine maps horizontal transport into depth space for diagnostic output. More... | |
subroutine, public | find_overlap (e, Z_top, Z_bot, k_max, k_start, k_top, k_bot, wt, z1, z2) |
This subroutine 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. More... | |
subroutine, public | find_limited_slope (val, e, slope, k) |
This subroutine determines a limited slope for val to be advected with a piecewise limited scheme. More... | |
subroutine, public | calc_zint_diags (h, in_ptrs, ids, num_diags, G, GV, CS) |
subroutine, public | register_z_tracer (tr_ptr, name, long_name, units, Time, G, CS, standard_name, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name) |
This subroutine registers a tracer to be output in depth space. More... | |
subroutine | register_z_tracer_low (tr_ptr, name, long_name, units, standard_name, Time, G, CS) |
This subroutine registers a tracer to be output in depth space. More... | |
subroutine, public | mom_diag_to_z_init (Time, G, GV, param_file, diag, CS) |
subroutine | get_z_depths (depth_file, int_depth_name, int_depth, cell_depth_name, z_axis_index, edge_index, nk_out) |
This subroutine reads the depths of the interfaces bounding the intended layers from a NetCDF file. If no appropriate file is found, -1 is returned as the number of layers in the output file. Also, a diag_manager axis is set up with the same information as this axis. More... | |
subroutine, public | mom_diag_to_z_end (CS) |
integer function, public | ocean_register_diag_with_z (tr_ptr, vardesc_tr, G, Time, CS) |
This subroutine registers a tracer to be output in depth space. More... | |
integer function | register_z_diag (var_desc, CS, day, missing) |
integer function, public | register_zint_diag (var_desc, CS, day) |
Variables | |
integer, parameter | no_zspace = -1 |
subroutine, public mom_diag_to_z::calc_zint_diags | ( | real, dimension(szi_(g),szj_(g),szk_(g)), intent(in) | h, |
type(p3d), dimension(:), intent(in) | in_ptrs, | ||
integer, dimension(:), intent(in) | ids, | ||
integer, intent(in) | num_diags, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(diag_to_z_cs), pointer | CS | ||
) |
[in] | g | The ocean's grid structure. |
[in] | h | Layer thicknesses, in H (usually m or kg m-2). |
[in] | gv | The ocean's vertical grid structure. |
Definition at line 811 of file MOM_diag_to_Z.F90.
References mom_error_handler::mom_error().
Referenced by mom_diabatic_driver::diabatic(), and mom_set_diffusivity::set_diffusivity().
subroutine, public mom_diag_to_z::calculate_z_diag_fields | ( | real, dimension(szib_(g),szj_(g),szk_(g)), intent(in) | u, |
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in) | v, | ||
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in) | h, | ||
real, dimension(szi_(g),szj_(g)), intent(in) | ssh_in, | ||
real, dimension(:,:), pointer | frac_shelf_h, | ||
type(ocean_grid_type), intent(inout) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(diag_to_z_cs), pointer | CS | ||
) |
This subroutine maps tracers and velocities into depth space for diagnostics.
[in,out] | g | The ocean's grid structure. |
[in] | gv | The ocean's vertical grid structure. |
[in] | u | The zonal velocity, in m s-1. |
[in] | v | The meridional velocity, in m s-1. |
[in] | h | Layer thicknesses, in H (usually m or kg m-2). |
[in] | ssh_in | Sea surface height (meter or kg/m2). |
cs | Control structure returned by previous call to diagnostics_init. |
Definition at line 167 of file MOM_diag_to_Z.F90.
References find_limited_slope(), find_overlap(), global_z_mean(), and mom_error_handler::mom_error().
subroutine, public mom_diag_to_z::calculate_z_transport | ( | real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in) | uh_int, |
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in) | vh_int, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in) | h, | ||
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(inout) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(diag_to_z_cs), pointer | CS | ||
) |
This subroutine maps horizontal transport into depth space for diagnostic output.
[in,out] | g | The ocean's grid structure. |
[in] | gv | The ocean's vertical grid structure. |
[in] | uh_int | Time integrated zonal transport (m3 or kg). |
[in] | vh_int | Time integrated meridional transport (m3 or kg). |
[in] | h | Layer thicknesses, in H (usually m or kg m-2). |
[in] | dt | The time difference in s since the last call to this subroutine. |
cs | Control structure returned by previous call to diagnostics_init. |
Definition at line 528 of file MOM_diag_to_Z.F90.
References mom_error_handler::mom_error().
subroutine, public mom_diag_to_z::find_limited_slope | ( | real, dimension(:), intent(in) | val, |
real, dimension(:), intent(in) | e, | ||
real, intent(out) | slope, | ||
integer, intent(in) | k | ||
) |
This subroutine determines a limited slope for val to be advected with a piecewise limited scheme.
[in] | val | A column of values that are being interpolated. |
[in] | e | Column interface heights (meter or kg/m2). |
[out] | slope | Normalized slope in the intracell distribution of val. |
[in] | k | Layer whose slope is being determined. |
Definition at line 777 of file MOM_diag_to_Z.F90.
Referenced by calculate_z_diag_fields(), and mom_tracer_z_init::tracer_z_init().
subroutine, public mom_diag_to_z::find_overlap | ( | real, dimension(:), intent(in) | e, |
real, intent(in) | Z_top, | ||
real, intent(in) | Z_bot, | ||
integer, intent(in) | k_max, | ||
integer, intent(in) | k_start, | ||
integer, intent(inout) | k_top, | ||
integer, intent(inout) | k_bot, | ||
real, dimension(:), intent(out) | wt, | ||
real, dimension(:), intent(out) | z1, | ||
real, dimension(:), intent(out) | z2 | ||
) |
This subroutine 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.
[in] | e | Column interface heights (meter or kg/m2). |
[in] | z_top | Top of range being mapped to (meter or kg/m2). |
[in] | z_bot | Bottom of range being mapped to (meter or kg/m2). |
[in] | k_max | Number of valid layers. |
[in] | k_start | Layer at which to start searching. |
[in,out] | k_top | Indices of top layers that overlap with the depth range. |
[in,out] | k_bot | Indices of bottom layers that overlap with the depth range. |
[out] | wt | Relative weights of each layer from k_top to k_bot. |
[out] | z2 | Depths of the top and bottom 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. |
Definition at line 701 of file MOM_diag_to_Z.F90.
Referenced by calculate_z_diag_fields(), and mom_tracer_z_init::tracer_z_init().
|
private |
This subroutine reads the depths of the interfaces bounding the intended layers from a NetCDF file. If no appropriate file is found, -1 is returned as the number of layers in the output file. Also, a diag_manager axis is set up with the same information as this axis.
Definition at line 1155 of file MOM_diag_to_Z.F90.
References mom_error_handler::mom_error().
Referenced by mom_diag_to_z_init().
|
private |
[in] | g | The ocean's grid structure |
Definition at line 119 of file MOM_diag_to_Z.F90.
Referenced by calculate_z_diag_fields().
subroutine, public mom_diag_to_z::mom_diag_to_z_end | ( | type(diag_to_z_cs), pointer | CS | ) |
Definition at line 1275 of file MOM_diag_to_Z.F90.
subroutine, public mom_diag_to_z::mom_diag_to_z_init | ( | type(time_type), intent(in) | Time, |
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(param_file_type), intent(in) | param_file, | ||
type(diag_ctrl), intent(inout), target | diag, | ||
type(diag_to_z_cs), pointer | CS | ||
) |
[in] | time | Current model time. |
[in] | g | The ocean's grid structure. |
[in] | gv | The ocean's vertical grid structure. |
[in] | param_file | A structure to parse for run-time parameters. |
[in,out] | diag | Struct to regulate diagnostic output. |
cs | Pointer to point to control structure for this module. |
Definition at line 1047 of file MOM_diag_to_Z.F90.
References get_z_depths(), and mom_error_handler::mom_error().
integer function, public mom_diag_to_z::ocean_register_diag_with_z | ( | real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in), target | tr_ptr, |
type(vardesc), intent(in) | vardesc_tr, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(time_type), intent(in) | Time, | ||
type(diag_to_z_cs), pointer | CS | ||
) |
This subroutine registers a tracer to be output in depth space.
[in] | g | The ocean's grid structure. |
[in] | tr_ptr | Tracer for translation to Z-space. |
[in] | vardesc_tr | Variable descriptor. |
[in] | time | Current model time. |
cs | Control struct returned by a previous call to diagnostics_init. |
Definition at line 1289 of file MOM_diag_to_Z.F90.
References mom_io::modify_vardesc(), mom_error_handler::mom_error(), no_zspace, mom_diag_mediator::ocean_register_diag(), mom_io::query_vardesc(), and register_z_diag().
|
private |
Definition at line 1354 of file MOM_diag_to_Z.F90.
References mom_error_handler::mom_error(), and mom_io::query_vardesc().
Referenced by ocean_register_diag_with_z().
subroutine, public mom_diag_to_z::register_z_tracer | ( | real, dimension(szi_(g),szj_(g),szk_(g)), intent(in), target | tr_ptr, |
character(len=*), intent(in) | name, | ||
character(len=*), intent(in) | long_name, | ||
character(len=*), intent(in) | units, | ||
type(time_type), intent(in) | Time, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(diag_to_z_cs), pointer | CS, | ||
character(len=*), intent(in), optional | standard_name, | ||
character(len=*), intent(in), optional | cmor_field_name, | ||
character(len=*), intent(in), optional | cmor_long_name, | ||
character(len=*), intent(in), optional | cmor_units, | ||
character(len=*), intent(in), optional | cmor_standard_name | ||
) |
This subroutine registers a tracer to be output in depth space.
[in] | g | The ocean's grid structure. |
[in] | tr_ptr | Tracer for translation to Z-space. |
[in] | name | name for the output tracer. |
[in] | long_name | Long name for the output tracer. |
[in] | units | Units of output tracer. |
[in] | time | Current model time. |
cs | Control struct returned by previous call to diagnostics_init. | |
[in] | cmor_field_name | cmor name of a field. |
[in] | cmor_long_name | cmor long name of a field. |
[in] | cmor_units | cmor units of a field. |
[in] | cmor_standard_name | cmor standardized name associated with a field. |
Definition at line 910 of file MOM_diag_to_Z.F90.
References register_z_tracer_low().
Referenced by advection_test_tracer::initialize_advection_test_tracer(), boundary_impulse_tracer::initialize_boundary_impulse_tracer(), dome_tracer::initialize_dome_tracer(), regional_dyes::initialize_dye_tracer(), ideal_age_example::initialize_ideal_age_tracer(), isomip_tracer::initialize_isomip_tracer(), mom_ocmip2_cfc::initialize_ocmip2_cfc(), oil_tracer::initialize_oil_tracer(), pseudo_salt_tracer::initialize_pseudo_salt_tracer(), and user_tracer_example::user_initialize_tracer().
|
private |
This subroutine registers a tracer to be output in depth space.
[in] | g | The ocean's grid structure. |
[in] | tr_ptr | Tracer for translation to Z-space. |
[in] | name | Name for the output tracer. |
[in] | long_name | Long name for output tracer. |
[in] | units | Units of output tracer. |
[in] | time | Current model time. |
cs | Control struct returned by previous call to diagnostics_init. |
Definition at line 981 of file MOM_diag_to_Z.F90.
References mom_error_handler::mom_error().
Referenced by register_z_tracer().
integer function, public mom_diag_to_z::register_zint_diag | ( | type(vardesc), intent(in) | var_desc, |
type(diag_to_z_cs), pointer | CS, | ||
type(time_type), intent(in) | day | ||
) |
Definition at line 1407 of file MOM_diag_to_Z.F90.
References mom_error_handler::mom_error(), and mom_io::query_vardesc().
Referenced by mom_diabatic_driver::diabatic_driver_init().
integer, parameter mom_diag_to_z::no_zspace = -1 |
Definition at line 114 of file MOM_diag_to_Z.F90.
Referenced by ocean_register_diag_with_z().