\hypertarget{MOM__unit__scaling_8F90_source}{}\section{M\+O\+M\+\_\+unit\+\_\+scaling.\+F90}
\label{MOM__unit__scaling_8F90_source}\index{/home/cermak/src/\+M\+O\+M6/src/framework/\+M\+O\+M\+\_\+unit\+\_\+scaling.\+F90@{/home/cermak/src/\+M\+O\+M6/src/framework/\+M\+O\+M\+\_\+unit\+\_\+scaling.\+F90}}

\begin{DoxyCode}
00001 \textcolor{comment}{!> Provides a transparent unit rescaling type to facilitate dimensional consistency testing}
\Hypertarget{MOM__unit__scaling_8F90_source_l00002}\hyperlink{namespacemom__unit__scaling}{00002} \textcolor{keyword}{module} \hyperlink{namespacemom__unit__scaling}{mom\_unit\_scaling}
00003 
00004 \textcolor{comment}{! This file is part of MOM6. See LICENSE.md for the license.}
00005 
00006 \textcolor{keywordtype}{use }mom\_error\_handler\textcolor{keywordtype}{, only} : mom\_error, mom\_mesg, fatal
00007 \textcolor{keywordtype}{use }mom\_file\_parser\textcolor{keywordtype}{, only} : get\_param, log\_param, log\_version, param\_file\_type
00008 
00009 \textcolor{keywordtype}{implicit none} ; \textcolor{keywordtype}{private}
00010 
00011 \textcolor{keywordtype}{public} \hyperlink{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}{unit\_scaling\_init}, \hyperlink{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}{unit\_scaling\_end}, 
      \hyperlink{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}{fix\_restart\_unit\_scaling}
00012 \textcolor{comment}{}
00013 \textcolor{comment}{!> Describes various unit conversion factors}
\Hypertarget{MOM__unit__scaling_8F90_source_l00014}\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{00014} \textcolor{keyword}{type}, \textcolor{keyword}{public} :: \hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\_scale\_type}
\Hypertarget{MOM__unit__scaling_8F90_source_l00015}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a9e9bcd35c60ffe1ad0ca65c39a0c9ee4}{00015}   \textcolor{keywordtype}{real} :: m\_to\_z\textcolor{comment}{ !< A constant that translates distances in meters to the units of depth.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00016}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_aa6683003710c410cb29766ca59e67386}{00016}   \textcolor{keywordtype}{real} :: z\_to\_m\textcolor{comment}{ !< A constant that translates distances in the units of depth to meters.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00017}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a1a734d24b7a23915805fee4745c56302}{00017}   \textcolor{keywordtype}{real} :: m\_to\_l\textcolor{comment}{ !< A constant that translates lengths in meters to the units of horizontal lengths.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00018}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_aa9c07449f18ad5b33d8c5ff99a9583b9}{00018}   \textcolor{keywordtype}{real} :: l\_to\_m\textcolor{comment}{ !< A constant that translates lengths in the units of horizontal lengths to meters.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00019}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_ad7a5c7cfc6e5aadec5404748ceeea8c4}{00019}   \textcolor{keywordtype}{real} :: s\_to\_t\textcolor{comment}{ !< A constant that translates time intervals in seconds to the units of time.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00020}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a64fe95d7ab8988296ef331d200286512}{00020}   \textcolor{keywordtype}{real} :: t\_to\_s\textcolor{comment}{ !< A constant that translates the units of time to seconds.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00021}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a8bf357de04042f736da22f050d39c029}{00021}   \textcolor{keywordtype}{real} :: r\_to\_kg\_m3\textcolor{comment}{ !< A constant that translates the units of density to kilograms per meter cubed.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00022}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_aeb294a6431c4125780ac7a774c5e5358}{00022}   \textcolor{keywordtype}{real} :: kg\_m3\_to\_r\textcolor{comment}{ !< A constant that translates kilograms per meter cubed to the units of density.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00023}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a9e1309c254ade9a943b94f2d08e307f5}{00023}   \textcolor{keywordtype}{real} :: q\_to\_j\_kg\textcolor{comment}{  !< A constant that translates the units of enthalpy to Joules per kilogram.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00024}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a61e20c873414ced5ce42ee1e6340a74a}{00024}   \textcolor{keywordtype}{real} :: j\_kg\_to\_q\textcolor{comment}{  !< A constant that translates Joules per kilogram to the units of enthalpy.}
00025 
00026   \textcolor{comment}{! These are useful combinations of the fundamental scale conversion factors above.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00027}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a933c52cfe6d06a15560df5c268c7236e}{00027}   \textcolor{keywordtype}{real} :: z\_to\_l\textcolor{comment}{          !< Convert vertical distances to lateral lengths}
\Hypertarget{MOM__unit__scaling_8F90_source_l00028}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a5e6fdba434e46632df141b3e59ed1fa3}{00028}   \textcolor{keywordtype}{real} :: l\_to\_z\textcolor{comment}{          !< Convert lateral lengths to vertical distances}
\Hypertarget{MOM__unit__scaling_8F90_source_l00029}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_aa3773a01be375b3141ee04f31970ea62}{00029}   \textcolor{keywordtype}{real} :: l\_t\_to\_m\_s\textcolor{comment}{      !< Convert lateral velocities from L T-1 to m s-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00030}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a13075ccb70a9a0b5126aa22dad3db3f7}{00030}   \textcolor{keywordtype}{real} :: m\_s\_to\_l\_t\textcolor{comment}{      !< Convert lateral velocities from m s-1 to L T-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00031}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a4d885c1acf4a3469143214c079f3c966}{00031}   \textcolor{keywordtype}{real} :: l\_t2\_to\_m\_s2\textcolor{comment}{    !< Convert lateral accelerations from L T-2 to m s-2.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00032}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a5ac4c9a5f24988048b87b8c8be00d473}{00032}   \textcolor{keywordtype}{real} :: z2\_t\_to\_m2\_s\textcolor{comment}{    !< Convert vertical diffusivities from Z2 T-1 to m2 s-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00033}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a5f746535d28fa4264acab5008788e8ac}{00033}   \textcolor{keywordtype}{real} :: m2\_s\_to\_z2\_t\textcolor{comment}{    !< Convert vertical diffusivities from m2 s-1 to Z2 T-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00034}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a51868180eba3d6ce17ecde6aaab0b4e7}{00034}   \textcolor{keywordtype}{real} :: w\_m2\_to\_qrz\_t\textcolor{comment}{   !< Convert heat fluxes from W m-2 to Q R Z T-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00035}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a096733eb3a3ffd970392a756671d5884}{00035}   \textcolor{keywordtype}{real} :: qrz\_t\_to\_w\_m2\textcolor{comment}{   !< Convert heat fluxes from Q R Z T-1 to W m-2.}
00036   \textcolor{comment}{! Not used enough:  real :: kg\_m2\_to\_RZ   !< Convert mass loads from kg m-2 to R Z.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00037}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a3da58884b89b3b8b57c87a8bea3c1318}{00037}   \textcolor{keywordtype}{real} :: rz\_to\_kg\_m2\textcolor{comment}{     !< Convert mass loads from R Z to kg m-2.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00038}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a0a7140cfea2ef9d5f7ef9ab7ce419567}{00038}   \textcolor{keywordtype}{real} :: kg\_m2s\_to\_rz\_t\textcolor{comment}{  !< Convert mass fluxes from kg m-2 s-1 to R Z T-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00039}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_af64fdb0914704e583b9f6813a0090118}{00039}   \textcolor{keywordtype}{real} :: rz\_t\_to\_kg\_m2s\textcolor{comment}{  !< Convert mass fluxes from R Z T-1 to kg m-2 s-1.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00040}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a50f4d02009d7bc90055f3b7f09342b15}{00040}   \textcolor{keywordtype}{real} :: rz3\_t3\_to\_w\_m2\textcolor{comment}{  !< Convert turbulent kinetic energy fluxes from R Z3 T-3 to W m-2.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00041}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a6281d21581543c235553106557446205}{00041}   \textcolor{keywordtype}{real} :: w\_m2\_to\_rz3\_t3\textcolor{comment}{  !< Convert turbulent kinetic energy fluxes from W m-2 to R Z3 T-3.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00042}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_af5b174ed064f32776f76e986a20c9033}{00042}   \textcolor{keywordtype}{real} :: rl2\_t2\_to\_pa\textcolor{comment}{    !< Convert pressures from R L2 T-2 to Pa.}
00043   \textcolor{comment}{! Not used enough:  real :: Pa\_to\_RL2\_T2    !< Convert pressures from Pa to R L2 T-2.}
00044 
00045   \textcolor{comment}{! These are used for changing scaling across restarts.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00046}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a75c2ecdf2cb841803dbcd5b1f6416d3e}{00046}   \textcolor{keywordtype}{real} :: m\_to\_z\_restart = 0.0\textcolor{comment}{ !< A copy of the m\_to\_Z that is used in restart files.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00047}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_acd8d43a5d748a42106918d86ec2a1945}{00047}   \textcolor{keywordtype}{real} :: m\_to\_l\_restart = 0.0\textcolor{comment}{ !< A copy of the m\_to\_L that is used in restart files.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00048}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_ae65d0572a709312821f1baa6e02349a7}{00048}   \textcolor{keywordtype}{real} :: s\_to\_t\_restart = 0.0\textcolor{comment}{ !< A copy of the s\_to\_T that is used in restart files.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00049}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a489d12f8352d574af7a20c377c794780}{00049}   \textcolor{keywordtype}{real} :: kg\_m3\_to\_r\_restart = 0.0\textcolor{comment}{ !< A copy of the kg\_m3\_to\_R that is used in restart files.}
\Hypertarget{MOM__unit__scaling_8F90_source_l00050}\hyperlink{structmom__unit__scaling_1_1unit__scale__type_a1f2f1512df2efd36aba57c0bf4e333a1}{00050}   \textcolor{keywordtype}{real} :: j\_kg\_to\_q\_restart = 0.0\textcolor{comment}{ !< A copy of the J\_kg\_to\_Q that is used in restart files.}
00051 \textcolor{keyword}{end type }\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\_scale\_type}
00052 
00053 \textcolor{keyword}{contains}
00054 \textcolor{comment}{}
00055 \textcolor{comment}{!> Allocates and initializes the ocean model unit scaling type}
00056 \textcolor{keyword}{subroutine }\hyperlink{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}{unit\_scaling\_init}( param\_file, US )
\Hypertarget{MOM__unit__scaling_8F90_source_l00057}\hyperlink{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}{00057}   \textcolor{keywordtype}{type}(param\_file\_type), \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)} :: param\_file\textcolor{comment}{ !< Parameter file handle/type}
00058   \textcolor{keywordtype}{type}(\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\_scale\_type}), \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{pointer}    :: US\textcolor{comment}{         !< A dimensional unit scaling
       type}
00059 
00060   \textcolor{comment}{! This routine initializes a unit\_scale\_type structure (US).}
00061 
00062   \textcolor{comment}{! Local variables}
00063   \textcolor{keywordtype}{integer} :: Z\_power, L\_power, T\_power, R\_power, Q\_power
00064   \textcolor{keywordtype}{real}    :: Z\_rescale\_factor, L\_rescale\_factor, T\_rescale\_factor, R\_rescale\_factor, Q\_rescale\_factor
00065   \textcolor{comment}{! This include declares and sets the variable "version".}
00066 \textcolor{preprocessor}{# include "version\_variable.h"}
00067 \textcolor{preprocessor}{}  \textcolor{keywordtype}{character(len=16)} :: mdl = \textcolor{stringliteral}{"MOM\_unit\_scaling"}
00068 
00069   \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{present}(us)) \textcolor{keywordflow}{return}
00070 
00071   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(us)) \textcolor{keyword}{call }mom\_error(fatal, &
00072      \textcolor{stringliteral}{'unit\_scaling\_init: called with an associated US pointer.'})
00073   \textcolor{keyword}{allocate}(us)
00074 
00075   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(param\_file)) \textcolor{keywordflow}{then}
00076     \textcolor{comment}{! Read all relevant parameters and write them to the model log.}
00077     \textcolor{keyword}{call }log\_version(param\_file, mdl, version, &
00078                  \textcolor{stringliteral}{"Parameters for doing unit scaling of variables."}, debugging=.true.)
00079     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"Z\_RESCALE\_POWER"}, z\_power, &
00080                  \textcolor{stringliteral}{"An integer power of 2 that is used to rescale the model's "}//&
00081                  \textcolor{stringliteral}{"internal units of depths and heights.  Valid values range from -300 to 300."}, &
00082                  units=\textcolor{stringliteral}{"nondim"}, default=0, debuggingparam=.true.)
00083     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"L\_RESCALE\_POWER"}, l\_power, &
00084                  \textcolor{stringliteral}{"An integer power of 2 that is used to rescale the model's "}//&
00085                  \textcolor{stringliteral}{"internal units of lateral distances.  Valid values range from -300 to 300."}, &
00086                  units=\textcolor{stringliteral}{"nondim"}, default=0, debuggingparam=.true.)
00087     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"T\_RESCALE\_POWER"}, t\_power, &
00088                  \textcolor{stringliteral}{"An integer power of 2 that is used to rescale the model's "}//&
00089                  \textcolor{stringliteral}{"internal units of time.  Valid values range from -300 to 300."}, &
00090                  units=\textcolor{stringliteral}{"nondim"}, default=0, debuggingparam=.true.)
00091     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"R\_RESCALE\_POWER"}, r\_power, &
00092                  \textcolor{stringliteral}{"An integer power of 2 that is used to rescale the model's "}//&
00093                  \textcolor{stringliteral}{"internal units of density.  Valid values range from -300 to 300."}, &
00094                  units=\textcolor{stringliteral}{"nondim"}, default=0, debuggingparam=.true.)
00095     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"Q\_RESCALE\_POWER"}, q\_power, &
00096                  \textcolor{stringliteral}{"An integer power of 2 that is used to rescale the model's "}//&
00097                    \textcolor{stringliteral}{"internal units of heat content.  Valid values range from -300 to 300."}, &
00098                    units=\textcolor{stringliteral}{"nondim"}, default=0, debuggingparam=.true.)
00099   \textcolor{keywordflow}{else}
00100     z\_power = 0 ; l\_power = 0 ; t\_power = 0 ; r\_power = 0 ; q\_power = 0
00101 \textcolor{keywordflow}{  endif}
00102 
00103   \textcolor{keywordflow}{if} (abs(z\_power) > 300) \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"unit\_scaling\_init: "}//&
00104                  \textcolor{stringliteral}{"Z\_RESCALE\_POWER is outside of the valid range of -300 to 300."})
00105   \textcolor{keywordflow}{if} (abs(l\_power) > 300) \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"unit\_scaling\_init: "}//&
00106                  \textcolor{stringliteral}{"L\_RESCALE\_POWER is outside of the valid range of -300 to 300."})
00107   \textcolor{keywordflow}{if} (abs(t\_power) > 300) \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"unit\_scaling\_init: "}//&
00108                  \textcolor{stringliteral}{"T\_RESCALE\_POWER is outside of the valid range of -300 to 300."})
00109   \textcolor{keywordflow}{if} (abs(r\_power) > 300) \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"unit\_scaling\_init: "}//&
00110                  \textcolor{stringliteral}{"R\_RESCALE\_POWER is outside of the valid range of -300 to 300."})
00111   \textcolor{keywordflow}{if} (abs(q\_power) > 300) \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"unit\_scaling\_init: "}//&
00112                  \textcolor{stringliteral}{"Q\_RESCALE\_POWER is outside of the valid range of -300 to 300."})
00113 
00114   z\_rescale\_factor = 1.0
00115   \textcolor{keywordflow}{if} (z\_power /= 0) z\_rescale\_factor = 2.0**z\_power
00116   us%Z\_to\_m = 1.0 * z\_rescale\_factor
00117   us%m\_to\_Z = 1.0 / z\_rescale\_factor
00118 
00119   l\_rescale\_factor = 1.0
00120   \textcolor{keywordflow}{if} (l\_power /= 0) l\_rescale\_factor = 2.0**l\_power
00121   us%L\_to\_m = 1.0 * l\_rescale\_factor
00122   us%m\_to\_L = 1.0 / l\_rescale\_factor
00123 
00124   t\_rescale\_factor = 1.0
00125   \textcolor{keywordflow}{if} (t\_power /= 0) t\_rescale\_factor = 2.0**t\_power
00126   us%T\_to\_s = 1.0 * t\_rescale\_factor
00127   us%s\_to\_T = 1.0 / t\_rescale\_factor
00128 
00129   r\_rescale\_factor = 1.0
00130   \textcolor{keywordflow}{if} (r\_power /= 0) r\_rescale\_factor = 2.0**r\_power
00131   us%R\_to\_kg\_m3 = 1.0 * r\_rescale\_factor
00132   us%kg\_m3\_to\_R = 1.0 / r\_rescale\_factor
00133 
00134   q\_rescale\_factor = 1.0
00135   \textcolor{keywordflow}{if} (q\_power /= 0) q\_rescale\_factor = 2.0**q\_power
00136   us%Q\_to\_J\_kg = 1.0 * q\_rescale\_factor
00137   us%J\_kg\_to\_Q = 1.0 / q\_rescale\_factor
00138 
00139   \textcolor{comment}{! These are useful combinations of the fundamental scale conversion factors set above.}
00140   us%Z\_to\_L = us%Z\_to\_m * us%m\_to\_L
00141   us%L\_to\_Z = us%L\_to\_m * us%m\_to\_Z
00142   \textcolor{comment}{! Horizontal velocities:}
00143   us%L\_T\_to\_m\_s = us%L\_to\_m * us%s\_to\_T
00144   us%m\_s\_to\_L\_T = us%m\_to\_L * us%T\_to\_s
00145   \textcolor{comment}{! Horizontal accelerations:}
00146   us%L\_T2\_to\_m\_s2 = us%L\_to\_m * us%s\_to\_T**2
00147     \textcolor{comment}{! It does not look like US%m\_s2\_to\_L\_T2 would be used, so it does not exist.}
00148   \textcolor{comment}{! Vertical diffusivities and viscosities:}
00149   us%Z2\_T\_to\_m2\_s = us%Z\_to\_m**2 * us%s\_to\_T
00150   us%m2\_s\_to\_Z2\_T = us%m\_to\_Z**2 * us%T\_to\_s
00151   \textcolor{comment}{! Column mass loads:}
00152   us%RZ\_to\_kg\_m2  = us%R\_to\_kg\_m3 * us%Z\_to\_m
00153     \textcolor{comment}{! It does not seem like US%kg\_m2\_to\_RZ would be used enough in MOM6 to justify its existence.}
00154   \textcolor{comment}{! Vertical mass fluxes:}
00155   us%kg\_m2s\_to\_RZ\_T = us%kg\_m3\_to\_R * us%m\_to\_Z * us%T\_to\_s
00156   us%RZ\_T\_to\_kg\_m2s = us%R\_to\_kg\_m3 * us%Z\_to\_m * us%s\_to\_T
00157   \textcolor{comment}{! Turbulent kinetic energy vertical fluxes:}
00158   us%RZ3\_T3\_to\_W\_m2 = us%R\_to\_kg\_m3 * us%Z\_to\_m**3 * us%s\_to\_T**3
00159   us%W\_m2\_to\_RZ3\_T3 = us%kg\_m3\_to\_R * us%m\_to\_Z**3 * us%T\_to\_s**3
00160   \textcolor{comment}{! Vertical heat fluxes:}
00161   us%W\_m2\_to\_QRZ\_T = us%J\_kg\_to\_Q * us%kg\_m3\_to\_R * us%m\_to\_Z * us%T\_to\_s
00162   us%QRZ\_T\_to\_W\_m2 = us%Q\_to\_J\_kg * us%R\_to\_kg\_m3 * us%Z\_to\_m * us%s\_to\_T
00163   \textcolor{comment}{! Pressures:}
00164   us%RL2\_T2\_to\_Pa = us%R\_to\_kg\_m3 * us%L\_T\_to\_m\_s**2
00165     \textcolor{comment}{! It does not seem like US%Pa\_to\_RL2\_T2 would be used enough in MOM6 to justify its existence.}
00166   \textcolor{comment}{! US%Pa\_to\_RL2\_T2 = US%kg\_m3\_to\_R * US%m\_s\_to\_L\_T**2}
00167 
00168 \textcolor{keyword}{end subroutine }\hyperlink{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}{unit\_scaling\_init}
00169 \textcolor{comment}{}
00170 \textcolor{comment}{!> Set the unit scaling factors for output to restart files to the unit scaling}
00171 \textcolor{comment}{!! factors for this run.}
00172 \textcolor{keyword}{subroutine }\hyperlink{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}{fix\_restart\_unit\_scaling}(US)
\Hypertarget{MOM__unit__scaling_8F90_source_l00173}\hyperlink{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}{00173}   \textcolor{keywordtype}{type}(\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\_scale\_type}), \textcolor{keywordtype}{intent(inout)} :: US\textcolor{comment}{ !< A dimensional unit scaling type}
00174 
00175   us%m\_to\_Z\_restart = us%m\_to\_Z
00176   us%m\_to\_L\_restart = us%m\_to\_L
00177   us%s\_to\_T\_restart = us%s\_to\_T
00178   us%kg\_m3\_to\_R\_restart = us%kg\_m3\_to\_R
00179   us%J\_kg\_to\_Q\_restart = us%J\_kg\_to\_Q
00180 
00181 \textcolor{keyword}{end subroutine }\hyperlink{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}{fix\_restart\_unit\_scaling}
00182 \textcolor{comment}{}
00183 \textcolor{comment}{!> Deallocates a unit scaling structure.}
00184 \textcolor{keyword}{subroutine }\hyperlink{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}{unit\_scaling\_end}( US )
\Hypertarget{MOM__unit__scaling_8F90_source_l00185}\hyperlink{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}{00185}   \textcolor{keywordtype}{type}(\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\_scale\_type}), \textcolor{keywordtype}{pointer} :: US\textcolor{comment}{ !< A dimensional unit scaling type}
00186 
00187   \textcolor{keyword}{deallocate}( us )
00188 
00189 \textcolor{keyword}{end subroutine }\hyperlink{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}{unit\_scaling\_end}
00190 
00191 \textcolor{keyword}{end module }\hyperlink{namespacemom__unit__scaling}{mom\_unit\_scaling}
\end{DoxyCode}
