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

109  type(unit_scale_type), intent(inout) :: US !< A dimensional unit scaling type
110 
111  us%m_to_Z_restart = us%m_to_Z
112  us%m_to_L_restart = us%m_to_L
113  us%s_to_T_restart = us%s_to_T
114 

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

119  type(unit_scale_type), pointer :: US !< A dimensional unit scaling type
120 
121  deallocate( us )
122 

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

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