MOM6
mom_unit_scaling Module Reference

Detailed Description

Provides a transparent unit rescaling type to facilitate dimensional consistency testing.

Data Types

type  unit_scale_type
 Describes various unit conversion factors. More...
 

Functions/Subroutines

subroutine, public unit_scaling_init (param_file, US)
 Allocates and initializes the ocean model unit scaling type. More...
 
subroutine, public fix_restart_unit_scaling (US)
 Set the unit scaling factors for output to restart files to the unit scaling factors for this run. More...
 
subroutine, public unit_scaling_end (US)
 Deallocates a unit scaling structure. More...
 

Function/Subroutine Documentation

◆ fix_restart_unit_scaling()

subroutine, public mom_unit_scaling::fix_restart_unit_scaling ( type(unit_scale_type), intent(inout)  US)

Set the unit scaling factors for output to restart files to the unit scaling factors for this run.

Parameters
[in,out]usA dimensional unit scaling type

Definition at line 123 of file MOM_unit_scaling.F90.

123  type(unit_scale_type), intent(inout) :: US !< A dimensional unit scaling type
124 
125  us%m_to_Z_restart = us%m_to_Z
126  us%m_to_L_restart = us%m_to_L
127  us%s_to_T_restart = us%s_to_T
128  us%kg_m3_to_R_restart = us%kg_m3_to_R
129 

Referenced by mom::finish_mom_initialization(), and mom_ice_shelf::initialize_ice_shelf().

Here is the caller graph for this function:

◆ unit_scaling_end()

subroutine, public mom_unit_scaling::unit_scaling_end ( type(unit_scale_type), pointer  US)

Deallocates a unit scaling structure.

Parameters
usA dimensional unit scaling type

Definition at line 134 of file MOM_unit_scaling.F90.

134  type(unit_scale_type), pointer :: US !< A dimensional unit scaling type
135 
136  deallocate( us )
137 

Referenced by mom::mom_end().

Here is the caller graph for this function:

◆ unit_scaling_init()

subroutine, public mom_unit_scaling::unit_scaling_init ( type(param_file_type), intent(in)  param_file,
type(unit_scale_type), pointer  US 
)

Allocates and initializes the ocean model unit scaling type.

Parameters
[in]param_fileParameter file handle/type
usA dimensional unit scaling type

Definition at line 44 of file MOM_unit_scaling.F90.

44  type(param_file_type), intent(in) :: param_file !< Parameter file handle/type
45  type(unit_scale_type), pointer :: US !< A dimensional unit scaling type
46 
47  ! This routine initializes a unit_scale_type structure (US).
48 
49  ! Local variables
50  integer :: Z_power, L_power, T_power, R_power
51  real :: Z_rescale_factor, L_rescale_factor, T_rescale_factor, R_rescale_factor
52  ! This include declares and sets the variable "version".
53 # include "version_variable.h"
54  character(len=16) :: mdl = "MOM_unit_scaling"
55 
56  if (associated(us)) call mom_error(fatal, &
57  'unit_scaling_init: called with an associated US pointer.')
58  allocate(us)
59 
60  ! Read all relevant parameters and write them to the model log.
61  call log_version(param_file, mdl, version, &
62  "Parameters for doing unit scaling of variables.")
63  call get_param(param_file, mdl, "Z_RESCALE_POWER", z_power, &
64  "An integer power of 2 that is used to rescale the model's "//&
65  "intenal units of depths and heights. Valid values range from -300 to 300.", &
66  units="nondim", default=0, debuggingparam=.true.)
67  call get_param(param_file, mdl, "L_RESCALE_POWER", l_power, &
68  "An integer power of 2 that is used to rescale the model's "//&
69  "intenal units of lateral distances. Valid values range from -300 to 300.", &
70  units="nondim", default=0, debuggingparam=.true.)
71  call get_param(param_file, mdl, "T_RESCALE_POWER", t_power, &
72  "An integer power of 2 that is used to rescale the model's "//&
73  "intenal units of time. Valid values range from -300 to 300.", &
74  units="nondim", default=0, debuggingparam=.true.)
75  call get_param(param_file, mdl, "R_RESCALE_POWER", r_power, &
76  "An integer power of 2 that is used to rescale the model's "//&
77  "intenal units of density. Valid values range from -300 to 300.", &
78  units="nondim", default=0, debuggingparam=.true.)
79  if (abs(z_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
80  "Z_RESCALE_POWER is outside of the valid range of -300 to 300.")
81  if (abs(l_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
82  "L_RESCALE_POWER is outside of the valid range of -300 to 300.")
83  if (abs(t_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
84  "T_RESCALE_POWER is outside of the valid range of -300 to 300.")
85  if (abs(r_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
86  "R_RESCALE_POWER is outside of the valid range of -300 to 300.")
87 
88  z_rescale_factor = 1.0
89  if (z_power /= 0) z_rescale_factor = 2.0**z_power
90  us%Z_to_m = 1.0 * z_rescale_factor
91  us%m_to_Z = 1.0 / z_rescale_factor
92 
93  l_rescale_factor = 1.0
94  if (l_power /= 0) l_rescale_factor = 2.0**l_power
95  us%L_to_m = 1.0 * l_rescale_factor
96  us%m_to_L = 1.0 / l_rescale_factor
97 
98  t_rescale_factor = 1.0
99  if (t_power /= 0) t_rescale_factor = 2.0**t_power
100  us%T_to_s = 1.0 * t_rescale_factor
101  us%s_to_T = 1.0 / t_rescale_factor
102 
103  r_rescale_factor = 1.0
104  if (r_power /= 0) r_rescale_factor = 2.0**r_power
105  us%R_to_kg_m3 = 1.0 * r_rescale_factor
106  us%kg_m3_to_R = 1.0 / r_rescale_factor
107 
108  ! These are useful combinations of the fundamental scale conversion factors set above.
109  us%Z_to_L = us%Z_to_m * us%m_to_L
110  us%L_to_Z = us%L_to_m * us%m_to_Z
111  us%L_T_to_m_s = us%L_to_m * us%s_to_T
112  us%m_s_to_L_T = us%m_to_L * us%T_to_s
113  us%L_T2_to_m_s2 = us%L_to_m * us%s_to_T**2
114  ! It does not look like US%m_s2_to_L_T2 would be used, so it does not exist.
115  us%Z2_T_to_m2_s = us%Z_to_m**2 * us%s_to_T
116  us%m2_s_to_Z2_T = us%m_to_Z**2 * us%T_to_s
117 

References mom_error_handler::mom_error().

Referenced by mom_oda_driver_mod::init_oda(), and mom_ice_shelf::initialize_ice_shelf().

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