\hypertarget{namespacemom__coriolisadv}{}\section{mom\+\_\+coriolisadv Module Reference}
\label{namespacemom__coriolisadv}\index{mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}}


\subsection{Detailed Description}
Accelerations due to the Coriolis force and momentum advection. 

This file contains the subroutine that calculates the time derivatives of the velocities due to Coriolis acceleration and momentum advection. This subroutine uses either a vorticity advection scheme from Arakawa and Hsu, Mon. Wea. Rev. 1990, or Sadourny\textquotesingle{}s (J\+AS 1975) energy conserving scheme. Both have been modified to use general orthogonal coordinates as described in Arakawa and Lamb, Mon. Wea. Rev. 1981. Both schemes are second order accurate, and allow for vanishingly small layer thicknesses. The Arakawa and Hsu scheme globally conserves both total energy and potential enstrophy in the limit of nondivergent flow. Sadourny\textquotesingle{}s energy conserving scheme conserves energy if the flow is nondivergent or centered difference thickness fluxes are used.

A small fragment of the grid is shown below\+: \begin{DoxyVerb}    j+1  x ^ x ^ x   At x:  q, CoriolisBu
    j+1  > o > o >   At ^:  v, CAv, vh
    j    x ^ x ^ x   At >:  u, CAu, uh, a, b, c, d
    j    > o > o >   At o:  h, KE
    j-1  x ^ x ^ x
        i-1  i  i+1  At x & ^:
           i  i+1    At > & o:\end{DoxyVerb}


The boundaries always run through q grid points (x). \subsection*{Data Types}
\begin{DoxyCompactItemize}
\item 
type \hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}
\begin{DoxyCompactList}\small\item\em Control structure for \hyperlink{namespacemom__coriolisadv}{mom\+\_\+coriolisadv}. \end{DoxyCompactList}\end{DoxyCompactItemize}
\subsection*{Functions/\+Subroutines}
\begin{DoxyCompactItemize}
\item 
subroutine, public \hyperlink{namespacemom__coriolisadv_ac677e9d644c881b7e8ce6413aa5450cd}{coradcalc} (u, v, h, uh, vh, C\+Au, C\+Av, O\+BC, AD, G, GV, US, CS)
\begin{DoxyCompactList}\small\item\em Calculates the Coriolis and momentum advection contributions to the acceleration. \end{DoxyCompactList}\item 
subroutine \hyperlink{namespacemom__coriolisadv_a87e4a437552052fa238260442af19868}{gradke} (u, v, h, KE, K\+Ex, K\+Ey, k, O\+BC, G, US, CS)
\begin{DoxyCompactList}\small\item\em Calculates the acceleration due to the gradient of kinetic energy. \end{DoxyCompactList}\item 
subroutine, public \hyperlink{namespacemom__coriolisadv_ae021ac8de3b3510ca4552314ec9e1a9a}{coriolisadv\+\_\+init} (Time, G, GV, US, param\+\_\+file, diag, AD, CS)
\begin{DoxyCompactList}\small\item\em Initializes the control structure for \hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}. \end{DoxyCompactList}\item 
subroutine, public \hyperlink{namespacemom__coriolisadv_a6252eaea90947c83b5a1900d31191b96}{coriolisadv\+\_\+end} (CS)
\begin{DoxyCompactList}\small\item\em Destructor for \hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}. \end{DoxyCompactList}\end{DoxyCompactItemize}
\subsection*{Variables}
\textbf{ }\par
\begin{DoxyCompactItemize}
\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a147b51b8c90695f067981c088fd97a7d}\label{namespacemom__coriolisadv_a147b51b8c90695f067981c088fd97a7d}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_a147b51b8c90695f067981c088fd97a7d}{sadourny75\+\_\+energy} = 1
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_ac6ce86a064b1931b11b1f0cab8e98ae9}\label{namespacemom__coriolisadv_ac6ce86a064b1931b11b1f0cab8e98ae9}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_ac6ce86a064b1931b11b1f0cab8e98ae9}{arakawa\+\_\+hsu90} = 2
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a7e3de4f6d1439feb52cad07e9321d531}\label{namespacemom__coriolisadv_a7e3de4f6d1439feb52cad07e9321d531}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_a7e3de4f6d1439feb52cad07e9321d531}{robust\+\_\+enstro} = 3
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_aa666114e7e3d763580cd2fb9e98a327e}\label{namespacemom__coriolisadv_aa666114e7e3d763580cd2fb9e98a327e}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_aa666114e7e3d763580cd2fb9e98a327e}{sadourny75\+\_\+enstro} = 4
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a668bd901852bbdee753d57f5a89418d2}\label{namespacemom__coriolisadv_a668bd901852bbdee753d57f5a89418d2}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_a668bd901852bbdee753d57f5a89418d2}{arakawa\+\_\+lamb81} = 5
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_ad7bda9fd606377eaef6bd65cd3739746}\label{namespacemom__coriolisadv_ad7bda9fd606377eaef6bd65cd3739746}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_ad7bda9fd606377eaef6bd65cd3739746}{al\+\_\+blend} = 6
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a1cc8656f19fd69ba7d46a1ef9f14dffe}\label{namespacemom__coriolisadv_a1cc8656f19fd69ba7d46a1ef9f14dffe}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a1cc8656f19fd69ba7d46a1ef9f14dffe}{sadourny75\+\_\+energy\+\_\+string} = \char`\"{}S\+A\+D\+O\+U\+R\+N\+Y75\+\_\+\+E\+N\+E\+R\+GY\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a83331235df2dc1d04ce4dcacc98ceac8}\label{namespacemom__coriolisadv_a83331235df2dc1d04ce4dcacc98ceac8}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a83331235df2dc1d04ce4dcacc98ceac8}{arakawa\+\_\+hsu\+\_\+string} = \char`\"{}A\+R\+A\+K\+A\+W\+A\+\_\+\+H\+S\+U90\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a9d8edb5a8f791d163be42e83e64318f0}\label{namespacemom__coriolisadv_a9d8edb5a8f791d163be42e83e64318f0}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a9d8edb5a8f791d163be42e83e64318f0}{robust\+\_\+enstro\+\_\+string} = \char`\"{}R\+O\+B\+U\+S\+T\+\_\+\+E\+N\+S\+T\+RO\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_ac5c369c1180e5f42d6e72abf129be7b3}\label{namespacemom__coriolisadv_ac5c369c1180e5f42d6e72abf129be7b3}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_ac5c369c1180e5f42d6e72abf129be7b3}{sadourny75\+\_\+enstro\+\_\+string} = \char`\"{}S\+A\+D\+O\+U\+R\+N\+Y75\+\_\+\+E\+N\+S\+T\+RO\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a72266cbde19f201a04aee9e657794e32}\label{namespacemom__coriolisadv_a72266cbde19f201a04aee9e657794e32}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a72266cbde19f201a04aee9e657794e32}{arakawa\+\_\+lamb\+\_\+string} = \char`\"{}A\+R\+A\+K\+A\+W\+A\+\_\+\+L\+A\+M\+B81\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a136855c10a40b48669032a7e46718302}\label{namespacemom__coriolisadv_a136855c10a40b48669032a7e46718302}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a136855c10a40b48669032a7e46718302}{al\+\_\+blend\+\_\+string} = \char`\"{}A\+R\+A\+K\+A\+W\+A\+\_\+\+L\+A\+M\+B\+\_\+\+B\+L\+E\+ND\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for Coriolis\+\_\+\+Scheme. \end{DoxyCompactList}\end{DoxyCompactItemize}

\textbf{ }\par
\begin{DoxyCompactItemize}
\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_ad4cf1614860611458fcbbf6f460c0b86}\label{namespacemom__coriolisadv_ad4cf1614860611458fcbbf6f460c0b86}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_ad4cf1614860611458fcbbf6f460c0b86}{ke\+\_\+arakawa} = 10
\begin{DoxyCompactList}\small\item\em Enumeration values for K\+E\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a7303f7b6bb505f800e828f6df8d31180}\label{namespacemom__coriolisadv_a7303f7b6bb505f800e828f6df8d31180}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_a7303f7b6bb505f800e828f6df8d31180}{ke\+\_\+simple\+\_\+gudonov} = 11
\begin{DoxyCompactList}\small\item\em Enumeration values for K\+E\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a3f4f30b7016d655cb6c247040b6fdb29}\label{namespacemom__coriolisadv_a3f4f30b7016d655cb6c247040b6fdb29}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_a3f4f30b7016d655cb6c247040b6fdb29}{ke\+\_\+gudonov} = 12
\begin{DoxyCompactList}\small\item\em Enumeration values for K\+E\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a90cbcf11d9b2bf9f2693b98bd061e959}\label{namespacemom__coriolisadv_a90cbcf11d9b2bf9f2693b98bd061e959}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a90cbcf11d9b2bf9f2693b98bd061e959}{ke\+\_\+arakawa\+\_\+string} = \char`\"{}K\+E\+\_\+\+A\+R\+A\+K\+A\+WA\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for K\+E\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_add83990bb15763a8db3134e73a999f02}\label{namespacemom__coriolisadv_add83990bb15763a8db3134e73a999f02}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_add83990bb15763a8db3134e73a999f02}{ke\+\_\+simple\+\_\+gudonov\+\_\+string} = \char`\"{}K\+E\+\_\+\+S\+I\+M\+P\+L\+E\+\_\+\+G\+U\+D\+O\+N\+OV\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for K\+E\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_acc382668558ee525d20e0e0deae4e1c2}\label{namespacemom__coriolisadv_acc382668558ee525d20e0e0deae4e1c2}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_acc382668558ee525d20e0e0deae4e1c2}{ke\+\_\+gudonov\+\_\+string} = \char`\"{}K\+E\+\_\+\+G\+U\+D\+O\+N\+OV\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for K\+E\+\_\+\+Scheme. \end{DoxyCompactList}\end{DoxyCompactItemize}

\textbf{ }\par
\begin{DoxyCompactItemize}
\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_ae3a6b2ec6bbd73efc51f89dfdccbe375}\label{namespacemom__coriolisadv_ae3a6b2ec6bbd73efc51f89dfdccbe375}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_ae3a6b2ec6bbd73efc51f89dfdccbe375}{pv\+\_\+adv\+\_\+centered} = 21
\begin{DoxyCompactList}\small\item\em Enumeration values for P\+V\+\_\+\+Adv\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a5bf435f8a66d139afe2446ef13e7fdc5}\label{namespacemom__coriolisadv_a5bf435f8a66d139afe2446ef13e7fdc5}} 
integer, parameter \hyperlink{namespacemom__coriolisadv_a5bf435f8a66d139afe2446ef13e7fdc5}{pv\+\_\+adv\+\_\+upwind1} = 22
\begin{DoxyCompactList}\small\item\em Enumeration values for P\+V\+\_\+\+Adv\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_a5f70330ffdb45608e567a6bcc748c0e5}\label{namespacemom__coriolisadv_a5f70330ffdb45608e567a6bcc748c0e5}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_a5f70330ffdb45608e567a6bcc748c0e5}{pv\+\_\+adv\+\_\+centered\+\_\+string} = \char`\"{}P\+V\+\_\+\+A\+D\+V\+\_\+\+C\+E\+N\+T\+E\+R\+ED\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for P\+V\+\_\+\+Adv\+\_\+\+Scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__coriolisadv_ae07c9a121db5179bd72fbcd940a09c59}\label{namespacemom__coriolisadv_ae07c9a121db5179bd72fbcd940a09c59}} 
character $\ast$(20), parameter \hyperlink{namespacemom__coriolisadv_ae07c9a121db5179bd72fbcd940a09c59}{pv\+\_\+adv\+\_\+upwind1\+\_\+string} = \char`\"{}P\+V\+\_\+\+A\+D\+V\+\_\+\+U\+P\+W\+I\+N\+D1\char`\"{}
\begin{DoxyCompactList}\small\item\em Enumeration values for P\+V\+\_\+\+Adv\+\_\+\+Scheme. \end{DoxyCompactList}\end{DoxyCompactItemize}



\subsection{Function/\+Subroutine Documentation}
\mbox{\Hypertarget{namespacemom__coriolisadv_ac677e9d644c881b7e8ce6413aa5450cd}\label{namespacemom__coriolisadv_ac677e9d644c881b7e8ce6413aa5450cd}} 
\index{mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}!coradcalc@{coradcalc}}
\index{coradcalc@{coradcalc}!mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}}
\subsubsection{\texorpdfstring{coradcalc()}{coradcalc()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+coriolisadv\+::coradcalc (\begin{DoxyParamCaption}\item[{real, dimension(szib\+\_\+(g),szj\+\_\+(g),szk\+\_\+(g)), intent(in)}]{u,  }\item[{real, dimension(szi\+\_\+(g),szjb\+\_\+(g),szk\+\_\+(g)), intent(in)}]{v,  }\item[{real, dimension(szi\+\_\+(g),szj\+\_\+(g),szk\+\_\+(g)), intent(in)}]{h,  }\item[{real, dimension(szib\+\_\+(g),szj\+\_\+(g),szk\+\_\+(g)), intent(in)}]{uh,  }\item[{real, dimension(szi\+\_\+(g),szjb\+\_\+(g),szk\+\_\+(g)), intent(in)}]{vh,  }\item[{real, dimension(szib\+\_\+(g),szj\+\_\+(g),szk\+\_\+(g)), intent(out)}]{C\+Au,  }\item[{real, dimension(szi\+\_\+(g),szjb\+\_\+(g),szk\+\_\+(g)), intent(out)}]{C\+Av,  }\item[{type(ocean\+\_\+obc\+\_\+type), pointer}]{O\+BC,  }\item[{type(accel\+\_\+diag\+\_\+ptrs), intent(inout)}]{AD,  }\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{type(\hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}), pointer}]{CS }\end{DoxyParamCaption})}



Calculates the Coriolis and momentum advection contributions to the acceleration. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em g} & Ocen grid structure\\
\hline
\mbox{\tt in}  & {\em gv} & Vertical grid structure\\
\hline
\mbox{\tt in}  & {\em u} & Zonal velocity \mbox{[}L T-\/1 $\sim$$>$ m s-\/1\mbox{]}\\
\hline
\mbox{\tt in}  & {\em v} & Meridional velocity \mbox{[}L T-\/1 $\sim$$>$ m s-\/1\mbox{]}\\
\hline
\mbox{\tt in}  & {\em h} & Layer thickness \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em uh} & Zonal transport u$\ast$h$\ast$dy \mbox{[}H L2 T-\/1 $\sim$$>$ m3 s-\/1 or kg s-\/1\mbox{]}\\
\hline
\mbox{\tt in}  & {\em vh} & Meridional transport v$\ast$h$\ast$dx \mbox{[}H L2 T-\/1 $\sim$$>$ m3 s-\/1 or kg s-\/1\mbox{]}\\
\hline
\mbox{\tt out}  & {\em cau} & Zonal acceleration due to Coriolis and momentum advection \mbox{[}L T-\/2 $\sim$$>$ m s-\/2\mbox{]}.\\
\hline
\mbox{\tt out}  & {\em cav} & Meridional acceleration due to Coriolis and momentum advection \mbox{[}L T-\/2 $\sim$$>$ m s-\/2\mbox{]}.\\
\hline
 & {\em obc} & Open boundary control structure\\
\hline
\mbox{\tt in,out}  & {\em ad} & Storage for acceleration diagnostics\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
 & {\em cs} & Control structure for M\+O\+M\+\_\+\+Coriolis\+Adv \\
\hline
\end{DoxyParams}


Definition at line 117 of file M\+O\+M\+\_\+\+Coriolis\+Adv.\+F90.


\begin{DoxyCode}
117   \textcolor{keywordtype}{type}(ocean\_grid\_type),                     \textcolor{keywordtype}{intent(in)}    :: g\textcolor{comment}{  !< Ocen grid structure}
118   \textcolor{keywordtype}{type}(verticalgrid\_type),                   \textcolor{keywordtype}{intent(in)}    :: gv\textcolor{comment}{ !< Vertical grid structure}
119   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G),SZK\_(G))}, \textcolor{keywordtype}{intent(in)}    :: u\textcolor{comment}{  !< Zonal velocity [L T-1 ~> m s-1]}
120   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJB\_(G),SZK\_(G))}, \textcolor{keywordtype}{intent(in)}    :: v\textcolor{comment}{  !< Meridional velocity [L T-1 ~> m s-1]}
121   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G))},  \textcolor{keywordtype}{intent(in)}    :: h\textcolor{comment}{  !< Layer thickness [H ~> m or kg m-2]}
122   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G),SZK\_(G))}, \textcolor{keywordtype}{intent(in)}    :: uh\textcolor{comment}{ !< Zonal transport u*h*dy}
123 \textcolor{comment}{                                                                 !! [H L2 T-1 ~> m3 s-1 or kg s-1]}
124   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJB\_(G),SZK\_(G))}, \textcolor{keywordtype}{intent(in)}    :: vh\textcolor{comment}{ !< Meridional transport v*h*dx}
125 \textcolor{comment}{                                                                 !! [H L2 T-1 ~> m3 s-1 or kg s-1]}
126   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G),SZK\_(G))}, \textcolor{keywordtype}{intent(out)}   :: cau\textcolor{comment}{ !< Zonal acceleration due to Coriolis}
127 \textcolor{comment}{                                                                  !! and momentum advection [L T-2 ~> m
       s-2].}
128   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJB\_(G),SZK\_(G))}, \textcolor{keywordtype}{intent(out)}   :: cav\textcolor{comment}{ !< Meridional acceleration due to
       Coriolis}
129 \textcolor{comment}{                                                                  !! and momentum advection [L T-2 ~> m
       s-2].}
130   \textcolor{keywordtype}{type}(ocean\_obc\_type),                      \textcolor{keywordtype}{pointer}       :: obc\textcolor{comment}{ !< Open boundary control structure}
131   \textcolor{keywordtype}{type}(accel\_diag\_ptrs),                     \textcolor{keywordtype}{intent(inout)} :: ad\textcolor{comment}{  !< Storage for acceleration diagnostics}
132   \textcolor{keywordtype}{type}(unit\_scale\_type),                     \textcolor{keywordtype}{intent(in)}    :: us\textcolor{comment}{  !< A dimensional unit scaling type}
133   \textcolor{keywordtype}{type}(coriolisadv\_cs),                      \textcolor{keywordtype}{pointer}       :: cs\textcolor{comment}{  !< Control structure for MOM\_CoriolisAdv}
134 
135   \textcolor{comment}{! Local variables}
136   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJB\_(G))} :: &
137     q, &        \textcolor{comment}{! Layer potential vorticity [H-1 T-1 ~> m-1 s-1 or m2 kg-1 s-1].}
138     ih\_q, &     \textcolor{comment}{! The inverse of thickness interpolated to q points [H-1 ~> m-1 or m2 kg-1].}
139     area\_q      \textcolor{comment}{! The sum of the ocean areas at the 4 adjacent thickness points [L2 ~> m2].}
140 
141   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G))} :: &
142     a, b, c, d  \textcolor{comment}{! a, b, c, & d are combinations of the potential vorticities}
143                 \textcolor{comment}{! surrounding an h grid point.  At small scales, a = q/4,}
144                 \textcolor{comment}{! b = q/4, etc.  All are in [H-1 T-1 ~> m-1 s-1 or m2 kg-1 s-1],}
145                 \textcolor{comment}{! and use the indexing of the corresponding u point.}
146 
147   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))} :: &
148     area\_h, &   \textcolor{comment}{! The ocean area at h points [L2 ~> m2].  Area\_h is used to find the}
149                 \textcolor{comment}{! average thickness in the denominator of q.  0 for land points.}
150     ke          \textcolor{comment}{! Kinetic energy per unit mass [L2 T-2 ~> m2 s-2], KE = (u^2 + v^2)/2.}
151   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G))} :: &
152     harea\_u, &  \textcolor{comment}{! The cell area weighted thickness interpolated to u points}
153                 \textcolor{comment}{! times the effective areas [H L2 ~> m3 or kg].}
154     kex, &      \textcolor{comment}{! The zonal gradient of Kinetic energy per unit mass [L T-2 ~> m s-2],}
155                 \textcolor{comment}{! KEx = d/dx KE.}
156     uh\_center   \textcolor{comment}{! Transport based on arithmetic mean h at u-points [H L2 T-1 ~> m3 s-1 or kg s-1]}
157   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJB\_(G))} :: &
158     harea\_v, &  \textcolor{comment}{! The cell area weighted thickness interpolated to v points}
159                 \textcolor{comment}{! times the effective areas [H L2 ~> m3 or kg].}
160     key, &      \textcolor{comment}{! The meridonal gradient of Kinetic energy per unit mass [L T-2 ~> m s-2],}
161                 \textcolor{comment}{! KEy = d/dy KE.}
162     vh\_center   \textcolor{comment}{! Transport based on arithmetic mean h at v-points [H L2 T-1 ~> m3 s-1 or kg s-1]}
163   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))} :: &
164     uh\_min, uh\_max, &   \textcolor{comment}{! The smallest and largest estimates of the volume}
165     vh\_min, vh\_max, &   \textcolor{comment}{! fluxes through the faces (i.e. u*h*dy & v*h*dx)}
166                         \textcolor{comment}{! [H L2 T-1 ~> m3 s-1 or kg s-1].}
167     ep\_u, ep\_v  \textcolor{comment}{! Additional pseudo-Coriolis terms in the Arakawa and Lamb}
168                 \textcolor{comment}{! discretization [H-1 s-1 ~> m-1 s-1 or m2 kg-1 s-1].}
169   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJB\_(G))} :: &
170     dvdx, dudy, & \textcolor{comment}{! Contributions to the circulation around q-points [L2 T-1 ~> m2 s-1]}
171     abs\_vort, & \textcolor{comment}{! Absolute vorticity at q-points [T-1 ~> s-1].}
172     q2, &       \textcolor{comment}{! Relative vorticity over thickness [H-1 T-1 ~> m-1 s-1 or m2 kg-1 s-1].}
173     max\_fvq, &  \textcolor{comment}{! The maximum of the adjacent values of (-u) times absolute vorticity [L T-2 ~> m s-2].}
174     min\_fvq, &  \textcolor{comment}{! The minimum of the adjacent values of (-u) times absolute vorticity [L T-2 ~> m s-2].}
175     max\_fuq, &  \textcolor{comment}{! The maximum of the adjacent values of u times absolute vorticity [L T-2 ~> m s-2].}
176     min\_fuq     \textcolor{comment}{! The minimum of the adjacent values of u times absolute vorticity [L T-2 ~> m s-2].}
177   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJB\_(G),SZK\_(G))} :: &
178     pv, &       \textcolor{comment}{! A diagnostic array of the potential vorticities [H-1 T-1 ~> m-1 s-1 or m2 kg-1 s-1].}
179     rv          \textcolor{comment}{! A diagnostic array of the relative vorticities [T-1 ~> s-1].}
180   \textcolor{keywordtype}{real} :: fv1, fv2, fu1, fu2   \textcolor{comment}{! (f+rv)*v or (f+rv)*u [L T-2 ~> m s-2].}
181   \textcolor{keywordtype}{real} :: max\_fv, max\_fu       \textcolor{comment}{! The maximum or minimum of the neighboring Coriolis}
182   \textcolor{keywordtype}{real} :: min\_fv, min\_fu       \textcolor{comment}{! accelerations [L T-2 ~> m s-2], i.e. max(min)\_fu(v)q.}
183 
184   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{parameter} :: c1\_12=1.0/12.0 \textcolor{comment}{! C1\_12 = 1/12}
185   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{parameter} :: c1\_24=1.0/24.0 \textcolor{comment}{! C1\_24 = 1/24}
186   \textcolor{keywordtype}{real} :: absolute\_vorticity     \textcolor{comment}{! Absolute vorticity [T-1 ~> s-1].}
187   \textcolor{keywordtype}{real} :: relative\_vorticity     \textcolor{comment}{! Relative vorticity [T-1 ~> s-1].}
188   \textcolor{keywordtype}{real} :: ih                     \textcolor{comment}{! Inverse of thickness [H-1 ~> m-1 or m2 kg-1].}
189   \textcolor{keywordtype}{real} :: max\_ihq, min\_ihq       \textcolor{comment}{! The maximum and minimum of the nearby Ihq [H-1 ~> m-1 or m2 kg-1].}
190   \textcolor{keywordtype}{real} :: harea\_q                \textcolor{comment}{! The sum of area times thickness of the cells}
191                                  \textcolor{comment}{! surrounding a q point [H L2 ~> m3 or kg].}
192   \textcolor{keywordtype}{real} :: h\_neglect              \textcolor{comment}{! A thickness that is so small it is usually}
193                                  \textcolor{comment}{! lost in roundoff and can be neglected [H ~> m or kg m-2].}
194   \textcolor{keywordtype}{real} :: temp1, temp2           \textcolor{comment}{! Temporary variables [L2 T-2 ~> m2 s-2].}
195   \textcolor{keywordtype}{real} :: eps\_vel                \textcolor{comment}{! A tiny, positive velocity [L T-1 ~> m s-1].}
196 
197   \textcolor{keywordtype}{real} :: uhc, vhc               \textcolor{comment}{! Centered estimates of uh and vh [H L2 T-1 ~> m3 s-1 or kg s-1].}
198   \textcolor{keywordtype}{real} :: uhm, vhm               \textcolor{comment}{! The input estimates of uh and vh [H L2 T-1 ~> m3 s-1 or kg s-1].}
199   \textcolor{keywordtype}{real} :: c1, c2, c3, slope      \textcolor{comment}{! Nondimensional parameters for the Coriolis limiter scheme.}
200 
201   \textcolor{keywordtype}{real} :: fe\_m2         \textcolor{comment}{! Nondimensional temporary variables asssociated with}
202   \textcolor{keywordtype}{real} :: rat\_lin       \textcolor{comment}{! the ARAKAWA\_LAMB\_BLEND scheme.}
203   \textcolor{keywordtype}{real} :: rat\_m1        \textcolor{comment}{! The ratio of the maximum neighboring inverse thickness}
204                         \textcolor{comment}{! to the minimum inverse thickness minus 1. rat\_m1 >= 0.}
205   \textcolor{keywordtype}{real} :: al\_wt         \textcolor{comment}{! The relative weight of the Arakawa & Lamb scheme to the}
206                         \textcolor{comment}{! Arakawa & Hsu scheme, nondimensional between 0 and 1.}
207   \textcolor{keywordtype}{real} :: sad\_wt        \textcolor{comment}{! The relative weight of the Sadourny energy scheme to}
208                         \textcolor{comment}{! the other two with the ARAKAWA\_LAMB\_BLEND scheme,}
209                         \textcolor{comment}{! nondimensional between 0 and 1.}
210 
211   \textcolor{keywordtype}{real} :: heff1, heff2  \textcolor{comment}{! Temporary effective H at U or V points [H ~> m or kg m-2].}
212   \textcolor{keywordtype}{real} :: heff3, heff4  \textcolor{comment}{! Temporary effective H at U or V points [H ~> m or kg m-2].}
213   \textcolor{keywordtype}{real} :: h\_tiny        \textcolor{comment}{! A very small thickness [H ~> m or kg m-2].}
214   \textcolor{keywordtype}{real} :: uheff, vheff  \textcolor{comment}{! More temporary variables [H L2 T-1 ~> m3 s-1 or kg s-1].}
215   \textcolor{keywordtype}{real} :: quheff,qvheff \textcolor{comment}{! More temporary variables [H L2 T-1 s-1 ~> m3 s-2 or kg s-2].}
216   \textcolor{keywordtype}{integer} :: i, j, k, n, is, ie, js, je, isq, ieq, jsq, jeq, nz
217 
218 \textcolor{comment}{! Diagnostics for fractional thickness-weighted terms}
219   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{allocatable}, \textcolor{keywordtype}{dimension(:,:)} :: &
220     hf\_gkeu\_2d, hf\_gkev\_2d, & \textcolor{comment}{! Depth sum of hf\_gKEu, hf\_gKEv [L T-2 ~> m s-2].}
221     hf\_rvxu\_2d, hf\_rvxv\_2d    \textcolor{comment}{! Depth sum of hf\_rvxu, hf\_rvxv [L T-2 ~> m s-2].}
222   \textcolor{comment}{!real, allocatable, dimension(:,:,:) :: &}
223   \textcolor{comment}{!  hf\_gKEu, hf\_gKEv, & ! accel. due to KE gradient x fract. thickness  [L T-2 ~> m s-2].}
224   \textcolor{comment}{!  hf\_rvxu, hf\_rvxv    ! accel. due to RV x fract. thickness [L T-2 ~> m s-2].}
225   \textcolor{comment}{! 3D diagnostics hf\_gKEu etc. are commented because there is no clarity on proper remapping grid option.}
226   \textcolor{comment}{! The code is retained for degugging purposes in the future.}
227 
228 \textcolor{comment}{! To work, the following fields must be set outside of the usual}
229 \textcolor{comment}{! is to ie range before this subroutine is called:}
230 \textcolor{comment}{!   v(is-1:ie+2,js-1:je+1), u(is-1:ie+1,js-1:je+2), h(is-1:ie+2,js-1:je+2),}
231 \textcolor{comment}{!   uh(is-1,ie,js:je+1) and vh(is:ie+1,js-1:je).}
232 
233   \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(cs)) \textcolor{keyword}{call }mom\_error(fatal, &
234          \textcolor{stringliteral}{"MOM\_CoriolisAdv: Module must be initialized before it is used."})
235   is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
236   isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nz = g%ke
237   h\_neglect = gv%H\_subroundoff
238   eps\_vel = 1.0e-10*us%m\_s\_to\_L\_T
239   h\_tiny = gv%Angstrom\_H  \textcolor{comment}{! Perhaps this should be set to h\_neglect instead.}
240 
241   \textcolor{comment}{!$OMP parallel do default(private) shared(Isq,Ieq,Jsq,Jeq,G,Area\_h)}
242   \textcolor{keywordflow}{do} j=jsq-1,jeq+2 ; \textcolor{keywordflow}{do} i=isq-1,ieq+2
243     area\_h(i,j) = g%mask2dT(i,j) * g%areaT(i,j)
244 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo}
245   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(obc)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} n=1,obc%number\_of\_segments
246     \textcolor{keywordflow}{if} (.not. obc%segment(n)%on\_pe) cycle
247     i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
248     \textcolor{keywordflow}{if} (obc%segment(n)%is\_N\_or\_S .and. (j >= jsq-1) .and. (j <= jeq+1)) \textcolor{keywordflow}{then}
249       \textcolor{keywordflow}{do} i = max(isq-1,obc%segment(n)%HI%isd), min(ieq+2,obc%segment(n)%HI%ied)
250         \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_n) \textcolor{keywordflow}{then}
251           area\_h(i,j+1) = area\_h(i,j)
252         \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_S)}
253           area\_h(i,j) = area\_h(i,j+1)
254 \textcolor{keywordflow}{        endif}
255 \textcolor{keywordflow}{      enddo}
256     \textcolor{keywordflow}{elseif} (obc%segment(n)%is\_E\_or\_W .and. (i >= isq-1) .and. (i <= ieq+1)) \textcolor{keywordflow}{then}
257       \textcolor{keywordflow}{do} j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+2,obc%segment(n)%HI%jed)
258         \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_e) \textcolor{keywordflow}{then}
259           area\_h(i+1,j) = area\_h(i,j)
260         \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_W)}
261           area\_h(i,j) = area\_h(i+1,j)
262 \textcolor{keywordflow}{        endif}
263 \textcolor{keywordflow}{      enddo}
264 \textcolor{keywordflow}{    endif}
265 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ endif}
266   \textcolor{comment}{!$OMP parallel do default(private) shared(Isq,Ieq,Jsq,Jeq,G,Area\_h,Area\_q)}
267   \textcolor{keywordflow}{do} j=jsq-1,jeq+1 ; \textcolor{keywordflow}{do} i=isq-1,ieq+1
268     area\_q(i,j) = (area\_h(i,j) + area\_h(i+1,j+1)) + &
269                   (area\_h(i+1,j) + area\_h(i,j+1))
270 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo}
271 
272   \textcolor{comment}{!$OMP parallel do default(private) shared(u,v,h,uh,vh,CAu,CAv,G,CS,AD,Area\_h,Area\_q,&}
273   \textcolor{comment}{!$OMP                        RV,PV,is,ie,js,je,Isq,Ieq,Jsq,Jeq,nz,h\_neglect,h\_tiny,OBC,eps\_vel)}
274   \textcolor{keywordflow}{do} k=1,nz
275 
276     \textcolor{comment}{! Here the second order accurate layer potential vorticities, q,}
277     \textcolor{comment}{! are calculated.  hq is  second order accurate in space.  Relative}
278     \textcolor{comment}{! vorticity is second order accurate everywhere with free slip b.c.s,}
279     \textcolor{comment}{! but only first order accurate at boundaries with no slip b.c.s.}
280     \textcolor{comment}{! First calculate the contributions to the circulation around the q-point.}
281     \textcolor{keywordflow}{do} j=jsq-1,jeq+1 ; \textcolor{keywordflow}{do} i=isq-1,ieq+1
282       dvdx(i,j) = (v(i+1,j,k)*g%dyCv(i+1,j) - v(i,j,k)*g%dyCv(i,j))
283       dudy(i,j) = (u(i,j+1,k)*g%dxCu(i,j+1) - u(i,j,k)*g%dxCu(i,j))
284 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
285     \textcolor{keywordflow}{do} j=jsq-1,jeq+1 ; \textcolor{keywordflow}{do} i=isq-1,ieq+2
286       harea\_v(i,j) = 0.5*(area\_h(i,j) * h(i,j,k) + area\_h(i,j+1) * h(i,j+1,k))
287 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
288     \textcolor{keywordflow}{do} j=jsq-1,jeq+2 ; \textcolor{keywordflow}{do} i=isq-1,ieq+1
289       harea\_u(i,j) = 0.5*(area\_h(i,j) * h(i,j,k) + area\_h(i+1,j) * h(i+1,j,k))
290 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
291     \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis) \textcolor{keywordflow}{then}
292       \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=is-1,ie
293         uh\_center(i,j) = 0.5 * (g%dy\_Cu(i,j) * u(i,j,k)) * (h(i,j,k) + h(i+1,j,k))
294 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
295       \textcolor{keywordflow}{do} j=js-1,je ; \textcolor{keywordflow}{do} i=isq,ieq+1
296         vh\_center(i,j) = 0.5 * (g%dx\_Cv(i,j) * v(i,j,k)) * (h(i,j,k) + h(i,j+1,k))
297 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
298 \textcolor{keywordflow}{    endif}
299 
300     \textcolor{comment}{! Adjust circulation components to relative vorticity and thickness projected onto}
301     \textcolor{comment}{! velocity points on open boundaries.}
302     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(obc)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} n=1,obc%number\_of\_segments
303       \textcolor{keywordflow}{if} (.not. obc%segment(n)%on\_pe) cycle
304       i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
305       \textcolor{keywordflow}{if} (obc%segment(n)%is\_N\_or\_S .and. (j >= jsq-1) .and. (j <= jeq+1)) \textcolor{keywordflow}{then}
306         \textcolor{keywordflow}{if} (obc%zero\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
307           dvdx(i,j) = 0. ; dudy(i,j) = 0.
308 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
309         \textcolor{keywordflow}{if} (obc%freeslip\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
310           dudy(i,j) = 0.
311 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
312         \textcolor{keywordflow}{if} (obc%computed\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
313           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_n) \textcolor{keywordflow}{then}
314             dudy(i,j) = 2.0*(obc%segment(n)%tangential\_vel(i,j,k) - u(i,j,k))*g%dxCu(i,j)
315           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_S)}
316             dudy(i,j) = 2.0*(u(i,j+1,k) - obc%segment(n)%tangential\_vel(i,j,k))*g%dxCu(i,j+1)
317 \textcolor{keywordflow}{          endif}
318 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
319         \textcolor{keywordflow}{if} (obc%specified\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
320           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_n) \textcolor{keywordflow}{then}
321             dudy(i,j) = obc%segment(n)%tangential\_grad(i,j,k)*g%dxCu(i,j)*g%dyBu(i,j)
322           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_S)}
323             dudy(i,j) = obc%segment(n)%tangential\_grad(i,j,k)*g%dxCu(i,j+1)*g%dyBu(i,j)
324 \textcolor{keywordflow}{          endif}
325 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
326 
327         \textcolor{comment}{! Project thicknesses across OBC points with a no-gradient condition.}
328         \textcolor{keywordflow}{do} i = max(isq-1,obc%segment(n)%HI%isd), min(ieq+2,obc%segment(n)%HI%ied)
329           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_n) \textcolor{keywordflow}{then}
330             harea\_v(i,j) = 0.5 * (area\_h(i,j) + area\_h(i,j+1)) * h(i,j,k)
331           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_S)}
332             harea\_v(i,j) = 0.5 * (area\_h(i,j) + area\_h(i,j+1)) * h(i,j+1,k)
333 \textcolor{keywordflow}{          endif}
334 \textcolor{keywordflow}{        enddo}
335 
336         \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis) \textcolor{keywordflow}{then}
337           \textcolor{keywordflow}{do} i = max(isq-1,obc%segment(n)%HI%isd), min(ieq+2,obc%segment(n)%HI%ied)
338             \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_n) \textcolor{keywordflow}{then}
339               vh\_center(i,j) = g%dx\_Cv(i,j) * v(i,j,k) * h(i,j,k)
340             \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_S)}
341               vh\_center(i,j) = g%dx\_Cv(i,j) * v(i,j,k) * h(i,j+1,k)
342 \textcolor{keywordflow}{            endif}
343 \textcolor{keywordflow}{          enddo}
344 \textcolor{keywordflow}{        endif}
345       \textcolor{keywordflow}{elseif} (obc%segment(n)%is\_E\_or\_W .and. (i >= isq-1) .and. (i <= ieq+1)) \textcolor{keywordflow}{then}
346         \textcolor{keywordflow}{if} (obc%zero\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
347           dvdx(i,j) = 0. ; dudy(i,j) = 0.
348 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
349         \textcolor{keywordflow}{if} (obc%freeslip\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
350           dvdx(i,j) = 0.
351 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
352         \textcolor{keywordflow}{if} (obc%computed\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
353           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_e) \textcolor{keywordflow}{then}
354             dvdx(i,j) = 2.0*(obc%segment(n)%tangential\_vel(i,j,k) - v(i,j,k))*g%dyCv(i,j)
355           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_W)}
356             dvdx(i,j) = 2.0*(v(i+1,j,k) - obc%segment(n)%tangential\_vel(i,j,k))*g%dyCv(i+1,j)
357 \textcolor{keywordflow}{          endif}
358 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
359         \textcolor{keywordflow}{if} (obc%specified\_vorticity) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
360           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_e) \textcolor{keywordflow}{then}
361             dvdx(i,j) = obc%segment(n)%tangential\_grad(i,j,k)*g%dyCv(i,j)*g%dxBu(i,j)
362           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_W)}
363             dvdx(i,j) = obc%segment(n)%tangential\_grad(i,j,k)*g%dyCv(i+1,j)*g%dxBu(i,j)
364 \textcolor{keywordflow}{          endif}
365 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ endif}
366 
367         \textcolor{comment}{! Project thicknesses across OBC points with a no-gradient condition.}
368         \textcolor{keywordflow}{do} j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+2,obc%segment(n)%HI%jed)
369           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_e) \textcolor{keywordflow}{then}
370             harea\_u(i,j) = 0.5*(area\_h(i,j) + area\_h(i+1,j)) * h(i,j,k)
371           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_W)}
372             harea\_u(i,j) = 0.5*(area\_h(i,j) + area\_h(i+1,j)) * h(i+1,j,k)
373 \textcolor{keywordflow}{          endif}
374 \textcolor{keywordflow}{        enddo}
375         \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis) \textcolor{keywordflow}{then}
376           \textcolor{keywordflow}{do} j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+2,obc%segment(n)%HI%jed)
377             \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_e) \textcolor{keywordflow}{then}
378               uh\_center(i,j) = g%dy\_Cu(i,j) * u(i,j,k) * h(i,j,k)
379             \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_W)}
380               uh\_center(i,j) = g%dy\_Cu(i,j) * u(i,j,k) * h(i+1,j,k)
381 \textcolor{keywordflow}{            endif}
382 \textcolor{keywordflow}{          enddo}
383 \textcolor{keywordflow}{        endif}
384 \textcolor{keywordflow}{      endif}
385 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ endif}
386 
387     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(obc)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} n=1,obc%number\_of\_segments
388       \textcolor{keywordflow}{if} (.not. obc%segment(n)%on\_pe) cycle
389       \textcolor{comment}{! Now project thicknesses across cell-corner points in the OBCs.  The two}
390       \textcolor{comment}{! projections have to occur in sequence and can not be combined easily.}
391       i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
392       \textcolor{keywordflow}{if} (obc%segment(n)%is\_N\_or\_S .and. (j >= jsq-1) .and. (j <= jeq+1)) \textcolor{keywordflow}{then}
393         \textcolor{keywordflow}{do} i = max(isq-1,obc%segment(n)%HI%IsdB), min(ieq+1,obc%segment(n)%HI%IedB)
394           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_n) \textcolor{keywordflow}{then}
395             \textcolor{keywordflow}{if} (area\_h(i,j) + area\_h(i+1,j) > 0.0) \textcolor{keywordflow}{then}
396               harea\_u(i,j+1) = harea\_u(i,j) * ((area\_h(i,j+1) + area\_h(i+1,j+1)) / &
397                                                (area\_h(i,j) + area\_h(i+1,j)))
398             \textcolor{keywordflow}{else} ; harea\_u(i,j+1) = 0.0 ;\textcolor{keywordflow}{ endif}
399           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_S)}
400             \textcolor{keywordflow}{if} (area\_h(i,j+1) + area\_h(i+1,j+1) > 0.0) \textcolor{keywordflow}{then}
401               harea\_u(i,j) = harea\_u(i,j+1) * ((area\_h(i,j) + area\_h(i+1,j)) / &
402                                                (area\_h(i,j+1) + area\_h(i+1,j+1)))
403             \textcolor{keywordflow}{else} ; harea\_u(i,j) = 0.0 ;\textcolor{keywordflow}{ endif}
404 \textcolor{keywordflow}{          endif}
405 \textcolor{keywordflow}{        enddo}
406       \textcolor{keywordflow}{elseif} (obc%segment(n)%is\_E\_or\_W .and. (i >= isq-1) .and. (i <= ieq+1)) \textcolor{keywordflow}{then}
407         \textcolor{keywordflow}{do} j = max(jsq-1,obc%segment(n)%HI%JsdB), min(jeq+1,obc%segment(n)%HI%JedB)
408           \textcolor{keywordflow}{if} (obc%segment(n)%direction == obc\_direction\_e) \textcolor{keywordflow}{then}
409             \textcolor{keywordflow}{if} (area\_h(i,j) + area\_h(i,j+1) > 0.0) \textcolor{keywordflow}{then}
410               harea\_v(i+1,j) = harea\_v(i,j) * ((area\_h(i+1,j) + area\_h(i+1,j+1)) / &
411                                                (area\_h(i,j) + area\_h(i,j+1)))
412             \textcolor{keywordflow}{else} ; harea\_v(i+1,j) = 0.0 ;\textcolor{keywordflow}{ endif}
413           \textcolor{keywordflow}{else} \textcolor{comment}{! (OBC%segment(n)%direction == OBC\_DIRECTION\_W)}
414             harea\_v(i,j) = 0.5 * (area\_h(i,j) + area\_h(i,j+1)) * h(i,j+1,k)
415             \textcolor{keywordflow}{if} (area\_h(i+1,j) + area\_h(i+1,j+1) > 0.0) \textcolor{keywordflow}{then}
416               harea\_v(i,j) = harea\_v(i+1,j) * ((area\_h(i,j) + area\_h(i,j+1)) / &
417                                                (area\_h(i+1,j) + area\_h(i+1,j+1)))
418             \textcolor{keywordflow}{else} ; harea\_v(i,j) = 0.0 ;\textcolor{keywordflow}{ endif}
419 \textcolor{keywordflow}{          endif}
420 \textcolor{keywordflow}{        enddo}
421 \textcolor{keywordflow}{      endif}
422 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ endif}
423 
424     \textcolor{keywordflow}{do} j=jsq-1,jeq+1 ; \textcolor{keywordflow}{do} i=isq-1,ieq+1
425       \textcolor{keywordflow}{if} (cs%no\_slip ) \textcolor{keywordflow}{then}
426         relative\_vorticity = (2.0-g%mask2dBu(i,j)) * (dvdx(i,j) - dudy(i,j)) * g%IareaBu(i,j)
427       \textcolor{keywordflow}{else}
428         relative\_vorticity = g%mask2dBu(i,j) * (dvdx(i,j) - dudy(i,j)) * g%IareaBu(i,j)
429 \textcolor{keywordflow}{      endif}
430       absolute\_vorticity = g%CoriolisBu(i,j) + relative\_vorticity
431       ih = 0.0
432       \textcolor{keywordflow}{if} (area\_q(i,j) > 0.0) \textcolor{keywordflow}{then}
433         harea\_q = (harea\_u(i,j) + harea\_u(i,j+1)) + (harea\_v(i,j) + harea\_v(i+1,j))
434         ih = area\_q(i,j) / (harea\_q + h\_neglect*area\_q(i,j))
435 \textcolor{keywordflow}{      endif}
436       q(i,j) = absolute\_vorticity * ih
437       abs\_vort(i,j) = absolute\_vorticity
438       ih\_q(i,j) = ih
439 
440       \textcolor{keywordflow}{if} (cs%bound\_Coriolis) \textcolor{keywordflow}{then}
441         fv1 = absolute\_vorticity * v(i+1,j,k)
442         fv2 = absolute\_vorticity * v(i,j,k)
443         fu1 = -absolute\_vorticity * u(i,j+1,k)
444         fu2 = -absolute\_vorticity * u(i,j,k)
445         \textcolor{keywordflow}{if} (fv1 > fv2) \textcolor{keywordflow}{then}
446           max\_fvq(i,j) = fv1 ; min\_fvq(i,j) = fv2
447         \textcolor{keywordflow}{else}
448           max\_fvq(i,j) = fv2 ; min\_fvq(i,j) = fv1
449 \textcolor{keywordflow}{        endif}
450         \textcolor{keywordflow}{if} (fu1 > fu2) \textcolor{keywordflow}{then}
451           max\_fuq(i,j) = fu1 ; min\_fuq(i,j) = fu2
452         \textcolor{keywordflow}{else}
453           max\_fuq(i,j) = fu2 ; min\_fuq(i,j) = fu1
454 \textcolor{keywordflow}{        endif}
455 \textcolor{keywordflow}{      endif}
456 
457       \textcolor{keywordflow}{if} (cs%id\_rv > 0) rv(i,j,k) = relative\_vorticity
458       \textcolor{keywordflow}{if} (cs%id\_PV > 0) pv(i,j,k) = q(i,j)
459       \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%rv\_x\_v) .or. \textcolor{keyword}{associated}(ad%rv\_x\_u)) &
460         q2(i,j) = relative\_vorticity * ih
461 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
462 
463     \textcolor{comment}{!   a, b, c, and d are combinations of neighboring potential}
464     \textcolor{comment}{! vorticities which form the Arakawa and Hsu vorticity advection}
465     \textcolor{comment}{! scheme.  All are defined at u grid points.}
466 
467     \textcolor{keywordflow}{if} (cs%Coriolis\_Scheme == arakawa\_hsu90) \textcolor{keywordflow}{then}
468       \textcolor{keywordflow}{do} j=jsq,jeq+1
469         \textcolor{keywordflow}{do} i=is-1,ieq
470           a(i,j) = (q(i,j) + (q(i+1,j) + q(i,j-1))) * c1\_12
471           d(i,j) = ((q(i,j) + q(i+1,j-1)) + q(i,j-1)) * c1\_12
472 \textcolor{keywordflow}{        enddo}
473         \textcolor{keywordflow}{do} i=isq,ieq
474           b(i,j) = (q(i,j) + (q(i-1,j) + q(i,j-1))) * c1\_12
475           c(i,j) = ((q(i,j) + q(i-1,j-1)) + q(i,j-1)) * c1\_12
476 \textcolor{keywordflow}{        enddo}
477 \textcolor{keywordflow}{      enddo}
478     \textcolor{keywordflow}{elseif} (cs%Coriolis\_Scheme == arakawa\_lamb81) \textcolor{keywordflow}{then}
479       \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=isq,ieq+1
480         a(i-1,j) = (2.0*(q(i,j) + q(i-1,j-1)) + (q(i-1,j) + q(i,j-1))) * c1\_24
481         d(i-1,j) = ((q(i,j) + q(i-1,j-1)) + 2.0*(q(i-1,j) + q(i,j-1))) * c1\_24
482         b(i,j) =   ((q(i,j) + q(i-1,j-1)) + 2.0*(q(i-1,j) + q(i,j-1))) * c1\_24
483         c(i,j) =   (2.0*(q(i,j) + q(i-1,j-1)) + (q(i-1,j) + q(i,j-1))) * c1\_24
484         ep\_u(i,j) = ((q(i,j) - q(i-1,j-1)) + (q(i-1,j) - q(i,j-1))) * c1\_24
485         ep\_v(i,j) = (-(q(i,j) - q(i-1,j-1)) + (q(i-1,j) - q(i,j-1))) * c1\_24
486 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
487     \textcolor{keywordflow}{elseif} (cs%Coriolis\_Scheme == al\_blend) \textcolor{keywordflow}{then}
488       fe\_m2 = cs%F\_eff\_max\_blend - 2.0
489       rat\_lin = 1.5 * fe\_m2 / max(cs%wt\_lin\_blend, 1.0e-16)
490 
491       \textcolor{comment}{! This allows the code to always give Sadourny Energy}
492       \textcolor{keywordflow}{if} (cs%F\_eff\_max\_blend <= 2.0) \textcolor{keywordflow}{then} ; fe\_m2 = -1. ; rat\_lin = -1.0 ;\textcolor{keywordflow}{ endif}
493 
494       \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=isq,ieq+1
495         min\_ihq = min(ih\_q(i-1,j-1), ih\_q(i,j-1), ih\_q(i-1,j), ih\_q(i,j))
496         max\_ihq = max(ih\_q(i-1,j-1), ih\_q(i,j-1), ih\_q(i-1,j), ih\_q(i,j))
497         rat\_m1 = 1.0e15
498         \textcolor{keywordflow}{if} (max\_ihq < 1.0e15*min\_ihq) rat\_m1 = max\_ihq / min\_ihq - 1.0
499         \textcolor{comment}{! The weights used here are designed to keep the effective Coriolis}
500         \textcolor{comment}{! acceleration from any one point on its neighbors within a factor}
501         \textcolor{comment}{! of F\_eff\_max.  The minimum permitted value is 2 (the factor for}
502         \textcolor{comment}{! Sadourny's energy conserving scheme).}
503 
504         \textcolor{comment}{! Determine the relative weights of Arakawa & Lamb vs. Arakawa and Hsu.}
505         \textcolor{keywordflow}{if} (rat\_m1 <= fe\_m2) \textcolor{keywordflow}{then} ; al\_wt = 1.0
506         \textcolor{keywordflow}{elseif} (rat\_m1 < 1.5*fe\_m2) \textcolor{keywordflow}{then} ; al\_wt = 3.0*fe\_m2 / rat\_m1 - 2.0
507         \textcolor{keywordflow}{else} ; al\_wt = 0.0 ;\textcolor{keywordflow}{ endif}
508 
509         \textcolor{comment}{! Determine the relative weights of Sadourny Energy vs. the other two.}
510         \textcolor{keywordflow}{if} (rat\_m1 <= 1.5*fe\_m2) \textcolor{keywordflow}{then} ; sad\_wt = 0.0
511         \textcolor{keywordflow}{elseif} (rat\_m1 <= rat\_lin) \textcolor{keywordflow}{then}
512           sad\_wt = 1.0 - (1.5*fe\_m2) / rat\_m1
513         \textcolor{keywordflow}{elseif} (rat\_m1 < 2.0*rat\_lin) \textcolor{keywordflow}{then}
514           sad\_wt = 1.0 - (cs%wt\_lin\_blend / rat\_lin) * (rat\_m1 - 2.0*rat\_lin)
515         \textcolor{keywordflow}{else} ; sad\_wt = 1.0 ;\textcolor{keywordflow}{ endif}
516 
517         a(i-1,j) = sad\_wt * 0.25 * q(i-1,j) + (1.0 - sad\_wt) * &
518                    ( ((2.0-al\_wt)* q(i-1,j) + al\_wt*q(i,j-1)) + &
519                       2.0 * (q(i,j) + q(i-1,j-1)) ) * c1\_24
520         d(i-1,j) = sad\_wt * 0.25 * q(i-1,j-1) + (1.0 - sad\_wt) * &
521                    ( ((2.0-al\_wt)* q(i-1,j-1) + al\_wt*q(i,j)) + &
522                       2.0 * (q(i-1,j) + q(i,j-1)) ) * c1\_24
523         b(i,j) =   sad\_wt * 0.25 * q(i,j) + (1.0 - sad\_wt) * &
524                    ( ((2.0-al\_wt)* q(i,j) + al\_wt*q(i-1,j-1)) + &
525                       2.0 * (q(i-1,j) + q(i,j-1)) ) * c1\_24
526         c(i,j) =   sad\_wt * 0.25 * q(i,j-1) + (1.0 - sad\_wt) * &
527                    ( ((2.0-al\_wt)* q(i,j-1) + al\_wt*q(i-1,j)) + &
528                       2.0 * (q(i,j) + q(i-1,j-1)) ) * c1\_24
529         ep\_u(i,j) = al\_wt  * ((q(i,j) - q(i-1,j-1)) + (q(i-1,j) - q(i,j-1))) * c1\_24
530         ep\_v(i,j) = al\_wt * (-(q(i,j) - q(i-1,j-1)) + (q(i-1,j) - q(i,j-1))) * c1\_24
531 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
532 \textcolor{keywordflow}{    endif}
533 
534     \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis) \textcolor{keywordflow}{then}
535     \textcolor{comment}{!  c1 = 1.0-1.5*RANGE ; c2 = 1.0-RANGE ; c3 = 2.0 ; slope = 0.5}
536       c1 = 1.0-1.5*0.5 ; c2 = 1.0-0.5 ; c3 = 2.0 ; slope = 0.5
537 
538       \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=is-1,ie
539         uhc = uh\_center(i,j)
540         uhm = uh(i,j,k)
541         \textcolor{comment}{! This sometimes matters with some types of open boundary conditions.}
542         \textcolor{keywordflow}{if} (g%dy\_Cu(i,j) == 0.0) uhc = uhm
543 
544         \textcolor{keywordflow}{if} (abs(uhc) < 0.1*abs(uhm)) \textcolor{keywordflow}{then}
545           uhm = 10.0*uhc
546         \textcolor{keywordflow}{elseif} (abs(uhc) > c1*abs(uhm)) \textcolor{keywordflow}{then}
547           \textcolor{keywordflow}{if} (abs(uhc) < c2*abs(uhm)) \textcolor{keywordflow}{then} ; uhc = (3.0*uhc+(1.0-c2*3.0)*uhm)
548           \textcolor{keywordflow}{elseif} (abs(uhc) <= c3*abs(uhm)) \textcolor{keywordflow}{then} ; uhc = uhm
549           \textcolor{keywordflow}{else} ; uhc = slope*uhc+(1.0-c3*slope)*uhm
550 \textcolor{keywordflow}{          endif}
551 \textcolor{keywordflow}{        endif}
552 
553         \textcolor{keywordflow}{if} (uhc > uhm) \textcolor{keywordflow}{then}
554           uh\_min(i,j) = uhm ; uh\_max(i,j) = uhc
555         \textcolor{keywordflow}{else}
556           uh\_max(i,j) = uhm ; uh\_min(i,j) = uhc
557 \textcolor{keywordflow}{        endif}
558 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
559       \textcolor{keywordflow}{do} j=js-1,je ; \textcolor{keywordflow}{do} i=isq,ieq+1
560         vhc = vh\_center(i,j)
561         vhm = vh(i,j,k)
562         \textcolor{comment}{! This sometimes matters with some types of open boundary conditions.}
563         \textcolor{keywordflow}{if} (g%dx\_Cv(i,j) == 0.0) vhc = vhm
564 
565         \textcolor{keywordflow}{if} (abs(vhc) < 0.1*abs(vhm)) \textcolor{keywordflow}{then}
566           vhm = 10.0*vhc
567         \textcolor{keywordflow}{elseif} (abs(vhc) > c1*abs(vhm)) \textcolor{keywordflow}{then}
568           \textcolor{keywordflow}{if} (abs(vhc) < c2*abs(vhm)) \textcolor{keywordflow}{then} ; vhc = (3.0*vhc+(1.0-c2*3.0)*vhm)
569           \textcolor{keywordflow}{elseif} (abs(vhc) <= c3*abs(vhm)) \textcolor{keywordflow}{then} ; vhc = vhm
570           \textcolor{keywordflow}{else} ; vhc = slope*vhc+(1.0-c3*slope)*vhm
571 \textcolor{keywordflow}{          endif}
572 \textcolor{keywordflow}{        endif}
573 
574         \textcolor{keywordflow}{if} (vhc > vhm) \textcolor{keywordflow}{then}
575           vh\_min(i,j) = vhm ; vh\_max(i,j) = vhc
576         \textcolor{keywordflow}{else}
577           vh\_max(i,j) = vhm ; vh\_min(i,j) = vhc
578 \textcolor{keywordflow}{        endif}
579 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
580 \textcolor{keywordflow}{    endif}
581 
582     \textcolor{comment}{! Calculate KE and the gradient of KE}
583     \textcolor{keyword}{call }gradke(u, v, h, ke, kex, key, k, obc, g, us, cs)
584 
585     \textcolor{comment}{! Calculate the tendencies of zonal velocity due to the Coriolis}
586     \textcolor{comment}{! force and momentum advection.  On a Cartesian grid, this is}
587     \textcolor{comment}{!     CAu =  q * vh - d(KE)/dx.}
588     \textcolor{keywordflow}{if} (cs%Coriolis\_Scheme == sadourny75\_energy) \textcolor{keywordflow}{then}
589       \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis) \textcolor{keywordflow}{then}
590         \textcolor{comment}{! Energy dissipating biased scheme, Hallberg 200x}
591         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
592           \textcolor{keywordflow}{if} (q(i,j)*u(i,j,k) == 0.0) \textcolor{keywordflow}{then}
593             temp1 = q(i,j) * ( (vh\_max(i,j)+vh\_max(i+1,j)) &
594                              + (vh\_min(i,j)+vh\_min(i+1,j)) )*0.5
595           \textcolor{keywordflow}{elseif} (q(i,j)*u(i,j,k) < 0.0) \textcolor{keywordflow}{then}
596             temp1 = q(i,j) * (vh\_max(i,j)+vh\_max(i+1,j))
597           \textcolor{keywordflow}{else}
598             temp1 = q(i,j) * (vh\_min(i,j)+vh\_min(i+1,j))
599 \textcolor{keywordflow}{          endif}
600           \textcolor{keywordflow}{if} (q(i,j-1)*u(i,j,k) == 0.0) \textcolor{keywordflow}{then}
601             temp2 = q(i,j-1) * ( (vh\_max(i,j-1)+vh\_max(i+1,j-1)) &
602                                + (vh\_min(i,j-1)+vh\_min(i+1,j-1)) )*0.5
603           \textcolor{keywordflow}{elseif} (q(i,j-1)*u(i,j,k) < 0.0) \textcolor{keywordflow}{then}
604             temp2 = q(i,j-1) * (vh\_max(i,j-1)+vh\_max(i+1,j-1))
605           \textcolor{keywordflow}{else}
606             temp2 = q(i,j-1) * (vh\_min(i,j-1)+vh\_min(i+1,j-1))
607 \textcolor{keywordflow}{          endif}
608           cau(i,j,k) = 0.25 * g%IdxCu(i,j) * (temp1 + temp2)
609 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
610       \textcolor{keywordflow}{else}
611         \textcolor{comment}{! Energy conserving scheme, Sadourny 1975}
612         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
613           cau(i,j,k) = 0.25 * &
614             (q(i,j) * (vh(i+1,j,k) + vh(i,j,k)) + &
615              q(i,j-1) * (vh(i,j-1,k) + vh(i+1,j-1,k))) * g%IdxCu(i,j)
616 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
617 \textcolor{keywordflow}{      endif}
618     \textcolor{keywordflow}{elseif} (cs%Coriolis\_Scheme == sadourny75\_enstro) \textcolor{keywordflow}{then}
619       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
620         cau(i,j,k) = 0.125 * (g%IdxCu(i,j) * (q(i,j) + q(i,j-1))) * &
621                      ((vh(i+1,j,k) + vh(i,j,k)) + (vh(i,j-1,k) + vh(i+1,j-1,k)))
622 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
623     \textcolor{keywordflow}{elseif} ((cs%Coriolis\_Scheme == arakawa\_hsu90) .or. &
624             (cs%Coriolis\_Scheme == arakawa\_lamb81) .or. &
625             (cs%Coriolis\_Scheme == al\_blend)) \textcolor{keywordflow}{then}
626       \textcolor{comment}{! (Global) Energy and (Local) Enstrophy conserving, Arakawa & Hsu 1990}
627       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
628         cau(i,j,k) = ((a(i,j) * vh(i+1,j,k) +  c(i,j) * vh(i,j-1,k))  + &
629                       (b(i,j) * vh(i,j,k) +  d(i,j) * vh(i+1,j-1,k))) * g%IdxCu(i,j)
630 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
631     \textcolor{keywordflow}{elseif} (cs%Coriolis\_Scheme == robust\_enstro) \textcolor{keywordflow}{then}
632       \textcolor{comment}{! An enstrophy conserving scheme robust to vanishing layers}
633       \textcolor{comment}{! Note: Heffs are in lieu of h\_at\_v that should be returned by the}
634       \textcolor{comment}{!       continuity solver. AJA}
635       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
636         heff1 = abs(vh(i,j,k) * g%IdxCv(i,j)) / (eps\_vel+abs(v(i,j,k)))
637         heff1 = max(heff1, min(h(i,j,k),h(i,j+1,k)))
638         heff1 = min(heff1, max(h(i,j,k),h(i,j+1,k)))
639         heff2 = abs(vh(i,j-1,k) * g%IdxCv(i,j-1)) / (eps\_vel+abs(v(i,j-1,k)))
640         heff2 = max(heff2, min(h(i,j-1,k),h(i,j,k)))
641         heff2 = min(heff2, max(h(i,j-1,k),h(i,j,k)))
642         heff3 = abs(vh(i+1,j,k) * g%IdxCv(i+1,j)) / (eps\_vel+abs(v(i+1,j,k)))
643         heff3 = max(heff3, min(h(i+1,j,k),h(i+1,j+1,k)))
644         heff3 = min(heff3, max(h(i+1,j,k),h(i+1,j+1,k)))
645         heff4 = abs(vh(i+1,j-1,k) * g%IdxCv(i+1,j-1)) / (eps\_vel+abs(v(i+1,j-1,k)))
646         heff4 = max(heff4, min(h(i+1,j-1,k),h(i+1,j,k)))
647         heff4 = min(heff4, max(h(i+1,j-1,k),h(i+1,j,k)))
648         \textcolor{keywordflow}{if} (cs%PV\_Adv\_Scheme == pv\_adv\_centered) \textcolor{keywordflow}{then}
649           cau(i,j,k) = 0.5*(abs\_vort(i,j)+abs\_vort(i,j-1)) * &
650                        ((vh(i,j,k) + vh(i+1,j-1,k)) + (vh(i,j-1,k) + vh(i+1,j,k)) ) /  &
651                        (h\_tiny + ((heff1+heff4) + (heff2+heff3)) ) * g%IdxCu(i,j)
652         \textcolor{keywordflow}{elseif} (cs%PV\_Adv\_Scheme == pv\_adv\_upwind1) \textcolor{keywordflow}{then}
653           vheff = ((vh(i,j,k) + vh(i+1,j-1,k)) + (vh(i,j-1,k) + vh(i+1,j,k)) )
654           qvheff = 0.5*( (abs\_vort(i,j)+abs\_vort(i,j-1))*vheff &
655                         -(abs\_vort(i,j)-abs\_vort(i,j-1))*abs(vheff) )
656           cau(i,j,k) = (qvheff / ( h\_tiny + ((heff1+heff4) + (heff2+heff3)) ) ) * g%IdxCu(i,j)
657 \textcolor{keywordflow}{        endif}
658 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
659 \textcolor{keywordflow}{    endif}
660     \textcolor{comment}{! Add in the additonal terms with Arakawa & Lamb.}
661     \textcolor{keywordflow}{if} ((cs%Coriolis\_Scheme == arakawa\_lamb81) .or. &
662         (cs%Coriolis\_Scheme == al\_blend)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
663       cau(i,j,k) = cau(i,j,k) + &
664             (ep\_u(i,j)*uh(i-1,j,k) - ep\_u(i+1,j)*uh(i+1,j,k)) * g%IdxCu(i,j)
665 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ endif}
666 
667 
668     \textcolor{keywordflow}{if} (cs%bound\_Coriolis) \textcolor{keywordflow}{then}
669       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
670         max\_fv = max(max\_fvq(i,j), max\_fvq(i,j-1))
671         min\_fv = min(min\_fvq(i,j), min\_fvq(i,j-1))
672        \textcolor{comment}{! CAu(I,j,k) = min( CAu(I,j,k), max\_fv )}
673        \textcolor{comment}{! CAu(I,j,k) = max( CAu(I,j,k), min\_fv )}
674         \textcolor{keywordflow}{if} (cau(i,j,k) > max\_fv) \textcolor{keywordflow}{then}
675             cau(i,j,k) = max\_fv
676         \textcolor{keywordflow}{else}
677           \textcolor{keywordflow}{if} (cau(i,j,k) < min\_fv) cau(i,j,k) = min\_fv
678 \textcolor{keywordflow}{        endif}
679 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
680 \textcolor{keywordflow}{    endif}
681 
682     \textcolor{comment}{! Term - d(KE)/dx.}
683     \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
684       cau(i,j,k) = cau(i,j,k) - kex(i,j)
685       \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%gradKEu)) ad%gradKEu(i,j,k) = -kex(i,j)
686 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
687 
688 
689     \textcolor{comment}{! Calculate the tendencies of meridional velocity due to the Coriolis}
690     \textcolor{comment}{! force and momentum advection.  On a Cartesian grid, this is}
691     \textcolor{comment}{!     CAv = - q * uh - d(KE)/dy.}
692     \textcolor{keywordflow}{if} (cs%Coriolis\_Scheme == sadourny75\_energy) \textcolor{keywordflow}{then}
693       \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis) \textcolor{keywordflow}{then}
694         \textcolor{comment}{! Energy dissipating biased scheme, Hallberg 200x}
695         \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
696           \textcolor{keywordflow}{if} (q(i-1,j)*v(i,j,k) == 0.0) \textcolor{keywordflow}{then}
697             temp1 = q(i-1,j) * ( (uh\_max(i-1,j)+uh\_max(i-1,j+1)) &
698                                + (uh\_min(i-1,j)+uh\_min(i-1,j+1)) )*0.5
699           \textcolor{keywordflow}{elseif} (q(i-1,j)*v(i,j,k) > 0.0) \textcolor{keywordflow}{then}
700             temp1 = q(i-1,j) * (uh\_max(i-1,j)+uh\_max(i-1,j+1))
701           \textcolor{keywordflow}{else}
702             temp1 = q(i-1,j) * (uh\_min(i-1,j)+uh\_min(i-1,j+1))
703 \textcolor{keywordflow}{          endif}
704           \textcolor{keywordflow}{if} (q(i,j)*v(i,j,k) == 0.0) \textcolor{keywordflow}{then}
705             temp2 = q(i,j) * ( (uh\_max(i,j)+uh\_max(i,j+1)) &
706                              + (uh\_min(i,j)+uh\_min(i,j+1)) )*0.5
707           \textcolor{keywordflow}{elseif} (q(i,j)*v(i,j,k) > 0.0) \textcolor{keywordflow}{then}
708             temp2 = q(i,j) * (uh\_max(i,j)+uh\_max(i,j+1))
709           \textcolor{keywordflow}{else}
710             temp2 = q(i,j) * (uh\_min(i,j)+uh\_min(i,j+1))
711 \textcolor{keywordflow}{          endif}
712           cav(i,j,k) = -0.25 * g%IdyCv(i,j) * (temp1 + temp2)
713 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
714       \textcolor{keywordflow}{else}
715         \textcolor{comment}{! Energy conserving scheme, Sadourny 1975}
716         \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
717           cav(i,j,k) = - 0.25* &
718               (q(i-1,j)*(uh(i-1,j,k) + uh(i-1,j+1,k)) + &
719                q(i,j)*(uh(i,j,k) + uh(i,j+1,k))) * g%IdyCv(i,j)
720 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
721 \textcolor{keywordflow}{      endif}
722     \textcolor{keywordflow}{elseif} (cs%Coriolis\_Scheme == sadourny75\_enstro) \textcolor{keywordflow}{then}
723       \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
724         cav(i,j,k) = -0.125 * (g%IdyCv(i,j) * (q(i-1,j) + q(i,j))) * &
725                      ((uh(i-1,j,k) + uh(i-1,j+1,k)) + (uh(i,j,k) + uh(i,j+1,k)))
726 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
727     \textcolor{keywordflow}{elseif} ((cs%Coriolis\_Scheme == arakawa\_hsu90) .or. &
728             (cs%Coriolis\_Scheme == arakawa\_lamb81) .or. &
729             (cs%Coriolis\_Scheme == al\_blend)) \textcolor{keywordflow}{then}
730       \textcolor{comment}{! (Global) Energy and (Local) Enstrophy conserving, Arakawa & Hsu 1990}
731       \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
732         cav(i,j,k) = - ((a(i-1,j)   * uh(i-1,j,k) + &
733                          c(i,j+1)   * uh(i,j+1,k))  &
734                       + (b(i,j)     * uh(i,j,k) +   &
735                          d(i-1,j+1) * uh(i-1,j+1,k))) * g%IdyCv(i,j)
736 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
737     \textcolor{keywordflow}{elseif} (cs%Coriolis\_Scheme == robust\_enstro) \textcolor{keywordflow}{then}
738       \textcolor{comment}{! An enstrophy conserving scheme robust to vanishing layers}
739       \textcolor{comment}{! Note: Heffs are in lieu of h\_at\_u that should be returned by the}
740       \textcolor{comment}{!       continuity solver. AJA}
741       \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
742         heff1 = abs(uh(i,j,k) * g%IdyCu(i,j)) / (eps\_vel+abs(u(i,j,k)))
743         heff1 = max(heff1, min(h(i,j,k),h(i+1,j,k)))
744         heff1 = min(heff1, max(h(i,j,k),h(i+1,j,k)))
745         heff2 = abs(uh(i-1,j,k) * g%IdyCu(i-1,j)) / (eps\_vel+abs(u(i-1,j,k)))
746         heff2 = max(heff2, min(h(i-1,j,k),h(i,j,k)))
747         heff2 = min(heff2, max(h(i-1,j,k),h(i,j,k)))
748         heff3 = abs(uh(i,j+1,k) * g%IdyCu(i,j+1)) / (eps\_vel+abs(u(i,j+1,k)))
749         heff3 = max(heff3, min(h(i,j+1,k),h(i+1,j+1,k)))
750         heff3 = min(heff3, max(h(i,j+1,k),h(i+1,j+1,k)))
751         heff4 = abs(uh(i-1,j+1,k) * g%IdyCu(i-1,j+1)) / (eps\_vel+abs(u(i-1,j+1,k)))
752         heff4 = max(heff4, min(h(i-1,j+1,k),h(i,j+1,k)))
753         heff4 = min(heff4, max(h(i-1,j+1,k),h(i,j+1,k)))
754         \textcolor{keywordflow}{if} (cs%PV\_Adv\_Scheme == pv\_adv\_centered) \textcolor{keywordflow}{then}
755           cav(i,j,k) = - 0.5*(abs\_vort(i,j)+abs\_vort(i-1,j)) * &
756                          ((uh(i  ,j  ,k)+uh(i-1,j+1,k)) +      &
757                           (uh(i-1,j  ,k)+uh(i  ,j+1,k)) ) /    &
758                       (h\_tiny + ((heff1+heff4) +(heff2+heff3)) ) * g%IdyCv(i,j)
759         \textcolor{keywordflow}{elseif} (cs%PV\_Adv\_Scheme == pv\_adv\_upwind1) \textcolor{keywordflow}{then}
760           uheff = ((uh(i  ,j  ,k)+uh(i-1,j+1,k)) +      &
761                    (uh(i-1,j  ,k)+uh(i  ,j+1,k)) )
762           quheff = 0.5*( (abs\_vort(i,j)+abs\_vort(i-1,j))*uheff &
763                         -(abs\_vort(i,j)-abs\_vort(i-1,j))*abs(uheff) )
764           cav(i,j,k) = - quheff / &
765                        (h\_tiny + ((heff1+heff4) +(heff2+heff3)) ) * g%IdyCv(i,j)
766 \textcolor{keywordflow}{        endif}
767 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
768 \textcolor{keywordflow}{    endif}
769     \textcolor{comment}{! Add in the additonal terms with Arakawa & Lamb.}
770     \textcolor{keywordflow}{if} ((cs%Coriolis\_Scheme == arakawa\_lamb81) .or. &
771         (cs%Coriolis\_Scheme == al\_blend)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
772       cav(i,j,k) = cav(i,j,k) + &
773             (ep\_v(i,j)*vh(i,j-1,k) - ep\_v(i,j+1)*vh(i,j+1,k)) * g%IdyCv(i,j)
774 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ endif}
775 
776     \textcolor{keywordflow}{if} (cs%bound\_Coriolis) \textcolor{keywordflow}{then}
777       \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
778         max\_fu = max(max\_fuq(i,j),max\_fuq(i-1,j))
779         min\_fu = min(min\_fuq(i,j),min\_fuq(i-1,j))
780         \textcolor{keywordflow}{if} (cav(i,j,k) > max\_fu) \textcolor{keywordflow}{then}
781             cav(i,j,k) = max\_fu
782         \textcolor{keywordflow}{else}
783           \textcolor{keywordflow}{if} (cav(i,j,k) < min\_fu) cav(i,j,k) = min\_fu
784 \textcolor{keywordflow}{        endif}
785 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
786 \textcolor{keywordflow}{    endif}
787 
788     \textcolor{comment}{! Term - d(KE)/dy.}
789     \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
790       cav(i,j,k) = cav(i,j,k) - key(i,j)
791       \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%gradKEv)) ad%gradKEv(i,j,k) = -key(i,j)
792 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
793 
794     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%rv\_x\_u) .or. \textcolor{keyword}{associated}(ad%rv\_x\_v)) \textcolor{keywordflow}{then}
795       \textcolor{comment}{! Calculate the Coriolis-like acceleration due to relative vorticity.}
796       \textcolor{keywordflow}{if} (cs%Coriolis\_Scheme == sadourny75\_energy) \textcolor{keywordflow}{then}
797         \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%rv\_x\_u)) \textcolor{keywordflow}{then}
798           \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
799             ad%rv\_x\_u(i,j,k) = - 0.25* &
800               (q2(i-1,j)*(uh(i-1,j,k) + uh(i-1,j+1,k)) + &
801                q2(i,j)*(uh(i,j,k) + uh(i,j+1,k))) * g%IdyCv(i,j)
802 \textcolor{keywordflow}{          enddo} ;\textcolor{keywordflow}{ enddo}
803 \textcolor{keywordflow}{        endif}
804 
805         \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%rv\_x\_v)) \textcolor{keywordflow}{then}
806           \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
807             ad%rv\_x\_v(i,j,k) = 0.25 * &
808               (q2(i,j) * (vh(i+1,j,k) + vh(i,j,k)) + &
809                q2(i,j-1) * (vh(i,j-1,k) + vh(i+1,j-1,k))) * g%IdxCu(i,j)
810 \textcolor{keywordflow}{          enddo} ;\textcolor{keywordflow}{ enddo}
811 \textcolor{keywordflow}{        endif}
812       \textcolor{keywordflow}{else}
813         \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%rv\_x\_u)) \textcolor{keywordflow}{then}
814           \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
815             ad%rv\_x\_u(i,j,k) = -g%IdyCv(i,j) * c1\_12 * &
816               ((q2(i,j) + q2(i-1,j) + q2(i-1,j-1)) * uh(i-1,j,k) + &
817                (q2(i-1,j) + q2(i,j) + q2(i,j-1)) * uh(i,j,k) + &
818                (q2(i-1,j) + q2(i,j+1) + q2(i,j)) * uh(i,j+1,k) + &
819                (q2(i,j) + q2(i-1,j+1) + q2(i-1,j)) * uh(i-1,j+1,k))
820 \textcolor{keywordflow}{          enddo} ;\textcolor{keywordflow}{ enddo}
821 \textcolor{keywordflow}{        endif}
822 
823         \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(ad%rv\_x\_v)) \textcolor{keywordflow}{then}
824           \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
825             ad%rv\_x\_v(i,j,k) = g%IdxCu(i,j) * c1\_12 * &
826               ((q2(i+1,j) + q2(i,j) + q2(i,j-1)) * vh(i+1,j,k) + &
827                (q2(i-1,j) + q2(i,j) + q2(i,j-1)) * vh(i,j,k) + &
828                (q2(i-1,j-1) + q2(i,j) + q2(i,j-1)) * vh(i,j-1,k) + &
829                (q2(i+1,j-1) + q2(i,j) + q2(i,j-1)) * vh(i+1,j-1,k))
830 \textcolor{keywordflow}{          enddo} ;\textcolor{keywordflow}{ enddo}
831 \textcolor{keywordflow}{        endif}
832 \textcolor{keywordflow}{      endif}
833 \textcolor{keywordflow}{    endif}
834 
835 \textcolor{keywordflow}{  enddo} \textcolor{comment}{! k-loop.}
836 
837   \textcolor{comment}{! Here the various Coriolis-related derived quantities are offered for averaging.}
838   \textcolor{keywordflow}{if} (query\_averaging\_enabled(cs%diag)) \textcolor{keywordflow}{then}
839     \textcolor{keywordflow}{if} (cs%id\_rv > 0) \textcolor{keyword}{call }post\_data(cs%id\_rv, rv, cs%diag)
840     \textcolor{keywordflow}{if} (cs%id\_PV > 0) \textcolor{keyword}{call }post\_data(cs%id\_PV, pv, cs%diag)
841     \textcolor{keywordflow}{if} (cs%id\_gKEu>0) \textcolor{keyword}{call }post\_data(cs%id\_gKEu, ad%gradKEu, cs%diag)
842     \textcolor{keywordflow}{if} (cs%id\_gKEv>0) \textcolor{keyword}{call }post\_data(cs%id\_gKEv, ad%gradKEv, cs%diag)
843     \textcolor{keywordflow}{if} (cs%id\_rvxu > 0) \textcolor{keyword}{call }post\_data(cs%id\_rvxu, ad%rv\_x\_u, cs%diag)
844     \textcolor{keywordflow}{if} (cs%id\_rvxv > 0) \textcolor{keyword}{call }post\_data(cs%id\_rvxv, ad%rv\_x\_v, cs%diag)
845 
846     \textcolor{comment}{! Diagnostics for terms multiplied by fractional thicknesses}
847 
848     \textcolor{comment}{! 3D diagnostics hf\_gKEu etc. are commented because there is no clarity on proper remapping grid
       option.}
849     \textcolor{comment}{! The code is retained for degugging purposes in the future.}
850     \textcolor{comment}{!if (CS%id\_hf\_gKEu > 0) then}
851     \textcolor{comment}{!  allocate(hf\_gKEu(G%IsdB:G%IedB,G%jsd:G%jed,G%ke))}
852     \textcolor{comment}{!  do k=1,nz ; do j=js,je ; do I=Isq,Ieq}
853     \textcolor{comment}{!    hf\_gKEu(I,j,k) = AD%gradKEu(I,j,k) * AD%diag\_hfrac\_u(I,j,k)}
854     \textcolor{comment}{!  enddo ; enddo ; enddo}
855     \textcolor{comment}{!  call post\_data(CS%id\_hf\_gKEu, hf\_gKEu, CS%diag)}
856     \textcolor{comment}{!endif}
857 
858     \textcolor{comment}{!if (CS%id\_hf\_gKEv > 0) then}
859     \textcolor{comment}{!  allocate(hf\_gKEv(G%isd:G%ied,G%JsdB:G%JedB,G%ke))}
860     \textcolor{comment}{!  do k=1,nz ; do J=Jsq,Jeq ; do i=is,ie}
861     \textcolor{comment}{!    hf\_gKEv(i,J,k) = AD%gradKEv(i,J,k) * AD%diag\_hfrac\_v(i,J,k)}
862     \textcolor{comment}{!  enddo ; enddo ; enddo}
863     \textcolor{comment}{!  call post\_data(CS%id\_hf\_gKEv, hf\_gKEv, CS%diag)}
864     \textcolor{comment}{!endif}
865 
866     \textcolor{keywordflow}{if} (cs%id\_hf\_gKEu\_2d > 0) \textcolor{keywordflow}{then}
867       \textcolor{keyword}{allocate}(hf\_gkeu\_2d(g%IsdB:g%IedB,g%jsd:g%jed))
868       hf\_gkeu\_2d(:,:) = 0.0
869       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
870         hf\_gkeu\_2d(i,j) = hf\_gkeu\_2d(i,j) + ad%gradKEu(i,j,k) * ad%diag\_hfrac\_u(i,j,k)
871 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
872       \textcolor{keyword}{call }post\_data(cs%id\_hf\_gKEu\_2d, hf\_gkeu\_2d, cs%diag)
873       \textcolor{keyword}{deallocate}(hf\_gkeu\_2d)
874 \textcolor{keywordflow}{    endif}
875 
876     \textcolor{keywordflow}{if} (cs%id\_hf\_gKEv\_2d > 0) \textcolor{keywordflow}{then}
877       \textcolor{keyword}{allocate}(hf\_gkev\_2d(g%isd:g%ied,g%JsdB:g%JedB))
878       hf\_gkev\_2d(:,:) = 0.0
879       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
880         hf\_gkev\_2d(i,j) = hf\_gkev\_2d(i,j) + ad%gradKEv(i,j,k) * ad%diag\_hfrac\_v(i,j,k)
881 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
882       \textcolor{keyword}{call }post\_data(cs%id\_hf\_gKEv\_2d, hf\_gkev\_2d, cs%diag)
883       \textcolor{keyword}{deallocate}(hf\_gkev\_2d)
884 \textcolor{keywordflow}{    endif}
885 
886     \textcolor{comment}{!if (CS%id\_hf\_rvxv > 0) then}
887     \textcolor{comment}{!  allocate(hf\_rvxv(G%IsdB:G%IedB,G%jsd:G%jed,G%ke))}
888     \textcolor{comment}{!  do k=1,nz ; do j=js,je ; do I=Isq,Ieq}
889     \textcolor{comment}{!    hf\_rvxv(I,j,k) = AD%rv\_x\_v(I,j,k) * AD%diag\_hfrac\_u(I,j,k)}
890     \textcolor{comment}{!  enddo ; enddo ; enddo}
891     \textcolor{comment}{!  call post\_data(CS%id\_hf\_rvxv, hf\_rvxv, CS%diag)}
892     \textcolor{comment}{!endif}
893 
894     \textcolor{comment}{!if (CS%id\_hf\_rvxu > 0) then}
895     \textcolor{comment}{!  allocate(hf\_rvxu(G%isd:G%ied,G%JsdB:G%JedB,G%ke))}
896     \textcolor{comment}{!  do k=1,nz ; do J=Jsq,Jeq ; do i=is,ie}
897     \textcolor{comment}{!    hf\_rvxu(i,J,k) = AD%rv\_x\_u(i,J,k) * AD%diag\_hfrac\_v(i,J,k)}
898     \textcolor{comment}{!  enddo ; enddo ; enddo}
899     \textcolor{comment}{!  call post\_data(CS%id\_hf\_rvxu, hf\_rvxu, CS%diag)}
900     \textcolor{comment}{!endif}
901 
902     \textcolor{keywordflow}{if} (cs%id\_hf\_rvxv\_2d > 0) \textcolor{keywordflow}{then}
903       \textcolor{keyword}{allocate}(hf\_rvxv\_2d(g%IsdB:g%IedB,g%jsd:g%jed))
904       hf\_rvxv\_2d(:,:) = 0.0
905       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
906         hf\_rvxv\_2d(i,j) = hf\_rvxv\_2d(i,j) + ad%rv\_x\_v(i,j,k) * ad%diag\_hfrac\_u(i,j,k)
907 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
908       \textcolor{keyword}{call }post\_data(cs%id\_hf\_rvxv\_2d, hf\_rvxv\_2d, cs%diag)
909       \textcolor{keyword}{deallocate}(hf\_rvxv\_2d)
910 \textcolor{keywordflow}{    endif}
911 
912     \textcolor{keywordflow}{if} (cs%id\_hf\_rvxu\_2d > 0) \textcolor{keywordflow}{then}
913       \textcolor{keyword}{allocate}(hf\_rvxu\_2d(g%isd:g%ied,g%JsdB:g%JedB))
914       hf\_rvxu\_2d(:,:) = 0.0
915       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
916         hf\_rvxu\_2d(i,j) = hf\_rvxu\_2d(i,j) + ad%rv\_x\_u(i,j,k) * ad%diag\_hfrac\_v(i,j,k)
917 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
918       \textcolor{keyword}{call }post\_data(cs%id\_hf\_rvxu\_2d, hf\_rvxu\_2d, cs%diag)
919       \textcolor{keyword}{deallocate}(hf\_rvxu\_2d)
920 \textcolor{keywordflow}{    endif}
921 \textcolor{keywordflow}{  endif}
922 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__coriolisadv_a6252eaea90947c83b5a1900d31191b96}\label{namespacemom__coriolisadv_a6252eaea90947c83b5a1900d31191b96}} 
\index{mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}!coriolisadv\+\_\+end@{coriolisadv\+\_\+end}}
\index{coriolisadv\+\_\+end@{coriolisadv\+\_\+end}!mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}}
\subsubsection{\texorpdfstring{coriolisadv\+\_\+end()}{coriolisadv\_end()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+coriolisadv\+::coriolisadv\+\_\+end (\begin{DoxyParamCaption}\item[{type(\hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}), pointer}]{CS }\end{DoxyParamCaption})}



Destructor for \hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}. 


\begin{DoxyParams}{Parameters}
{\em cs} & Control structure fro M\+O\+M\+\_\+\+Coriolis\+Adv \\
\hline
\end{DoxyParams}


Definition at line 1248 of file M\+O\+M\+\_\+\+Coriolis\+Adv.\+F90.


\begin{DoxyCode}
1248   \textcolor{keywordtype}{type}(coriolisadv\_cs), \textcolor{keywordtype}{pointer} :: cs\textcolor{comment}{ !< Control structure fro MOM\_CoriolisAdv}
1249   \textcolor{keyword}{deallocate}(cs)
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__coriolisadv_ae021ac8de3b3510ca4552314ec9e1a9a}\label{namespacemom__coriolisadv_ae021ac8de3b3510ca4552314ec9e1a9a}} 
\index{mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}!coriolisadv\+\_\+init@{coriolisadv\+\_\+init}}
\index{coriolisadv\+\_\+init@{coriolisadv\+\_\+init}!mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}}
\subsubsection{\texorpdfstring{coriolisadv\+\_\+init()}{coriolisadv\_init()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+coriolisadv\+::coriolisadv\+\_\+init (\begin{DoxyParamCaption}\item[{type(time\+\_\+type), intent(in), target}]{Time,  }\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{type(param\+\_\+file\+\_\+type), intent(in)}]{param\+\_\+file,  }\item[{type(diag\+\_\+ctrl), intent(inout), target}]{diag,  }\item[{type(accel\+\_\+diag\+\_\+ptrs), intent(inout), target}]{AD,  }\item[{type(\hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}), pointer}]{CS }\end{DoxyParamCaption})}



Initializes the control structure for \hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em time} & Current model time\\
\hline
\mbox{\tt in}  & {\em g} & Ocean grid structure\\
\hline
\mbox{\tt in}  & {\em gv} & Vertical grid structure\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
\mbox{\tt in}  & {\em param\+\_\+file} & Runtime parameter handles\\
\hline
\mbox{\tt in,out}  & {\em diag} & Diagnostics control structure\\
\hline
\mbox{\tt in,out}  & {\em ad} & Strorage for acceleration diagnostics\\
\hline
 & {\em cs} & Control structure fro M\+O\+M\+\_\+\+Coriolis\+Adv \\
\hline
\end{DoxyParams}


Definition at line 1012 of file M\+O\+M\+\_\+\+Coriolis\+Adv.\+F90.


\begin{DoxyCode}
1012   \textcolor{keywordtype}{type}(time\_type), \textcolor{keywordtype}{target}, \textcolor{keywordtype}{intent(in)}    :: time\textcolor{comment}{ !< Current model time}
1013   \textcolor{keywordtype}{type}(ocean\_grid\_type),   \textcolor{keywordtype}{intent(in)}    :: g\textcolor{comment}{  !< Ocean grid structure}
1014   \textcolor{keywordtype}{type}(verticalgrid\_type), \textcolor{keywordtype}{intent(in)}    :: gv\textcolor{comment}{ !< Vertical grid structure}
1015   \textcolor{keywordtype}{type}(unit\_scale\_type),   \textcolor{keywordtype}{intent(in)}    :: us\textcolor{comment}{  !< A dimensional unit scaling type}
1016   \textcolor{keywordtype}{type}(param\_file\_type),   \textcolor{keywordtype}{intent(in)}    :: param\_file\textcolor{comment}{ !< Runtime parameter handles}
1017   \textcolor{keywordtype}{type}(diag\_ctrl), \textcolor{keywordtype}{target}, \textcolor{keywordtype}{intent(inout)} :: diag\textcolor{comment}{ !< Diagnostics control structure}
1018   \textcolor{keywordtype}{type}(accel\_diag\_ptrs),   \textcolor{keywordtype}{target}, \textcolor{keywordtype}{intent(inout)} :: ad\textcolor{comment}{ !< Strorage for acceleration diagnostics}
1019   \textcolor{keywordtype}{type}(coriolisadv\_cs),    \textcolor{keywordtype}{pointer}       :: cs\textcolor{comment}{ !< Control structure fro MOM\_CoriolisAdv}
1020   \textcolor{comment}{! Local variables}
1021 \textcolor{comment}{! This include declares and sets the variable "version".}
1022 \textcolor{preprocessor}{#include "version\_variable.h"}
1023 \textcolor{preprocessor}{}  \textcolor{keywordtype}{character(len=40)}  :: mdl = \textcolor{stringliteral}{"MOM\_CoriolisAdv"} \textcolor{comment}{! This module's name.}
1024   \textcolor{keywordtype}{character(len=20)}  :: tmpstr
1025   \textcolor{keywordtype}{character(len=400)} :: mesg
1026   \textcolor{keywordtype}{integer} :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb, nz
1027 
1028   isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = gv%ke
1029   isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1030 
1031   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(cs)) \textcolor{keywordflow}{then}
1032     \textcolor{keyword}{call }mom\_error(warning, \textcolor{stringliteral}{"CoriolisAdv\_init called with associated control structure."})
1033     \textcolor{keywordflow}{return}
1034 \textcolor{keywordflow}{  endif}
1035   \textcolor{keyword}{allocate}(cs)
1036 
1037   cs%diag => diag ; cs%Time => time
1038 
1039   \textcolor{comment}{! Read all relevant parameters and write them to the model log.}
1040   \textcolor{keyword}{call }log\_version(param\_file, mdl, version, \textcolor{stringliteral}{""})
1041   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"NOSLIP"}, cs%no\_slip, &
1042                  \textcolor{stringliteral}{"If true, no slip boundary conditions are used; otherwise "}//&
1043                  \textcolor{stringliteral}{"free slip boundary conditions are assumed. The "}//&
1044                  \textcolor{stringliteral}{"implementation of the free slip BCs on a C-grid is much "}//&
1045                  \textcolor{stringliteral}{"cleaner than the no slip BCs. The use of free slip BCs "}//&
1046                  \textcolor{stringliteral}{"is strongly encouraged, and no slip BCs are not used with "}//&
1047                  \textcolor{stringliteral}{"the biharmonic viscosity."}, default=.false.)
1048 
1049   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"CORIOLIS\_EN\_DIS"}, cs%Coriolis\_En\_Dis, &
1050                  \textcolor{stringliteral}{"If true, two estimates of the thickness fluxes are used "}//&
1051                  \textcolor{stringliteral}{"to estimate the Coriolis term, and the one that "}//&
1052                  \textcolor{stringliteral}{"dissipates energy relative to the other one is used."}, &
1053                  default=.false.)
1054 
1055   \textcolor{comment}{! Set %Coriolis\_Scheme}
1056   \textcolor{comment}{! (Select the baseline discretization for the Coriolis term)}
1057   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"CORIOLIS\_SCHEME"}, tmpstr, &
1058                  \textcolor{stringliteral}{"CORIOLIS\_SCHEME selects the discretization for the "}//&
1059                  \textcolor{stringliteral}{"Coriolis terms. Valid values are: \(\backslash\)n"}//&
1060                  \textcolor{stringliteral}{"\(\backslash\)t SADOURNY75\_ENERGY - Sadourny, 1975; energy cons. \(\backslash\)n"}//&
1061                  \textcolor{stringliteral}{"\(\backslash\)t ARAKAWA\_HSU90     - Arakawa & Hsu, 1990 \(\backslash\)n"}//&
1062                  \textcolor{stringliteral}{"\(\backslash\)t SADOURNY75\_ENSTRO - Sadourny, 1975; enstrophy cons. \(\backslash\)n"}//&
1063                  \textcolor{stringliteral}{"\(\backslash\)t ARAKAWA\_LAMB81    - Arakawa & Lamb, 1981; En. + Enst.\(\backslash\)n"}//&
1064                  \textcolor{stringliteral}{"\(\backslash\)t ARAKAWA\_LAMB\_BLEND - A blend of Arakawa & Lamb with \(\backslash\)n"}//&
1065                  \textcolor{stringliteral}{"\(\backslash\)t                      Arakawa & Hsu and Sadourny energy"}, &
1066                  default=sadourny75\_energy\_string)
1067   tmpstr = uppercase(tmpstr)
1068   \textcolor{keywordflow}{select case} (tmpstr)
1069     \textcolor{keywordflow}{case} (sadourny75\_energy\_string)
1070       cs%Coriolis\_Scheme = sadourny75\_energy
1071     \textcolor{keywordflow}{case} (arakawa\_hsu\_string)
1072       cs%Coriolis\_Scheme = arakawa\_hsu90
1073     \textcolor{keywordflow}{case} (sadourny75\_enstro\_string)
1074       cs%Coriolis\_Scheme = sadourny75\_enstro
1075     \textcolor{keywordflow}{case} (arakawa\_lamb\_string)
1076       cs%Coriolis\_Scheme = arakawa\_lamb81
1077     \textcolor{keywordflow}{case} (al\_blend\_string)
1078       cs%Coriolis\_Scheme = al\_blend
1079     \textcolor{keywordflow}{case} (robust\_enstro\_string)
1080       cs%Coriolis\_Scheme = robust\_enstro
1081       cs%Coriolis\_En\_Dis = .false.
1082 \textcolor{keywordflow}{    case default}
1083       \textcolor{keyword}{call }mom\_mesg(\textcolor{stringliteral}{'CoriolisAdv\_init: Coriolis\_Scheme ="'}//trim(tmpstr)//\textcolor{stringliteral}{'"'}, 0)
1084       \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"CoriolisAdv\_init: Unrecognized setting "}// &
1085             \textcolor{stringliteral}{"#define CORIOLIS\_SCHEME "}//trim(tmpstr)//\textcolor{stringliteral}{" found in input file."})
1086 \textcolor{keywordflow}{  end select}
1087   \textcolor{keywordflow}{if} (cs%Coriolis\_Scheme == al\_blend) \textcolor{keywordflow}{then}
1088     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"CORIOLIS\_BLEND\_WT\_LIN"}, cs%wt\_lin\_blend, &
1089                  \textcolor{stringliteral}{"A weighting value for the ratio of inverse thicknesses, "}//&
1090                  \textcolor{stringliteral}{"beyond which the blending between Sadourny Energy and "}//&
1091                  \textcolor{stringliteral}{"Arakawa & Hsu goes linearly to 0 when CORIOLIS\_SCHEME "}//&
1092                  \textcolor{stringliteral}{"is ARAWAKA\_LAMB\_BLEND. This must be between 1 and 1e-16."}, &
1093                  units=\textcolor{stringliteral}{"nondim"}, default=0.125)
1094     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"CORIOLIS\_BLEND\_F\_EFF\_MAX"}, cs%F\_eff\_max\_blend, &
1095                  \textcolor{stringliteral}{"The factor by which the maximum effective Coriolis "}//&
1096                  \textcolor{stringliteral}{"acceleration from any point can be increased when "}//&
1097                  \textcolor{stringliteral}{"blending different discretizations with the "}//&
1098                  \textcolor{stringliteral}{"ARAKAWA\_LAMB\_BLEND Coriolis scheme.  This must be "}//&
1099                  \textcolor{stringliteral}{"greater than 2.0 (the max value for Sadourny energy)."}, &
1100                  units=\textcolor{stringliteral}{"nondim"}, default=4.0)
1101     cs%wt\_lin\_blend = min(1.0, max(cs%wt\_lin\_blend,1e-16))
1102     \textcolor{keywordflow}{if} (cs%F\_eff\_max\_blend < 2.0) \textcolor{keyword}{call }mom\_error(warning, \textcolor{stringliteral}{"CoriolisAdv\_init: "}//&
1103            \textcolor{stringliteral}{"CORIOLIS\_BLEND\_F\_EFF\_MAX should be at least 2."})
1104 \textcolor{keywordflow}{  endif}
1105 
1106   mesg = \textcolor{stringliteral}{"If true, the Coriolis terms at u-points are bounded by "}//&
1107          \textcolor{stringliteral}{"the four estimates of (f+rv)v from the four neighboring "}//&
1108          \textcolor{stringliteral}{"v-points, and similarly at v-points."}
1109   \textcolor{keywordflow}{if} (cs%Coriolis\_En\_Dis .and. (cs%Coriolis\_Scheme == sadourny75\_energy)) \textcolor{keywordflow}{then}
1110     mesg = trim(mesg)//\textcolor{stringliteral}{"  This option is "}//&
1111                  \textcolor{stringliteral}{"always effectively false with CORIOLIS\_EN\_DIS defined and "}//&
1112                  \textcolor{stringliteral}{"CORIOLIS\_SCHEME set to "}//trim(sadourny75\_energy\_string)//\textcolor{stringliteral}{"."}
1113   \textcolor{keywordflow}{else}
1114     mesg = trim(mesg)//\textcolor{stringliteral}{"  This option would "}//&
1115                  \textcolor{stringliteral}{"have no effect on the SADOURNY Coriolis scheme if it "}//&
1116                  \textcolor{stringliteral}{"were possible to use centered difference thickness fluxes."}
1117 \textcolor{keywordflow}{  endif}
1118   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"BOUND\_CORIOLIS"}, cs%bound\_Coriolis, mesg, &
1119                  default=.false.)
1120   \textcolor{keywordflow}{if} ((cs%Coriolis\_En\_Dis .and. (cs%Coriolis\_Scheme == sadourny75\_energy)) .or. &
1121       (cs%Coriolis\_Scheme == robust\_enstro)) cs%bound\_Coriolis = .false.
1122 
1123   \textcolor{comment}{! Set KE\_Scheme (selects discretization of KE)}
1124   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"KE\_SCHEME"}, tmpstr, &
1125                  \textcolor{stringliteral}{"KE\_SCHEME selects the discretization for acceleration "}//&
1126                  \textcolor{stringliteral}{"due to the kinetic energy gradient. Valid values are: \(\backslash\)n"}//&
1127                  \textcolor{stringliteral}{"\(\backslash\)t KE\_ARAKAWA, KE\_SIMPLE\_GUDONOV, KE\_GUDONOV"}, &
1128                  default=ke\_arakawa\_string)
1129   tmpstr = uppercase(tmpstr)
1130   \textcolor{keywordflow}{select case} (tmpstr)
1131     \textcolor{keywordflow}{case} (ke\_arakawa\_string); cs%KE\_Scheme = ke\_arakawa
1132     \textcolor{keywordflow}{case} (ke\_simple\_gudonov\_string); cs%KE\_Scheme = ke\_simple\_gudonov
1133     \textcolor{keywordflow}{case} (ke\_gudonov\_string); cs%KE\_Scheme = ke\_gudonov
1134 \textcolor{keywordflow}{    case default}
1135       \textcolor{keyword}{call }mom\_mesg(\textcolor{stringliteral}{'CoriolisAdv\_init: KE\_Scheme ="'}//trim(tmpstr)//\textcolor{stringliteral}{'"'}, 0)
1136       \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"CoriolisAdv\_init: "}// &
1137                \textcolor{stringliteral}{"#define KE\_SCHEME "}//trim(tmpstr)//\textcolor{stringliteral}{" in input file is invalid."})
1138 \textcolor{keywordflow}{  end select}
1139 
1140   \textcolor{comment}{! Set PV\_Adv\_Scheme (selects discretization of PV advection)}
1141   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PV\_ADV\_SCHEME"}, tmpstr, &
1142                  \textcolor{stringliteral}{"PV\_ADV\_SCHEME selects the discretization for PV "}//&
1143                  \textcolor{stringliteral}{"advection. Valid values are: \(\backslash\)n"}//&
1144                  \textcolor{stringliteral}{"\(\backslash\)t PV\_ADV\_CENTERED - centered (aka Sadourny, 75) \(\backslash\)n"}//&
1145                  \textcolor{stringliteral}{"\(\backslash\)t PV\_ADV\_UPWIND1  - upwind, first order"}, &
1146                  default=pv\_adv\_centered\_string)
1147   \textcolor{keywordflow}{select case} (uppercase(tmpstr))
1148     \textcolor{keywordflow}{case} (pv\_adv\_centered\_string)
1149       cs%PV\_Adv\_Scheme = pv\_adv\_centered
1150     \textcolor{keywordflow}{case} (pv\_adv\_upwind1\_string)
1151       cs%PV\_Adv\_Scheme = pv\_adv\_upwind1
1152 \textcolor{keywordflow}{    case default}
1153       \textcolor{keyword}{call }mom\_mesg(\textcolor{stringliteral}{'CoriolisAdv\_init: PV\_Adv\_Scheme ="'}//trim(tmpstr)//\textcolor{stringliteral}{'"'}, 0)
1154       \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"CoriolisAdv\_init: "}// &
1155                      \textcolor{stringliteral}{"#DEFINE PV\_ADV\_SCHEME in input file is invalid."})
1156 \textcolor{keywordflow}{  end select}
1157 
1158   cs%id\_rv = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'RV'}, diag%axesBL, time, &
1159      \textcolor{stringliteral}{'Relative Vorticity'}, \textcolor{stringliteral}{'s-1'}, conversion=us%s\_to\_T)
1160 
1161   cs%id\_PV = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'PV'}, diag%axesBL, time, &
1162      \textcolor{stringliteral}{'Potential Vorticity'}, \textcolor{stringliteral}{'m-1 s-1'}, conversion=gv%m\_to\_H*us%s\_to\_T)
1163 
1164   cs%id\_gKEu = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'gKEu'}, diag%axesCuL, time, &
1165      \textcolor{stringliteral}{'Zonal Acceleration from Grad. Kinetic Energy'}, \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1166   \textcolor{keywordflow}{if} (cs%id\_gKEu > 0) \textcolor{keyword}{call }safe\_alloc\_ptr(ad%gradKEu,isdb,iedb,jsd,jed,nz)
1167 
1168   cs%id\_gKEv = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'gKEv'}, diag%axesCvL, time, &
1169      \textcolor{stringliteral}{'Meridional Acceleration from Grad. Kinetic Energy'}, \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1170   \textcolor{keywordflow}{if} (cs%id\_gKEv > 0) \textcolor{keyword}{call }safe\_alloc\_ptr(ad%gradKEv,isd,ied,jsdb,jedb,nz)
1171 
1172   cs%id\_rvxu = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'rvxu'}, diag%axesCvL, time, &
1173      \textcolor{stringliteral}{'Meridional Acceleration from Relative Vorticity'}, \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1174   \textcolor{keywordflow}{if} (cs%id\_rvxu > 0) \textcolor{keyword}{call }safe\_alloc\_ptr(ad%rv\_x\_u,isd,ied,jsdb,jedb,nz)
1175 
1176   cs%id\_rvxv = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'rvxv'}, diag%axesCuL, time, &
1177      \textcolor{stringliteral}{'Zonal Acceleration from Relative Vorticity'}, \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1178   \textcolor{keywordflow}{if} (cs%id\_rvxv > 0) \textcolor{keyword}{call }safe\_alloc\_ptr(ad%rv\_x\_v,isdb,iedb,jsd,jed,nz)
1179 
1180   \textcolor{comment}{!CS%id\_hf\_gKEu = register\_diag\_field('ocean\_model', 'hf\_gKEu', diag%axesCuL, Time, &}
1181   \textcolor{comment}{!   'Fractional Thickness-weighted Zonal Acceleration from Grad. Kinetic Energy', &}
1182   \textcolor{comment}{!   'm-1 s-2', v\_extensive=.true., conversion=US%L\_T2\_to\_m\_s2)}
1183   \textcolor{comment}{!if (CS%id\_hf\_gKEu > 0) then}
1184   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%gradKEu,IsdB,IedB,jsd,jed,nz)}
1185   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%diag\_hfrac\_u,IsdB,IedB,jsd,jed,nz)}
1186   \textcolor{comment}{!endif}
1187 
1188   \textcolor{comment}{!CS%id\_hf\_gKEv = register\_diag\_field('ocean\_model', 'hf\_gKEv', diag%axesCvL, Time, &}
1189   \textcolor{comment}{!   'Fractional Thickness-weighted Meridional Acceleration from Grad. Kinetic Energy', &}
1190   \textcolor{comment}{!   'm-1 s-2', v\_extensive=.true., conversion=US%L\_T2\_to\_m\_s2)}
1191   \textcolor{comment}{!if (CS%id\_hf\_gKEv > 0) then}
1192   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%gradKEv,isd,ied,JsdB,JedB,nz)}
1193   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%diag\_hfrac\_v,isd,ied,Jsd,JedB,nz)}
1194   \textcolor{comment}{!endif}
1195 
1196   cs%id\_hf\_gKEu\_2d = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'hf\_gKEu\_2d'}, diag%axesCu1, time, &
1197      \textcolor{stringliteral}{'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Grad. Kinetic Energy'}, &
1198      \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1199   \textcolor{keywordflow}{if} (cs%id\_hf\_gKEu\_2d > 0) \textcolor{keywordflow}{then}
1200     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%gradKEu,isdb,iedb,jsd,jed,nz)
1201     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%diag\_hfrac\_u,isdb,iedb,jsd,jed,nz)
1202 \textcolor{keywordflow}{  endif}
1203 
1204   cs%id\_hf\_gKEv\_2d = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'hf\_gKEv\_2d'}, diag%axesCv1, time, &
1205      \textcolor{stringliteral}{'Depth-sum Fractional Thickness-weighted Meridional Acceleration from Grad. Kinetic Energy'}, &
1206      \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1207   \textcolor{keywordflow}{if} (cs%id\_hf\_gKEv\_2d > 0) \textcolor{keywordflow}{then}
1208     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%gradKEv,isd,ied,jsdb,jedb,nz)
1209     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%diag\_hfrac\_v,isd,ied,jsd,jedb,nz)
1210 \textcolor{keywordflow}{  endif}
1211 
1212   \textcolor{comment}{!CS%id\_hf\_rvxu = register\_diag\_field('ocean\_model', 'hf\_rvxu', diag%axesCvL, Time, &}
1213   \textcolor{comment}{!   'Fractional Thickness-weighted Meridional Acceleration from Relative Vorticity', &}
1214   \textcolor{comment}{!   'm-1 s-2', v\_extensive=.true., conversion=US%L\_T2\_to\_m\_s2)}
1215   \textcolor{comment}{!if (CS%id\_hf\_rvxu > 0) then}
1216   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%rv\_x\_u,isd,ied,JsdB,JedB,nz)}
1217   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%diag\_hfrac\_v,isd,ied,Jsd,JedB,nz)}
1218   \textcolor{comment}{!endif}
1219 
1220   \textcolor{comment}{!CS%id\_hf\_rvxv = register\_diag\_field('ocean\_model', 'hf\_rvxv', diag%axesCuL, Time, &}
1221   \textcolor{comment}{!   'Fractional Thickness-weighted Zonal Acceleration from Relative Vorticity', &}
1222   \textcolor{comment}{!   'm-1 s-2', v\_extensive=.true., conversion=US%L\_T2\_to\_m\_s2)}
1223   \textcolor{comment}{!if (CS%id\_hf\_rvxv > 0) then}
1224   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%rv\_x\_v,IsdB,IedB,jsd,jed,nz)}
1225   \textcolor{comment}{!  call safe\_alloc\_ptr(AD%diag\_hfrac\_u,IsdB,IedB,jsd,jed,nz)}
1226   \textcolor{comment}{!endif}
1227 
1228   cs%id\_hf\_rvxu\_2d = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'hf\_rvxu\_2d'}, diag%axesCv1, time, &
1229      \textcolor{stringliteral}{'Fractional Thickness-weighted Meridional Acceleration from Relative Vorticity'}, &
1230      \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1231   \textcolor{keywordflow}{if} (cs%id\_hf\_rvxu\_2d > 0) \textcolor{keywordflow}{then}
1232     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%rv\_x\_u,isd,ied,jsdb,jedb,nz)
1233     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%diag\_hfrac\_v,isd,ied,jsd,jedb,nz)
1234 \textcolor{keywordflow}{  endif}
1235 
1236   cs%id\_hf\_rvxv\_2d = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'hf\_rvxv\_2d'}, diag%axesCu1, time, &
1237      \textcolor{stringliteral}{'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Relative Vorticity'}, &
1238      \textcolor{stringliteral}{'m-1 s-2'}, conversion=us%L\_T2\_to\_m\_s2)
1239   \textcolor{keywordflow}{if} (cs%id\_hf\_rvxv\_2d > 0) \textcolor{keywordflow}{then}
1240     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%rv\_x\_v,isdb,iedb,jsd,jed,nz)
1241     \textcolor{keyword}{call }safe\_alloc\_ptr(ad%diag\_hfrac\_u,isdb,iedb,jsd,jed,nz)
1242 \textcolor{keywordflow}{  endif}
1243 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__coriolisadv_a87e4a437552052fa238260442af19868}\label{namespacemom__coriolisadv_a87e4a437552052fa238260442af19868}} 
\index{mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}!gradke@{gradke}}
\index{gradke@{gradke}!mom\+\_\+coriolisadv@{mom\+\_\+coriolisadv}}
\subsubsection{\texorpdfstring{gradke()}{gradke()}}
{\footnotesize\ttfamily subroutine mom\+\_\+coriolisadv\+::gradke (\begin{DoxyParamCaption}\item[{real, dimension( g \%isdb\+: g \%iedb, g \%jsd\+: g \%jed, g \%ke), intent(in)}]{u,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsdb\+: g \%jedb, g \%ke), intent(in)}]{v,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed, g \%ke), intent(in)}]{h,  }\item[{real, dimension( g \%isd\+: g \%ied , g \%jsd\+: g \%jed ), intent(out)}]{KE,  }\item[{real, dimension( g \%isdb\+: g \%iedb, g \%jsd\+: g \%jed ), intent(out)}]{K\+Ex,  }\item[{real, dimension( g \%isd\+: g \%ied , g \%jsdb\+: g \%jedb), intent(out)}]{K\+Ey,  }\item[{integer, intent(in)}]{k,  }\item[{type(ocean\+\_\+obc\+\_\+type), pointer}]{O\+BC,  }\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{type(\hyperlink{structmom__coriolisadv_1_1coriolisadv__cs}{coriolisadv\+\_\+cs}), pointer}]{CS }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



Calculates the acceleration due to the gradient of kinetic energy. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em g} & Ocen grid structure\\
\hline
\mbox{\tt in}  & {\em u} & Zonal velocity \mbox{[}L T-\/1 $\sim$$>$ m s-\/1\mbox{]}\\
\hline
\mbox{\tt in}  & {\em v} & Meridional velocity \mbox{[}L T-\/1 $\sim$$>$ m s-\/1\mbox{]}\\
\hline
\mbox{\tt in}  & {\em h} & Layer thickness \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}\\
\hline
\mbox{\tt out}  & {\em ke} & Kinetic energy per unit mass \mbox{[}L2 T-\/2 $\sim$$>$ m2 s-\/2\mbox{]}\\
\hline
\mbox{\tt out}  & {\em kex} & Zonal acceleration due to kinetic energy gradient \mbox{[}L T-\/2 $\sim$$>$ m s-\/2\mbox{]}\\
\hline
\mbox{\tt out}  & {\em key} & Meridional acceleration due to kinetic energy gradient \mbox{[}L T-\/2 $\sim$$>$ m s-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em k} & Layer number to calculate for\\
\hline
 & {\em obc} & Open boundary control structure\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
 & {\em cs} & Control structure for M\+O\+M\+\_\+\+Coriolis\+Adv \\
\hline
\end{DoxyParams}


Definition at line 928 of file M\+O\+M\+\_\+\+Coriolis\+Adv.\+F90.


\begin{DoxyCode}
928   \textcolor{keywordtype}{type}(ocean\_grid\_type),                      \textcolor{keywordtype}{intent(in)}  :: g\textcolor{comment}{ !< Ocen grid structure}
929   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G),SZK\_(G))},  \textcolor{keywordtype}{intent(in)}  :: u\textcolor{comment}{ !< Zonal velocity [L T-1 ~> m s-1]}
930   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJB\_(G),SZK\_(G))},  \textcolor{keywordtype}{intent(in)}  :: v\textcolor{comment}{ !< Meridional velocity [L T-1 ~> m s-1]}
931   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(G))},   \textcolor{keywordtype}{intent(in)}  :: h\textcolor{comment}{ !< Layer thickness [H ~> m or kg m-2]}
932   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G) ,SZJ\_(G) )},         \textcolor{keywordtype}{intent(out)} :: ke\textcolor{comment}{ !< Kinetic energy per unit mass [L2 T-2 ~>
       m2 s-2]}
933   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZIB\_(G),SZJ\_(G) )},         \textcolor{keywordtype}{intent(out)} :: kex\textcolor{comment}{ !< Zonal acceleration due to kinetic}
934 \textcolor{comment}{                                                                 !! energy gradient [L T-2 ~> m s-2]}
935   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G) ,SZJB\_(G))},         \textcolor{keywordtype}{intent(out)} :: key\textcolor{comment}{ !< Meridional acceleration due to kinetic}
936 \textcolor{comment}{                                                                 !! energy gradient [L T-2 ~> m s-2]}
937   \textcolor{keywordtype}{integer},                                    \textcolor{keywordtype}{intent(in)}  :: k\textcolor{comment}{ !< Layer number to calculate for}
938   \textcolor{keywordtype}{type}(ocean\_obc\_type),                       \textcolor{keywordtype}{pointer}     :: obc\textcolor{comment}{ !< Open boundary control structure}
939   \textcolor{keywordtype}{type}(unit\_scale\_type),                      \textcolor{keywordtype}{intent(in)}  :: us\textcolor{comment}{  !< A dimensional unit scaling type}
940   \textcolor{keywordtype}{type}(coriolisadv\_cs),                       \textcolor{keywordtype}{pointer}     :: cs\textcolor{comment}{ !< Control structure for MOM\_CoriolisAdv}
941   \textcolor{comment}{! Local variables}
942   \textcolor{keywordtype}{real} :: um, up, vm, vp         \textcolor{comment}{! Temporary variables [L T-1 ~> m s-1].}
943   \textcolor{keywordtype}{real} :: um2, up2, vm2, vp2     \textcolor{comment}{! Temporary variables [L2 T-2 ~> m2 s-2].}
944   \textcolor{keywordtype}{real} :: um2a, up2a, vm2a, vp2a \textcolor{comment}{! Temporary variables [L4 T-2 ~> m4 s-2].}
945   \textcolor{keywordtype}{integer} :: i, j, is, ie, js, je, isq, ieq, jsq, jeq, nz, n
946 
947   is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
948   isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
949 
950 
951   \textcolor{comment}{! Calculate KE (Kinetic energy for use in the -grad(KE) acceleration term).}
952   \textcolor{keywordflow}{if} (cs%KE\_Scheme == ke\_arakawa) \textcolor{keywordflow}{then}
953     \textcolor{comment}{! The following calculation of Kinetic energy includes the metric terms}
954     \textcolor{comment}{! identified in Arakawa & Lamb 1982 as important for KE conservation.  It}
955     \textcolor{comment}{! also includes the possibility of partially-blocked tracer cell faces.}
956     \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=isq,ieq+1
957       ke(i,j) = ( ( g%areaCu( i ,j)*(u( i ,j,k)*u( i ,j,k)) + &
958                     g%areaCu(i-1,j)*(u(i-1,j,k)*u(i-1,j,k)) ) + &
959                   ( g%areaCv(i, j )*(v(i, j ,k)*v(i, j ,k)) + &
960                     g%areaCv(i,j-1)*(v(i,j-1,k)*v(i,j-1,k)) ) )*0.25*g%IareaT(i,j)
961 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
962   \textcolor{keywordflow}{elseif} (cs%KE\_Scheme == ke\_simple\_gudonov) \textcolor{keywordflow}{then}
963     \textcolor{comment}{! The following discretization of KE is based on the one-dimensinal Gudonov}
964     \textcolor{comment}{! scheme which does not take into account any geometric factors}
965     \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=isq,ieq+1
966       up = 0.5*( u(i-1,j,k) + abs( u(i-1,j,k) ) ) ; up2 = up*up
967       um = 0.5*( u( i ,j,k) - abs( u( i ,j,k) ) ) ; um2 = um*um
968       vp = 0.5*( v(i,j-1,k) + abs( v(i,j-1,k) ) ) ; vp2 = vp*vp
969       vm = 0.5*( v(i, j ,k) - abs( v(i, j ,k) ) ) ; vm2 = vm*vm
970       ke(i,j) = ( max(up2,um2) + max(vp2,vm2) ) *0.5
971 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
972   \textcolor{keywordflow}{elseif} (cs%KE\_Scheme == ke\_gudonov) \textcolor{keywordflow}{then}
973     \textcolor{comment}{! The following discretization of KE is based on the one-dimensinal Gudonov}
974     \textcolor{comment}{! scheme but has been adapted to take horizontal grid factors into account}
975     \textcolor{keywordflow}{do} j=jsq,jeq+1 ; \textcolor{keywordflow}{do} i=isq,ieq+1
976       up = 0.5*( u(i-1,j,k) + abs( u(i-1,j,k) ) ) ; up2a = up*up*g%areaCu(i-1,j)
977       um = 0.5*( u( i ,j,k) - abs( u( i ,j,k) ) ) ; um2a = um*um*g%areaCu( i ,j)
978       vp = 0.5*( v(i,j-1,k) + abs( v(i,j-1,k) ) ) ; vp2a = vp*vp*g%areaCv(i,j-1)
979       vm = 0.5*( v(i, j ,k) - abs( v(i, j ,k) ) ) ; vm2a = vm*vm*g%areaCv(i, j )
980       ke(i,j) = ( max(um2a,up2a) + max(vm2a,vp2a) )*0.5*g%IareaT(i,j)
981 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
982 \textcolor{keywordflow}{  endif}
983 
984   \textcolor{comment}{! Term - d(KE)/dx.}
985   \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=isq,ieq
986     kex(i,j) = (ke(i+1,j) - ke(i,j)) * g%IdxCu(i,j)
987 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo}
988 
989   \textcolor{comment}{! Term - d(KE)/dy.}
990   \textcolor{keywordflow}{do} j=jsq,jeq ; \textcolor{keywordflow}{do} i=is,ie
991     key(i,j) = (ke(i,j+1) - ke(i,j)) * g%IdyCv(i,j)
992 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo}
993 
994   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(obc)) \textcolor{keywordflow}{then}
995     \textcolor{keywordflow}{do} n=1,obc%number\_of\_segments
996       \textcolor{keywordflow}{if} (obc%segment(n)%is\_N\_or\_S) \textcolor{keywordflow}{then}
997         \textcolor{keywordflow}{do} i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
998           key(i,obc%segment(n)%HI%JsdB) = 0.
999 \textcolor{keywordflow}{        enddo}
1000       \textcolor{keywordflow}{elseif} (obc%segment(n)%is\_E\_or\_W) \textcolor{keywordflow}{then}
1001         \textcolor{keywordflow}{do} j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
1002           kex(obc%segment(n)%HI%IsdB,j) = 0.
1003 \textcolor{keywordflow}{        enddo}
1004 \textcolor{keywordflow}{      endif}
1005 \textcolor{keywordflow}{    enddo}
1006 \textcolor{keywordflow}{  endif}
1007 
\end{DoxyCode}
