\hypertarget{namespacemom__interface__heights}{}\section{mom\+\_\+interface\+\_\+heights Module Reference}
\label{namespacemom__interface__heights}\index{mom\+\_\+interface\+\_\+heights@{mom\+\_\+interface\+\_\+heights}}


\subsection{Detailed Description}
Functions for calculating interface heights, including free surface height. \subsection*{Data Types}
\begin{DoxyCompactItemize}
\item 
interface \mbox{\hyperlink{interfacemom__interface__heights_1_1find__eta}{find\+\_\+eta}}
\begin{DoxyCompactList}\small\item\em Calculates the heights of sruface or all interfaces from layer thicknesses. \end{DoxyCompactList}\end{DoxyCompactItemize}
\subsection*{Functions/\+Subroutines}
\begin{DoxyCompactItemize}
\item 
subroutine \mbox{\hyperlink{namespacemom__interface__heights_a534ce6d17993d55a7e2b108198bfcfd4}{find\+\_\+eta\+\_\+3d}} (h, tv, G, GV, US, eta, eta\+\_\+bt, halo\+\_\+size, eta\+\_\+to\+\_\+m)
\begin{DoxyCompactList}\small\item\em Calculates the heights of all interfaces between layers, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, these height may be dilated for consistency with the corresponding time-\/average quantity from the barotropic calculation. \end{DoxyCompactList}\item 
subroutine \mbox{\hyperlink{namespacemom__interface__heights_a0262e4e5e2d7bd849a2117f91b644ea1}{find\+\_\+eta\+\_\+2d}} (h, tv, G, GV, US, eta, eta\+\_\+bt, halo\+\_\+size, eta\+\_\+to\+\_\+m)
\begin{DoxyCompactList}\small\item\em Calculates the free surface height, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, the sea surface height may be adjusted for consistency with the corresponding time-\/average quantity from the barotropic calculation. \end{DoxyCompactList}\end{DoxyCompactItemize}


\subsection{Function/\+Subroutine Documentation}
\mbox{\Hypertarget{namespacemom__interface__heights_a0262e4e5e2d7bd849a2117f91b644ea1}\label{namespacemom__interface__heights_a0262e4e5e2d7bd849a2117f91b644ea1}} 
\index{mom\+\_\+interface\+\_\+heights@{mom\+\_\+interface\+\_\+heights}!find\+\_\+eta\+\_\+2d@{find\+\_\+eta\+\_\+2d}}
\index{find\+\_\+eta\+\_\+2d@{find\+\_\+eta\+\_\+2d}!mom\+\_\+interface\+\_\+heights@{mom\+\_\+interface\+\_\+heights}}
\subsubsection{\texorpdfstring{find\+\_\+eta\+\_\+2d()}{find\_eta\_2d()}}
{\footnotesize\ttfamily subroutine mom\+\_\+interface\+\_\+heights\+::find\+\_\+eta\+\_\+2d (\begin{DoxyParamCaption}\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed, g \%ke), intent(in)}]{h,  }\item[{type(thermo\+\_\+var\+\_\+ptrs), intent(in)}]{tv,  }\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed), intent(out)}]{eta,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed), intent(in), optional}]{eta\+\_\+bt,  }\item[{integer, intent(in), optional}]{halo\+\_\+size,  }\item[{real, intent(in), optional}]{eta\+\_\+to\+\_\+m }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



Calculates the free surface height, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, the sea surface height may be adjusted for consistency with the corresponding time-\/average quantity from the barotropic calculation. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
\mbox{\tt in}  & {\em h} & Layer thicknesses \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em tv} & A structure pointing to various thermodynamic variables.\\
\hline
\mbox{\tt out}  & {\em eta} & free surface height relative to mean sea level (z=0) often \mbox{[}Z $\sim$$>$ m\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em eta\+\_\+bt} & optional barotropic variable that gives the \char`\"{}correct\char`\"{} free surface height (Boussinesq) or total water column mass per unit area (non-\/\+Boussinesq) \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em halo\+\_\+size} & width of halo points on which to calculate eta.\\
\hline
\mbox{\tt in}  & {\em eta\+\_\+to\+\_\+m} & The conversion factor from the units of eta to m; by default this is USZ\+\_\+to\+\_\+m. \\
\hline
\end{DoxyParams}


Definition at line 149 of file M\+O\+M\+\_\+interface\+\_\+heights.\+F90.


\begin{DoxyCode}
149   \textcolor{keywordtype}{type}(ocean\_grid\_type),                      \textcolor{keywordtype}{intent(in)}  :: G\textcolor{comment}{   !< The ocean's grid structure.}
150   \textcolor{keywordtype}{type}(verticalGrid\_type),                    \textcolor{keywordtype}{intent(in)}  :: GV\textcolor{comment}{  !< The ocean's vertical grid structure.}
151   \textcolor{keywordtype}{type}(unit\_scale\_type),                      \textcolor{keywordtype}{intent(in)}  :: US\textcolor{comment}{  !< A dimensional unit scaling type}
152   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G))},   \textcolor{keywordtype}{intent(in)}  :: h\textcolor{comment}{   !< Layer thicknesses [H ~> m or kg m-2]}
153   \textcolor{keywordtype}{type}(thermo\_var\_ptrs),                      \textcolor{keywordtype}{intent(in)}  :: tv\textcolor{comment}{  !< A structure pointing to various}
154 \textcolor{comment}{                                                                 !! thermodynamic variables.}
155   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))},           \textcolor{keywordtype}{intent(out)} :: eta\textcolor{comment}{ !< free surface height relative to}
156 \textcolor{comment}{                                                                 !! mean sea level (z=0) often [Z ~> m].}
157   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: eta\_bt\textcolor{comment}{ !< optional barotropic}
158 \textcolor{comment}{                   !! variable that gives the "correct" free surface height (Boussinesq) or total}
159 \textcolor{comment}{                   !! water column mass per unit area (non-Boussinesq) [H ~> m or kg m-2].}
160   \textcolor{keywordtype}{integer},                          \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: halo\_size\textcolor{comment}{ !< width of halo points on}
161 \textcolor{comment}{                                                                 !! which to calculate eta.}
162   \textcolor{keywordtype}{real},                             \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: eta\_to\_m\textcolor{comment}{  !< The conversion factor from}
163 \textcolor{comment}{             !! the units of eta to m; by default this is US%Z\_to\_m.}
164   \textcolor{comment}{! Local variables}
165   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G)+1)} :: &
166     p          \textcolor{comment}{! Hydrostatic pressure at each interface [R L2 T-2 ~> Pa]}
167   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G))} :: &
168     dz\_geo     \textcolor{comment}{! The change in geopotential height across a layer [L2 T-2 ~> m2 s-2].}
169   \textcolor{keywordtype}{real} :: htot(SZI\_(G))  \textcolor{comment}{! The sum of all layers' thicknesses [H ~> m or kg m-2].}
170   \textcolor{keywordtype}{real} :: I\_gEarth          \textcolor{comment}{! The inverse of the gravitational acceleration times the}
171                             \textcolor{comment}{! rescaling factor derived from eta\_to\_m [T2 Z L-2 ~> s2 m-1]}
172   \textcolor{keywordtype}{real} :: Z\_to\_eta, H\_to\_eta, H\_to\_rho\_eta \textcolor{comment}{! Unit conversion factors with obvious names.}
173   \textcolor{keywordtype}{integer} i, j, k, is, ie, js, je, nz, halo
174 
175   halo = 0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(halo\_size)) halo = max(0,halo\_size)
176   is = g%isc-halo ; ie = g%iec+halo ; js = g%jsc-halo ; je = g%jec+halo
177   nz = g%ke
178 
179   z\_to\_eta = 1.0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eta\_to\_m)) z\_to\_eta = us%Z\_to\_m / eta\_to\_m
180   h\_to\_eta = gv%H\_to\_Z * z\_to\_eta
181   h\_to\_rho\_eta =  gv%H\_to\_RZ * z\_to\_eta
182   i\_gearth = z\_to\_eta / gv%g\_Earth
183 
184 \textcolor{comment}{!$OMP parallel default(shared) private(htot)}
185 \textcolor{comment}{!$OMP do}
186   \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; eta(i,j) = -z\_to\_eta*g%bathyT(i,j) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
187 
188   \textcolor{keywordflow}{if} (gv%Boussinesq) \textcolor{keywordflow}{then}
189     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eta\_bt)) \textcolor{keywordflow}{then}
190 \textcolor{comment}{!$OMP do}
191       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
192         eta(i,j) = h\_to\_eta*eta\_bt(i,j)
193 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
194     \textcolor{keywordflow}{else}
195 \textcolor{comment}{!$OMP do}
196       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie
197         eta(i,j) = eta(i,j) + h(i,j,k)*h\_to\_eta
198 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
199 \textcolor{keywordflow}{    endif}
200   \textcolor{keywordflow}{else}
201     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(tv%eqn\_of\_state)) \textcolor{keywordflow}{then}
202 \textcolor{comment}{!$OMP do}
203       \textcolor{keywordflow}{do} j=js,je
204         \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(tv%p\_surf)) \textcolor{keywordflow}{then}
205           \textcolor{keywordflow}{do} i=is,ie ; p(i,j,1) = tv%p\_surf(i,j) ;\textcolor{keywordflow}{ enddo}
206         \textcolor{keywordflow}{else}
207           \textcolor{keywordflow}{do} i=is,ie ; p(i,j,1) = 0.0 ;\textcolor{keywordflow}{ enddo}
208 \textcolor{keywordflow}{        endif}
209 
210         \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie
211           p(i,j,k+1) = p(i,j,k) + gv%g\_Earth*gv%H\_to\_RZ*h(i,j,k)
212 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
213 \textcolor{keywordflow}{      enddo}
214 \textcolor{comment}{!$OMP do}
215       \textcolor{keywordflow}{do} k = 1, nz
216         \textcolor{keyword}{call }int\_specific\_vol\_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), 0.0, &
217                                  g%HI, tv%eqn\_of\_state, us, dz\_geo(:,:,k), halo\_size=halo)
218 \textcolor{keywordflow}{      enddo}
219 \textcolor{comment}{!$OMP do}
220       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie
221           eta(i,j) = eta(i,j) + i\_gearth * dz\_geo(i,j,k)
222 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
223     \textcolor{keywordflow}{else}
224 \textcolor{comment}{!$OMP do}
225       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie
226         eta(i,j) = eta(i,j) + h\_to\_rho\_eta*h(i,j,k) / gv%Rlay(k)
227 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
228 \textcolor{keywordflow}{    endif}
229     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eta\_bt)) \textcolor{keywordflow}{then}
230       \textcolor{comment}{!   Dilate the water column to agree with the time-averaged column}
231       \textcolor{comment}{! mass from the barotropic solution.}
232 \textcolor{comment}{!$OMP do}
233       \textcolor{keywordflow}{do} j=js,je
234         \textcolor{keywordflow}{do} i=is,ie ; htot(i) = gv%H\_subroundoff ;\textcolor{keywordflow}{ enddo}
235         \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie ; htot(i) = htot(i) + h(i,j,k) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
236         \textcolor{keywordflow}{do} i=is,ie
237           eta(i,j) = (eta\_bt(i,j) / htot(i)) * (eta(i,j) + z\_to\_eta*g%bathyT(i,j)) - &
238                      z\_to\_eta*g%bathyT(i,j)
239 \textcolor{keywordflow}{        enddo}
240 \textcolor{keywordflow}{      enddo}
241 \textcolor{keywordflow}{    endif}
242 \textcolor{keywordflow}{  endif}
243 \textcolor{comment}{!$OMP end parallel}
244 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__interface__heights_a534ce6d17993d55a7e2b108198bfcfd4}\label{namespacemom__interface__heights_a534ce6d17993d55a7e2b108198bfcfd4}} 
\index{mom\+\_\+interface\+\_\+heights@{mom\+\_\+interface\+\_\+heights}!find\+\_\+eta\+\_\+3d@{find\+\_\+eta\+\_\+3d}}
\index{find\+\_\+eta\+\_\+3d@{find\+\_\+eta\+\_\+3d}!mom\+\_\+interface\+\_\+heights@{mom\+\_\+interface\+\_\+heights}}
\subsubsection{\texorpdfstring{find\+\_\+eta\+\_\+3d()}{find\_eta\_3d()}}
{\footnotesize\ttfamily subroutine mom\+\_\+interface\+\_\+heights\+::find\+\_\+eta\+\_\+3d (\begin{DoxyParamCaption}\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed, g \%ke), intent(in)}]{h,  }\item[{type(thermo\+\_\+var\+\_\+ptrs), intent(in)}]{tv,  }\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed, g \%ke+1), intent(out)}]{eta,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed), intent(in), optional}]{eta\+\_\+bt,  }\item[{integer, intent(in), optional}]{halo\+\_\+size,  }\item[{real, intent(in), optional}]{eta\+\_\+to\+\_\+m }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



Calculates the heights of all interfaces between layers, using the appropriate form for consistency with the calculation of the pressure gradient forces. Additionally, these height may be dilated for consistency with the corresponding time-\/average quantity from the barotropic calculation. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
\mbox{\tt in}  & {\em h} & Layer thicknesses \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em tv} & A structure pointing to various thermodynamic variables.\\
\hline
\mbox{\tt out}  & {\em eta} & layer interface heights \mbox{[}Z $\sim$$>$ m\mbox{]} or 1/eta\+\_\+to\+\_\+m m).\\
\hline
\mbox{\tt in}  & {\em eta\+\_\+bt} & optional barotropic variable that gives the \char`\"{}correct\char`\"{} free surface height (Boussinesq) or total water column mass per unit area (non-\/\+Boussinesq). This is used to dilate the layer. thicknesses when calculating interfaceheights \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em halo\+\_\+size} & width of halo points on which to calculate eta.\\
\hline
\mbox{\tt in}  & {\em eta\+\_\+to\+\_\+m} & The conversion factor from the units of eta to m; by default this is USZ\+\_\+to\+\_\+m. \\
\hline
\end{DoxyParams}


Definition at line 32 of file M\+O\+M\+\_\+interface\+\_\+heights.\+F90.


\begin{DoxyCode}
32   \textcolor{keywordtype}{type}(ocean\_grid\_type),                      \textcolor{keywordtype}{intent(in)}  :: G\textcolor{comment}{   !< The ocean's grid structure.}
33   \textcolor{keywordtype}{type}(verticalGrid\_type),                    \textcolor{keywordtype}{intent(in)}  :: GV\textcolor{comment}{  !< The ocean's vertical grid structure.}
34   \textcolor{keywordtype}{type}(unit\_scale\_type),                      \textcolor{keywordtype}{intent(in)}  :: US\textcolor{comment}{  !< A dimensional unit scaling type}
35   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G))},   \textcolor{keywordtype}{intent(in)}  :: h\textcolor{comment}{   !< Layer thicknesses [H ~> m or kg m-2]}
36   \textcolor{keywordtype}{type}(thermo\_var\_ptrs),                      \textcolor{keywordtype}{intent(in)}  :: tv\textcolor{comment}{  !< A structure pointing to various}
37 \textcolor{comment}{                                                                 !! thermodynamic variables.}
38   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G)+1)}, \textcolor{keywordtype}{intent(out)} :: eta\textcolor{comment}{ !< layer interface heights}
39 \textcolor{comment}{                                                                 !! [Z ~> m] or 1/eta\_to\_m m).}
40   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: eta\_bt\textcolor{comment}{ !< optional barotropic}
41 \textcolor{comment}{             !! variable that gives the "correct" free surface height (Boussinesq) or total water}
42 \textcolor{comment}{             !! column mass per unit area (non-Boussinesq).  This is used to dilate the layer.}
43 \textcolor{comment}{             !! thicknesses when calculating interfaceheights [H ~> m or kg m-2].}
44   \textcolor{keywordtype}{integer},                          \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: halo\_size\textcolor{comment}{ !< width of halo points on}
45 \textcolor{comment}{                                                                 !! which to calculate eta.}
46   \textcolor{keywordtype}{real},                             \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: eta\_to\_m\textcolor{comment}{  !< The conversion factor from}
47 \textcolor{comment}{             !! the units of eta to m; by default this is US%Z\_to\_m.}
48 
49   \textcolor{comment}{! Local variables}
50   \textcolor{keywordtype}{real} :: p(SZI\_(G),SZJ\_(G),SZK\_(G)+1)    \textcolor{comment}{! Hydrostatic pressure at each interface [R L2 T-2 ~> Pa]}
51   \textcolor{keywordtype}{real} :: dz\_geo(SZI\_(G),SZJ\_(G),SZK\_(G)) \textcolor{comment}{! The change in geopotential height}
52                                           \textcolor{comment}{! across a layer [L2 T-2 ~> m2 s-2].}
53   \textcolor{keywordtype}{real} :: dilate(SZI\_(G))                 \textcolor{comment}{! non-dimensional dilation factor}
54   \textcolor{keywordtype}{real} :: htot(SZI\_(G))                   \textcolor{comment}{! total thickness [H ~> m or kg m-2]}
55   \textcolor{keywordtype}{real} :: I\_gEarth          \textcolor{comment}{! The inverse of the gravitational acceleration times the}
56                             \textcolor{comment}{! rescaling factor derived from eta\_to\_m [T2 Z L-2 ~> s2 m-1]}
57   \textcolor{keywordtype}{real} :: Z\_to\_eta, H\_to\_eta, H\_to\_rho\_eta \textcolor{comment}{! Unit conversion factors with obvious names.}
58   \textcolor{keywordtype}{integer} i, j, k, isv, iev, jsv, jev, nz, halo
59 
60   halo = 0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(halo\_size)) halo = max(0,halo\_size)
61 
62   isv = g%isc-halo ; iev = g%iec+halo ; jsv = g%jsc-halo ; jev = g%jec+halo
63   nz = g%ke
64 
65   \textcolor{keywordflow}{if} ((isv<g%isd) .or. (iev>g%ied) .or. (jsv<g%jsd) .or. (jev>g%jed)) &
66     \textcolor{keyword}{call }mom\_error(fatal,\textcolor{stringliteral}{"find\_eta called with an overly large halo\_size."})
67 
68   z\_to\_eta = 1.0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eta\_to\_m)) z\_to\_eta = us%Z\_to\_m / eta\_to\_m
69   h\_to\_eta = gv%H\_to\_Z * z\_to\_eta
70   h\_to\_rho\_eta =  gv%H\_to\_RZ * z\_to\_eta
71   i\_gearth = z\_to\_eta / gv%g\_Earth
72 
73 \textcolor{comment}{!$OMP parallel default(shared) private(dilate,htot)}
74 \textcolor{comment}{!$OMP do}
75   \textcolor{keywordflow}{do} j=jsv,jev ; \textcolor{keywordflow}{do} i=isv,iev ; eta(i,j,nz+1) = -z\_to\_eta*g%bathyT(i,j) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
76 
77   \textcolor{keywordflow}{if} (gv%Boussinesq) \textcolor{keywordflow}{then}
78 \textcolor{comment}{!$OMP do}
79     \textcolor{keywordflow}{do} j=jsv,jev ; \textcolor{keywordflow}{do} k=nz,1,-1; \textcolor{keywordflow}{do} i=isv,iev
80       eta(i,j,k) = eta(i,j,k+1) + h(i,j,k)*h\_to\_eta
81 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
82     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eta\_bt)) \textcolor{keywordflow}{then}
83       \textcolor{comment}{! Dilate the water column to agree with the free surface height}
84       \textcolor{comment}{! that is used for the dynamics.}
85 \textcolor{comment}{!$OMP do}
86       \textcolor{keywordflow}{do} j=jsv,jev
87         \textcolor{keywordflow}{do} i=isv,iev
88           dilate(i) = (eta\_bt(i,j)*h\_to\_eta + z\_to\_eta*g%bathyT(i,j)) / &
89                       (eta(i,j,1) + z\_to\_eta*g%bathyT(i,j))
90 \textcolor{keywordflow}{        enddo}
91         \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=isv,iev
92           eta(i,j,k) = dilate(i) * (eta(i,j,k) + z\_to\_eta*g%bathyT(i,j)) - z\_to\_eta*g%bathyT(i,j)
93 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
94 \textcolor{keywordflow}{      enddo}
95 \textcolor{keywordflow}{    endif}
96   \textcolor{keywordflow}{else}
97     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(tv%eqn\_of\_state)) \textcolor{keywordflow}{then}
98 \textcolor{comment}{!$OMP do}
99       \textcolor{keywordflow}{do} j=jsv,jev
100         \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(tv%p\_surf)) \textcolor{keywordflow}{then}
101           \textcolor{keywordflow}{do} i=isv,iev ; p(i,j,1) = tv%p\_surf(i,j) ;\textcolor{keywordflow}{ enddo}
102         \textcolor{keywordflow}{else}
103           \textcolor{keywordflow}{do} i=isv,iev ; p(i,j,1) = 0.0 ;\textcolor{keywordflow}{ enddo}
104 \textcolor{keywordflow}{        endif}
105         \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=isv,iev
106           p(i,j,k+1) = p(i,j,k) + gv%g\_Earth*gv%H\_to\_RZ*h(i,j,k)
107 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
108 \textcolor{keywordflow}{      enddo}
109 \textcolor{comment}{!$OMP do}
110       \textcolor{keywordflow}{do} k=1,nz
111         \textcolor{keyword}{call }int\_specific\_vol\_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), &
112                                  0.0, g%HI, tv%eqn\_of\_state, us, dz\_geo(:,:,k), halo\_size=halo)
113 \textcolor{keywordflow}{      enddo}
114 \textcolor{comment}{!$OMP do}
115       \textcolor{keywordflow}{do} j=jsv,jev
116         \textcolor{keywordflow}{do} k=nz,1,-1 ; \textcolor{keywordflow}{do} i=isv,iev
117           eta(i,j,k) = eta(i,j,k+1) + i\_gearth * dz\_geo(i,j,k)
118 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
119 \textcolor{keywordflow}{      enddo}
120     \textcolor{keywordflow}{else}
121 \textcolor{comment}{!$OMP do}
122       \textcolor{keywordflow}{do} j=jsv,jev ;  \textcolor{keywordflow}{do} k=nz,1,-1; \textcolor{keywordflow}{do} i=isv,iev
123         eta(i,j,k) = eta(i,j,k+1) + h\_to\_rho\_eta*h(i,j,k) / gv%Rlay(k)
124 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
125 \textcolor{keywordflow}{    endif}
126     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eta\_bt)) \textcolor{keywordflow}{then}
127       \textcolor{comment}{! Dilate the water column to agree with the free surface height}
128       \textcolor{comment}{! from the time-averaged barotropic solution.}
129 \textcolor{comment}{!$OMP do}
130       \textcolor{keywordflow}{do} j=jsv,jev
131         \textcolor{keywordflow}{do} i=isv,iev ; htot(i) = gv%H\_subroundoff ;\textcolor{keywordflow}{ enddo}
132         \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=isv,iev ; htot(i) = htot(i) + h(i,j,k) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
133         \textcolor{keywordflow}{do} i=isv,iev ; dilate(i) = eta\_bt(i,j) / htot(i) ;\textcolor{keywordflow}{ enddo}
134         \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=isv,iev
135           eta(i,j,k) = dilate(i) * (eta(i,j,k) + z\_to\_eta*g%bathyT(i,j)) - z\_to\_eta*g%bathyT(i,j)
136 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
137 \textcolor{keywordflow}{      enddo}
138 \textcolor{keywordflow}{    endif}
139 \textcolor{keywordflow}{  endif}
140 \textcolor{comment}{!$OMP end parallel}
141 
\end{DoxyCode}
