\hypertarget{namespacemom__opacity}{}\section{mom\+\_\+opacity Module Reference}
\label{namespacemom__opacity}\index{mom\+\_\+opacity@{mom\+\_\+opacity}}


\subsection{Detailed Description}
Routines used to calculate the opacity of the ocean. 

opacity\+\_\+from\+\_\+chl\+: In this routine, the Morel (modified) or Manizza (modified) schemes use the \char`\"{}blue\char`\"{} band in the paramterizations to determine the e-\/folding depth of the incoming shortwave attenuation. The red portion is lumped into the net heating at the surface.

Morel, A., 1988\+: Optical modeling of the upper ocean in relation to its biogenous matter content (case-\/i waters)., J. Geo. Res., 93, 10,749-\/10,768.

Manizza, M., C. Le\+Quere, A. J. Watson, and E. T. Buitenhuis, 2005\+: Bio-\/optical feedbacks among phytoplankton, upper ocean physics and sea-\/ice in a global model, Geophys. Res. Let., 32, L05603, doi\+:10.\+1029/2004\+G\+L020778. \subsection*{Data Types}
\begin{DoxyCompactItemize}
\item 
type \mbox{\hyperlink{structmom__opacity_1_1opacity__cs}{opacity\+\_\+cs}}
\begin{DoxyCompactList}\small\item\em The control structure with paramters for the M\+O\+M\+\_\+opacity module. \end{DoxyCompactList}\item 
type \mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}
\begin{DoxyCompactList}\small\item\em This type is used to store information about ocean optical properties. \end{DoxyCompactList}\end{DoxyCompactItemize}
\subsection*{Functions/\+Subroutines}
\begin{DoxyCompactItemize}
\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_a05ef9c5d86adff869fad832f3083bba4}{set\+\_\+opacity}} (optics, sw\+\_\+total, sw\+\_\+vis\+\_\+dir, sw\+\_\+vis\+\_\+dif, sw\+\_\+nir\+\_\+dir, sw\+\_\+nir\+\_\+dif, G, GV, US, CS, chl\+\_\+2d, chl\+\_\+3d)
\begin{DoxyCompactList}\small\item\em This sets the opacity of sea water based based on one of several different schemes. \end{DoxyCompactList}\item 
subroutine \mbox{\hyperlink{namespacemom__opacity_aa14b4bfe45861d86740c24e474ce682c}{opacity\+\_\+from\+\_\+chl}} (optics, sw\+\_\+total, sw\+\_\+vis\+\_\+dir, sw\+\_\+vis\+\_\+dif, sw\+\_\+nir\+\_\+dir, sw\+\_\+nir\+\_\+dif, G, GV, US, CS, chl\+\_\+2d, chl\+\_\+3d)
\begin{DoxyCompactList}\small\item\em This sets the \char`\"{}blue\char`\"{} band opacity based on chloophyll A concencentrations The red portion is lumped into the net heating at the surface. \end{DoxyCompactList}\item 
real function, public \mbox{\hyperlink{namespacemom__opacity_a4498b4bb6fcf1b7d849f89aa87c0332e}{opacity\+\_\+morel}} (chl\+\_\+data)
\begin{DoxyCompactList}\small\item\em This sets the blue-\/wavelength opacity according to the scheme proposed by Morel and Antoine (1994). \end{DoxyCompactList}\item 
real function \mbox{\hyperlink{namespacemom__opacity_a0017241c03e4536115674fc5fc9608bf}{sw\+\_\+pen\+\_\+frac\+\_\+morel}} (chl\+\_\+data)
\begin{DoxyCompactList}\small\item\em This sets the penetrating shortwave fraction according to the scheme proposed by Morel and Antoine (1994). \end{DoxyCompactList}\item 
real function, public \mbox{\hyperlink{namespacemom__opacity_a27873e3942c2ff7980cf516944f78f03}{opacity\+\_\+manizza}} (chl\+\_\+data)
\begin{DoxyCompactList}\small\item\em This sets the blue-\/wavelength opacity according to the scheme proposed by Manizza, M. et al, 2005. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_a4c1942f798619a9ad854d1152ebcab63}{extract\+\_\+optics\+\_\+slice}} (optics, j, G, GV, opacity, opacity\+\_\+scale, pen\+S\+W\+\_\+top, pen\+S\+W\+\_\+scale)
\begin{DoxyCompactList}\small\item\em This subroutine returns a 2-\/d slice at constant j of fields from an \mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}, with the potential for rescaling these fields. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_aefdfc303272a4f4fc07f84a8aea2b0f1}{extract\+\_\+optics\+\_\+fields}} (optics, nbands)
\begin{DoxyCompactList}\small\item\em Set arguments to fields from the optics type. \end{DoxyCompactList}\item 
integer function, public \mbox{\hyperlink{namespacemom__opacity_a349c6934f113d238e4e2ef229b931a0c}{optics\+\_\+nbands}} (optics)
\begin{DoxyCompactList}\small\item\em Return the number of bands of penetrating shortwave radiation. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_a21db9da24cea8b875040ba1e7e8b2e9b}{absorbremainingsw}} (G, GV, US, h, opacity\+\_\+band, nsw, optics, j, dt, H\+\_\+limit\+\_\+fluxes, adjust\+Absorption\+Profile, absorb\+All\+SW, T, Pen\+\_\+\+S\+W\+\_\+bnd, eps, ksort, htot, Ttot, T\+KE, d\+S\+V\+\_\+dT)
\begin{DoxyCompactList}\small\item\em Apply shortwave heating below the boundary layer (when running with the bulk mixed layer inhereted from G\+O\+LD) or throughout the water column. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_ad27db4bd0d010d98a3f5a54902c7a05e}{sumswoverbands}} (G, GV, US, h, nsw, optics, j, dt, H\+\_\+limit\+\_\+fluxes, absorb\+All\+SW, i\+Pen\+\_\+\+S\+W\+\_\+bnd, net\+Pen)
\begin{DoxyCompactList}\small\item\em This subroutine calculates the total shortwave heat flux integrated over bands as a function of depth. This routine is only called for computing buoyancy fluxes for use in K\+PP. This routine does not updat e the state. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_a39fce7bd33a469e3e9fe7cfeb51825b5}{opacity\+\_\+init}} (Time, G, GV, US, param\+\_\+file, diag, CS, optics)
\begin{DoxyCompactList}\small\item\em This routine initalizes the opacity module, including an \mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__opacity_ab5c0caabbf8a806a95bcc416da673841}{opacity\+\_\+end}} (CS, optics)
\end{DoxyCompactItemize}
\subsection*{Variables}
\begin{DoxyCompactItemize}
\item 
\mbox{\Hypertarget{namespacemom__opacity_aecb5a81b64d808f474cb4e5c722dc76b}\label{namespacemom__opacity_aecb5a81b64d808f474cb4e5c722dc76b}} 
character $\ast$(10), parameter \mbox{\hyperlink{namespacemom__opacity_aecb5a81b64d808f474cb4e5c722dc76b}{manizza\+\_\+05\+\_\+string}} = \char`\"{}M\+A\+N\+I\+Z\+Z\+A\+\_\+05\char`\"{}
\begin{DoxyCompactList}\small\item\em String to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_a5e5e40ad5acbf12058294f2ea17a7368}\label{namespacemom__opacity_a5e5e40ad5acbf12058294f2ea17a7368}} 
character $\ast$(10), parameter \mbox{\hyperlink{namespacemom__opacity_a5e5e40ad5acbf12058294f2ea17a7368}{morel\+\_\+88\+\_\+string}} = \char`\"{}M\+O\+R\+E\+L\+\_\+88\char`\"{}
\begin{DoxyCompactList}\small\item\em String to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_ad430c567c09952e6ecfd07d6c05c8f69}\label{namespacemom__opacity_ad430c567c09952e6ecfd07d6c05c8f69}} 
character $\ast$(10), parameter \mbox{\hyperlink{namespacemom__opacity_ad430c567c09952e6ecfd07d6c05c8f69}{single\+\_\+exp\+\_\+string}} = \char`\"{}S\+I\+N\+G\+L\+E\+\_\+\+E\+XP\char`\"{}
\begin{DoxyCompactList}\small\item\em String to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_a5a3d62dfb71c22ac0b13a312c8d23c70}\label{namespacemom__opacity_a5a3d62dfb71c22ac0b13a312c8d23c70}} 
character $\ast$(10), parameter \mbox{\hyperlink{namespacemom__opacity_a5a3d62dfb71c22ac0b13a312c8d23c70}{double\+\_\+exp\+\_\+string}} = \char`\"{}D\+O\+U\+B\+L\+E\+\_\+\+E\+XP\char`\"{}
\begin{DoxyCompactList}\small\item\em String to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_a9e2e33c15dd65f100263da9efc0934a6}\label{namespacemom__opacity_a9e2e33c15dd65f100263da9efc0934a6}} 
real, parameter \mbox{\hyperlink{namespacemom__opacity_a9e2e33c15dd65f100263da9efc0934a6}{op\+\_\+diag\+\_\+len}} = 1e-\/10
\begin{DoxyCompactList}\small\item\em Lengthscale L used to remap opacity from op to 1/L $\ast$ tanh(op $\ast$ L) \end{DoxyCompactList}\end{DoxyCompactItemize}
\textbf{ }\par
\begin{DoxyCompactItemize}
\item 
\mbox{\Hypertarget{namespacemom__opacity_a948b6b8f52bd40409385ee13f68c9626}\label{namespacemom__opacity_a948b6b8f52bd40409385ee13f68c9626}} 
integer, parameter \mbox{\hyperlink{namespacemom__opacity_a948b6b8f52bd40409385ee13f68c9626}{no\+\_\+scheme}} = 0
\begin{DoxyCompactList}\small\item\em Coded integers to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_a28a3fe4f0c5e8c81882449450f6ba785}\label{namespacemom__opacity_a28a3fe4f0c5e8c81882449450f6ba785}} 
integer, parameter \mbox{\hyperlink{namespacemom__opacity_a28a3fe4f0c5e8c81882449450f6ba785}{manizza\+\_\+05}} = 1
\begin{DoxyCompactList}\small\item\em Coded integers to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_aa47b1866014e588f0306199f7ea7e73d}\label{namespacemom__opacity_aa47b1866014e588f0306199f7ea7e73d}} 
integer, parameter \mbox{\hyperlink{namespacemom__opacity_aa47b1866014e588f0306199f7ea7e73d}{morel\+\_\+88}} = 2
\begin{DoxyCompactList}\small\item\em Coded integers to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_a77646ea3d3da72b4e00694b7a121a771}\label{namespacemom__opacity_a77646ea3d3da72b4e00694b7a121a771}} 
integer, parameter \mbox{\hyperlink{namespacemom__opacity_a77646ea3d3da72b4e00694b7a121a771}{single\+\_\+exp}} = 3
\begin{DoxyCompactList}\small\item\em Coded integers to specify the opacity scheme. \end{DoxyCompactList}\item 
\mbox{\Hypertarget{namespacemom__opacity_aba8c34bcba1b2f2d80a48cd7bce261ec}\label{namespacemom__opacity_aba8c34bcba1b2f2d80a48cd7bce261ec}} 
integer, parameter \mbox{\hyperlink{namespacemom__opacity_aba8c34bcba1b2f2d80a48cd7bce261ec}{double\+\_\+exp}} = 4
\begin{DoxyCompactList}\small\item\em Coded integers to specify the opacity scheme. \end{DoxyCompactList}\end{DoxyCompactItemize}



\subsection{Function/\+Subroutine Documentation}
\mbox{\Hypertarget{namespacemom__opacity_a21db9da24cea8b875040ba1e7e8b2e9b}\label{namespacemom__opacity_a21db9da24cea8b875040ba1e7e8b2e9b}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!absorbremainingsw@{absorbremainingsw}}
\index{absorbremainingsw@{absorbremainingsw}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{absorbremainingsw()}{absorbremainingsw()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::absorbremainingsw (\begin{DoxyParamCaption}\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke), intent(in)}]{h,  }\item[{real, dimension(max(1,nsw), g \%isd\+: g \%ied, gv \%ke), intent(in)}]{opacity\+\_\+band,  }\item[{integer, intent(in)}]{nsw,  }\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(in)}]{optics,  }\item[{integer, intent(in)}]{j,  }\item[{real, intent(in)}]{dt,  }\item[{real, intent(in)}]{H\+\_\+limit\+\_\+fluxes,  }\item[{logical, intent(in)}]{adjust\+Absorption\+Profile,  }\item[{logical, intent(in)}]{absorb\+All\+SW,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke), intent(inout)}]{T,  }\item[{real, dimension(max(1,nsw), g \%isd\+: g \%ied), intent(inout)}]{Pen\+\_\+\+S\+W\+\_\+bnd,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke), intent(in), optional}]{eps,  }\item[{integer, dimension( g \%isd\+: g \%ied, gv \%ke), intent(in), optional}]{ksort,  }\item[{real, dimension( g \%isd\+: g \%ied), intent(in), optional}]{htot,  }\item[{real, dimension( g \%isd\+: g \%ied), intent(inout), optional}]{Ttot,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke), intent(inout), optional}]{T\+KE,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke), intent(in), optional}]{d\+S\+V\+\_\+dT }\end{DoxyParamCaption})}



Apply shortwave heating below the boundary layer (when running with the bulk mixed layer inhereted from G\+O\+LD) or throughout the water column. 

In addition, it causes all of the remaining SW radiation to be absorbed, provided that the total water column thickness is greater than H\+\_\+limit\+\_\+fluxes. For thinner water columns, the heating is scaled down proportionately, the assumption being that the remaining heating (which is left in Pen\+\_\+\+SW) should go into an (absent for now) ocean bottom sediment layer.


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
\mbox{\tt in}  & {\em nsw} & Number of bands of penetrating shortwave radiation.\\
\hline
\mbox{\tt in}  & {\em h} & Layer thicknesses \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em opacity\+\_\+band} & Opacity in each band of penetrating shortwave radiation \mbox{[}H-\/1 $\sim$$>$ m-\/1 or m2 kg-\/1\mbox{]}. The indicies are band, i, k.\\
\hline
\mbox{\tt in}  & {\em optics} & An optics structure that has values of opacities and shortwave fluxes.\\
\hline
\mbox{\tt in}  & {\em j} & j-\/index to work on.\\
\hline
\mbox{\tt in}  & {\em dt} & Time step \mbox{[}T $\sim$$>$ s\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em h\+\_\+limit\+\_\+fluxes} & If the total ocean depth is less than this, they are scaled away to avoid numerical instabilities \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}. This would not be necessary if a finite heat capacity mud-\/layer were added.\\
\hline
\mbox{\tt in}  & {\em adjustabsorptionprofile} & If true, apply heating above the layers in which it should have occurred to get the correct mean depth (and potential energy change) of the shortwave that should be absorbed by each layer.\\
\hline
\mbox{\tt in}  & {\em absorballsw} & If true, apply heating above the layers in which it should have occurred to get the correct mean depth (and potential energy change) of the shortwave that should be absorbed by each layer.\\
\hline
\mbox{\tt in,out}  & {\em t} & Layer potential/conservative temperatures \mbox{[}degC\mbox{]}\\
\hline
\mbox{\tt in,out}  & {\em pen\+\_\+sw\+\_\+bnd} & Penetrating shortwave heating in each band that hits the bottom and will will be redistributed through the water column \mbox{[}degC H $\sim$$>$ degC m or degC kg m-\/2\mbox{]}, size nsw x G isd\+: G ied.\\
\hline
\mbox{\tt in}  & {\em eps} & Small thickness that must remain in each layer, and which will not be subject to heating \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em ksort} & Density-\/sorted k-\/indicies.\\
\hline
\mbox{\tt in}  & {\em htot} & Total mixed layer thickness \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}.\\
\hline
\mbox{\tt in,out}  & {\em ttot} & Depth integrated mixed layer temperature \mbox{[}degC H $\sim$$>$ degC m or degC kg m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em dsv\+\_\+dt} & The partial derivative of specific volume with temperature \mbox{[}R-\/1 deg\+C-\/1\mbox{]}.\\
\hline
\mbox{\tt in,out}  & {\em tke} & The T\+KE sink from mixing the heating throughout a layer \mbox{[}R Z3 T-\/2 $\sim$$>$ J m-\/2\mbox{]}. \\
\hline
\end{DoxyParams}


Definition at line 511 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
511 
512   \textcolor{keywordtype}{type}(ocean\_grid\_type),             \textcolor{keywordtype}{intent(in)}    :: G\textcolor{comment}{    !< The ocean's grid structure.}
513   \textcolor{keywordtype}{type}(verticalGrid\_type),           \textcolor{keywordtype}{intent(in)}    :: GV\textcolor{comment}{   !< The ocean's vertical grid structure.}
514   \textcolor{keywordtype}{type}(unit\_scale\_type),             \textcolor{keywordtype}{intent(in)}    :: US\textcolor{comment}{   !< A dimensional unit scaling type}
515   \textcolor{keywordtype}{integer},                           \textcolor{keywordtype}{intent(in)}    :: nsw\textcolor{comment}{  !< Number of bands of penetrating}
516 \textcolor{comment}{                                                           !! shortwave radiation.}
517   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{intent(in)}    :: h\textcolor{comment}{    !< Layer thicknesses [H ~> m or kg m-2].}
518   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(max(1,nsw),SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{intent(in)} :: opacity\_band\textcolor{comment}{ !< Opacity in each band of
       penetrating}
519 \textcolor{comment}{                                                           !! shortwave radiation [H-1 ~> m-1 or m2 kg-1].}
520 \textcolor{comment}{                                                           !! The indicies are band, i, k.}
521   \textcolor{keywordtype}{type}(optics\_type),                 \textcolor{keywordtype}{intent(in)}    :: optics\textcolor{comment}{ !< An optics structure that has values of}
522 \textcolor{comment}{                                                           !! opacities and shortwave fluxes.}
523   \textcolor{keywordtype}{integer},                           \textcolor{keywordtype}{intent(in)}    :: j\textcolor{comment}{    !< j-index to work on.}
524   \textcolor{keywordtype}{real},                              \textcolor{keywordtype}{intent(in)}    :: dt\textcolor{comment}{   !< Time step [T ~> s].}
525   \textcolor{keywordtype}{real},                              \textcolor{keywordtype}{intent(in)}    :: H\_limit\_fluxes\textcolor{comment}{ !< If the total ocean depth is}
526 \textcolor{comment}{                                                           !! less than this, they are scaled away}
527 \textcolor{comment}{                                                           !! to avoid numerical instabilities}
528 \textcolor{comment}{                                                           !! [H ~> m or kg m-2]. This would}
529 \textcolor{comment}{                                                           !! not be necessary if a finite heat}
530 \textcolor{comment}{                                                           !! capacity mud-layer were added.}
531   \textcolor{keywordtype}{logical},                          \textcolor{keywordtype}{intent(in)}    :: adjustAbsorptionProfile\textcolor{comment}{ !< If true, apply}
532 \textcolor{comment}{                                                           !! heating above the layers in which it}
533 \textcolor{comment}{                                                           !! should have occurred to get the}
534 \textcolor{comment}{                                                           !! correct mean depth (and potential}
535 \textcolor{comment}{                                                           !! energy change) of the shortwave that}
536 \textcolor{comment}{                                                           !! should be absorbed by each layer.}
537   \textcolor{keywordtype}{logical},                          \textcolor{keywordtype}{intent(in)}    :: absorbAllSW\textcolor{comment}{ !< If true, apply heating above the}
538 \textcolor{comment}{                                                           !! layers in which it should have occurred}
539 \textcolor{comment}{                                                           !! to get the correct mean depth (and}
540 \textcolor{comment}{                                                           !! potential energy change) of the}
541 \textcolor{comment}{                                                           !! shortwave that should be absorbed by}
542 \textcolor{comment}{                                                           !! each layer.}
543   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{intent(inout)} :: T\textcolor{comment}{    !< Layer potential/conservative}
544 \textcolor{comment}{                                                           !! temperatures [degC]}
545   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(max(1,nsw),SZI\_(G))}, \textcolor{keywordtype}{intent(inout)} :: Pen\_SW\_bnd\textcolor{comment}{ !< Penetrating shortwave heating in}
546 \textcolor{comment}{                                                           !! each band that hits the bottom and will}
547 \textcolor{comment}{                                                           !! will be redistributed through the water}
548 \textcolor{comment}{                                                           !! column [degC H ~> degC m or degC kg m-2],}
549 \textcolor{comment}{                                                           !! size nsw x SZI\_(G).}
550   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)} :: eps\textcolor{comment}{ !< Small thickness that must remain in}
551 \textcolor{comment}{                                                           !! each layer, and which will not be}
552 \textcolor{comment}{                                                           !! subject to heating [H ~> m or kg m-2]}
553   \textcolor{keywordtype}{integer}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)} :: ksort\textcolor{comment}{ !< Density-sorted k-indicies.}
554   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}    :: htot\textcolor{comment}{ !< Total mixed layer thickness [H ~> m or kg
       m-2].}
555   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(inout)} :: Ttot\textcolor{comment}{ !< Depth integrated mixed layer}
556 \textcolor{comment}{                                                           !! temperature [degC H ~> degC m or degC kg m-2]}
557   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)} :: dSV\_dT\textcolor{comment}{ !< The partial derivative of specific}
558 \textcolor{comment}{                                                           !! volume with temperature [R-1 degC-1].}
559   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(inout)} :: TKE\textcolor{comment}{ !< The TKE sink from mixing the heating}
560 \textcolor{comment}{                                                           !! throughout a layer [R Z3 T-2 ~> J m-2].}
561 
562   \textcolor{comment}{! Local variables}
563   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))} :: &
564     T\_chg\_above    \textcolor{comment}{! A temperature change that will be applied to all the thick}
565                    \textcolor{comment}{! layers above a given layer [degC].  This is only nonzero if}
566                    \textcolor{comment}{! adjustAbsorptionProfile is true, in which case the net}
567                    \textcolor{comment}{! change in the temperature of a layer is the sum of the}
568                    \textcolor{comment}{! direct heating of that layer plus T\_chg\_above from all of}
569                    \textcolor{comment}{! the layers below, plus any contribution from absorbing}
570                    \textcolor{comment}{! radiation that hits the bottom.}
571   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G))} :: &
572     h\_heat, &      \textcolor{comment}{! The thickness of the water column that will be heated by}
573                    \textcolor{comment}{! any remaining shortwave radiation [H ~> m or kg m-2].}
574     t\_chg, &       \textcolor{comment}{! The temperature change of thick layers due to the remaining}
575                    \textcolor{comment}{! shortwave radiation and contributions from T\_chg\_above [degC].}
576     pen\_sw\_rem     \textcolor{comment}{! The sum across all wavelength bands of the penetrating shortwave}
577                    \textcolor{comment}{! heating that hits the bottom and will be redistributed through}
578                    \textcolor{comment}{! the water column [degC H ~> degC m or degC kg m-2]}
579   \textcolor{keywordtype}{real} :: SW\_trans          \textcolor{comment}{! fraction of shortwave radiation that is not}
580                             \textcolor{comment}{! absorbed in a layer [nondim]}
581   \textcolor{keywordtype}{real} :: unabsorbed        \textcolor{comment}{! fraction of the shortwave radiation that}
582                             \textcolor{comment}{! is not absorbed because the layers are too thin}
583   \textcolor{keywordtype}{real} :: Ih\_limit          \textcolor{comment}{! inverse of the total depth at which the}
584                             \textcolor{comment}{! surface fluxes start to be limited [H-1 ~> m-1 or m2 kg-1]}
585   \textcolor{keywordtype}{real} :: h\_min\_heat        \textcolor{comment}{! minimum thickness layer that should get heated [H ~> m or kg m-2]}
586   \textcolor{keywordtype}{real} :: opt\_depth         \textcolor{comment}{! optical depth of a layer [nondim]}
587   \textcolor{keywordtype}{real} :: exp\_OD            \textcolor{comment}{! exp(-opt\_depth) [nondim]}
588   \textcolor{keywordtype}{real} :: heat\_bnd          \textcolor{comment}{! heating due to absorption in the current}
589                             \textcolor{comment}{! layer by the current band, including any piece that}
590                             \textcolor{comment}{! is moved upward [degC H ~> degC m or degC kg m-2]}
591   \textcolor{keywordtype}{real} :: SWa               \textcolor{comment}{! fraction of the absorbed shortwave that is}
592                             \textcolor{comment}{! moved to layers above with adjustAbsorptionProfile [nondim]}
593   \textcolor{keywordtype}{real} :: coSWa\_frac        \textcolor{comment}{! The fraction of SWa that is actually moved upward.}
594   \textcolor{keywordtype}{real} :: min\_SW\_heat       \textcolor{comment}{! A minimum remaining shortwave heating within a timestep that will be simply}
595                             \textcolor{comment}{! absorbed in the next layer for computational efficiency, instead of}
596                             \textcolor{comment}{! continuing to penetrate [degC H ~> degC m or degC kg m-2].}
597   \textcolor{keywordtype}{real} :: I\_Habs            \textcolor{comment}{! The inverse of the absorption length for a minimal flux [H-1 ~> m-1 or m2
       kg-1]}
598   \textcolor{keywordtype}{real} :: epsilon           \textcolor{comment}{! A small thickness that must remain in each}
599                             \textcolor{comment}{! layer, and which will not be subject to heating [H ~> m or kg m-2]}
600   \textcolor{keywordtype}{real} :: g\_Hconv2          \textcolor{comment}{! A conversion factor for use in the TKE calculation}
601                             \textcolor{comment}{! in units of [Z3 R2 T-2 H-2 ~> kg2 m-5 s-2 or m s-2].}
602   \textcolor{keywordtype}{logical} :: SW\_Remains     \textcolor{comment}{! If true, some column has shortwave radiation that}
603                             \textcolor{comment}{! was not entirely absorbed.}
604   \textcolor{keywordtype}{logical} :: TKE\_calc       \textcolor{comment}{! If true, calculate the implications to the}
605                             \textcolor{comment}{! TKE budget of the shortwave heating.}
606   \textcolor{keywordtype}{real} :: C1\_6, C1\_60
607   \textcolor{keywordtype}{integer} :: is, ie, nz, i, k, ks, n
608   sw\_remains = .false.
609 
610   min\_sw\_heat = optics%PenSW\_flux\_absorb * dt
611   i\_habs = optics%PenSW\_absorb\_Invlen
612 
613   h\_min\_heat = 2.0*gv%Angstrom\_H + gv%H\_subroundoff
614   is = g%isc ; ie = g%iec ; nz = g%ke
615   c1\_6 = 1.0 / 6.0 ; c1\_60 = 1.0 / 60.0
616 
617   tke\_calc = (\textcolor{keyword}{present}(tke) .and. \textcolor{keyword}{present}(dsv\_dt))
618 
619   \textcolor{keywordflow}{if} (optics%answers\_2018) \textcolor{keywordflow}{then}
620     g\_hconv2 = (us%L\_to\_Z**2*gv%g\_Earth * gv%H\_to\_RZ) * gv%H\_to\_RZ
621   \textcolor{keywordflow}{else}
622     g\_hconv2 = us%L\_to\_Z**2*gv%g\_Earth * gv%H\_to\_RZ**2
623 \textcolor{keywordflow}{  endif}
624 
625   h\_heat(:) = 0.0
626   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(htot)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} i=is,ie ; h\_heat(i) = htot(i) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ endif}
627 
628   \textcolor{comment}{! Apply penetrating SW radiation to remaining parts of layers.}
629   \textcolor{comment}{! Excessively thin layers are not heated to avoid runaway temps.}
630   \textcolor{keywordflow}{do} ks=1,nz ; \textcolor{keywordflow}{do} i=is,ie
631     k = ks
632     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(ksort)) \textcolor{keywordflow}{then}
633       \textcolor{keywordflow}{if} (ksort(i,ks) <= 0) cycle
634       k = ksort(i,ks)
635 \textcolor{keywordflow}{    endif}
636     epsilon = 0.0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(eps)) epsilon = eps(i,k)
637 
638     t\_chg\_above(i,k) = 0.0
639 
640     \textcolor{keywordflow}{if} (h(i,k) > 1.5*epsilon) \textcolor{keywordflow}{then}
641       \textcolor{keywordflow}{do} n=1,nsw ; \textcolor{keywordflow}{if} (pen\_sw\_bnd(n,i) > 0.0) \textcolor{keywordflow}{then}
642         \textcolor{comment}{! SW\_trans is the SW that is transmitted THROUGH the layer}
643         opt\_depth = h(i,k) * opacity\_band(n,i,k)
644         exp\_od = exp(-opt\_depth)
645         sw\_trans = exp\_od
646 
647         \textcolor{comment}{! Heating at a very small rate can be absorbed by a sufficiently thick layer or several}
648         \textcolor{comment}{! thin layers without further penetration.}
649         \textcolor{keywordflow}{if} (optics%answers\_2018) \textcolor{keywordflow}{then}
650           \textcolor{keywordflow}{if} (nsw*pen\_sw\_bnd(n,i)*sw\_trans < min\_sw\_heat*min(1.0, i\_habs*h(i,k)) ) sw\_trans = 0.0
651         \textcolor{keywordflow}{elseif} ((nsw*pen\_sw\_bnd(n,i)*sw\_trans < min\_sw\_heat) .and. (h(i,k) > h\_min\_heat)) \textcolor{keywordflow}{then}
652           \textcolor{keywordflow}{if} (nsw*pen\_sw\_bnd(n,i) <= min\_sw\_heat * (i\_habs*(h(i,k) - h\_min\_heat))) \textcolor{keywordflow}{then}
653             sw\_trans = 0.0
654           \textcolor{keywordflow}{else}
655             sw\_trans = min(sw\_trans, &
656                            1.0 - (min\_sw\_heat*(i\_habs*(h(i,k) - h\_min\_heat))) / (nsw*pen\_sw\_bnd(n,i)))
657 \textcolor{keywordflow}{          endif}
658 \textcolor{keywordflow}{        endif}
659 
660         heat\_bnd = pen\_sw\_bnd(n,i) * (1.0 - sw\_trans)
661         \textcolor{keywordflow}{if} (adjustabsorptionprofile .and. (h\_heat(i) > 0.0)) \textcolor{keywordflow}{then}
662           \textcolor{comment}{!   In this case, a fraction of the heating is applied to the}
663           \textcolor{comment}{! overlying water so that the mean pressure at which the shortwave}
664           \textcolor{comment}{! heating occurs is exactly what it would have been with a careful}
665           \textcolor{comment}{! pressure-weighted averaging of the exponential heating profile,}
666           \textcolor{comment}{! hence there should be no TKE budget requirements due to this}
667           \textcolor{comment}{! layer.  Very clever, but this is also limited so that the}
668           \textcolor{comment}{! water above is not heated at a faster rate than the layer}
669           \textcolor{comment}{! actually being heated, i.e., SWA <= h\_heat / (h\_heat + h(i,k))}
670           \textcolor{comment}{! and takes the energetics of the rest of the heating into account.}
671           \textcolor{comment}{! (-RWH, ~7 years later.)}
672           \textcolor{keywordflow}{if} (opt\_depth > 1e-5) \textcolor{keywordflow}{then}
673             swa = ((opt\_depth + (opt\_depth + 2.0)*exp\_od) - 2.0) / &
674               ((opt\_depth + opacity\_band(n,i,k) * h\_heat(i)) * &
675                (1.0 - exp\_od))
676           \textcolor{keywordflow}{else}
677             \textcolor{comment}{! Use Taylor series expansion of the expression above for a}
678             \textcolor{comment}{! more accurate form with very small layer optical depths.}
679             swa = h(i,k) * (opt\_depth * (1.0 - opt\_depth)) / &
680               ((h\_heat(i) + h(i,k)) * (6.0 - 3.0*opt\_depth))
681 \textcolor{keywordflow}{          endif}
682           coswa\_frac = 0.0
683           \textcolor{keywordflow}{if} (swa*(h\_heat(i) + h(i,k)) > h\_heat(i)) \textcolor{keywordflow}{then}
684             coswa\_frac = (swa*(h\_heat(i) + h(i,k)) - h\_heat(i) ) / &
685                          (swa*(h\_heat(i) + h(i,k)))
686             swa = h\_heat(i) / (h\_heat(i) + h(i,k))
687 \textcolor{keywordflow}{          endif}
688 
689           t\_chg\_above(i,k) = t\_chg\_above(i,k) + (swa * heat\_bnd) / h\_heat(i)
690           t(i,k) = t(i,k) + ((1.0 - swa) * heat\_bnd) / h(i,k)
691         \textcolor{keywordflow}{else}
692           coswa\_frac = 1.0
693           t(i,k) = t(i,k) + pen\_sw\_bnd(n,i) * (1.0 - sw\_trans) / h(i,k)
694 \textcolor{keywordflow}{        endif}
695 
696         \textcolor{keywordflow}{if} (tke\_calc) \textcolor{keywordflow}{then}
697           \textcolor{keywordflow}{if} (opt\_depth > 1e-2) \textcolor{keywordflow}{then}
698             tke(i,k) = tke(i,k) - coswa\_frac*heat\_bnd*dsv\_dt(i,k)* &
699                (0.5*h(i,k)*g\_hconv2) * &
700                (opt\_depth*(1.0+exp\_od) - 2.0*(1.0-exp\_od)) / (opt\_depth*(1.0-exp\_od))
701           \textcolor{keywordflow}{else}
702             \textcolor{comment}{! Use Taylor series-derived approximation to the above expression}
703             \textcolor{comment}{! that is well behaved and more accurate when opt\_depth is small.}
704             tke(i,k) = tke(i,k) - coswa\_frac*heat\_bnd*dsv\_dt(i,k)* &
705                (0.5*h(i,k)*g\_hconv2) * &
706                (c1\_6*opt\_depth * (1.0 - c1\_60*opt\_depth**2))
707 \textcolor{keywordflow}{          endif}
708 \textcolor{keywordflow}{        endif}
709 
710         pen\_sw\_bnd(n,i) = pen\_sw\_bnd(n,i) * sw\_trans
711 \textcolor{keywordflow}{      endif} ;\textcolor{keywordflow}{ enddo}
712 \textcolor{keywordflow}{    endif}
713 
714     \textcolor{comment}{! Add to the accumulated thickness above that could be heated.}
715     \textcolor{comment}{! Only layers greater than h\_min\_heat thick should get heated.}
716     \textcolor{keywordflow}{if} (h(i,k) >= 2.0*h\_min\_heat) \textcolor{keywordflow}{then}
717       h\_heat(i) = h\_heat(i) + h(i,k)
718     \textcolor{keywordflow}{elseif} (h(i,k) > h\_min\_heat) \textcolor{keywordflow}{then}
719       h\_heat(i) = h\_heat(i) + (2.0*h(i,k) - 2.0*h\_min\_heat)
720 \textcolor{keywordflow}{    endif}
721 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo} \textcolor{comment}{! i & k loops}
722 
723 \textcolor{comment}{! if (.not.absorbAllSW .and. .not.adjustAbsorptionProfile) return}
724 
725   \textcolor{comment}{! Unless modified, there is no temperature change due to fluxes from the bottom.}
726   \textcolor{keywordflow}{do} i=is,ie ; t\_chg(i) = 0.0 ;\textcolor{keywordflow}{ enddo}
727 
728   \textcolor{keywordflow}{if} (absorballsw) \textcolor{keywordflow}{then}
729     \textcolor{comment}{! If there is still shortwave radiation at this point, it could go into}
730     \textcolor{comment}{! the bottom (with a bottom mud model), or it could be redistributed back}
731     \textcolor{comment}{! through the water column.}
732     \textcolor{keywordflow}{do} i=is,ie
733       pen\_sw\_rem(i) = pen\_sw\_bnd(1,i)
734       \textcolor{keywordflow}{do} n=2,nsw ; pen\_sw\_rem(i) = pen\_sw\_rem(i) + pen\_sw\_bnd(n,i) ;\textcolor{keywordflow}{ enddo}
735 \textcolor{keywordflow}{    enddo}
736     \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{if} (pen\_sw\_rem(i) > 0.0) sw\_remains = .true. ;\textcolor{keywordflow}{ enddo}
737 
738     ih\_limit = 1.0 / h\_limit\_fluxes
739     \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{if} ((pen\_sw\_rem(i) > 0.0) .and. (h\_heat(i) > 0.0)) \textcolor{keywordflow}{then}
740       \textcolor{keywordflow}{if} (h\_heat(i)*ih\_limit >= 1.0) \textcolor{keywordflow}{then}
741         t\_chg(i) = pen\_sw\_rem(i) / h\_heat(i) ; unabsorbed = 0.0
742       \textcolor{keywordflow}{else}
743         t\_chg(i) = pen\_sw\_rem(i) * ih\_limit
744         unabsorbed = 1.0 - h\_heat(i)*ih\_limit
745 \textcolor{keywordflow}{      endif}
746       \textcolor{keywordflow}{do} n=1,nsw ; pen\_sw\_bnd(n,i) = unabsorbed * pen\_sw\_bnd(n,i) ;\textcolor{keywordflow}{ enddo}
747 \textcolor{keywordflow}{    endif} ;\textcolor{keywordflow}{ enddo}
748 \textcolor{keywordflow}{  endif} \textcolor{comment}{! absorbAllSW}
749 
750   \textcolor{keywordflow}{if} (absorballsw .or. adjustabsorptionprofile) \textcolor{keywordflow}{then}
751     \textcolor{keywordflow}{do} ks=nz,1,-1 ; \textcolor{keywordflow}{do} i=is,ie
752       k = ks
753       \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(ksort)) \textcolor{keywordflow}{then}
754         \textcolor{keywordflow}{if} (ksort(i,ks) <= 0) cycle
755         k = ksort(i,ks)
756 \textcolor{keywordflow}{      endif}
757 
758       \textcolor{keywordflow}{if} (t\_chg(i) > 0.0) \textcolor{keywordflow}{then}
759         \textcolor{comment}{! Only layers greater than h\_min\_heat thick should get heated.}
760         \textcolor{keywordflow}{if} (h(i,k) >= 2.0*h\_min\_heat) \textcolor{keywordflow}{then} ; t(i,k) = t(i,k) + t\_chg(i)
761         \textcolor{keywordflow}{elseif} (h(i,k) > h\_min\_heat) \textcolor{keywordflow}{then}
762           t(i,k) = t(i,k) + t\_chg(i) * (2.0 - 2.0*h\_min\_heat/h(i,k))
763 \textcolor{keywordflow}{        endif}
764 \textcolor{keywordflow}{      endif}
765       \textcolor{comment}{! Increase the heating for layers above.}
766       t\_chg(i) = t\_chg(i) + t\_chg\_above(i,k)
767 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
768     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(htot) .and. \textcolor{keyword}{present}(ttot)) \textcolor{keywordflow}{then}
769       \textcolor{keywordflow}{do} i=is,ie ; ttot(i) = ttot(i) + t\_chg(i) * htot(i) ;\textcolor{keywordflow}{ enddo}
770 \textcolor{keywordflow}{    endif}
771 \textcolor{keywordflow}{  endif} \textcolor{comment}{! absorbAllSW .or. adjustAbsorptionProfile}
772 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_aefdfc303272a4f4fc07f84a8aea2b0f1}\label{namespacemom__opacity_aefdfc303272a4f4fc07f84a8aea2b0f1}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!extract\+\_\+optics\+\_\+fields@{extract\+\_\+optics\+\_\+fields}}
\index{extract\+\_\+optics\+\_\+fields@{extract\+\_\+optics\+\_\+fields}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{extract\+\_\+optics\+\_\+fields()}{extract\_optics\_fields()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::extract\+\_\+optics\+\_\+fields (\begin{DoxyParamCaption}\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(in)}]{optics,  }\item[{integer, intent(out), optional}]{nbands }\end{DoxyParamCaption})}



Set arguments to fields from the optics type. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em optics} & An optics structure that has values of opacities and shortwave fluxes.\\
\hline
\mbox{\tt out}  & {\em nbands} & The number of penetrating bands of SW radiation \\
\hline
\end{DoxyParams}


Definition at line 484 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
484   \textcolor{keywordtype}{type}(optics\_type),       \textcolor{keywordtype}{intent(in)}  :: optics\textcolor{comment}{ !< An optics structure that has values of opacities}
485 \textcolor{comment}{                                                 !! and shortwave fluxes.}
486   \textcolor{keywordtype}{integer}, \textcolor{keywordtype}{optional},       \textcolor{keywordtype}{intent(out)} :: nbands\textcolor{comment}{ !< The number of penetrating bands of SW radiation}
487 
488   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(nbands)) nbands = optics%nbands
489 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a4c1942f798619a9ad854d1152ebcab63}\label{namespacemom__opacity_a4c1942f798619a9ad854d1152ebcab63}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!extract\+\_\+optics\+\_\+slice@{extract\+\_\+optics\+\_\+slice}}
\index{extract\+\_\+optics\+\_\+slice@{extract\+\_\+optics\+\_\+slice}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{extract\+\_\+optics\+\_\+slice()}{extract\_optics\_slice()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::extract\+\_\+optics\+\_\+slice (\begin{DoxyParamCaption}\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(in)}]{optics,  }\item[{integer, intent(in)}]{j,  }\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{real, dimension(max(optics\%nbands,1), g \%isd\+: g \%ied, gv \%ke), intent(out), optional}]{opacity,  }\item[{real, intent(in), optional}]{opacity\+\_\+scale,  }\item[{real, dimension(max(optics\%nbands,1), g \%isd\+: g \%ied), intent(out), optional}]{pen\+S\+W\+\_\+top,  }\item[{real, intent(in), optional}]{pen\+S\+W\+\_\+scale }\end{DoxyParamCaption})}



This subroutine returns a 2-\/d slice at constant j of fields from an \mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}, with the potential for rescaling these fields. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em optics} & An optics structure that has values of opacities and shortwave fluxes.\\
\hline
\mbox{\tt in}  & {\em j} & j-\/index to extract\\
\hline
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt out}  & {\em opacity} & The opacity in each band, i-\/point, and layer\\
\hline
\mbox{\tt in}  & {\em opacity\+\_\+scale} & A factor by which to rescale the opacity.\\
\hline
\mbox{\tt out}  & {\em pensw\+\_\+top} & The shortwave radiation \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em pensw\+\_\+scale} & A factor by which to rescale the shortwave flux. \\
\hline
\end{DoxyParams}


Definition at line 446 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
446   \textcolor{keywordtype}{type}(optics\_type),       \textcolor{keywordtype}{intent(in)}  :: optics\textcolor{comment}{ !< An optics structure that has values of opacities}
447 \textcolor{comment}{                                                 !! and shortwave fluxes.}
448   \textcolor{keywordtype}{integer},                 \textcolor{keywordtype}{intent(in)}  :: j\textcolor{comment}{      !< j-index to extract}
449   \textcolor{keywordtype}{type}(ocean\_grid\_type),   \textcolor{keywordtype}{intent(in)}  :: G\textcolor{comment}{      !< The ocean's grid structure.}
450   \textcolor{keywordtype}{type}(verticalGrid\_type), \textcolor{keywordtype}{intent(in)}  :: GV\textcolor{comment}{     !< The ocean's vertical grid structure.}
451   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(max(optics%nbands,1),SZI\_(G),SZK\_(GV))}, &
452                  \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(out)} :: opacity\textcolor{comment}{   !< The opacity in each band, i-point, and layer}
453   \textcolor{keywordtype}{real},          \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: opacity\_scale\textcolor{comment}{ !< A factor by which to rescale the opacity.}
454   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(max(optics%nbands,1),SZI\_(G))}, &
455                  \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(out)} :: penSW\_top\textcolor{comment}{ !< The shortwave radiation [Q R Z T-1 ~> W m-2]}
456 \textcolor{comment}{                                                    !! at the surface in each of the nbands bands}
457 \textcolor{comment}{                                                    !! that penetrates beyond the surface skin layer.}
458   \textcolor{keywordtype}{real},          \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}  :: penSW\_scale\textcolor{comment}{ !< A factor by which to rescale the shortwave flux.}
459 
460   \textcolor{comment}{! Local variables}
461   \textcolor{keywordtype}{real} :: scale\_opacity, scale\_penSW \textcolor{comment}{! Rescaling factors}
462   \textcolor{keywordtype}{integer} :: i, is, ie, k, nz, n
463   is = g%isc ; ie = g%iec ; nz = g%ke
464 
465   scale\_opacity = 1.0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(opacity\_scale)) scale\_opacity = opacity\_scale
466   scale\_pensw = 1.0 ; \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(pensw\_scale)) scale\_pensw = pensw\_scale
467 
468   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(opacity)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie
469     \textcolor{keywordflow}{do} n=1,optics%nbands
470       opacity(n,i,k) = scale\_opacity * optics%opacity\_band(n,i,j,k)
471 \textcolor{keywordflow}{    enddo}
472 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ endif}
473 
474   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(pensw\_top)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} i=is,ie
475     \textcolor{keywordflow}{do} n=1,optics%nbands
476       pensw\_top(n,i) = scale\_pensw * optics%sw\_pen\_band(n,i,j)
477 \textcolor{keywordflow}{    enddo}
478 \textcolor{keywordflow}{  enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ endif}
479 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_ab5c0caabbf8a806a95bcc416da673841}\label{namespacemom__opacity_ab5c0caabbf8a806a95bcc416da673841}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!opacity\+\_\+end@{opacity\+\_\+end}}
\index{opacity\+\_\+end@{opacity\+\_\+end}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{opacity\+\_\+end()}{opacity\_end()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::opacity\+\_\+end (\begin{DoxyParamCaption}\item[{type(\mbox{\hyperlink{structmom__opacity_1_1opacity__cs}{opacity\+\_\+cs}}), pointer}]{CS,  }\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), optional, pointer}]{optics }\end{DoxyParamCaption})}


\begin{DoxyParams}{Parameters}
{\em cs} & An opacity control structure that should be deallocated.\\
\hline
{\em optics} & An optics type structure that should be deallocated. \\
\hline
\end{DoxyParams}


Definition at line 1119 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
1119   \textcolor{keywordtype}{type}(opacity\_CS),  \textcolor{keywordtype}{pointer}           :: CS\textcolor{comment}{ !< An opacity control structure that should be deallocated.}
1120   \textcolor{keywordtype}{type}(optics\_type), \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{pointer} :: optics\textcolor{comment}{ !< An optics type structure that should be deallocated.}
1121 
1122   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(cs%id\_opacity)) \textcolor{keyword}{deallocate}(cs%id\_opacity)
1123   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(cs)) \textcolor{keyword}{deallocate}(cs)
1124 
1125   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(optics)) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(optics)) \textcolor{keywordflow}{then}
1126     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(optics%opacity\_band)) \textcolor{keyword}{deallocate}(optics%opacity\_band)
1127     \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(optics%sw\_pen\_band)) \textcolor{keyword}{deallocate}(optics%sw\_pen\_band)
1128 \textcolor{keywordflow}{  endif} ;\textcolor{keywordflow}{ endif}
1129 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_aa14b4bfe45861d86740c24e474ce682c}\label{namespacemom__opacity_aa14b4bfe45861d86740c24e474ce682c}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!opacity\+\_\+from\+\_\+chl@{opacity\+\_\+from\+\_\+chl}}
\index{opacity\+\_\+from\+\_\+chl@{opacity\+\_\+from\+\_\+chl}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{opacity\+\_\+from\+\_\+chl()}{opacity\_from\_chl()}}
{\footnotesize\ttfamily subroutine mom\+\_\+opacity\+::opacity\+\_\+from\+\_\+chl (\begin{DoxyParamCaption}\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(inout)}]{optics,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+total,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+vis\+\_\+dir,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+vis\+\_\+dif,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+nir\+\_\+dir,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+nir\+\_\+dif,  }\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(\mbox{\hyperlink{structmom__opacity_1_1opacity__cs}{opacity\+\_\+cs}}), pointer}]{CS,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed), intent(in), optional}]{chl\+\_\+2d,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed, gv \%ke), intent(in), optional}]{chl\+\_\+3d }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



This sets the \char`\"{}blue\char`\"{} band opacity based on chloophyll A concencentrations The red portion is lumped into the net heating at the surface. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in,out}  & {\em optics} & An optics structure that has values set based on the opacities.\\
\hline
 & {\em sw\+\_\+total} & Total shortwave flux into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+vis\+\_\+dir} & Visible, direct shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+vis\+\_\+dif} & Visible, diffuse shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+nir\+\_\+dir} & Near-\/\+IR, direct shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+nir\+\_\+dif} & Near-\/\+IR, diffuse shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
 & {\em cs} & The control structure.\\
\hline
\mbox{\tt in}  & {\em chl\+\_\+2d} & Vertically uniform chlorophyll-\/A concentractions \mbox{[}mg m-\/3\mbox{]}\\
\hline
\mbox{\tt in}  & {\em chl\+\_\+3d} & A 3-\/d field of chlorophyll-\/A concentractions \mbox{[}mg m-\/3\mbox{]} \\
\hline
\end{DoxyParams}


Definition at line 222 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
222   \textcolor{keywordtype}{type}(optics\_type),       \textcolor{keywordtype}{intent(inout)} :: optics\textcolor{comment}{ !< An optics structure that has values}
223 \textcolor{comment}{                                                 !! set based on the opacities.}
224   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_total\textcolor{comment}{ !< Total shortwave flux into the ocean [Q R Z T-1 ~> W
       m-2]}
225   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_vis\_dir\textcolor{comment}{ !< Visible, direct shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
226   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_vis\_dif\textcolor{comment}{ !< Visible, diffuse shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
227   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_nir\_dir\textcolor{comment}{ !< Near-IR, direct shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
228   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_nir\_dif\textcolor{comment}{ !< Near-IR, diffuse shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
229   \textcolor{keywordtype}{type}(ocean\_grid\_type),   \textcolor{keywordtype}{intent(in)}    :: G\textcolor{comment}{      !< The ocean's grid structure.}
230   \textcolor{keywordtype}{type}(verticalGrid\_type), \textcolor{keywordtype}{intent(in)}    :: GV\textcolor{comment}{     !< The ocean's vertical grid structure.}
231   \textcolor{keywordtype}{type}(unit\_scale\_type),   \textcolor{keywordtype}{intent(in)}    :: US\textcolor{comment}{     !< A dimensional unit scaling type}
232   \textcolor{keywordtype}{type}(opacity\_CS),        \textcolor{keywordtype}{pointer}       :: CS\textcolor{comment}{     !< The control structure.}
233   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))}, &
234                  \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}    :: chl\_2d\textcolor{comment}{ !< Vertically uniform chlorophyll-A concentractions [mg
       m-3]}
235   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(GV))}, &
236                  \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}    :: chl\_3d\textcolor{comment}{ !< A 3-d field of chlorophyll-A concentractions [mg m-3]}
237 
238   \textcolor{keywordtype}{real} :: chl\_data(SZI\_(G),SZJ\_(G)) \textcolor{comment}{! The chlorophyll A concentrations in a layer [mg m-3].}
239   \textcolor{keywordtype}{real} :: Inv\_nbands        \textcolor{comment}{! The inverse of the number of bands of penetrating}
240                             \textcolor{comment}{! shortwave radiation.}
241   \textcolor{keywordtype}{real} :: Inv\_nbands\_nir    \textcolor{comment}{! The inverse of the number of bands of penetrating}
242                             \textcolor{comment}{! near-infrafed radiation.}
243   \textcolor{keywordtype}{real} :: SW\_pen\_tot        \textcolor{comment}{! The sum across the bands of the penetrating}
244                             \textcolor{comment}{! shortwave radiation [Q R Z T-1 ~> W m-2].}
245   \textcolor{keywordtype}{real} :: SW\_vis\_tot        \textcolor{comment}{! The sum across the visible bands of shortwave}
246                             \textcolor{comment}{! radiation [Q R Z T-1 ~> W m-2].}
247   \textcolor{keywordtype}{real} :: SW\_nir\_tot        \textcolor{comment}{! The sum across the near infrared bands of shortwave}
248                             \textcolor{comment}{! radiation [Q R Z T-1 ~> W m-2].}
249   \textcolor{keywordtype}{type}(time\_type) :: day
250   \textcolor{keywordtype}{character(len=128)} :: mesg
251   \textcolor{keywordtype}{integer} :: i, j, k, n, is, ie, js, je, nz, nbands
252   \textcolor{keywordtype}{logical} :: multiband\_vis\_input, multiband\_nir\_input
253 
254   is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
255 
256 \textcolor{comment}{!   In this model, the Morel (modified) and Manizza (modified) schemes}
257 \textcolor{comment}{! use the "blue" band in the parameterizations to determine the e-folding}
258 \textcolor{comment}{! depth of the incoming shortwave attenuation. The red portion is lumped}
259 \textcolor{comment}{! into the net heating at the surface.}
260 \textcolor{comment}{!}
261 \textcolor{comment}{! Morel, A., Optical modeling of the upper ocean in relation to its biogenous}
262 \textcolor{comment}{!   matter content (case-i waters).,J. Geo. Res., \{93\}, 10,749--10,768, 1988.}
263 \textcolor{comment}{!}
264 \textcolor{comment}{! Manizza, M., C.~L. Quere, A.~Watson, and E.~T. Buitenhuis, Bio-optical}
265 \textcolor{comment}{!   feedbacks among phytoplankton, upper ocean physics and sea-ice in a}
266 \textcolor{comment}{!   global model, Geophys. Res. Let., , L05,603, 2005.}
267 
268   nbands = optics%nbands
269 
270   \textcolor{keywordflow}{if} (nbands <= 1) \textcolor{keywordflow}{then} ; inv\_nbands = 1.0
271   \textcolor{keywordflow}{else} ; inv\_nbands = 1.0 / \textcolor{keywordtype}{real(nbands)} ; endif
272 
273   \textcolor{keywordflow}{if} (nbands <= 2) \textcolor{keywordflow}{then} ; inv\_nbands\_nir = 0.0
274   \textcolor{keywordflow}{else} ; inv\_nbands\_nir = 1.0 / \textcolor{keywordtype}{real(nbands - 2.0)} ; endif
275 
276   multiband\_vis\_input = (\textcolor{keyword}{associated}(sw\_vis\_dir) .and. &
277                          \textcolor{keyword}{associated}(sw\_vis\_dif))
278   multiband\_nir\_input = (\textcolor{keyword}{associated}(sw\_nir\_dir) .and. &
279                          \textcolor{keyword}{associated}(sw\_nir\_dif))
280 
281   chl\_data(:,:) = 0.0
282   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(chl\_3d)) \textcolor{keywordflow}{then}
283     \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; chl\_data(i,j) = chl\_3d(i,j,1) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
284     \textcolor{keywordflow}{do} k=1,nz; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
285       \textcolor{keywordflow}{if} ((g%mask2dT(i,j) > 0.5) .and. (chl\_3d(i,j,k) < 0.0)) \textcolor{keywordflow}{then}
286         \textcolor{keyword}{write}(mesg,\textcolor{stringliteral}{'(" Negative chl\_3d of ",(1pe12.4)," found at i,j,k = ", &}
287 \textcolor{stringliteral}{}\textcolor{stringliteral}{                  & 3(1x,i3), " lon/lat = ",(1pe12.4)," E ", (1pe12.4), " N.")'}) &
288                    chl\_3d(i,j,k), i, j, k, g%geoLonT(i,j), g%geoLatT(i,j)
289         \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"MOM\_opacity opacity\_from\_chl: "}//trim(mesg))
290 \textcolor{keywordflow}{      endif}
291 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
292   \textcolor{keywordflow}{elseif} (\textcolor{keyword}{present}(chl\_2d)) \textcolor{keywordflow}{then}
293     \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; chl\_data(i,j) = chl\_2d(i,j) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
294     \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
295       \textcolor{keywordflow}{if} ((g%mask2dT(i,j) > 0.5) .and. (chl\_2d(i,j) < 0.0)) \textcolor{keywordflow}{then}
296         \textcolor{keyword}{write}(mesg,\textcolor{stringliteral}{'(" Negative chl\_2d of ",(1pe12.4)," at i,j = ", &}
297 \textcolor{stringliteral}{}\textcolor{stringliteral}{                  & 2(i3), "lon/lat = ",(1pe12.4)," E ", (1pe12.4), " N.")'}) &
298                    chl\_data(i,j), i, j, g%geoLonT(i,j), g%geoLatT(i,j)
299         \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"MOM\_opacity opacity\_from\_chl: "}//trim(mesg))
300 \textcolor{keywordflow}{      endif}
301 \textcolor{keywordflow}{    enddo} ;\textcolor{keywordflow}{ enddo}
302   \textcolor{keywordflow}{else}
303     \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"Either chl\_2d or chl\_3d must be preesnt in a call to opacity\_form\_chl."})
304 \textcolor{keywordflow}{  endif}
305 
306   \textcolor{keywordflow}{select case} (cs%opacity\_scheme)
307     \textcolor{keywordflow}{case} (manizza\_05)
308       \textcolor{comment}{!$OMP parallel do default(shared) private(SW\_vis\_tot,SW\_nir\_tot)}
309       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
310         sw\_vis\_tot = 0.0 ; sw\_nir\_tot = 0.0
311         \textcolor{keywordflow}{if} (g%mask2dT(i,j) > 0.5) \textcolor{keywordflow}{then}
312           \textcolor{keywordflow}{if} (multiband\_vis\_input) \textcolor{keywordflow}{then}
313             sw\_vis\_tot = sw\_vis\_dir(i,j) + sw\_vis\_dif(i,j)
314           \textcolor{keywordflow}{else}  \textcolor{comment}{! Follow Manizza 05 in assuming that 42% of SW is visible.}
315             sw\_vis\_tot = 0.42 * sw\_total(i,j)
316 \textcolor{keywordflow}{          endif}
317           \textcolor{keywordflow}{if} (multiband\_nir\_input) \textcolor{keywordflow}{then}
318             sw\_nir\_tot = sw\_nir\_dir(i,j) + sw\_nir\_dif(i,j)
319           \textcolor{keywordflow}{else}
320             sw\_nir\_tot = sw\_total(i,j) - sw\_vis\_tot
321 \textcolor{keywordflow}{          endif}
322 \textcolor{keywordflow}{        endif}
323 
324         \textcolor{comment}{! Band 1 is Manizza blue.}
325         optics%sw\_pen\_band(1,i,j) = cs%blue\_frac*sw\_vis\_tot
326         \textcolor{comment}{! Band 2 (if used) is Manizza red.}
327         \textcolor{keywordflow}{if} (nbands > 1) &
328           optics%sw\_pen\_band(2,i,j) = (1.0-cs%blue\_frac)*sw\_vis\_tot
329         \textcolor{comment}{! All remaining bands are NIR, for lack of something better to do.}
330         \textcolor{keywordflow}{do} n=3,nbands
331           optics%sw\_pen\_band(n,i,j) = inv\_nbands\_nir * sw\_nir\_tot
332 \textcolor{keywordflow}{        enddo}
333 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
334     \textcolor{keywordflow}{case} (morel\_88)
335       \textcolor{comment}{!$OMP parallel do default(shared) private(SW\_pen\_tot)}
336       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
337         sw\_pen\_tot = 0.0
338         \textcolor{keywordflow}{if} (g%mask2dT(i,j) > 0.5) \textcolor{keywordflow}{then} ; \textcolor{keywordflow}{if} (multiband\_vis\_input) \textcolor{keywordflow}{then}
339             sw\_pen\_tot = sw\_pen\_frac\_morel(chl\_data(i,j)) * (sw\_vis\_dir(i,j) + sw\_vis\_dif(i,j))
340           \textcolor{keywordflow}{else}
341             sw\_pen\_tot = sw\_pen\_frac\_morel(chl\_data(i,j)) * 0.5*sw\_total(i,j)
342 \textcolor{keywordflow}{        endif} ;\textcolor{keywordflow}{ endif}
343 
344         \textcolor{keywordflow}{do} n=1,nbands
345           optics%sw\_pen\_band(n,i,j) = inv\_nbands*sw\_pen\_tot
346 \textcolor{keywordflow}{        enddo}
347 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
348 \textcolor{keywordflow}{    case default}
349       \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"opacity\_from\_chl: CS%opacity\_scheme is not valid."})
350 \textcolor{keywordflow}{  end select}
351 
352   \textcolor{comment}{!$OMP parallel do default(shared) firstprivate(chl\_data)}
353   \textcolor{keywordflow}{do} k=1,nz
354     \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(chl\_3d)) \textcolor{keywordflow}{then}
355       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; chl\_data(i,j) = chl\_3d(i,j,k) ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
356 \textcolor{keywordflow}{    endif}
357 
358     \textcolor{keywordflow}{select case} (cs%opacity\_scheme)
359       \textcolor{keywordflow}{case} (manizza\_05)
360         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
361           \textcolor{keywordflow}{if} (g%mask2dT(i,j) <= 0.5) \textcolor{keywordflow}{then}
362             \textcolor{keywordflow}{do} n=1,optics%nbands
363               optics%opacity\_band(n,i,j,k) = cs%opacity\_land\_value
364 \textcolor{keywordflow}{            enddo}
365           \textcolor{keywordflow}{else}
366             \textcolor{comment}{! Band 1 is Manizza blue.}
367             optics%opacity\_band(1,i,j,k) = 0.0232 + 0.074*chl\_data(i,j)**0.674
368             \textcolor{keywordflow}{if} (nbands >= 2) &  \textcolor{comment}{!  Band 2 is Manizza red.}
369               optics%opacity\_band(2,i,j,k) = 0.225 + 0.037*chl\_data(i,j)**0.629
370             \textcolor{comment}{! All remaining bands are NIR, for lack of something better to do.}
371             \textcolor{keywordflow}{do} n=3,nbands ; optics%opacity\_band(n,i,j,k) = 2.86 ;\textcolor{keywordflow}{ enddo}
372 \textcolor{keywordflow}{          endif}
373 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
374       \textcolor{keywordflow}{case} (morel\_88)
375         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
376           optics%opacity\_band(1,i,j,k) = cs%opacity\_land\_value
377           \textcolor{keywordflow}{if} (g%mask2dT(i,j) > 0.5) &
378             optics%opacity\_band(1,i,j,k) = opacity\_morel(chl\_data(i,j))
379 
380           \textcolor{keywordflow}{do} n=2,optics%nbands
381             optics%opacity\_band(n,i,j,k) = optics%opacity\_band(1,i,j,k)
382 \textcolor{keywordflow}{          enddo}
383 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
384 
385 \textcolor{keywordflow}{      case default}
386         \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"opacity\_from\_chl: CS%opacity\_scheme is not valid."})
387 \textcolor{keywordflow}{    end select}
388 \textcolor{keywordflow}{  enddo}
389 
390 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a39fce7bd33a469e3e9fe7cfeb51825b5}\label{namespacemom__opacity_a39fce7bd33a469e3e9fe7cfeb51825b5}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!opacity\+\_\+init@{opacity\+\_\+init}}
\index{opacity\+\_\+init@{opacity\+\_\+init}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{opacity\+\_\+init()}{opacity\_init()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::opacity\+\_\+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(\mbox{\hyperlink{structmom__opacity_1_1opacity__cs}{opacity\+\_\+cs}}), pointer}]{CS,  }\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), pointer}]{optics }\end{DoxyParamCaption})}



This routine initalizes the opacity module, including an \mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em time} & The current model time.\\
\hline
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & model vertical grid structure\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
\mbox{\tt in}  & {\em param\+\_\+file} & A structure to parse for run-\/time parameters.\\
\hline
\mbox{\tt in,out}  & {\em diag} & A structure that is used to regulate diagnostic output.\\
\hline
 & {\em cs} & A pointer that is set to point to the control structure for this module.\\
\hline
 & {\em optics} & An optics structure that has parameters set and arrays allocated here. \\
\hline
\end{DoxyParams}


Definition at line 920 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
920   \textcolor{keywordtype}{type}(time\_type), \textcolor{keywordtype}{target}, \textcolor{keywordtype}{intent(in)}    :: Time\textcolor{comment}{ !< The current model time.}
921   \textcolor{keywordtype}{type}(ocean\_grid\_type),   \textcolor{keywordtype}{intent(in)}    :: G\textcolor{comment}{    !< The ocean's grid structure.}
922   \textcolor{keywordtype}{type}(verticalGrid\_type), \textcolor{keywordtype}{intent(in)}    :: GV\textcolor{comment}{   !< model vertical grid structure}
923   \textcolor{keywordtype}{type}(unit\_scale\_type),   \textcolor{keywordtype}{intent(in)}    :: US\textcolor{comment}{   !< A dimensional unit scaling type}
924   \textcolor{keywordtype}{type}(param\_file\_type),   \textcolor{keywordtype}{intent(in)}    :: param\_file\textcolor{comment}{ !< A structure to parse for run-time}
925 \textcolor{comment}{                                                 !! parameters.}
926   \textcolor{keywordtype}{type}(diag\_ctrl), \textcolor{keywordtype}{target}, \textcolor{keywordtype}{intent(inout)} :: diag\textcolor{comment}{ !< A structure that is used to regulate diagnostic}
927 \textcolor{comment}{                                                 !! output.}
928   \textcolor{keywordtype}{type}(opacity\_CS),        \textcolor{keywordtype}{pointer}       :: CS\textcolor{comment}{   !< A pointer that is set to point to the control}
929 \textcolor{comment}{                                                 !! structure for this module.}
930   \textcolor{keywordtype}{type}(optics\_type),       \textcolor{keywordtype}{pointer}       :: optics\textcolor{comment}{ !< An optics structure that has parameters}
931 \textcolor{comment}{                                                 !! set and arrays allocated here.}
932   \textcolor{comment}{! Local variables}
933   \textcolor{keywordtype}{character(len=200)} :: tmpstr
934   \textcolor{keywordtype}{character(len=40)}  :: mdl = \textcolor{stringliteral}{"MOM\_opacity"}
935   \textcolor{keywordtype}{character(len=40)}  :: bandnum, shortname
936   \textcolor{keywordtype}{character(len=200)} :: longname
937   \textcolor{keywordtype}{character(len=40)}  :: scheme\_string
938   \textcolor{comment}{! This include declares and sets the variable "version".}
939 \textcolor{preprocessor}{# include "version\_variable.h"}
940 \textcolor{preprocessor}{}  \textcolor{keywordtype}{real} :: PenSW\_absorb\_minthick \textcolor{comment}{! A thickness that is used to absorb the remaining shortwave heat}
941                                 \textcolor{comment}{! flux when that flux drops below PEN\_SW\_FLUX\_ABSORB [m].}
942   \textcolor{keywordtype}{real} :: PenSW\_minthick\_dflt \textcolor{comment}{! The default for PenSW\_absorb\_minthick [m]}
943   \textcolor{keywordtype}{logical} :: default\_2018\_answers
944   \textcolor{keywordtype}{logical} :: use\_scheme
945   \textcolor{keywordtype}{integer} :: isd, ied, jsd, jed, nz, n
946   isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed ; nz = g%ke
947 
948   \textcolor{keywordflow}{if} (\textcolor{keyword}{associated}(cs)) \textcolor{keywordflow}{then}
949     \textcolor{keyword}{call }mom\_error(warning, \textcolor{stringliteral}{"opacity\_init called with an associated"}// &
950                              \textcolor{stringliteral}{"associated control structure."})
951     \textcolor{keywordflow}{return}
952   \textcolor{keywordflow}{else} ; \textcolor{keyword}{allocate}(cs) ;\textcolor{keywordflow}{ endif}
953 
954   cs%diag => diag
955 
956   \textcolor{comment}{! Read all relevant parameters and write them to the model log.}
957   \textcolor{keyword}{call }log\_version(param\_file, mdl, version, \textcolor{stringliteral}{''})
958 
959 \textcolor{comment}{! parameters for CHL\_A routines}
960   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"VAR\_PEN\_SW"}, cs%var\_pen\_sw, &
961                  \textcolor{stringliteral}{"If true, use one of the CHL\_A schemes specified by "}//&
962                  \textcolor{stringliteral}{"OPACITY\_SCHEME to determine the e-folding depth of "}//&
963                  \textcolor{stringliteral}{"incoming short wave radiation."}, default=.false.)
964 
965   cs%opacity\_scheme = no\_scheme ; scheme\_string = \textcolor{stringliteral}{''}
966   \textcolor{keywordflow}{if} (cs%var\_pen\_sw) \textcolor{keywordflow}{then}
967     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"OPACITY\_SCHEME"}, tmpstr, &
968                  \textcolor{stringliteral}{"This character string specifies how chlorophyll "}//&
969                  \textcolor{stringliteral}{"concentrations are translated into opacities. Currently "}//&
970                  \textcolor{stringliteral}{"valid options include:\(\backslash\)n"}//&
971                  \textcolor{stringliteral}{" \(\backslash\)t\(\backslash\)t  MANIZZA\_05 - Use Manizza et al., GRL, 2005. \(\backslash\)n"}//&
972                  \textcolor{stringliteral}{" \(\backslash\)t\(\backslash\)t  MOREL\_88 - Use Morel, JGR, 1988."}, &
973                  default=manizza\_05\_string)
974     \textcolor{keywordflow}{if} (len\_trim(tmpstr) > 0) \textcolor{keywordflow}{then}
975       tmpstr = uppercase(tmpstr)
976       \textcolor{keywordflow}{select case} (tmpstr)
977         \textcolor{keywordflow}{case} (manizza\_05\_string)
978           cs%opacity\_scheme = manizza\_05 ; scheme\_string = manizza\_05\_string
979         \textcolor{keywordflow}{case} (morel\_88\_string)
980           cs%opacity\_scheme = morel\_88 ; scheme\_string = morel\_88\_string
981 \textcolor{keywordflow}{        case default}
982           \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"opacity\_init: #DEFINE OPACITY\_SCHEME "}//&
983                                   trim(tmpstr) // \textcolor{stringliteral}{"in input file is invalid."})
984 \textcolor{keywordflow}{      end select}
985       \textcolor{keyword}{call }mom\_mesg(\textcolor{stringliteral}{'opacity\_init: opacity scheme set to "'}//trim(tmpstr)//\textcolor{stringliteral}{'".'}, 5)
986 \textcolor{keywordflow}{    endif}
987     \textcolor{keywordflow}{if} (cs%opacity\_scheme == no\_scheme) \textcolor{keywordflow}{then}
988       \textcolor{keyword}{call }mom\_error(warning, \textcolor{stringliteral}{"opacity\_init: No scheme has successfully "}//&
989                \textcolor{stringliteral}{"been specified for the opacity.  Using the default MANIZZA\_05."})
990       cs%opacity\_scheme = manizza\_05 ; scheme\_string = manizza\_05\_string
991 \textcolor{keywordflow}{    endif}
992 
993     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"BLUE\_FRAC\_SW"}, cs%blue\_frac, &
994                  \textcolor{stringliteral}{"The fraction of the penetrating shortwave radiation "}//&
995                  \textcolor{stringliteral}{"that is in the blue band."}, default=0.5, units=\textcolor{stringliteral}{"nondim"})
996   \textcolor{keywordflow}{else}
997     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"EXP\_OPACITY\_SCHEME"}, tmpstr, &
998                  \textcolor{stringliteral}{"This character string specifies which exponential "}//&
999                  \textcolor{stringliteral}{"opacity scheme to utilize. Currently "}//&
1000                  \textcolor{stringliteral}{"valid options include:\(\backslash\)n"}//&
1001                  \textcolor{stringliteral}{" \(\backslash\)t\(\backslash\)t  SINGLE\_EXP - Single Exponent decay. \(\backslash\)n"}//&
1002                  \textcolor{stringliteral}{" \(\backslash\)t\(\backslash\)t  DOUBLE\_EXP - Double Exponent decay."}, &
1003                  default=single\_exp\_string)\textcolor{comment}{!New default for "else" above (non-Chl scheme)}
1004     \textcolor{keywordflow}{if} (len\_trim(tmpstr) > 0) \textcolor{keywordflow}{then}
1005       tmpstr = uppercase(tmpstr)
1006       \textcolor{keywordflow}{select case} (tmpstr)
1007         \textcolor{keywordflow}{case} (single\_exp\_string)
1008           cs%opacity\_scheme = single\_exp ; scheme\_string = single\_exp\_string
1009         \textcolor{keywordflow}{case} (double\_exp\_string)
1010           cs%opacity\_scheme = double\_exp ; scheme\_string = double\_exp\_string
1011 \textcolor{keywordflow}{      end select}
1012       \textcolor{keyword}{call }mom\_mesg(\textcolor{stringliteral}{'opacity\_init: opacity scheme set to "'}//trim(tmpstr)//\textcolor{stringliteral}{'".'}, 5)
1013 \textcolor{keywordflow}{    endif}
1014     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PEN\_SW\_SCALE"}, cs%pen\_sw\_scale, &
1015                  \textcolor{stringliteral}{"The vertical absorption e-folding depth of the "}//&
1016                  \textcolor{stringliteral}{"penetrating shortwave radiation."}, units=\textcolor{stringliteral}{"m"}, default=0.0)
1017     \textcolor{comment}{!BGR/ Added for opacity\_scheme==double\_exp read in 2nd exp-decay and fraction}
1018     \textcolor{keywordflow}{if} (cs%Opacity\_scheme == double\_exp ) \textcolor{keywordflow}{then}
1019       \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PEN\_SW\_SCALE\_2ND"}, cs%pen\_sw\_scale\_2nd, &
1020                  \textcolor{stringliteral}{"The (2nd) vertical absorption e-folding depth of the "}//&
1021                  \textcolor{stringliteral}{"penetrating shortwave radiation "}//&
1022                  \textcolor{stringliteral}{"(use if SW\_EXP\_MODE==double.)"},&
1023                  units=\textcolor{stringliteral}{"m"}, default=0.0)
1024       \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"SW\_1ST\_EXP\_RATIO"}, cs%sw\_1st\_exp\_ratio, &
1025                  \textcolor{stringliteral}{"The fraction of 1st vertical absorption e-folding depth "}//&
1026                  \textcolor{stringliteral}{"penetrating shortwave radiation if SW\_EXP\_MODE==double."},&
1027                   units=\textcolor{stringliteral}{"m"}, default=0.0)
1028     \textcolor{keywordflow}{elseif} (cs%OPACITY\_SCHEME == single\_exp) \textcolor{keywordflow}{then}
1029       \textcolor{comment}{!/Else disable 2nd\_exp scheme}
1030       cs%pen\_sw\_scale\_2nd = 0.0
1031       cs%sw\_1st\_exp\_ratio = 1.0
1032 \textcolor{keywordflow}{    endif}
1033     \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PEN\_SW\_FRAC"}, cs%pen\_sw\_frac, &
1034                  \textcolor{stringliteral}{"The fraction of the shortwave radiation that penetrates "}//&
1035                  \textcolor{stringliteral}{"below the surface."}, units=\textcolor{stringliteral}{"nondim"}, default=0.0)
1036 
1037 \textcolor{keywordflow}{  endif}
1038   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PEN\_SW\_NBANDS"}, optics%nbands, &
1039                  \textcolor{stringliteral}{"The number of bands of penetrating shortwave radiation."}, &
1040                  default=1)
1041   \textcolor{keywordflow}{if} (cs%Opacity\_scheme == double\_exp ) \textcolor{keywordflow}{then}
1042     \textcolor{keywordflow}{if} (optics%nbands /= 2) \textcolor{keyword}{call }mom\_error(fatal, &
1043         \textcolor{stringliteral}{"set\_opacity: \(\backslash\)Cannot use a double\_exp opacity scheme with nbands!=2."})
1044   \textcolor{keywordflow}{elseif} (cs%Opacity\_scheme == single\_exp ) \textcolor{keywordflow}{then}
1045     \textcolor{keywordflow}{if} (optics%nbands /= 1) \textcolor{keyword}{call }mom\_error(fatal, &
1046         \textcolor{stringliteral}{"set\_opacity: \(\backslash\)Cannot use a single\_exp opacity scheme with nbands!=1."})
1047 \textcolor{keywordflow}{  endif}
1048 
1049   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"DEFAULT\_2018\_ANSWERS"}, default\_2018\_answers, &
1050                  \textcolor{stringliteral}{"This sets the default value for the various \_2018\_ANSWERS parameters."}, &
1051                  default=.false.)
1052   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"OPTICS\_2018\_ANSWERS"}, optics%answers\_2018, &
1053                  \textcolor{stringliteral}{"If true, use the order of arithmetic and expressions that recover the "}//&
1054                  \textcolor{stringliteral}{"answers from the end of 2018.  Otherwise, use updated expressions for "}//&
1055                  \textcolor{stringliteral}{"handling the absorption of small remaining shortwave fluxes."}, &
1056                  default=default\_2018\_answers)
1057 
1058   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PEN\_SW\_FLUX\_ABSORB"}, optics%PenSW\_flux\_absorb, &
1059                  \textcolor{stringliteral}{"A minimum remaining shortwave heating rate that will be simply absorbed in "}//&
1060                  \textcolor{stringliteral}{"the next sufficiently thick layers for computational efficiency, instead of "}//&
1061                  \textcolor{stringliteral}{"continuing to penetrate.  The default, 2.5e-11 degC m s-1, is about 1e-4 W m-2 "}//&
1062                  \textcolor{stringliteral}{"or 0.08 degC m century-1, but 0 is also a valid value."}, &
1063                  default=2.5e-11, units=\textcolor{stringliteral}{"degC m s-1"}, scale=gv%m\_to\_H*us%T\_to\_s)
1064 
1065   \textcolor{keywordflow}{if} (optics%answers\_2018) \textcolor{keywordflow}{then} ; pensw\_minthick\_dflt = 0.001 ; \textcolor{keywordflow}{else} ; pensw\_minthick\_dflt = 1.0 ;\textcolor{keywordflow}{ endif}
1066   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"PEN\_SW\_ABSORB\_MINTHICK"}, pensw\_absorb\_minthick, &
1067                  \textcolor{stringliteral}{"A thickness that is used to absorb the remaining penetrating shortwave heat "}//&
1068                  \textcolor{stringliteral}{"flux when it drops below PEN\_SW\_FLUX\_ABSORB."}, &
1069                  default=pensw\_minthick\_dflt, units=\textcolor{stringliteral}{"m"}, scale=gv%m\_to\_H)
1070   optics%PenSW\_absorb\_Invlen = 1.0 / (pensw\_absorb\_minthick + gv%H\_subroundoff)
1071 
1072   \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(optics%min\_wavelength\_band)) &
1073     \textcolor{keyword}{allocate}(optics%min\_wavelength\_band(optics%nbands))
1074   \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(optics%max\_wavelength\_band)) &
1075     \textcolor{keyword}{allocate}(optics%max\_wavelength\_band(optics%nbands))
1076 
1077   \textcolor{keywordflow}{if} (cs%opacity\_scheme == manizza\_05) \textcolor{keywordflow}{then}
1078     optics%min\_wavelength\_band(1) =0
1079     optics%max\_wavelength\_band(1) =550
1080     \textcolor{keywordflow}{if} (optics%nbands >= 2) \textcolor{keywordflow}{then}
1081       optics%min\_wavelength\_band(2)=550
1082       optics%max\_wavelength\_band(2)=700
1083 \textcolor{keywordflow}{    endif}
1084     \textcolor{keywordflow}{if} (optics%nbands > 2) \textcolor{keywordflow}{then}
1085       \textcolor{keywordflow}{do} n=3,optics%nbands
1086         optics%min\_wavelength\_band(n) =700
1087         optics%max\_wavelength\_band(n) =2800
1088 \textcolor{keywordflow}{      enddo}
1089 \textcolor{keywordflow}{    endif}
1090 \textcolor{keywordflow}{  endif}
1091 
1092   \textcolor{keyword}{call }get\_param(param\_file, mdl, \textcolor{stringliteral}{"OPACITY\_LAND\_VALUE"}, cs%opacity\_land\_value, &
1093                  \textcolor{stringliteral}{"The value to use for opacity over land. The default is "}//&
1094                  \textcolor{stringliteral}{"10 m-1 - a value for muddy water."}, units=\textcolor{stringliteral}{"m-1"}, default=10.0)
1095 
1096   \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(optics%opacity\_band)) &
1097     \textcolor{keyword}{allocate}(optics%opacity\_band(optics%nbands,isd:ied,jsd:jed,nz))
1098   \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(optics%sw\_pen\_band)) &
1099     \textcolor{keyword}{allocate}(optics%sw\_pen\_band(optics%nbands,isd:ied,jsd:jed))
1100   \textcolor{keyword}{allocate}(cs%id\_opacity(optics%nbands)) ; cs%id\_opacity(:) = -1
1101 
1102   cs%id\_sw\_pen = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'SW\_pen'}, diag%axesT1, time, &
1103       \textcolor{stringliteral}{'Penetrating shortwave radiation flux into ocean'}, \textcolor{stringliteral}{'W m-2'}, conversion=us%QRZ\_T\_to\_W\_m2)
1104   cs%id\_sw\_vis\_pen = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, \textcolor{stringliteral}{'SW\_vis\_pen'}, diag%axesT1, time, &
1105       \textcolor{stringliteral}{'Visible penetrating shortwave radiation flux into ocean'}, \textcolor{stringliteral}{'W m-2'}, conversion=us%QRZ\_T\_to\_W\_m2)
1106   \textcolor{keywordflow}{do} n=1,optics%nbands
1107     \textcolor{keyword}{write}(bandnum,\textcolor{stringliteral}{'(i3)'}) n
1108     shortname = \textcolor{stringliteral}{'opac\_'}//trim(adjustl(bandnum))
1109     longname = \textcolor{stringliteral}{'Opacity for shortwave radiation in band '}//trim(adjustl(bandnum)) &
1110       // \textcolor{stringliteral}{', saved as L^-1 tanh(Opacity * L) for L = 10^-10 m'}
1111     cs%id\_opacity(n) = register\_diag\_field(\textcolor{stringliteral}{'ocean\_model'}, shortname, diag%axesTL, time, &
1112       longname, \textcolor{stringliteral}{'m-1'})
1113 \textcolor{keywordflow}{  enddo}
1114 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a27873e3942c2ff7980cf516944f78f03}\label{namespacemom__opacity_a27873e3942c2ff7980cf516944f78f03}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!opacity\+\_\+manizza@{opacity\+\_\+manizza}}
\index{opacity\+\_\+manizza@{opacity\+\_\+manizza}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{opacity\+\_\+manizza()}{opacity\_manizza()}}
{\footnotesize\ttfamily real function, public mom\+\_\+opacity\+::opacity\+\_\+manizza (\begin{DoxyParamCaption}\item[{real, intent(in)}]{chl\+\_\+data }\end{DoxyParamCaption})}



This sets the blue-\/wavelength opacity according to the scheme proposed by Manizza, M. et al, 2005. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em chl\+\_\+data} & The chlorophyll-\/A concentration in mg m-\/3.\\
\hline
\end{DoxyParams}
\begin{DoxyReturn}{Returns}
The returned opacity \mbox{[}m-\/1\mbox{]} 
\end{DoxyReturn}


Definition at line 436 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
436   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{intent(in)}  :: chl\_data\textcolor{comment}{ !< The chlorophyll-A concentration in mg m-3.}
437   \textcolor{keywordtype}{real} :: opacity\_manizza\textcolor{comment}{ !< The returned opacity [m-1]}
438 \textcolor{comment}{!   This sets the blue-wavelength opacity according to the scheme proposed by Manizza, M. et al, 2005.}
439 
440   opacity\_manizza = 0.0232 + 0.074*chl\_data**0.674
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a4498b4bb6fcf1b7d849f89aa87c0332e}\label{namespacemom__opacity_a4498b4bb6fcf1b7d849f89aa87c0332e}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!opacity\+\_\+morel@{opacity\+\_\+morel}}
\index{opacity\+\_\+morel@{opacity\+\_\+morel}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{opacity\+\_\+morel()}{opacity\_morel()}}
{\footnotesize\ttfamily real function, public mom\+\_\+opacity\+::opacity\+\_\+morel (\begin{DoxyParamCaption}\item[{real, intent(in)}]{chl\+\_\+data }\end{DoxyParamCaption})}



This sets the blue-\/wavelength opacity according to the scheme proposed by Morel and Antoine (1994). 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em chl\+\_\+data} & The chlorophyll-\/A concentration in mg m-\/3.\\
\hline
\end{DoxyParams}
\begin{DoxyReturn}{Returns}
The returned opacity \mbox{[}m-\/1\mbox{]} 
\end{DoxyReturn}


Definition at line 396 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
396   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{intent(in)}  :: chl\_data\textcolor{comment}{ !< The chlorophyll-A concentration in mg m-3.}
397   \textcolor{keywordtype}{real} :: opacity\_morel\textcolor{comment}{ !< The returned opacity [m-1]}
398 
399   \textcolor{comment}{!   The following are coefficients for the optical model taken from Morel and}
400   \textcolor{comment}{! Antoine (1994). These coeficients represent a non uniform distribution of}
401   \textcolor{comment}{! chlorophyll-a through the water column.  Other approaches may be more}
402   \textcolor{comment}{! appropriate when using an interactive ecosystem model that predicts}
403   \textcolor{comment}{! three-dimensional chl-a values.}
404   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(6)}, \textcolor{keywordtype}{parameter} :: &
405        Z2\_coef=(/7.925, -6.644, 3.662, -1.815, -0.218,  0.502/)
406   \textcolor{keywordtype}{real} :: Chl, Chl2 \textcolor{comment}{! The log10 of chl\_data (in mg m-3), and Chl^2.}
407 
408   chl = log10(min(max(chl\_data,0.02),60.0)) ; chl2 = chl*chl
409   opacity\_morel = 1.0 / ( (z2\_coef(1) + z2\_coef(2)*chl) + chl2 * &
410       ((z2\_coef(3) + chl*z2\_coef(4)) + chl2*(z2\_coef(5) + chl*z2\_coef(6))) )
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a349c6934f113d238e4e2ef229b931a0c}\label{namespacemom__opacity_a349c6934f113d238e4e2ef229b931a0c}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!optics\+\_\+nbands@{optics\+\_\+nbands}}
\index{optics\+\_\+nbands@{optics\+\_\+nbands}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{optics\+\_\+nbands()}{optics\_nbands()}}
{\footnotesize\ttfamily integer function, public mom\+\_\+opacity\+::optics\+\_\+nbands (\begin{DoxyParamCaption}\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(in)}]{optics }\end{DoxyParamCaption})}



Return the number of bands of penetrating shortwave radiation. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em optics} & An optics structure that has values of opacities and shortwave fluxes.\\
\hline
\end{DoxyParams}
\begin{DoxyReturn}{Returns}
The number of penetrating bands of SW radiation 
\end{DoxyReturn}


Definition at line 494 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
494   \textcolor{keywordtype}{type}(optics\_type),       \textcolor{keywordtype}{intent(in)}  :: optics\textcolor{comment}{ !< An optics structure that has values of opacities}
495 \textcolor{comment}{                                                 !! and shortwave fluxes.}
496   \textcolor{keywordtype}{integer} :: optics\_nbands\textcolor{comment}{ !< The number of penetrating bands of SW radiation}
497 
498   optics\_nbands = optics%nbands
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a05ef9c5d86adff869fad832f3083bba4}\label{namespacemom__opacity_a05ef9c5d86adff869fad832f3083bba4}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!set\+\_\+opacity@{set\+\_\+opacity}}
\index{set\+\_\+opacity@{set\+\_\+opacity}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{set\+\_\+opacity()}{set\_opacity()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::set\+\_\+opacity (\begin{DoxyParamCaption}\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(inout)}]{optics,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+total,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+vis\+\_\+dir,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+vis\+\_\+dif,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+nir\+\_\+dir,  }\item[{real, dimension(\+:,\+:), pointer}]{sw\+\_\+nir\+\_\+dif,  }\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(\mbox{\hyperlink{structmom__opacity_1_1opacity__cs}{opacity\+\_\+cs}}), pointer}]{CS,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed), intent(in), optional}]{chl\+\_\+2d,  }\item[{real, dimension( g \%isd\+: g \%ied, g \%jsd\+: g \%jed, gv \%ke), intent(in), optional}]{chl\+\_\+3d }\end{DoxyParamCaption})}



This sets the opacity of sea water based based on one of several different schemes. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in,out}  & {\em optics} & An optics structure that has values set based on the opacities.\\
\hline
 & {\em sw\+\_\+total} & Total shortwave flux into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+vis\+\_\+dir} & Visible, direct shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+vis\+\_\+dif} & Visible, diffuse shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+nir\+\_\+dir} & Near-\/\+IR, direct shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
 & {\em sw\+\_\+nir\+\_\+dif} & Near-\/\+IR, diffuse shortwave into the ocean \mbox{[}Q R Z T-\/1 $\sim$$>$ W m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
 & {\em cs} & The control structure earlier set up by opacity\+\_\+init.\\
\hline
\mbox{\tt in}  & {\em chl\+\_\+2d} & Vertically uniform chlorophyll-\/A concentractions \mbox{[}mg m-\/3\mbox{]}\\
\hline
\mbox{\tt in}  & {\em chl\+\_\+3d} & The chlorophyll-\/A concentractions of each layer \mbox{[}mg m-\/3\mbox{]} \\
\hline
\end{DoxyParams}


Definition at line 93 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
93   \textcolor{keywordtype}{type}(optics\_type),       \textcolor{keywordtype}{intent(inout)} :: optics\textcolor{comment}{ !< An optics structure that has values}
94 \textcolor{comment}{                                                   !! set based on the opacities.}
95   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_total\textcolor{comment}{ !< Total shortwave flux into the ocean [Q R Z T-1 ~> W
       m-2]}
96   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_vis\_dir\textcolor{comment}{ !< Visible, direct shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
97   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_vis\_dif\textcolor{comment}{ !< Visible, diffuse shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
98   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_nir\_dir\textcolor{comment}{ !< Near-IR, direct shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
99   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(:,:)},    \textcolor{keywordtype}{pointer}       :: sw\_nir\_dif\textcolor{comment}{ !< Near-IR, diffuse shortwave into the ocean [Q R Z
       T-1 ~> W m-2]}
100   \textcolor{keywordtype}{type}(ocean\_grid\_type),   \textcolor{keywordtype}{intent(in)}    :: G\textcolor{comment}{      !< The ocean's grid structure.}
101   \textcolor{keywordtype}{type}(verticalGrid\_type), \textcolor{keywordtype}{intent(in)}    :: GV\textcolor{comment}{     !< The ocean's vertical grid structure.}
102   \textcolor{keywordtype}{type}(unit\_scale\_type),   \textcolor{keywordtype}{intent(in)}    :: US\textcolor{comment}{     !< A dimensional unit scaling type}
103   \textcolor{keywordtype}{type}(opacity\_CS),        \textcolor{keywordtype}{pointer}       :: CS\textcolor{comment}{     !< The control structure earlier set up by opacity\_init.}
104   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G))}, &
105                  \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}    :: chl\_2d\textcolor{comment}{ !< Vertically uniform chlorophyll-A concentractions [mg
       m-3]}
106   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZJ\_(G),SZK\_(GV))}, &
107                  \textcolor{keywordtype}{optional}, \textcolor{keywordtype}{intent(in)}    :: chl\_3d\textcolor{comment}{ !< The chlorophyll-A concentractions of each layer [mg
       m-3]}
108 
109   \textcolor{comment}{! Local variables}
110   \textcolor{keywordtype}{integer} :: i, j, k, n, is, ie, js, je, nz
111   \textcolor{keywordtype}{real} :: inv\_sw\_pen\_scale  \textcolor{comment}{! The inverse of the e-folding scale [m-1].}
112   \textcolor{keywordtype}{real} :: Inv\_nbands        \textcolor{comment}{! The inverse of the number of bands of penetrating}
113                             \textcolor{comment}{! shortwave radiation.}
114   \textcolor{keywordtype}{logical} :: call\_for\_surface  \textcolor{comment}{! if horizontal slice is the surface layer}
115   \textcolor{keywordtype}{real} :: tmp(SZI\_(G),SZJ\_(G),SZK\_(GV)) \textcolor{comment}{! A 3-d temporary array.}
116   \textcolor{keywordtype}{real} :: chl(SZI\_(G),SZJ\_(G),SZK\_(GV)) \textcolor{comment}{! The concentration of chlorophyll-A [mg m-3].}
117   \textcolor{keywordtype}{real} :: Pen\_SW\_tot(SZI\_(G),SZJ\_(G))   \textcolor{comment}{! The penetrating shortwave radiation}
118                                         \textcolor{comment}{! summed across all bands [Q R Z T-1 ~> W m-2].}
119   is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
120 
121   \textcolor{keywordflow}{if} (.not. \textcolor{keyword}{associated}(cs)) \textcolor{keyword}{call }mom\_error(fatal, \textcolor{stringliteral}{"set\_opacity: "}// &
122          \textcolor{stringliteral}{"Module must be initialized via opacity\_init before it is used."})
123 
124   \textcolor{keywordflow}{if} (\textcolor{keyword}{present}(chl\_2d) .or. \textcolor{keyword}{present}(chl\_3d)) \textcolor{keywordflow}{then}
125     \textcolor{comment}{! The optical properties are based on cholophyll concentrations.}
126     \textcolor{keyword}{call }opacity\_from\_chl(optics, sw\_total, sw\_vis\_dir, sw\_vis\_dif, sw\_nir\_dir, sw\_nir\_dif, &
127                           g, gv, us, cs, chl\_2d, chl\_3d)
128   \textcolor{keywordflow}{else} \textcolor{comment}{! Use sw e-folding scale set by MOM\_input}
129     \textcolor{keywordflow}{if} (optics%nbands <= 1) \textcolor{keywordflow}{then} ; inv\_nbands = 1.0
130     \textcolor{keywordflow}{else} ; inv\_nbands = 1.0 / \textcolor{keywordtype}{real(optics%nbands)} ; endif
131 
132     \textcolor{comment}{! Make sure there is no division by 0.}
133     inv\_sw\_pen\_scale = 1.0 / max(cs%pen\_sw\_scale, 0.1*gv%Angstrom\_m, &
134                                  gv%H\_to\_m*gv%H\_subroundoff)
135     \textcolor{keywordflow}{if} ( cs%Opacity\_scheme == double\_exp ) \textcolor{keywordflow}{then}
136       \textcolor{comment}{!$OMP parallel do default(shared)}
137       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
138         optics%opacity\_band(1,i,j,k) = inv\_sw\_pen\_scale
139         optics%opacity\_band(2,i,j,k) = 1.0 / max(cs%pen\_sw\_scale\_2nd, &
140              0.1*gv%Angstrom\_m,gv%H\_to\_m*gv%H\_subroundoff)
141 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
142       \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(sw\_total) .or. (cs%pen\_SW\_scale <= 0.0)) \textcolor{keywordflow}{then}
143         \textcolor{comment}{!$OMP parallel do default(shared)}
144         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{do} n=1,optics%nbands
145           optics%sw\_pen\_band(n,i,j) = 0.0
146 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
147       \textcolor{keywordflow}{else}
148         \textcolor{comment}{!$OMP parallel do default(shared)}
149         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
150           optics%sw\_pen\_band(1,i,j) = (cs%SW\_1st\_EXP\_RATIO) * sw\_total(i,j)
151           optics%sw\_pen\_band(2,i,j) = (1.-cs%SW\_1st\_EXP\_RATIO) * sw\_total(i,j)
152 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
153 \textcolor{keywordflow}{      endif}
154     \textcolor{keywordflow}{else}
155       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie  ; \textcolor{keywordflow}{do} n=1,optics%nbands
156         optics%opacity\_band(n,i,j,k) = inv\_sw\_pen\_scale
157 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
158       \textcolor{keywordflow}{if} (.not.\textcolor{keyword}{associated}(sw\_total) .or. (cs%pen\_SW\_scale <= 0.0)) \textcolor{keywordflow}{then}
159         \textcolor{comment}{!$OMP parallel do default(shared)}
160         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{do} n=1,optics%nbands
161           optics%sw\_pen\_band(n,i,j) = 0.0
162 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
163       \textcolor{keywordflow}{else}
164         \textcolor{comment}{!$OMP parallel do default(shared)}
165         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{do} n=1,optics%nbands
166           optics%sw\_pen\_band(n,i,j) = cs%pen\_SW\_frac * inv\_nbands * sw\_total(i,j)
167 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
168 \textcolor{keywordflow}{      endif}
169 \textcolor{keywordflow}{    endif}
170 \textcolor{keywordflow}{  endif}
171 
172   \textcolor{keywordflow}{if} (query\_averaging\_enabled(cs%diag)) \textcolor{keywordflow}{then}
173     \textcolor{keywordflow}{if} (cs%id\_sw\_pen > 0) \textcolor{keywordflow}{then}
174       \textcolor{comment}{!$OMP parallel do default(shared)}
175       \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
176         pen\_sw\_tot(i,j) = 0.0
177         \textcolor{keywordflow}{do} n=1,optics%nbands
178           pen\_sw\_tot(i,j) = pen\_sw\_tot(i,j) + optics%sw\_pen\_band(n,i,j)
179 \textcolor{keywordflow}{        enddo}
180 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo}
181       \textcolor{keyword}{call }post\_data(cs%id\_sw\_pen, pen\_sw\_tot, cs%diag)
182 \textcolor{keywordflow}{    endif}
183     \textcolor{keywordflow}{if} (cs%id\_sw\_vis\_pen > 0) \textcolor{keywordflow}{then}
184       \textcolor{keywordflow}{if} (cs%opacity\_scheme == manizza\_05) \textcolor{keywordflow}{then}
185         \textcolor{comment}{!$OMP parallel do default(shared)}
186         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
187           pen\_sw\_tot(i,j) = 0.0
188           \textcolor{keywordflow}{do} n=1,min(optics%nbands,2)
189             pen\_sw\_tot(i,j) = pen\_sw\_tot(i,j) + optics%sw\_pen\_band(n,i,j)
190 \textcolor{keywordflow}{          enddo}
191 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
192       \textcolor{keywordflow}{else}
193         \textcolor{comment}{!$OMP parallel do default(shared)}
194         \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
195           pen\_sw\_tot(i,j) = 0.0
196           \textcolor{keywordflow}{do} n=1,optics%nbands
197             pen\_sw\_tot(i,j) = pen\_sw\_tot(i,j) + optics%sw\_pen\_band(n,i,j)
198 \textcolor{keywordflow}{          enddo}
199 \textcolor{keywordflow}{        enddo} ;\textcolor{keywordflow}{ enddo}
200 \textcolor{keywordflow}{      endif}
201       \textcolor{keyword}{call }post\_data(cs%id\_sw\_vis\_pen, pen\_sw\_tot, cs%diag)
202 \textcolor{keywordflow}{    endif}
203     \textcolor{keywordflow}{do} n=1,optics%nbands ; \textcolor{keywordflow}{if} (cs%id\_opacity(n) > 0) \textcolor{keywordflow}{then}
204       \textcolor{comment}{!$OMP parallel do default(shared)}
205       \textcolor{keywordflow}{do} k=1,nz ; \textcolor{keywordflow}{do} j=js,je ; \textcolor{keywordflow}{do} i=is,ie
206         \textcolor{comment}{! Remap opacity (op) to 1/L * tanh(op * L) where L is one Angstrom.}
207         \textcolor{comment}{! This gives a nearly identical value when op << 1/L but allows one to}
208         \textcolor{comment}{! store the values when opacity is divergent (i.e. opaque).}
209         tmp(i,j,k) = tanh(op\_diag\_len * optics%opacity\_band(n,i,j,k)) / op\_diag\_len
210 \textcolor{keywordflow}{      enddo} ;\textcolor{keywordflow}{ enddo} ;\textcolor{keywordflow}{ enddo}
211       \textcolor{keyword}{call }post\_data(cs%id\_opacity(n), tmp, cs%diag)
212 \textcolor{keywordflow}{    endif} ;\textcolor{keywordflow}{ enddo}
213 \textcolor{keywordflow}{  endif}
214 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_ad27db4bd0d010d98a3f5a54902c7a05e}\label{namespacemom__opacity_ad27db4bd0d010d98a3f5a54902c7a05e}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!sumswoverbands@{sumswoverbands}}
\index{sumswoverbands@{sumswoverbands}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{sumswoverbands()}{sumswoverbands()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+opacity\+::sumswoverbands (\begin{DoxyParamCaption}\item[{type(ocean\+\_\+grid\+\_\+type), intent(in)}]{G,  }\item[{type(verticalgrid\+\_\+type), intent(in)}]{GV,  }\item[{type(unit\+\_\+scale\+\_\+type), intent(in)}]{US,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke), intent(in)}]{h,  }\item[{integer, intent(in)}]{nsw,  }\item[{type(\mbox{\hyperlink{structmom__opacity_1_1optics__type}{optics\+\_\+type}}), intent(in)}]{optics,  }\item[{integer, intent(in)}]{j,  }\item[{real, intent(in)}]{dt,  }\item[{real, intent(in)}]{H\+\_\+limit\+\_\+fluxes,  }\item[{logical, intent(in)}]{absorb\+All\+SW,  }\item[{real, dimension(max(nsw,1), g \%isd\+: g \%ied), intent(in)}]{i\+Pen\+\_\+\+S\+W\+\_\+bnd,  }\item[{real, dimension( g \%isd\+: g \%ied, gv \%ke+1), intent(inout)}]{net\+Pen }\end{DoxyParamCaption})}



This subroutine calculates the total shortwave heat flux integrated over bands as a function of depth. This routine is only called for computing buoyancy fluxes for use in K\+PP. This routine does not updat e the state. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em g} & The ocean\textquotesingle{}s grid structure.\\
\hline
\mbox{\tt in}  & {\em gv} & The ocean\textquotesingle{}s vertical grid structure.\\
\hline
\mbox{\tt in}  & {\em us} & A dimensional unit scaling type\\
\hline
\mbox{\tt in}  & {\em h} & Layer thicknesses \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em nsw} & The number of bands of penetrating shortwave radiation, perhaps from optics\+\_\+nbands(optics),\\
\hline
\mbox{\tt in}  & {\em optics} & An optics structure that has values set based on the opacities.\\
\hline
\mbox{\tt in}  & {\em j} & j-\/index to work on.\\
\hline
\mbox{\tt in}  & {\em dt} & Time step \mbox{[}T $\sim$$>$ s\mbox{]}.\\
\hline
\mbox{\tt in}  & {\em h\+\_\+limit\+\_\+fluxes} & the total depth at which the surface fluxes start to be limited to avoid excessive heating of a thin ocean \mbox{[}H $\sim$$>$ m or kg m-\/2\mbox{]}\\
\hline
\mbox{\tt in}  & {\em absorballsw} & If true, ensure that all shortwave radiation is absorbed in the ocean water column.\\
\hline
\mbox{\tt in}  & {\em ipen\+\_\+sw\+\_\+bnd} & The incident penetrating shortwave heating in each band that hits the bottom and will be redistributed through the water column \mbox{[}degC H $\sim$$>$ degC m or degC kg m-\/2\mbox{]}; size nsw x G isd\+: G ied.\\
\hline
\mbox{\tt in,out}  & {\em netpen} & Net penetrating shortwave heat flux at each \\
\hline
\end{DoxyParams}


Definition at line 781 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
781   \textcolor{keywordtype}{type}(ocean\_grid\_type),    \textcolor{keywordtype}{intent(in)}    :: G\textcolor{comment}{   !< The ocean's grid structure.}
782   \textcolor{keywordtype}{type}(verticalGrid\_type),  \textcolor{keywordtype}{intent(in)}    :: GV\textcolor{comment}{  !< The ocean's vertical grid structure.}
783   \textcolor{keywordtype}{type}(unit\_scale\_type),    \textcolor{keywordtype}{intent(in)}    :: US\textcolor{comment}{    !< A dimensional unit scaling type}
784   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV))}, &
785                             \textcolor{keywordtype}{intent(in)}    :: h\textcolor{comment}{   !< Layer thicknesses [H ~> m or kg m-2].}
786   \textcolor{keywordtype}{integer},                  \textcolor{keywordtype}{intent(in)}    :: nsw\textcolor{comment}{ !< The number of bands of penetrating shortwave}
787 \textcolor{comment}{                                                 !! radiation, perhaps from optics\_nbands(optics),}
788   \textcolor{keywordtype}{type}(optics\_type),        \textcolor{keywordtype}{intent(in)}    :: optics\textcolor{comment}{ !< An optics structure that has values}
789 \textcolor{comment}{                                                   !! set based on the opacities.}
790   \textcolor{keywordtype}{integer},                  \textcolor{keywordtype}{intent(in)}    :: j\textcolor{comment}{   !< j-index to work on.}
791   \textcolor{keywordtype}{real},                     \textcolor{keywordtype}{intent(in)}    :: dt\textcolor{comment}{  !< Time step [T ~> s].}
792   \textcolor{keywordtype}{real},                     \textcolor{keywordtype}{intent(in)}    :: H\_limit\_fluxes\textcolor{comment}{ !< the total depth at which the}
793 \textcolor{comment}{                                                 !! surface fluxes start to be limited to avoid}
794 \textcolor{comment}{                                                 !! excessive heating of a thin ocean [H ~> m or kg m-2]}
795   \textcolor{keywordtype}{logical},                  \textcolor{keywordtype}{intent(in)}    :: absorbAllSW\textcolor{comment}{ !< If true, ensure that all shortwave}
796 \textcolor{comment}{                                                 !! radiation is absorbed in the ocean water column.}
797   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(max(nsw,1),SZI\_(G))}, \textcolor{keywordtype}{intent(in)} :: iPen\_SW\_bnd\textcolor{comment}{ !< The incident penetrating shortwave}
798 \textcolor{comment}{                                                 !! heating in each band that hits the bottom and}
799 \textcolor{comment}{                                                 !! will be redistributed through the water column}
800 \textcolor{comment}{                                                 !! [degC H ~> degC m or degC kg m-2]; size nsw x SZI\_(G).}
801   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(SZI\_(G),SZK\_(GV)+1)}, &
802                              \textcolor{keywordtype}{intent(inout)} :: netPen\textcolor{comment}{ !< Net penetrating shortwave heat flux at each}
803 \textcolor{comment}{                                                 !! interface, summed across all bands}
804 \textcolor{comment}{                                                 !! [degC H ~> degC m or degC kg m-2].}
805   \textcolor{comment}{! Local variables}
806   \textcolor{keywordtype}{real} :: h\_heat(SZI\_(G))     \textcolor{comment}{! thickness of the water column that receives}
807                               \textcolor{comment}{! remaining shortwave radiation [H ~> m or kg m-2].}
808   \textcolor{keywordtype}{real} :: Pen\_SW\_rem(SZI\_(G)) \textcolor{comment}{! sum across all wavelength bands of the}
809                               \textcolor{comment}{! penetrating shortwave heating that hits the bottom}
810                               \textcolor{comment}{! and will be redistributed through the water column}
811                               \textcolor{comment}{! [degC H ~> degC m or degC kg m-2]}
812 
813   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(max(nsw,1),SZI\_(G))} :: Pen\_SW\_bnd \textcolor{comment}{! The remaining penetrating shortwave radiation}
814                           \textcolor{comment}{! in each band, initially iPen\_SW\_bnd [degC H ~> degC m or degC kg m-2]}
815   \textcolor{keywordtype}{real} :: SW\_trans        \textcolor{comment}{! fraction of shortwave radiation not}
816                           \textcolor{comment}{! absorbed in a layer [nondim]}
817   \textcolor{keywordtype}{real} :: unabsorbed      \textcolor{comment}{! fraction of the shortwave radiation}
818                           \textcolor{comment}{! not absorbed because the layers are too thin.}
819   \textcolor{keywordtype}{real} :: Ih\_limit        \textcolor{comment}{! inverse of the total depth at which the}
820                           \textcolor{comment}{! surface fluxes start to be limited [H-1 ~> m-1 or m2 kg-1]}
821   \textcolor{keywordtype}{real} :: min\_SW\_heat     \textcolor{comment}{! A minimum remaining shortwave heating within a timestep that will be simply}
822                           \textcolor{comment}{! absorbed in the next layer for computational efficiency, instead of}
823                           \textcolor{comment}{! continuing to penetrate [degC H ~> degC m or degC kg m-2].}
824   \textcolor{keywordtype}{real} :: I\_Habs            \textcolor{comment}{! The inverse of the absorption length for a minimal flux [H-1 ~> m-1 or m2
       kg-1]}
825   \textcolor{keywordtype}{real} :: h\_min\_heat      \textcolor{comment}{! minimum thickness layer that should get heated [H ~> m or kg m-2]}
826   \textcolor{keywordtype}{real} :: opt\_depth       \textcolor{comment}{! optical depth of a layer [nondim]}
827   \textcolor{keywordtype}{real} :: exp\_OD          \textcolor{comment}{! exp(-opt\_depth) [nondim]}
828   \textcolor{keywordtype}{logical} :: SW\_Remains   \textcolor{comment}{! If true, some column has shortwave radiation that}
829                           \textcolor{comment}{! was not entirely absorbed.}
830 
831   \textcolor{keywordtype}{integer} :: is, ie, nz, i, k, ks, n
832   sw\_remains = .false.
833 
834   min\_sw\_heat = optics%PenSW\_flux\_absorb*dt \textcolor{comment}{! Default of 2.5e-11*US%T\_to\_s*GV%m\_to\_H}
835   i\_habs = 1e3*gv%H\_to\_m \textcolor{comment}{! optics%PenSW\_absorb\_Invlen}
836 
837   h\_min\_heat = 2.0*gv%Angstrom\_H + gv%H\_subroundoff
838   is = g%isc ; ie = g%iec ; nz = g%ke
839 
840   pen\_sw\_bnd(:,:) = ipen\_sw\_bnd(:,:)
841   \textcolor{keywordflow}{do} i=is,ie ; h\_heat(i) = 0.0 ;\textcolor{keywordflow}{ enddo}
842   \textcolor{keywordflow}{do} i=is,ie
843     netpen(i,1) = 0.
844     \textcolor{keywordflow}{do} n=1,max(nsw,1)
845       netpen(i,1) = netpen(i,1) + pen\_sw\_bnd(n,i)   \textcolor{comment}{! Surface interface}
846 \textcolor{keywordflow}{    enddo}
847 \textcolor{keywordflow}{  enddo}
848 
849   \textcolor{comment}{! Apply penetrating SW radiation to remaining parts of layers.}
850   \textcolor{comment}{! Excessively thin layers are not heated to avoid runaway temps.}
851   \textcolor{keywordflow}{do} k=1,nz
852 
853     \textcolor{keywordflow}{do} i=is,ie
854       netpen(i,k+1) = 0.
855 
856       \textcolor{keywordflow}{if} (h(i,k) > 0.0) \textcolor{keywordflow}{then}
857         \textcolor{keywordflow}{do} n=1,nsw ; \textcolor{keywordflow}{if} (pen\_sw\_bnd(n,i) > 0.0) \textcolor{keywordflow}{then}
858           \textcolor{comment}{! SW\_trans is the SW that is transmitted THROUGH the layer}
859           opt\_depth = h(i,k)*gv%H\_to\_m * optics%opacity\_band(n,i,j,k)
860           exp\_od = exp(-opt\_depth)
861           sw\_trans = exp\_od
862 
863           \textcolor{comment}{! Heating at a very small rate can be absorbed by a sufficiently thick layer or several}
864           \textcolor{comment}{! thin layers without further penetration.}
865           \textcolor{keywordflow}{if} (optics%answers\_2018) \textcolor{keywordflow}{then}
866             \textcolor{keywordflow}{if} (nsw*pen\_sw\_bnd(n,i)*sw\_trans < min\_sw\_heat*min(1.0, i\_habs*h(i,k)) ) sw\_trans = 0.0
867           \textcolor{keywordflow}{elseif} ((nsw*pen\_sw\_bnd(n,i)*sw\_trans < min\_sw\_heat) .and. (h(i,k) > h\_min\_heat)) \textcolor{keywordflow}{then}
868             \textcolor{keywordflow}{if} (nsw*pen\_sw\_bnd(n,i) <= min\_sw\_heat * (i\_habs*(h(i,k) - h\_min\_heat))) \textcolor{keywordflow}{then}
869               sw\_trans = 0.0
870             \textcolor{keywordflow}{else}
871               sw\_trans = min(sw\_trans, &
872                              1.0 - (min\_sw\_heat*(i\_habs*(h(i,k) - h\_min\_heat))) / (nsw*pen\_sw\_bnd(n,i)))
873 \textcolor{keywordflow}{            endif}
874 \textcolor{keywordflow}{          endif}
875 
876           pen\_sw\_bnd(n,i) = pen\_sw\_bnd(n,i) * sw\_trans
877           netpen(i,k+1)   = netpen(i,k+1) + pen\_sw\_bnd(n,i)
878 \textcolor{keywordflow}{        endif} ;\textcolor{keywordflow}{ enddo}
879 \textcolor{keywordflow}{      endif} \textcolor{comment}{! h(i,k) > 0.0}
880 
881       \textcolor{comment}{! Add to the accumulated thickness above that could be heated.}
882       \textcolor{comment}{! Only layers greater than h\_min\_heat thick should get heated.}
883       \textcolor{keywordflow}{if} (h(i,k) >= 2.0*h\_min\_heat) \textcolor{keywordflow}{then}
884         h\_heat(i) = h\_heat(i) + h(i,k)
885       \textcolor{keywordflow}{elseif} (h(i,k) > h\_min\_heat) \textcolor{keywordflow}{then}
886         h\_heat(i) = h\_heat(i) + (2.0*h(i,k) - 2.0*h\_min\_heat)
887 \textcolor{keywordflow}{      endif}
888 \textcolor{keywordflow}{    enddo} \textcolor{comment}{! i loop}
889 \textcolor{keywordflow}{  enddo} \textcolor{comment}{! k loop}
890 
891   \textcolor{keywordflow}{if} (absorballsw) \textcolor{keywordflow}{then}
892 
893     \textcolor{comment}{! If there is still shortwave radiation at this point, it could go into}
894     \textcolor{comment}{! the bottom (with a bottom mud model), or it could be redistributed back}
895     \textcolor{comment}{! through the water column.}
896     \textcolor{keywordflow}{do} i=is,ie
897       pen\_sw\_rem(i) = pen\_sw\_bnd(1,i)
898       \textcolor{keywordflow}{do} n=2,nsw ; pen\_sw\_rem(i) = pen\_sw\_rem(i) + pen\_sw\_bnd(n,i) ;\textcolor{keywordflow}{ enddo}
899 \textcolor{keywordflow}{    enddo}
900     \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{if} (pen\_sw\_rem(i) > 0.0) sw\_remains = .true. ;\textcolor{keywordflow}{ enddo}
901 
902     ih\_limit = 1.0 / h\_limit\_fluxes
903     \textcolor{keywordflow}{do} i=is,ie ; \textcolor{keywordflow}{if} ((pen\_sw\_rem(i) > 0.0) .and. (h\_heat(i) > 0.0)) \textcolor{keywordflow}{then}
904       \textcolor{keywordflow}{if} (h\_heat(i)*ih\_limit < 1.0) \textcolor{keywordflow}{then}
905         unabsorbed = 1.0 - h\_heat(i)*ih\_limit
906       \textcolor{keywordflow}{else}
907         unabsorbed = 0.0
908 \textcolor{keywordflow}{      endif}
909       \textcolor{keywordflow}{do} n=1,nsw ; pen\_sw\_bnd(n,i) = unabsorbed * pen\_sw\_bnd(n,i) ;\textcolor{keywordflow}{ enddo}
910 \textcolor{keywordflow}{    endif} ;\textcolor{keywordflow}{ enddo}
911 
912 \textcolor{keywordflow}{  endif} \textcolor{comment}{! absorbAllSW}
913 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__opacity_a0017241c03e4536115674fc5fc9608bf}\label{namespacemom__opacity_a0017241c03e4536115674fc5fc9608bf}} 
\index{mom\+\_\+opacity@{mom\+\_\+opacity}!sw\+\_\+pen\+\_\+frac\+\_\+morel@{sw\+\_\+pen\+\_\+frac\+\_\+morel}}
\index{sw\+\_\+pen\+\_\+frac\+\_\+morel@{sw\+\_\+pen\+\_\+frac\+\_\+morel}!mom\+\_\+opacity@{mom\+\_\+opacity}}
\subsubsection{\texorpdfstring{sw\+\_\+pen\+\_\+frac\+\_\+morel()}{sw\_pen\_frac\_morel()}}
{\footnotesize\ttfamily real function mom\+\_\+opacity\+::sw\+\_\+pen\+\_\+frac\+\_\+morel (\begin{DoxyParamCaption}\item[{real, intent(in)}]{chl\+\_\+data }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



This sets the penetrating shortwave fraction according to the scheme proposed by Morel and Antoine (1994). 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em chl\+\_\+data} & The chlorophyll-\/A concentration in mg m-\/3.\\
\hline
\end{DoxyParams}
\begin{DoxyReturn}{Returns}
The returned penetrating shortwave fraction \mbox{[}nondim\mbox{]} 
\end{DoxyReturn}


Definition at line 416 of file M\+O\+M\+\_\+opacity.\+F90.


\begin{DoxyCode}
416   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{intent(in)}  :: chl\_data\textcolor{comment}{ !< The chlorophyll-A concentration in mg m-3.}
417   \textcolor{keywordtype}{real} :: SW\_pen\_frac\_morel\textcolor{comment}{     !< The returned penetrating shortwave fraction [nondim]}
418 
419   \textcolor{comment}{!   The following are coefficients for the optical model taken from Morel and}
420   \textcolor{comment}{! Antoine (1994). These coeficients represent a non uniform distribution of}
421   \textcolor{comment}{! chlorophyll-a through the water column.  Other approaches may be more}
422   \textcolor{comment}{! appropriate when using an interactive ecosystem model that predicts}
423   \textcolor{comment}{! three-dimensional chl-a values.}
424   \textcolor{keywordtype}{real} :: Chl, Chl2         \textcolor{comment}{! The log10 of chl\_data in mg m-3, and Chl^2.}
425   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(6)}, \textcolor{keywordtype}{parameter} :: &
426        V1\_coef=(/0.321,  0.008, 0.132,  0.038, -0.017, -0.007/)
427 
428   chl = log10(min(max(chl\_data,0.02),60.0)) ; chl2 = chl*chl
429   sw\_pen\_frac\_morel = 1.0 - ( (v1\_coef(1) + v1\_coef(2)*chl) + chl2 * &
430        ((v1\_coef(3) + chl*v1\_coef(4)) + chl2*(v1\_coef(5) + chl*v1\_coef(6))) )
\end{DoxyCode}
