MOM6
mom_diag_manager_wrapper Module Reference

Detailed Description

A simple (very thin) wrapper for register_diag_field to avoid a compiler bug with PGI.

This module simply wraps register_diag_field() from FMS's diag_manager_mod. We used to be able to import register_diag_field and rename it to register_diag_field_fms with a simple "use, only : register_diag_field_fms => register_diag_field" but PGI 16.5 has a bug that refuses to compile this - earlier versions did work.

Data Types

interface  register_diag_field_fms
 A wrapper for register_diag_field_array() More...
 

Functions/Subroutines

integer function register_diag_field_array_fms (module_name, field_name, axes, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, area)
 An integer handle for a diagnostic array returned by register_diag_field() More...
 
integer function register_diag_field_scalar_fms (module_name, field_name, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, area)
 An integer handle for a diagnostic scalar array returned by register_diag_field() More...
 

Function/Subroutine Documentation

◆ register_diag_field_array_fms()

integer function mom_diag_manager_wrapper::register_diag_field_array_fms ( character(len=*), intent(in)  module_name,
character(len=*), intent(in)  field_name,
integer, dimension(:), intent(in)  axes,
type(time_type), intent(in)  init_time,
character(len=*), intent(in), optional  long_name,
character(len=*), intent(in), optional  units,
real, intent(in), optional  missing_value,
real, dimension(2), intent(in), optional  range,
logical, intent(in), optional  mask_variant,
character(len=*), intent(in), optional  standard_name,
logical, intent(in), optional  verbose,
logical, intent(in), optional  do_not_log,
character(len=*), intent(out), optional  err_msg,
character(len=*), intent(in), optional  interp_method,
integer, intent(in), optional  tile_count,
integer, intent(in), optional  area 
)

An integer handle for a diagnostic array returned by register_diag_field()

Parameters
[in]module_nameName of this module, usually "ocean_model" or "ice_shelf_model"
[in]field_nameName of the diagnostic field
[in]axesContainer w/ up to 3 integer handles that indicates axes for this field
[in]init_timeTime at which a field is first available?
[in]long_nameLong name of a field.
[in]unitsUnits of a field.
[in]standard_nameStandardized name associated with a field
[in]missing_valueA value that indicates missing values.
[in]rangeValid range of a variable (not used in MOM?)
[in]mask_variantIf true a logical mask must be provided with post_data calls (not used in MOM?)
[in]verboseIf true, FMS is verbose (not used in MOM?)
[in]do_not_logIf true, do not log something (not used in MOM?)
[out]err_msgString into which an error message might be placed (not used in MOM?)
[in]interp_methodIf 'none' indicates the field should not be interpolated as a scalar
[in]tile_countno clue (not used in MOM?)
[in]areaThe FMS id of cell area

Definition at line 22 of file MOM_diag_manager_wrapper.F90.

22  character(len=*), intent(in) :: module_name !< Name of this module, usually "ocean_model" or "ice_shelf_model"
23  character(len=*), intent(in) :: field_name !< Name of the diagnostic field
24  integer, intent(in) :: axes(:) !< Container w/ up to 3 integer handles that indicates axes for this field
25  type(time_type), intent(in) :: init_time !< Time at which a field is first available?
26  character(len=*), optional, intent(in) :: long_name !< Long name of a field.
27  character(len=*), optional, intent(in) :: units !< Units of a field.
28  character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
29  real, optional, intent(in) :: missing_value !< A value that indicates missing values.
30  real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
31  logical, optional, intent(in) :: mask_variant !< If true a logical mask must be provided with post_data calls (not used in MOM?)
32  logical, optional, intent(in) :: verbose !< If true, FMS is verbose (not used in MOM?)
33  logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
34  character(len=*), optional, intent(out):: err_msg !< String into which an error message might be placed (not used in MOM?)
35  character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should not be interpolated as a scalar
36  integer, optional, intent(in) :: tile_count !< no clue (not used in MOM?)
37  integer, optional, intent(in) :: area !< The FMS id of cell area
38  ! Local variables
39 
40  register_diag_field_array_fms = register_diag_field(module_name, field_name, axes, &
41  init_time, long_name=long_name, units=units, missing_value=missing_value, &
42  mask_variant=mask_variant, standard_name=standard_name, &
43  verbose=verbose, do_not_log=do_not_log, err_msg=err_msg, &
44  interp_method=interp_method)
45 

◆ register_diag_field_scalar_fms()

integer function mom_diag_manager_wrapper::register_diag_field_scalar_fms ( character(len=*), intent(in)  module_name,
character(len=*), intent(in)  field_name,
type(time_type), intent(in)  init_time,
character(len=*), intent(in), optional  long_name,
character(len=*), intent(in), optional  units,
real, intent(in), optional  missing_value,
real, dimension(2), intent(in), optional  range,
logical, intent(in), optional  mask_variant,
character(len=*), intent(in), optional  standard_name,
logical, intent(in), optional  verbose,
logical, intent(in), optional  do_not_log,
character(len=*), intent(out), optional  err_msg,
character(len=*), intent(in), optional  interp_method,
integer, intent(in), optional  tile_count,
integer, intent(in), optional  area 
)

An integer handle for a diagnostic scalar array returned by register_diag_field()

Parameters
[in]module_nameName of this module, usually "ocean_model" or "ice_shelf_model"
[in]field_nameName of the diagnostic field
[in]init_timeTime at which a field is first available?
[in]long_nameLong name of a field.
[in]unitsUnits of a field.
[in]standard_nameStandardized name associated with a field
[in]missing_valueA value that indicates missing values.
[in]rangeValid range of a variable (not used in MOM?)
[in]mask_variantIf true a logical mask must be provided with post_data calls (not used in MOM?)
[in]verboseIf true, FMS is verbose (not used in MOM?)
[in]do_not_logIf true, do not log something (not used in MOM?)
[out]err_msgString into which an error message might be placed (not used in MOM?)
[in]interp_methodIf 'none' indicates the field should not be interpolated as a scalar
[in]tile_countno clue (not used in MOM?)
[in]areaThe FMS id of cell area (not used for scalars)

Definition at line 52 of file MOM_diag_manager_wrapper.F90.

52  character(len=*), intent(in) :: module_name !< Name of this module, usually "ocean_model" or "ice_shelf_model"
53  character(len=*), intent(in) :: field_name !< Name of the diagnostic field
54  type(time_type), intent(in) :: init_time !< Time at which a field is first available?
55  character(len=*), optional, intent(in) :: long_name !< Long name of a field.
56  character(len=*), optional, intent(in) :: units !< Units of a field.
57  character(len=*), optional, intent(in) :: standard_name !< Standardized name associated with a field
58  real, optional, intent(in) :: missing_value !< A value that indicates missing values.
59  real, optional, intent(in) :: range(2) !< Valid range of a variable (not used in MOM?)
60  logical, optional, intent(in) :: mask_variant !< If true a logical mask must be provided with post_data calls (not used in MOM?)
61  logical, optional, intent(in) :: verbose !< If true, FMS is verbose (not used in MOM?)
62  logical, optional, intent(in) :: do_not_log !< If true, do not log something (not used in MOM?)
63  character(len=*), optional, intent(out):: err_msg !< String into which an error message might be placed (not used in MOM?)
64  character(len=*), optional, intent(in) :: interp_method !< If 'none' indicates the field should not be interpolated as a scalar
65  integer, optional, intent(in) :: tile_count !< no clue (not used in MOM?)
66  integer, optional, intent(in) :: area !< The FMS id of cell area (not used for scalars)
67  ! Local variables
68 
69  register_diag_field_scalar_fms = register_diag_field(module_name, field_name, &
70  init_time, long_name=long_name, units=units, missing_value=missing_value, &
71  standard_name=standard_name, do_not_log=do_not_log, err_msg=err_msg)
72