MOM6
MOM_unit_scaling.F90
1 !> Provides a transparent unit rescaling type to facilitate dimensional consistency testing
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, mom_mesg, fatal
8 
9 implicit none ; private
10 
11 public unit_scaling_init, unit_scaling_end, fix_restart_unit_scaling
12 
13 !> Describes various unit conversion factors
14 type, public :: unit_scale_type
15  real :: m_to_z !< A constant that translates distances in meters to the units of depth.
16  real :: z_to_m !< A constant that translates distances in the units of depth to meters.
17  real :: m_to_l !< A constant that translates lengths in meters to the units of horizontal lengths.
18  real :: l_to_m !< A constant that translates lengths in the units of horizontal lengths to meters.
19  real :: s_to_t !< A constant that translates time intervals in seconds to the units of time.
20  real :: t_to_s !< A constant that translates the units of time to seconds.
21  real :: r_to_kg_m3 !< A constant that translates the units of density to kilograms per meter cubed.
22  real :: kg_m3_to_r !< A constant that translates kilograms per meter cubed to the units of density.
23  real :: q_to_j_kg !< A constant that translates the units of enthalpy to Joules per kilogram.
24  real :: j_kg_to_q !< A constant that translates Joules per kilogram to the units of enthalpy.
25 
26  ! These are useful combinations of the fundamental scale conversion factors above.
27  real :: z_to_l !< Convert vertical distances to lateral lengths
28  real :: l_to_z !< Convert lateral lengths to vertical distances
29  real :: l_t_to_m_s !< Convert lateral velocities from L T-1 to m s-1.
30  real :: m_s_to_l_t !< Convert lateral velocities from m s-1 to L T-1.
31  real :: l_t2_to_m_s2 !< Convert lateral accelerations from L T-2 to m s-2.
32  real :: z2_t_to_m2_s !< Convert vertical diffusivities from Z2 T-1 to m2 s-1.
33  real :: m2_s_to_z2_t !< Convert vertical diffusivities from m2 s-1 to Z2 T-1.
34  real :: w_m2_to_qrz_t !< Convert heat fluxes from W m-2 to Q R Z T-1.
35  real :: qrz_t_to_w_m2 !< Convert heat fluxes from Q R Z T-1 to W m-2.
36  ! Not used enough: real :: kg_m2_to_RZ !< Convert mass loads from kg m-2 to R Z.
37  real :: rz_to_kg_m2 !< Convert mass loads from R Z to kg m-2.
38  real :: kg_m2s_to_rz_t !< Convert mass fluxes from kg m-2 s-1 to R Z T-1.
39  real :: rz_t_to_kg_m2s !< Convert mass fluxes from R Z T-1 to kg m-2 s-1.
40  real :: rz3_t3_to_w_m2 !< Convert turbulent kinetic energy fluxes from R Z3 T-3 to W m-2.
41  real :: w_m2_to_rz3_t3 !< Convert turbulent kinetic energy fluxes from W m-2 to R Z3 T-3.
42  real :: rl2_t2_to_pa !< Convert pressures from R L2 T-2 to Pa.
43  ! Not used enough: real :: Pa_to_RL2_T2 !< Convert pressures from Pa to R L2 T-2.
44 
45  ! These are used for changing scaling across restarts.
46  real :: m_to_z_restart = 0.0 !< A copy of the m_to_Z that is used in restart files.
47  real :: m_to_l_restart = 0.0 !< A copy of the m_to_L that is used in restart files.
48  real :: s_to_t_restart = 0.0 !< A copy of the s_to_T that is used in restart files.
49  real :: kg_m3_to_r_restart = 0.0 !< A copy of the kg_m3_to_R that is used in restart files.
50  real :: j_kg_to_q_restart = 0.0 !< A copy of the J_kg_to_Q that is used in restart files.
51 end type unit_scale_type
52 
53 contains
54 
55 !> Allocates and initializes the ocean model unit scaling type
56 subroutine unit_scaling_init( param_file, US )
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 
168 end subroutine unit_scaling_init
169 
170 !> Set the unit scaling factors for output to restart files to the unit scaling
171 !! factors for this run.
172 subroutine fix_restart_unit_scaling(US)
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 
181 end subroutine fix_restart_unit_scaling
182 
183 !> Deallocates a unit scaling structure.
184 subroutine unit_scaling_end( US )
185  type(unit_scale_type), pointer :: us !< A dimensional unit scaling type
186 
187  deallocate( us )
188 
189 end subroutine unit_scaling_end
190 
191 end module mom_unit_scaling
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2