\hypertarget{namespacemom__unit__scaling}{}\section{mom\+\_\+unit\+\_\+scaling Module Reference}
\label{namespacemom__unit__scaling}\index{mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}}


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


\subsection*{Data Types}
\begin{DoxyCompactItemize}
\item 
type \hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\+\_\+scale\+\_\+type}
\begin{DoxyCompactList}\small\item\em Describes various unit conversion factors. \end{DoxyCompactList}\end{DoxyCompactItemize}
\subsection*{Functions/\+Subroutines}
\begin{DoxyCompactItemize}
\item 
subroutine, public \hyperlink{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}{unit\+\_\+scaling\+\_\+init} (param\+\_\+file, US)
\begin{DoxyCompactList}\small\item\em Allocates and initializes the ocean model unit scaling type. \end{DoxyCompactList}\item 
subroutine, public \hyperlink{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}{fix\+\_\+restart\+\_\+unit\+\_\+scaling} (US)
\begin{DoxyCompactList}\small\item\em Set the unit scaling factors for output to restart files to the unit scaling factors for this run. \end{DoxyCompactList}\item 
subroutine, public \hyperlink{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}{unit\+\_\+scaling\+\_\+end} (US)
\begin{DoxyCompactList}\small\item\em Deallocates a unit scaling structure. \end{DoxyCompactList}\end{DoxyCompactItemize}


\subsection{Detailed Description}
Provides a transparent unit rescaling type to facilitate dimensional consistency testing. 

\subsection{Function/\+Subroutine Documentation}
\mbox{\Hypertarget{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}\label{namespacemom__unit__scaling_a0d99ae286970838e8f4cd534e3a2744c}} 
\index{mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}!fix\+\_\+restart\+\_\+unit\+\_\+scaling@{fix\+\_\+restart\+\_\+unit\+\_\+scaling}}
\index{fix\+\_\+restart\+\_\+unit\+\_\+scaling@{fix\+\_\+restart\+\_\+unit\+\_\+scaling}!mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}}
\subsubsection{\texorpdfstring{fix\+\_\+restart\+\_\+unit\+\_\+scaling()}{fix\_restart\_unit\_scaling()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+unit\+\_\+scaling\+::fix\+\_\+restart\+\_\+unit\+\_\+scaling (\begin{DoxyParamCaption}\item[{type(\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\+\_\+scale\+\_\+type}), intent(inout)}]{US }\end{DoxyParamCaption})}



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


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in,out}  & {\em us} & A dimensional unit scaling type \\
\hline
\end{DoxyParams}


Definition at line \hyperlink{MOM__unit__scaling_8F90_source_l00173}{173} of file \hyperlink{MOM__unit__scaling_8F90_source}{M\+O\+M\+\_\+unit\+\_\+scaling.\+F90}.



Referenced by \hyperlink{namespacemom__ice__shelf_a5990f9918493ff4984245eac74e5f4d9}{mom\+\_\+ice\+\_\+shelf\+::initialize\+\_\+ice\+\_\+shelf()}.


\begin{DoxyCode}
00173   \textcolor{keywordtype}{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 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}\label{namespacemom__unit__scaling_a6b58ce1b6a08d07a84da1257cd8e8694}} 
\index{mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}!unit\+\_\+scaling\+\_\+end@{unit\+\_\+scaling\+\_\+end}}
\index{unit\+\_\+scaling\+\_\+end@{unit\+\_\+scaling\+\_\+end}!mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}}
\subsubsection{\texorpdfstring{unit\+\_\+scaling\+\_\+end()}{unit\_scaling\_end()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+unit\+\_\+scaling\+::unit\+\_\+scaling\+\_\+end (\begin{DoxyParamCaption}\item[{type(\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\+\_\+scale\+\_\+type}), pointer}]{US }\end{DoxyParamCaption})}



Deallocates a unit scaling structure. 


\begin{DoxyParams}{Parameters}
{\em us} & A dimensional unit scaling type \\
\hline
\end{DoxyParams}


Definition at line \hyperlink{MOM__unit__scaling_8F90_source_l00185}{185} of file \hyperlink{MOM__unit__scaling_8F90_source}{M\+O\+M\+\_\+unit\+\_\+scaling.\+F90}.


\begin{DoxyCode}
00185   \textcolor{keywordtype}{type}(unit\_scale\_type), \textcolor{keywordtype}{pointer} :: us\textcolor{comment}{ !< A dimensional unit scaling type}
00186 
00187   \textcolor{keyword}{deallocate}( us )
00188 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}\label{namespacemom__unit__scaling_a74867ddf628f93dcee854980e08bbe21}} 
\index{mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}!unit\+\_\+scaling\+\_\+init@{unit\+\_\+scaling\+\_\+init}}
\index{unit\+\_\+scaling\+\_\+init@{unit\+\_\+scaling\+\_\+init}!mom\+\_\+unit\+\_\+scaling@{mom\+\_\+unit\+\_\+scaling}}
\subsubsection{\texorpdfstring{unit\+\_\+scaling\+\_\+init()}{unit\_scaling\_init()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+unit\+\_\+scaling\+::unit\+\_\+scaling\+\_\+init (\begin{DoxyParamCaption}\item[{type(param\+\_\+file\+\_\+type), intent(in), optional}]{param\+\_\+file,  }\item[{type(\hyperlink{structmom__unit__scaling_1_1unit__scale__type}{unit\+\_\+scale\+\_\+type}), optional, pointer}]{US }\end{DoxyParamCaption})}



Allocates and initializes the ocean model unit scaling type. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em param\+\_\+file} & Parameter file handle/type\\
\hline
 & {\em us} & A dimensional unit scaling type \\
\hline
\end{DoxyParams}


Definition at line \hyperlink{MOM__unit__scaling_8F90_source_l00057}{57} of file \hyperlink{MOM__unit__scaling_8F90_source}{M\+O\+M\+\_\+unit\+\_\+scaling.\+F90}.



Referenced by \hyperlink{namespacemom__ice__shelf_a5990f9918493ff4984245eac74e5f4d9}{mom\+\_\+ice\+\_\+shelf\+::initialize\+\_\+ice\+\_\+shelf()}.


\begin{DoxyCode}
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}(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 
\end{DoxyCode}
