\hypertarget{namespacemom__diag__vkernels}{}\section{mom\+\_\+diag\+\_\+vkernels Module Reference}
\label{namespacemom__diag__vkernels}\index{mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}}


\subsection{Detailed Description}
Provides kernels for single-\/column interpolation, re-\/integration (re-\/mapping of integrated quantities) and intensive-\/variable remapping in the vertical. \subsection*{Functions/\+Subroutines}
\begin{DoxyCompactItemize}
\item 
subroutine, public \mbox{\hyperlink{namespacemom__diag__vkernels_a4f42f472a725a147f8d97a68b2028c5b}{interpolate\+\_\+column}} (nsrc, h\+\_\+src, u\+\_\+src, ndest, h\+\_\+dest, missing\+\_\+value, u\+\_\+dest)
\begin{DoxyCompactList}\small\item\em Linearly interpolate interface data, u\+\_\+src, from grid h\+\_\+src to a grid h\+\_\+dest. \end{DoxyCompactList}\item 
subroutine, public \mbox{\hyperlink{namespacemom__diag__vkernels_a89369e4bf4b7795f6e384762d11f0d23}{reintegrate\+\_\+column}} (nsrc, h\+\_\+src, uh\+\_\+src, ndest, h\+\_\+dest, missing\+\_\+value, uh\+\_\+dest)
\begin{DoxyCompactList}\small\item\em Conservatively calculate integrated data, uh\+\_\+dest, on grid h\+\_\+dest, from layer-\/integrated data, uh\+\_\+src, on grid h\+\_\+src. \end{DoxyCompactList}\item 
logical function, public \mbox{\hyperlink{namespacemom__diag__vkernels_a6001aaa22610f7fd690106fc737775dd}{diag\+\_\+vkernels\+\_\+unit\+\_\+tests}} (verbose)
\begin{DoxyCompactList}\small\item\em Returns true if any unit tests for module M\+O\+M\+\_\+diag\+\_\+vkernels fail. \end{DoxyCompactList}\item 
logical function \mbox{\hyperlink{namespacemom__diag__vkernels_abeef457cda28b20c03a89ad402bcd434}{test\+\_\+interp}} (verbose, missing\+\_\+value, msg, nsrc, h\+\_\+src, u\+\_\+src, ndest, h\+\_\+dest, u\+\_\+true)
\begin{DoxyCompactList}\small\item\em Returns true if a test of \mbox{\hyperlink{namespacemom__diag__vkernels_a4f42f472a725a147f8d97a68b2028c5b}{interpolate\+\_\+column()}} produces the wrong answer. \end{DoxyCompactList}\item 
logical function \mbox{\hyperlink{namespacemom__diag__vkernels_aac7aad5875fcf17273b211c90ac0cbee}{test\+\_\+reintegrate}} (verbose, missing\+\_\+value, msg, nsrc, h\+\_\+src, uh\+\_\+src, ndest, h\+\_\+dest, uh\+\_\+true)
\begin{DoxyCompactList}\small\item\em Returns true if a test of \mbox{\hyperlink{namespacemom__diag__vkernels_a89369e4bf4b7795f6e384762d11f0d23}{reintegrate\+\_\+column()}} produces the wrong answer. \end{DoxyCompactList}\end{DoxyCompactItemize}


\subsection{Function/\+Subroutine Documentation}
\mbox{\Hypertarget{namespacemom__diag__vkernels_a6001aaa22610f7fd690106fc737775dd}\label{namespacemom__diag__vkernels_a6001aaa22610f7fd690106fc737775dd}} 
\index{mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}!diag\+\_\+vkernels\+\_\+unit\+\_\+tests@{diag\+\_\+vkernels\+\_\+unit\+\_\+tests}}
\index{diag\+\_\+vkernels\+\_\+unit\+\_\+tests@{diag\+\_\+vkernels\+\_\+unit\+\_\+tests}!mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}}
\subsubsection{\texorpdfstring{diag\+\_\+vkernels\+\_\+unit\+\_\+tests()}{diag\_vkernels\_unit\_tests()}}
{\footnotesize\ttfamily logical function, public mom\+\_\+diag\+\_\+vkernels\+::diag\+\_\+vkernels\+\_\+unit\+\_\+tests (\begin{DoxyParamCaption}\item[{logical, intent(in)}]{verbose }\end{DoxyParamCaption})}



Returns true if any unit tests for module M\+O\+M\+\_\+diag\+\_\+vkernels fail. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em verbose} & If true, write results to stdout \\
\hline
\end{DoxyParams}


Definition at line 171 of file M\+O\+M\+\_\+diag\+\_\+vkernels.\+F90.


\begin{DoxyCode}
171   \textcolor{keywordtype}{logical}, \textcolor{keywordtype}{intent(in)} :: verbose\textcolor{comment}{ !< If true, write results to stdout}
172   \textcolor{comment}{! Local variables}
173   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{parameter} :: mv=-9.999999999e9 \textcolor{comment}{! Value to use for vanished layers}
174   \textcolor{keywordtype}{logical} :: fail, v
175 
176   v = verbose
177 
178   \textcolor{keyword}{write}(stdout,*) \textcolor{stringliteral}{'==== MOM\_diag\_kernels: diag\_vkernels\_unit\_tests =========='}
179   \textcolor{keywordflow}{if} (v) \textcolor{keyword}{write}(stdout,*) \textcolor{stringliteral}{'- - - - - - - - - - interpolation tests  - - - - - - - - -'}
180 
181   fail = test\_interp(v,mv,\textcolor{stringliteral}{'Identity: 3 layer'}, &
182                      3, (/1.,2.,3./), (/1.,2.,3.,4./), &
183                      3, (/1.,2.,3./), (/1.,2.,3.,4./) )
184   diag\_vkernels\_unit\_tests = fail
185 
186   fail = test\_interp(v,mv,\textcolor{stringliteral}{'A: 3 layer to 2'}, &
187                      3, (/1.,1.,1./), (/1.,2.,3.,4./), &
188                      2, (/1.5,1.5/), (/1.,2.5,4./) )
189   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
190 
191   fail = test\_interp(v,mv,\textcolor{stringliteral}{'B: 2 layer to 3'}, &
192                      2, (/1.5,1.5/), (/1.,4.,7./), &
193                      3, (/1.,1.,1./), (/1.,3.,5.,7./) )
194   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
195 
196   fail = test\_interp(v,mv,\textcolor{stringliteral}{'C: 3 layer (vanished middle) to 2'}, &
197                      3, (/1.,0.,2./), (/1.,2.,2.,3./), &
198                      2, (/1.,2./), (/1.,2.,3./) )
199   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
200 
201   fail = test\_interp(v,mv,\textcolor{stringliteral}{'D: 3 layer (deep) to 3'}, &
202                      3, (/1.,2.,3./), (/1.,2.,4.,7./), &
203                      2, (/2.,2./), (/1.,3.,5./) )
204   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
205 
206   fail = test\_interp(v,mv,\textcolor{stringliteral}{'E: 3 layer to 3 (deep)'}, &
207                      3, (/1.,2.,4./), (/1.,2.,4.,8./), &
208                      3, (/2.,3.,4./), (/1.,3.,6.,8./) )
209   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
210 
211   fail = test\_interp(v,mv,\textcolor{stringliteral}{'F: 3 layer to 4 with vanished top/botton'}, &
212                      3, (/1.,2.,4./), (/1.,2.,4.,8./), &
213                      4, (/0.,2.,5.,0./), (/mv,1.,3.,8.,mv/) )
214   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
215 
216   fail = test\_interp(v,mv,\textcolor{stringliteral}{'Fs: 3 layer to 4 with vanished top/botton (shallow)'}, &
217                      3, (/1.,2.,4./), (/1.,2.,4.,8./), &
218                      4, (/0.,2.,4.,0./), (/mv,1.,3.,7.,mv/) )
219   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
220 
221   fail = test\_interp(v,mv,\textcolor{stringliteral}{'Fd: 3 layer to 4 with vanished top/botton (deep)'}, &
222                      3, (/1.,2.,4./), (/1.,2.,4.,8./), &
223                      4, (/0.,2.,6.,0./), (/mv,1.,3.,8.,mv/) )
224   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
225 
226   \textcolor{keywordflow}{if} (v) \textcolor{keyword}{write}(stdout,*) \textcolor{stringliteral}{'- - - - - - - - - - reintegration tests  - - - - - - - - -'}
227 
228   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'Identity: 3 layer'}, &
229                      3, (/1.,2.,3./), (/-5.,2.,1./), &
230                      3, (/1.,2.,3./), (/-5.,2.,1./) )
231   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
232 
233   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'A: 3 layer to 2'}, &
234                      3, (/2.,2.,2./), (/-5.,2.,1./), &
235                      2, (/3.,3./), (/-4.,2./) )
236   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
237 
238   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'A: 3 layer to 2 (deep)'}, &
239                      3, (/2.,2.,2./), (/-5.,2.,1./), &
240                      2, (/3.,4./), (/-4.,2./) )
241   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
242 
243   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'A: 3 layer to 2 (shallow)'}, &
244                      3, (/2.,2.,2./), (/-5.,2.,1./), &
245                      2, (/3.,2./), (/-4.,1.5/) )
246   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
247 
248   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'B: 3 layer to 4 with vanished top/bottom'}, &
249                      3, (/2.,2.,2./), (/-5.,2.,1./), &
250                      4, (/0.,3.,3.,0./), (/0.,-4.,2.,0./) )
251   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
252 
253   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'C: 3 layer to 4 with vanished top//middle/bottom'}, &
254                      3, (/2.,2.,2./), (/-5.,2.,1./), &
255                      5, (/0.,3.,0.,3.,0./), (/0.,-4.,0.,2.,0./) )
256   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
257 
258   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'D: 3 layer to 3 (vanished)'}, &
259                      3, (/2.,2.,2./), (/-5.,2.,1./), &
260                      3, (/0.,0.,0./), (/0.,0.,0./) )
261   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
262 
263   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'D: 3 layer (vanished) to 3'}, &
264                      3, (/0.,0.,0./), (/-5.,2.,1./), &
265                      3, (/2.,2.,2./), (/mv, mv, mv/) )
266   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
267 
268   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'D: 3 layer (vanished) to 3 (vanished)'}, &
269                      3, (/0.,0.,0./), (/-5.,2.,1./), &
270                      3, (/0.,0.,0./), (/mv, mv, mv/) )
271   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
272 
273   fail = test\_reintegrate(v,mv,\textcolor{stringliteral}{'D: 3 layer (vanished) to 3 (vanished)'}, &
274                      3, (/0.,0.,0./), (/0.,0.,0./), &
275                      3, (/0.,0.,0./), (/mv, mv, mv/) )
276   diag\_vkernels\_unit\_tests = diag\_vkernels\_unit\_tests .or. fail
277 
278   \textcolor{keywordflow}{if} (.not. fail) \textcolor{keyword}{write}(stdout,*) \textcolor{stringliteral}{'Pass'}
279 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__diag__vkernels_a4f42f472a725a147f8d97a68b2028c5b}\label{namespacemom__diag__vkernels_a4f42f472a725a147f8d97a68b2028c5b}} 
\index{mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}!interpolate\+\_\+column@{interpolate\+\_\+column}}
\index{interpolate\+\_\+column@{interpolate\+\_\+column}!mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}}
\subsubsection{\texorpdfstring{interpolate\+\_\+column()}{interpolate\_column()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+diag\+\_\+vkernels\+::interpolate\+\_\+column (\begin{DoxyParamCaption}\item[{integer, intent(in)}]{nsrc,  }\item[{real, dimension(nsrc), intent(in)}]{h\+\_\+src,  }\item[{real, dimension(nsrc+1), intent(in)}]{u\+\_\+src,  }\item[{integer, intent(in)}]{ndest,  }\item[{real, dimension(ndest), intent(in)}]{h\+\_\+dest,  }\item[{real, intent(in)}]{missing\+\_\+value,  }\item[{real, dimension(ndest+1), intent(inout)}]{u\+\_\+dest }\end{DoxyParamCaption})}



Linearly interpolate interface data, u\+\_\+src, from grid h\+\_\+src to a grid h\+\_\+dest. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em nsrc} & Number of source cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+src} & Thickness of source cells\\
\hline
\mbox{\tt in}  & {\em u\+\_\+src} & Values at source cell interfaces\\
\hline
\mbox{\tt in}  & {\em ndest} & Number of destination cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+dest} & Thickness of destination cells\\
\hline
\mbox{\tt in}  & {\em missing\+\_\+value} & Value to assign in vanished cells\\
\hline
\mbox{\tt in,out}  & {\em u\+\_\+dest} & Interpolated value at destination cell interfaces \\
\hline
\end{DoxyParams}


Definition at line 19 of file M\+O\+M\+\_\+diag\+\_\+vkernels.\+F90.


\begin{DoxyCode}
19   \textcolor{keywordtype}{integer},                  \textcolor{keywordtype}{intent(in)}    :: nsrc\textcolor{comment}{ !< Number of source cells}
20   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc)},    \textcolor{keywordtype}{intent(in)}    :: h\_src\textcolor{comment}{ !< Thickness of source cells}
21   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc+1)},  \textcolor{keywordtype}{intent(in)}    :: u\_src\textcolor{comment}{ !< Values at source cell interfaces}
22   \textcolor{keywordtype}{integer},                  \textcolor{keywordtype}{intent(in)}    :: ndest\textcolor{comment}{ !< Number of destination cells}
23   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)},   \textcolor{keywordtype}{intent(in)}    :: h\_dest\textcolor{comment}{ !< Thickness of destination cells}
24   \textcolor{keywordtype}{real},                     \textcolor{keywordtype}{intent(in)}    :: missing\_value\textcolor{comment}{ !< Value to assign in vanished cells}
25   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest+1)}, \textcolor{keywordtype}{intent(inout)} :: u\_dest\textcolor{comment}{ !< Interpolated value at destination cell interfaces}
26   \textcolor{comment}{! Local variables}
27   \textcolor{keywordtype}{real} :: x\_dest \textcolor{comment}{! Relative position of target interface}
28   \textcolor{keywordtype}{real} :: dh \textcolor{comment}{! Source cell thickness}
29   \textcolor{keywordtype}{real} :: u1, u2 \textcolor{comment}{! Values to interpolate between}
30   \textcolor{keywordtype}{real} :: weight\_a, weight\_b \textcolor{comment}{! Weights for interpolation}
31   \textcolor{keywordtype}{integer} :: k\_src, k\_dest \textcolor{comment}{! Index of cell in src and dest columns}
32   \textcolor{keywordtype}{logical} :: still\_vanished \textcolor{comment}{! Used for figuring out what to mask as missing}
33 
34   \textcolor{comment}{! Initial values for the loop}
35   still\_vanished = .true.
36 
37   \textcolor{comment}{! The following forces the "do while" loop to do one cycle that will set u1, u2, dh.}
38   k\_src = 0
39   dh = 0.
40   x\_dest = 0.
41 
42   \textcolor{keywordflow}{do} k\_dest=1, ndest+1
43     \textcolor{keywordflow}{do} \textcolor{keywordflow}{while} (dh<=x\_dest .and. k\_src<nsrc)
44       \textcolor{comment}{! Move positions pointers forward until the interval 0 .. dh spans x\_dest.}
45       x\_dest = x\_dest - dh
46       k\_src = k\_src + 1
47       dh = h\_src(k\_src) \textcolor{comment}{! Thickness of layer k\_src}
48 
49       \textcolor{comment}{! Values that will be used for the interpolation.}
50       u1 = u\_src(k\_src) \textcolor{comment}{! Value on left of source cell}
51       u2 = u\_src(k\_src+1) \textcolor{comment}{! Value on right of source cell}
52 \textcolor{keywordflow}{    enddo}
53 
54     \textcolor{keywordflow}{if} (dh>0.) \textcolor{keywordflow}{then}
55       weight\_a = max(0., ( dh - x\_dest ) / dh) \textcolor{comment}{! Weight of u1}
56       weight\_b = min(1., x\_dest / dh) \textcolor{comment}{! Weight of u2}
57       u\_dest(k\_dest) = weight\_a * u1 + weight\_b * u2 \textcolor{comment}{! Linear interpolation between u1 and u2}
58     \textcolor{keywordflow}{else}
59       u\_dest(k\_dest) = 0.5 * ( u1 + u2 ) \textcolor{comment}{! For a vanished layer we need to do something reasonable...}
60 \textcolor{keywordflow}{    endif}
61 
62     \textcolor{comment}{! Mask vanished layers at the surface which would be under an ice-shelf.}
63     \textcolor{comment}{! TODO: Need to figure out what to do for an isopycnal coordinate diagnostic that could}
64     \textcolor{comment}{!       also have vanished layers at the surface.}
65     \textcolor{keywordflow}{if} (k\_dest<=ndest) \textcolor{keywordflow}{then}
66       x\_dest = x\_dest + h\_dest(k\_dest) \textcolor{comment}{! Position of interface k\_dest+1}
67       \textcolor{keywordflow}{if} (still\_vanished .and. h\_dest(k\_dest)==0.) \textcolor{keywordflow}{then}
68         \textcolor{comment}{! When the layer k\_dest is vanished and all layers above are also vanished, the k\_dest}
69         \textcolor{comment}{! interface value should be missing.}
70         u\_dest(k\_dest) = missing\_value
71       \textcolor{keywordflow}{else}
72         still\_vanished = .false.
73 \textcolor{keywordflow}{      endif}
74 \textcolor{keywordflow}{    endif}
75 
76 \textcolor{keywordflow}{  enddo}
77 
78   \textcolor{comment}{! Mask vanished layers on topography}
79   still\_vanished = .true.
80   \textcolor{keywordflow}{do} k\_dest=ndest, 1, -1
81     \textcolor{keywordflow}{if} (still\_vanished .and. h\_dest(k\_dest)==0.) \textcolor{keywordflow}{then}
82       \textcolor{comment}{! When the layer k\_dest is vanished and all layers below are also vanished, the k\_dest+1}
83       \textcolor{comment}{! interface value should be missing.}
84       u\_dest(k\_dest+1) = missing\_value
85     \textcolor{keywordflow}{else}
86       \textcolor{keywordflow}{exit}
87 \textcolor{keywordflow}{    endif}
88 \textcolor{keywordflow}{  enddo}
89 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__diag__vkernels_a89369e4bf4b7795f6e384762d11f0d23}\label{namespacemom__diag__vkernels_a89369e4bf4b7795f6e384762d11f0d23}} 
\index{mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}!reintegrate\+\_\+column@{reintegrate\+\_\+column}}
\index{reintegrate\+\_\+column@{reintegrate\+\_\+column}!mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}}
\subsubsection{\texorpdfstring{reintegrate\+\_\+column()}{reintegrate\_column()}}
{\footnotesize\ttfamily subroutine, public mom\+\_\+diag\+\_\+vkernels\+::reintegrate\+\_\+column (\begin{DoxyParamCaption}\item[{integer, intent(in)}]{nsrc,  }\item[{real, dimension(nsrc), intent(in)}]{h\+\_\+src,  }\item[{real, dimension(nsrc), intent(in)}]{uh\+\_\+src,  }\item[{integer, intent(in)}]{ndest,  }\item[{real, dimension(ndest), intent(in)}]{h\+\_\+dest,  }\item[{real, intent(in)}]{missing\+\_\+value,  }\item[{real, dimension(ndest), intent(inout)}]{uh\+\_\+dest }\end{DoxyParamCaption})}



Conservatively calculate integrated data, uh\+\_\+dest, on grid h\+\_\+dest, from layer-\/integrated data, uh\+\_\+src, on grid h\+\_\+src. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em nsrc} & Number of source cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+src} & Thickness of source cells\\
\hline
\mbox{\tt in}  & {\em uh\+\_\+src} & Values at source cell interfaces\\
\hline
\mbox{\tt in}  & {\em ndest} & Number of destination cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+dest} & Thickness of destination cells\\
\hline
\mbox{\tt in}  & {\em missing\+\_\+value} & Value to assign in vanished cells\\
\hline
\mbox{\tt in,out}  & {\em uh\+\_\+dest} & Interpolated value at destination cell interfaces \\
\hline
\end{DoxyParams}


Definition at line 94 of file M\+O\+M\+\_\+diag\+\_\+vkernels.\+F90.


\begin{DoxyCode}
94   \textcolor{keywordtype}{integer},                \textcolor{keywordtype}{intent(in)}    :: nsrc\textcolor{comment}{ !< Number of source cells}
95   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc)},  \textcolor{keywordtype}{intent(in)}    :: h\_src\textcolor{comment}{ !< Thickness of source cells}
96   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc)},  \textcolor{keywordtype}{intent(in)}    :: uh\_src\textcolor{comment}{ !< Values at source cell interfaces}
97   \textcolor{keywordtype}{integer},                \textcolor{keywordtype}{intent(in)}    :: ndest\textcolor{comment}{ !< Number of destination cells}
98   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)}, \textcolor{keywordtype}{intent(in)}    :: h\_dest\textcolor{comment}{ !< Thickness of destination cells}
99   \textcolor{keywordtype}{real},                   \textcolor{keywordtype}{intent(in)}    :: missing\_value\textcolor{comment}{ !< Value to assign in vanished cells}
100   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)}, \textcolor{keywordtype}{intent(inout)} :: uh\_dest\textcolor{comment}{ !< Interpolated value at destination cell interfaces}
101   \textcolor{comment}{! Local variables}
102   \textcolor{keywordtype}{real} :: x\_dest \textcolor{comment}{! Relative position of target interface}
103   \textcolor{keywordtype}{real} :: h\_src\_rem, h\_dest\_rem, dh \textcolor{comment}{! Incremental thicknesses}
104   \textcolor{keywordtype}{real} :: uh\_src\_rem, uh\_dest\_rem, duh \textcolor{comment}{! Incremental amounts of stuff}
105   \textcolor{keywordtype}{integer} :: k\_src, k\_dest \textcolor{comment}{! Index of cell in src and dest columns}
106   \textcolor{keywordtype}{integer} :: iter
107   \textcolor{keywordtype}{logical} :: src\_ran\_out, src\_exists
108 
109   uh\_dest(:) = missing\_value
110 
111   k\_src = 0
112   k\_dest = 0
113   h\_dest\_rem = 0.
114   h\_src\_rem = 0.
115   src\_ran\_out = .false.
116   src\_exists = .false.
117 
118   \textcolor{keywordflow}{do} \textcolor{keywordflow}{while}(.true.)
119     \textcolor{keywordflow}{if} (h\_src\_rem==0. .and. k\_src<nsrc) \textcolor{keywordflow}{then}
120       \textcolor{comment}{! Supply is empty so move to the next source cell}
121       k\_src = k\_src + 1
122       h\_src\_rem = h\_src(k\_src)
123       uh\_src\_rem = uh\_src(k\_src)
124       \textcolor{keywordflow}{if} (h\_src\_rem==0.) cycle
125       src\_exists = .true. \textcolor{comment}{! This stops us masking out the entire column}
126 \textcolor{keywordflow}{    endif}
127     \textcolor{keywordflow}{if} (h\_dest\_rem==0. .and. k\_dest<ndest) \textcolor{keywordflow}{then}
128       \textcolor{comment}{! Sink has no capacity so move to the next destination cell}
129       k\_dest = k\_dest + 1
130       h\_dest\_rem = h\_dest(k\_dest)
131       uh\_dest(k\_dest) = 0.
132       \textcolor{keywordflow}{if} (h\_dest\_rem==0.) cycle
133 \textcolor{keywordflow}{    endif}
134     \textcolor{keywordflow}{if} (k\_src==nsrc .and. h\_src\_rem==0.) \textcolor{keywordflow}{then}
135       \textcolor{keywordflow}{if} (src\_ran\_out) \textcolor{keywordflow}{exit} \textcolor{comment}{! This is the second time implying there is no more src}
136       src\_ran\_out = .true.
137       cycle
138 \textcolor{keywordflow}{    endif}
139     duh = 0.
140     \textcolor{keywordflow}{if} (h\_src\_rem<h\_dest\_rem) \textcolor{keywordflow}{then}
141       \textcolor{comment}{! The source cell is fully within the destination cell}
142       dh = h\_src\_rem
143       \textcolor{keywordflow}{if} (dh>0.) duh = uh\_src\_rem
144       h\_src\_rem = 0.
145       uh\_src\_rem = 0.
146       h\_dest\_rem = max(0., h\_dest\_rem - dh)
147     \textcolor{keywordflow}{elseif} (h\_src\_rem>h\_dest\_rem) \textcolor{keywordflow}{then}
148       \textcolor{comment}{! Only part of the source cell can be used up}
149       dh = h\_dest\_rem
150       duh = (dh / h\_src\_rem) * uh\_src\_rem
151       h\_src\_rem = max(0., h\_src\_rem - dh)
152       uh\_src\_rem = uh\_src\_rem - duh
153       h\_dest\_rem = 0.
154     \textcolor{keywordflow}{else} \textcolor{comment}{! h\_src\_rem==h\_dest\_rem}
155       \textcolor{comment}{! The source cell exactly fits the destination cell}
156       duh = uh\_src\_rem
157       h\_src\_rem = 0.
158       uh\_src\_rem = 0.
159       h\_dest\_rem = 0.
160 \textcolor{keywordflow}{    endif}
161     uh\_dest(k\_dest) = uh\_dest(k\_dest) + duh
162     \textcolor{keywordflow}{if} (k\_dest==ndest .and. (k\_src==nsrc .or. h\_dest\_rem==0.)) \textcolor{keywordflow}{exit}
163 \textcolor{keywordflow}{  enddo}
164 
165   \textcolor{keywordflow}{if} (.not. src\_exists) uh\_dest(1:ndest) = missing\_value
166 
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__diag__vkernels_abeef457cda28b20c03a89ad402bcd434}\label{namespacemom__diag__vkernels_abeef457cda28b20c03a89ad402bcd434}} 
\index{mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}!test\+\_\+interp@{test\+\_\+interp}}
\index{test\+\_\+interp@{test\+\_\+interp}!mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}}
\subsubsection{\texorpdfstring{test\+\_\+interp()}{test\_interp()}}
{\footnotesize\ttfamily logical function mom\+\_\+diag\+\_\+vkernels\+::test\+\_\+interp (\begin{DoxyParamCaption}\item[{logical, intent(in)}]{verbose,  }\item[{real, intent(in)}]{missing\+\_\+value,  }\item[{character(len=$\ast$), intent(in)}]{msg,  }\item[{integer, intent(in)}]{nsrc,  }\item[{real, dimension(nsrc), intent(in)}]{h\+\_\+src,  }\item[{real, dimension(nsrc+1), intent(in)}]{u\+\_\+src,  }\item[{integer, intent(in)}]{ndest,  }\item[{real, dimension(ndest), intent(in)}]{h\+\_\+dest,  }\item[{real, dimension(ndest+1), intent(in)}]{u\+\_\+true }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



Returns true if a test of \mbox{\hyperlink{namespacemom__diag__vkernels_a4f42f472a725a147f8d97a68b2028c5b}{interpolate\+\_\+column()}} produces the wrong answer. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em verbose} & If true, write results to stdout\\
\hline
\mbox{\tt in}  & {\em missing\+\_\+value} & Value to indicate missing data\\
\hline
\mbox{\tt in}  & {\em msg} & Message to label test\\
\hline
\mbox{\tt in}  & {\em nsrc} & Number of source cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+src} & Thickness of source cells\\
\hline
\mbox{\tt in}  & {\em u\+\_\+src} & Values at source cell interfaces\\
\hline
\mbox{\tt in}  & {\em ndest} & Number of destination cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+dest} & Thickness of destination cells\\
\hline
\mbox{\tt in}  & {\em u\+\_\+true} & Correct value at destination cell interfaces \\
\hline
\end{DoxyParams}


Definition at line 284 of file M\+O\+M\+\_\+diag\+\_\+vkernels.\+F90.


\begin{DoxyCode}
284   \textcolor{keywordtype}{logical},                  \textcolor{keywordtype}{intent(in)} :: verbose\textcolor{comment}{ !< If true, write results to stdout}
285   \textcolor{keywordtype}{real},                     \textcolor{keywordtype}{intent(in)} :: missing\_value\textcolor{comment}{ !< Value to indicate missing data}
286   \textcolor{keywordtype}{character(len=*)},         \textcolor{keywordtype}{intent(in)} :: msg\textcolor{comment}{ !< Message to label test}
287   \textcolor{keywordtype}{integer},                  \textcolor{keywordtype}{intent(in)} :: nsrc\textcolor{comment}{ !< Number of source cells}
288   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc)},    \textcolor{keywordtype}{intent(in)} :: h\_src\textcolor{comment}{ !< Thickness of source cells}
289   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc+1)},  \textcolor{keywordtype}{intent(in)} :: u\_src\textcolor{comment}{ !< Values at source cell interfaces}
290   \textcolor{keywordtype}{integer},                  \textcolor{keywordtype}{intent(in)} :: ndest\textcolor{comment}{ !< Number of destination cells}
291   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)},   \textcolor{keywordtype}{intent(in)} :: h\_dest\textcolor{comment}{ !< Thickness of destination cells}
292   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest+1)}, \textcolor{keywordtype}{intent(in)} :: u\_true\textcolor{comment}{ !< Correct value at destination cell interfaces}
293   \textcolor{comment}{! Local variables}
294   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest+1)} :: u\_dest \textcolor{comment}{! Interpolated value at destination cell interfaces}
295   \textcolor{keywordtype}{integer} :: k
296   \textcolor{keywordtype}{real} :: error
297   \textcolor{keywordtype}{logical} :: print\_results
298 
299   \textcolor{comment}{! Interpolate from src to dest}
300   \textcolor{keyword}{call }interpolate\_column(nsrc, h\_src, u\_src, ndest, h\_dest, missing\_value, u\_dest)
301 
302   test\_interp = .false.
303   \textcolor{keywordflow}{do} k=1,ndest+1
304     \textcolor{keywordflow}{if} (u\_dest(k)/=u\_true(k)) test\_interp = .true.
305 \textcolor{keywordflow}{  enddo}
306   \textcolor{keywordflow}{if} (verbose .or. test\_interp) \textcolor{keywordflow}{then}
307     \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(2a)'}) \textcolor{stringliteral}{' Test: '},msg
308     \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(a3,3(a24))'}) \textcolor{stringliteral}{'k'},\textcolor{stringliteral}{'u\_result'},\textcolor{stringliteral}{'u\_true'},\textcolor{stringliteral}{'error'}
309     \textcolor{keywordflow}{do} k=1,ndest+1
310       error = u\_dest(k)-u\_true(k)
311       \textcolor{keywordflow}{if} (error==0.) \textcolor{keywordflow}{then}
312         \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(i3,3(1pe24.16))'}) k,u\_dest(k),u\_true(k),u\_dest(k)-u\_true(k)
313       \textcolor{keywordflow}{else}
314         \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(i3,3(1pe24.16),x,a)'}) k,u\_dest(k),u\_true(k),u\_dest(k)-u\_true(k),\textcolor{stringliteral}{'<--- WRONG!'}
315         \textcolor{keyword}{write}(stderr,\textcolor{stringliteral}{'(i3,3(1pe24.16),x,a)'}) k,u\_dest(k),u\_true(k),u\_dest(k)-u\_true(k),\textcolor{stringliteral}{'<--- WRONG!'}
316 \textcolor{keywordflow}{      endif}
317 \textcolor{keywordflow}{    enddo}
318 \textcolor{keywordflow}{  endif}
\end{DoxyCode}
\mbox{\Hypertarget{namespacemom__diag__vkernels_aac7aad5875fcf17273b211c90ac0cbee}\label{namespacemom__diag__vkernels_aac7aad5875fcf17273b211c90ac0cbee}} 
\index{mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}!test\+\_\+reintegrate@{test\+\_\+reintegrate}}
\index{test\+\_\+reintegrate@{test\+\_\+reintegrate}!mom\+\_\+diag\+\_\+vkernels@{mom\+\_\+diag\+\_\+vkernels}}
\subsubsection{\texorpdfstring{test\+\_\+reintegrate()}{test\_reintegrate()}}
{\footnotesize\ttfamily logical function mom\+\_\+diag\+\_\+vkernels\+::test\+\_\+reintegrate (\begin{DoxyParamCaption}\item[{logical, intent(in)}]{verbose,  }\item[{real, intent(in)}]{missing\+\_\+value,  }\item[{character(len=$\ast$), intent(in)}]{msg,  }\item[{integer, intent(in)}]{nsrc,  }\item[{real, dimension(nsrc), intent(in)}]{h\+\_\+src,  }\item[{real, dimension(nsrc), intent(in)}]{uh\+\_\+src,  }\item[{integer, intent(in)}]{ndest,  }\item[{real, dimension(ndest), intent(in)}]{h\+\_\+dest,  }\item[{real, dimension(ndest), intent(in)}]{uh\+\_\+true }\end{DoxyParamCaption})\hspace{0.3cm}{\ttfamily [private]}}



Returns true if a test of \mbox{\hyperlink{namespacemom__diag__vkernels_a89369e4bf4b7795f6e384762d11f0d23}{reintegrate\+\_\+column()}} produces the wrong answer. 


\begin{DoxyParams}[1]{Parameters}
\mbox{\tt in}  & {\em verbose} & If true, write results to stdout\\
\hline
\mbox{\tt in}  & {\em missing\+\_\+value} & Value to indicate missing data\\
\hline
\mbox{\tt in}  & {\em msg} & Message to label test\\
\hline
\mbox{\tt in}  & {\em nsrc} & Number of source cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+src} & Thickness of source cells\\
\hline
\mbox{\tt in}  & {\em uh\+\_\+src} & Values of source cell stuff\\
\hline
\mbox{\tt in}  & {\em ndest} & Number of destination cells\\
\hline
\mbox{\tt in}  & {\em h\+\_\+dest} & Thickness of destination cells\\
\hline
\mbox{\tt in}  & {\em uh\+\_\+true} & Correct value of destination cell stuff \\
\hline
\end{DoxyParams}


Definition at line 323 of file M\+O\+M\+\_\+diag\+\_\+vkernels.\+F90.


\begin{DoxyCode}
323   \textcolor{keywordtype}{logical},                \textcolor{keywordtype}{intent(in)} :: verbose\textcolor{comment}{ !< If true, write results to stdout}
324   \textcolor{keywordtype}{real},                   \textcolor{keywordtype}{intent(in)} :: missing\_value\textcolor{comment}{ !< Value to indicate missing data}
325   \textcolor{keywordtype}{character(len=*)},       \textcolor{keywordtype}{intent(in)} :: msg\textcolor{comment}{ !< Message to label test}
326   \textcolor{keywordtype}{integer},                \textcolor{keywordtype}{intent(in)} :: nsrc\textcolor{comment}{ !< Number of source cells}
327   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc)},  \textcolor{keywordtype}{intent(in)} :: h\_src\textcolor{comment}{ !< Thickness of source cells}
328   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(nsrc)},  \textcolor{keywordtype}{intent(in)} :: uh\_src\textcolor{comment}{ !< Values of source cell stuff}
329   \textcolor{keywordtype}{integer},                \textcolor{keywordtype}{intent(in)} :: ndest\textcolor{comment}{ !< Number of destination cells}
330   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)}, \textcolor{keywordtype}{intent(in)} :: h\_dest\textcolor{comment}{ !< Thickness of destination cells}
331   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)}, \textcolor{keywordtype}{intent(in)} :: uh\_true\textcolor{comment}{ !< Correct value of destination cell stuff}
332   \textcolor{comment}{! Local variables}
333   \textcolor{keywordtype}{real}, \textcolor{keywordtype}{dimension(ndest)} :: uh\_dest \textcolor{comment}{! Reintegrated value on destination cells}
334   \textcolor{keywordtype}{integer} :: k
335   \textcolor{keywordtype}{real} :: error
336   \textcolor{keywordtype}{logical} :: print\_results
337 
338   \textcolor{comment}{! Interpolate from src to dest}
339   \textcolor{keyword}{call }reintegrate\_column(nsrc, h\_src, uh\_src, ndest, h\_dest, missing\_value, uh\_dest)
340 
341   test\_reintegrate = .false.
342   \textcolor{keywordflow}{do} k=1,ndest
343     \textcolor{keywordflow}{if} (uh\_dest(k)/=uh\_true(k)) test\_reintegrate = .true.
344 \textcolor{keywordflow}{  enddo}
345   \textcolor{keywordflow}{if} (verbose .or. test\_reintegrate) \textcolor{keywordflow}{then}
346     \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(2a)'}) \textcolor{stringliteral}{' Test: '},msg
347     \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(a3,3(a24))'}) \textcolor{stringliteral}{'k'},\textcolor{stringliteral}{'uh\_result'},\textcolor{stringliteral}{'uh\_true'},\textcolor{stringliteral}{'error'}
348     \textcolor{keywordflow}{do} k=1,ndest
349       error = uh\_dest(k)-uh\_true(k)
350       \textcolor{keywordflow}{if} (error==0.) \textcolor{keywordflow}{then}
351         \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(i3,3(1pe24.16))'}) k,uh\_dest(k),uh\_true(k),uh\_dest(k)-uh\_true(k)
352       \textcolor{keywordflow}{else}
353         \textcolor{keyword}{write}(stdout,\textcolor{stringliteral}{'(i3,3(1pe24.16),x,a)'}) k,uh\_dest(k),uh\_true(k),uh\_dest(k)-uh\_true(k),\textcolor{stringliteral}{'<--- WRONG!'}
354         \textcolor{keyword}{write}(stderr,\textcolor{stringliteral}{'(i3,3(1pe24.16),x,a)'}) k,uh\_dest(k),uh\_true(k),uh\_dest(k)-uh\_true(k),\textcolor{stringliteral}{'<--- WRONG!'}
355 \textcolor{keywordflow}{      endif}
356 \textcolor{keywordflow}{    enddo}
357 \textcolor{keywordflow}{  endif}
\end{DoxyCode}
