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 172 of file MOM_unit_scaling.F90.

173  type(unit_scale_type), intent(inout) :: US !< A dimensional unit scaling type
174 
175  us%m_to_Z_restart = us%m_to_Z
176  us%m_to_L_restart = us%m_to_L
177  us%s_to_T_restart = us%s_to_T
178  us%kg_m3_to_R_restart = us%kg_m3_to_R
179  us%J_kg_to_Q_restart = us%J_kg_to_Q
180 

◆ 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 184 of file MOM_unit_scaling.F90.

185  type(unit_scale_type), pointer :: US !< A dimensional unit scaling type
186 
187  deallocate( us )
188 

◆ unit_scaling_init()

subroutine, public mom_unit_scaling::unit_scaling_init ( type(param_file_type), intent(in), optional  param_file,
type(unit_scale_type), optional, 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 56 of file MOM_unit_scaling.F90.

57  type(param_file_type), optional, intent(in) :: param_file !< Parameter file handle/type
58  type(unit_scale_type), optional, pointer :: US !< A dimensional unit scaling type
59 
60  ! This routine initializes a unit_scale_type structure (US).
61 
62  ! Local variables
63  integer :: Z_power, L_power, T_power, R_power, Q_power
64  real :: Z_rescale_factor, L_rescale_factor, T_rescale_factor, R_rescale_factor, Q_rescale_factor
65  ! This include declares and sets the variable "version".
66 # include "version_variable.h"
67  character(len=16) :: mdl = "MOM_unit_scaling"
68 
69  if (.not.present(us)) return
70 
71  if (associated(us)) call mom_error(fatal, &
72  'unit_scaling_init: called with an associated US pointer.')
73  allocate(us)
74 
75  if (present(param_file)) then
76  ! Read all relevant parameters and write them to the model log.
77  call log_version(param_file, mdl, version, &
78  "Parameters for doing unit scaling of variables.", debugging=.true.)
79  call get_param(param_file, mdl, "Z_RESCALE_POWER", z_power, &
80  "An integer power of 2 that is used to rescale the model's "//&
81  "internal units of depths and heights. Valid values range from -300 to 300.", &
82  units="nondim", default=0, debuggingparam=.true.)
83  call get_param(param_file, mdl, "L_RESCALE_POWER", l_power, &
84  "An integer power of 2 that is used to rescale the model's "//&
85  "internal units of lateral distances. Valid values range from -300 to 300.", &
86  units="nondim", default=0, debuggingparam=.true.)
87  call get_param(param_file, mdl, "T_RESCALE_POWER", t_power, &
88  "An integer power of 2 that is used to rescale the model's "//&
89  "internal units of time. Valid values range from -300 to 300.", &
90  units="nondim", default=0, debuggingparam=.true.)
91  call get_param(param_file, mdl, "R_RESCALE_POWER", r_power, &
92  "An integer power of 2 that is used to rescale the model's "//&
93  "internal units of density. Valid values range from -300 to 300.", &
94  units="nondim", default=0, debuggingparam=.true.)
95  call get_param(param_file, mdl, "Q_RESCALE_POWER", q_power, &
96  "An integer power of 2 that is used to rescale the model's "//&
97  "internal units of heat content. Valid values range from -300 to 300.", &
98  units="nondim", default=0, debuggingparam=.true.)
99  else
100  z_power = 0 ; l_power = 0 ; t_power = 0 ; r_power = 0 ; q_power = 0
101  endif
102 
103  if (abs(z_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
104  "Z_RESCALE_POWER is outside of the valid range of -300 to 300.")
105  if (abs(l_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
106  "L_RESCALE_POWER is outside of the valid range of -300 to 300.")
107  if (abs(t_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
108  "T_RESCALE_POWER is outside of the valid range of -300 to 300.")
109  if (abs(r_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
110  "R_RESCALE_POWER is outside of the valid range of -300 to 300.")
111  if (abs(q_power) > 300) call mom_error(fatal, "unit_scaling_init: "//&
112  "Q_RESCALE_POWER is outside of the valid range of -300 to 300.")
113 
114  z_rescale_factor = 1.0
115  if (z_power /= 0) z_rescale_factor = 2.0**z_power
116  us%Z_to_m = 1.0 * z_rescale_factor
117  us%m_to_Z = 1.0 / z_rescale_factor
118 
119  l_rescale_factor = 1.0
120  if (l_power /= 0) l_rescale_factor = 2.0**l_power
121  us%L_to_m = 1.0 * l_rescale_factor
122  us%m_to_L = 1.0 / l_rescale_factor
123 
124  t_rescale_factor = 1.0
125  if (t_power /= 0) t_rescale_factor = 2.0**t_power
126  us%T_to_s = 1.0 * t_rescale_factor
127  us%s_to_T = 1.0 / t_rescale_factor
128 
129  r_rescale_factor = 1.0
130  if (r_power /= 0) r_rescale_factor = 2.0**r_power
131  us%R_to_kg_m3 = 1.0 * r_rescale_factor
132  us%kg_m3_to_R = 1.0 / r_rescale_factor
133 
134  q_rescale_factor = 1.0
135  if (q_power /= 0) q_rescale_factor = 2.0**q_power
136  us%Q_to_J_kg = 1.0 * q_rescale_factor
137  us%J_kg_to_Q = 1.0 / q_rescale_factor
138 
139  ! These are useful combinations of the fundamental scale conversion factors set above.
140  us%Z_to_L = us%Z_to_m * us%m_to_L
141  us%L_to_Z = us%L_to_m * us%m_to_Z
142  ! Horizontal velocities:
143  us%L_T_to_m_s = us%L_to_m * us%s_to_T
144  us%m_s_to_L_T = us%m_to_L * us%T_to_s
145  ! Horizontal accelerations:
146  us%L_T2_to_m_s2 = us%L_to_m * us%s_to_T**2
147  ! It does not look like US%m_s2_to_L_T2 would be used, so it does not exist.
148  ! Vertical diffusivities and viscosities:
149  us%Z2_T_to_m2_s = us%Z_to_m**2 * us%s_to_T
150  us%m2_s_to_Z2_T = us%m_to_Z**2 * us%T_to_s
151  ! Column mass loads:
152  us%RZ_to_kg_m2 = us%R_to_kg_m3 * us%Z_to_m
153  ! It does not seem like US%kg_m2_to_RZ would be used enough in MOM6 to justify its existence.
154  ! Vertical mass fluxes:
155  us%kg_m2s_to_RZ_T = us%kg_m3_to_R * us%m_to_Z * us%T_to_s
156  us%RZ_T_to_kg_m2s = us%R_to_kg_m3 * us%Z_to_m * us%s_to_T
157  ! Turbulent kinetic energy vertical fluxes:
158  us%RZ3_T3_to_W_m2 = us%R_to_kg_m3 * us%Z_to_m**3 * us%s_to_T**3
159  us%W_m2_to_RZ3_T3 = us%kg_m3_to_R * us%m_to_Z**3 * us%T_to_s**3
160  ! Vertical heat fluxes:
161  us%W_m2_to_QRZ_T = us%J_kg_to_Q * us%kg_m3_to_R * us%m_to_Z * us%T_to_s
162  us%QRZ_T_to_W_m2 = us%Q_to_J_kg * us%R_to_kg_m3 * us%Z_to_m * us%s_to_T
163  ! Pressures:
164  us%RL2_T2_to_Pa = us%R_to_kg_m3 * us%L_T_to_m_s**2
165  ! It does not seem like US%Pa_to_RL2_T2 would be used enough in MOM6 to justify its existence.
166  ! US%Pa_to_RL2_T2 = US%kg_m3_to_R * US%m_s_to_L_T**2
167