MOM6
mom_obsolete_diagnostics Module Reference

Detailed Description

Provides a mechanism for recording diagnostic variables that are no longer valid, along with their replacement name if appropriate.

Functions/Subroutines

subroutine, public register_obsolete_diagnostics (param_file, diag)
 Scan through the diag_table searching for obsolete parameters and issue informational messages and optionallly a FATAL error. More...
 
logical function found_in_diagtable (diag, varName, newVarName)
 Fakes a register of a diagnostic to find out if an obsolete parameter appears in the diag_table. More...
 

Function/Subroutine Documentation

◆ found_in_diagtable()

logical function mom_obsolete_diagnostics::found_in_diagtable ( type(diag_ctrl), intent(in)  diag,
character(len=*), intent(in)  varName,
character(len=*), intent(in), optional  newVarName 
)
private

Fakes a register of a diagnostic to find out if an obsolete parameter appears in the diag_table.

Parameters
[in]diagA structure used to control diagnostics.
[in]varnameThe obsolete diagnostic name
[in]newvarnameThe valid name of this diagnostic

Definition at line 67 of file MOM_obsolete_diagnostics.F90.

67  type(diag_ctrl), intent(in) :: diag !< A structure used to control diagnostics.
68  character(len=*), intent(in) :: varName !< The obsolete diagnostic name
69  character(len=*), optional, intent(in) :: newVarName !< The valid name of this diagnostic
70  ! Local
71  integer :: handle ! Integer handle returned from diag_manager
72 
73  ! We use register_static_field_fms() instead of register_static_field() so
74  ! that the diagnostic does not appear in the available diagnostics list.
75  handle = register_static_field_fms('ocean_model', varname, &
76  diag%axesT1%handles, 'Obsolete parameter', 'N/A')
77 
78  found_in_diagtable = (handle>0)
79 
80  if (handle>0 .and. is_root_pe()) then
81  if (present(newvarname)) then
82  call mom_error(warning, 'MOM_obsolete_params: '// &
83  'diag_table entry "'//trim(varname)//'" found. Use '// &
84  '"'//trim(newvarname)//'" instead.' )
85  else
86  call mom_error(warning, 'MOM_obsolete_params: '// &
87  'diag_table entry "'//trim(varname)//'" is obsolete.' )
88  endif
89  endif
90 

References mom_error_handler::is_root_pe(), and mom_error_handler::mom_error().

Referenced by register_obsolete_diagnostics().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ register_obsolete_diagnostics()

subroutine, public mom_obsolete_diagnostics::register_obsolete_diagnostics ( type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(in)  diag 
)

Scan through the diag_table searching for obsolete parameters and issue informational messages and optionallly a FATAL error.

Parameters
[in]param_fileThe parameter file handle.
[in]diagA structure used to control diagnostics.

Definition at line 23 of file MOM_obsolete_diagnostics.F90.

23  type(param_file_type), intent(in) :: param_file !< The parameter file handle.
24  type(diag_ctrl), intent(in) :: diag !< A structure used to control diagnostics.
25 ! This include declares and sets the variable "version".
26 #include "version_variable.h"
27  ! Local variables
28  character(len=40) :: mdl = "MOM_obsolete_diagnostics" !< This module's name.
29  logical :: foundEntry, causeFatal
30  integer :: errType
31 
32  call log_version(param_file, mdl, version)
33  call get_param(param_file, mdl, "OBSOLETE_DIAGNOSTIC_IS_FATAL", causefatal, &
34  "If an obsolete diagnostic variable appears in the diag_table, "// &
35  "cause a FATAL error rather than issue a WARNING.", default=.true.)
36 
37  foundentry = .false.
38  ! Each obsolete entry, with replacement name is available.
39  if (found_in_diagtable(diag, 'Net_Heat', 'net_heat_surface or net_heat_coupler')) foundentry = .true.
40  if (found_in_diagtable(diag, 'PmE', 'PRCmE')) foundentry = .true.
41  if (found_in_diagtable(diag, 'froz_precip', 'fprec')) foundentry = .true.
42  if (found_in_diagtable(diag, 'liq_precip', 'lprec')) foundentry = .true.
43  if (found_in_diagtable(diag, 'virt_precip', 'vprec')) foundentry = .true.
44  if (found_in_diagtable(diag, 'froz_runoff', 'frunoff')) foundentry = .true.
45  if (found_in_diagtable(diag, 'liq_runoff', 'lrunoff')) foundentry = .true.
46  if (found_in_diagtable(diag, 'calving_heat_content', 'heat_content_frunoff')) foundentry = .true.
47  if (found_in_diagtable(diag, 'precip_heat_content', 'heat_content_lprec')) foundentry = .true.
48  if (found_in_diagtable(diag, 'evap_heat_content', 'heat_content_massout')) foundentry = .true.
49  if (found_in_diagtable(diag, 'runoff_heat_content', 'heat_content_lrunoff')) foundentry = .true.
50  if (found_in_diagtable(diag, 'latent_fprec')) foundentry = .true.
51  if (found_in_diagtable(diag, 'latent_calve')) foundentry = .true.
52  if (found_in_diagtable(diag, 'heat_rest', 'heat_restore')) foundentry = .true.
53  if (found_in_diagtable(diag, 'KPP_dTdt', 'KPP_NLT_dTdt')) foundentry = .true.
54  if (found_in_diagtable(diag, 'KPP_dSdt', 'KPP_NLT_dSdt')) foundentry = .true.
55 
56  if (causefatal) then; errtype = fatal
57  else ; errtype = warning ; endif
58  if (foundentry .and. is_root_pe()) &
59  call mom_error(errtype, 'MOM_obsolete_diagnostics: '//&
60  'Obsolete diagnostics found in diag_table')
61 

References found_in_diagtable(), mom_error_handler::is_root_pe(), and mom_error_handler::mom_error().

Referenced by mom::initialize_mom().

Here is the call graph for this function:
Here is the caller graph for this function: