1 |
jmc |
1.25 |
% $Header: /u/gcmpack/manual/s_algorithm/text/spatial-discrete.tex,v 1.24 2010/08/30 23:09:19 jmc Exp $ |
2 |
adcroft |
1.2 |
% $Name: $ |
3 |
adcroft |
1.1 |
|
4 |
|
|
\section{Spatial discretization of the dynamical equations} |
5 |
edhill |
1.17 |
\begin{rawhtml} |
6 |
|
|
<!-- CMIREDIR:spatial_discretization_of_dyn_eq: --> |
7 |
|
|
\end{rawhtml} |
8 |
adcroft |
1.1 |
|
9 |
adcroft |
1.2 |
Spatial discretization is carried out using the finite volume |
10 |
|
|
method. This amounts to a grid-point method (namely second-order |
11 |
|
|
centered finite difference) in the fluid interior but allows |
12 |
|
|
boundaries to intersect a regular grid allowing a more accurate |
13 |
|
|
representation of the position of the boundary. We treat the |
14 |
cnh |
1.10 |
horizontal and vertical directions as separable and differently. |
15 |
adcroft |
1.2 |
|
16 |
|
|
|
17 |
|
|
\subsection{The finite volume method: finite volumes versus finite difference} |
18 |
afe |
1.13 |
\begin{rawhtml} |
19 |
afe |
1.14 |
<!-- CMIREDIR:finite_volume: --> |
20 |
afe |
1.13 |
\end{rawhtml} |
21 |
|
|
|
22 |
|
|
|
23 |
adcroft |
1.2 |
|
24 |
|
|
The finite volume method is used to discretize the equations in |
25 |
|
|
space. The expression ``finite volume'' actually has two meanings; one |
26 |
adcroft |
1.8 |
is the method of embedded or intersecting boundaries (shaved or lopped |
27 |
|
|
cells in our terminology) and the other is non-linear interpolation |
28 |
|
|
methods that can deal with non-smooth solutions such as shocks |
29 |
|
|
(i.e. flux limiters for advection). Both make use of the integral form |
30 |
|
|
of the conservation laws to which the {\it weak solution} is a |
31 |
|
|
solution on each finite volume of (sub-domain). The weak solution can |
32 |
|
|
be constructed out of piece-wise constant elements or be |
33 |
adcroft |
1.3 |
differentiable. The differentiable equations can not be satisfied by |
34 |
|
|
piece-wise constant functions. |
35 |
|
|
|
36 |
|
|
As an example, the 1-D constant coefficient advection-diffusion |
37 |
|
|
equation: |
38 |
adcroft |
1.2 |
\begin{displaymath} |
39 |
|
|
\partial_t \theta + \partial_x ( u \theta - \kappa \partial_x \theta ) = 0 |
40 |
|
|
\end{displaymath} |
41 |
adcroft |
1.3 |
can be discretized by integrating over finite sub-domains, i.e. |
42 |
|
|
the lengths $\Delta x_i$: |
43 |
adcroft |
1.2 |
\begin{displaymath} |
44 |
|
|
\Delta x \partial_t \theta + \delta_i ( F ) = 0 |
45 |
|
|
\end{displaymath} |
46 |
cnh |
1.10 |
is exact if $\theta(x)$ is piece-wise constant over the interval |
47 |
adcroft |
1.3 |
$\Delta x_i$ or more generally if $\theta_i$ is defined as the average |
48 |
|
|
over the interval $\Delta x_i$. |
49 |
|
|
|
50 |
|
|
The flux, $F_{i-1/2}$, must be approximated: |
51 |
adcroft |
1.2 |
\begin{displaymath} |
52 |
|
|
F = u \overline{\theta} - \frac{\kappa}{\Delta x_c} \partial_i \theta |
53 |
|
|
\end{displaymath} |
54 |
adcroft |
1.3 |
and this is where truncation errors can enter the solution. The |
55 |
|
|
method for obtaining $\overline{\theta}$ is unspecified and a wide |
56 |
|
|
range of possibilities exist including centered and upwind |
57 |
|
|
interpolation, polynomial fits based on the the volume average |
58 |
|
|
definitions of quantities and non-linear interpolation such as |
59 |
|
|
flux-limiters. |
60 |
|
|
|
61 |
|
|
Choosing simple centered second-order interpolation and differencing |
62 |
|
|
recovers the same ODE's resulting from finite differencing for the |
63 |
|
|
interior of a fluid. Differences arise at boundaries where a boundary |
64 |
|
|
is not positioned on a regular or smoothly varying grid. This method |
65 |
|
|
is used to represent the topography using lopped cell, see |
66 |
adcroft |
1.11 |
\cite{adcroft:97}. Subtle difference also appear in more than one |
67 |
adcroft |
1.3 |
dimension away from boundaries. This happens because the each |
68 |
cnh |
1.10 |
direction is discretized independently in the finite difference method |
69 |
adcroft |
1.3 |
while the integrating over finite volume implicitly treats all |
70 |
|
|
directions simultaneously. Illustration of this is given in |
71 |
adcroft |
1.11 |
\cite{ac:02}. |
72 |
adcroft |
1.2 |
|
73 |
adcroft |
1.1 |
\subsection{C grid staggering of variables} |
74 |
|
|
|
75 |
|
|
\begin{figure} |
76 |
adcroft |
1.7 |
\begin{center} |
77 |
jmc |
1.23 |
\resizebox{!}{2in}{ \includegraphics{s_algorithm/figs/cgrid3d.eps}} |
78 |
adcroft |
1.7 |
\end{center} |
79 |
adcroft |
1.1 |
\caption{Three dimensional staggering of velocity components. This |
80 |
|
|
facilitates the natural discretization of the continuity and tracer |
81 |
|
|
equations. } |
82 |
adcroft |
1.2 |
\label{fig:cgrid3d} |
83 |
adcroft |
1.1 |
\end{figure} |
84 |
|
|
|
85 |
adcroft |
1.2 |
The basic algorithm employed for stepping forward the momentum |
86 |
|
|
equations is based on retaining non-divergence of the flow at all |
87 |
|
|
times. This is most naturally done if the components of flow are |
88 |
adcroft |
1.11 |
staggered in space in the form of an Arakawa C grid \cite{arakawa:77}. |
89 |
adcroft |
1.2 |
|
90 |
|
|
Fig. \ref{fig:cgrid3d} shows the components of flow ($u$,$v$,$w$) |
91 |
|
|
staggered in space such that the zonal component falls on the |
92 |
cnh |
1.10 |
interface between continuity cells in the zonal direction. Similarly |
93 |
|
|
for the meridional and vertical directions. The continuity cell is |
94 |
adcroft |
1.2 |
synonymous with tracer cells (they are one and the same). |
95 |
|
|
|
96 |
|
|
|
97 |
|
|
|
98 |
adcroft |
1.5 |
\subsection{Grid initialization and data} |
99 |
|
|
|
100 |
|
|
Initialization of grid data is controlled by subroutine {\em |
101 |
|
|
INI\_GRID} which in calls {\em INI\_VERTICAL\_GRID} to initialize the |
102 |
|
|
vertical grid, and then either of {\em INI\_CARTESIAN\_GRID}, {\em |
103 |
|
|
INI\_SPHERICAL\_POLAR\_GRID} or {\em INI\_CURV\-ILINEAR\_GRID} to |
104 |
|
|
initialize the horizontal grid for cartesian, spherical-polar or |
105 |
|
|
curvilinear coordinates respectively. |
106 |
|
|
|
107 |
|
|
The reciprocals of all grid quantities are pre-calculated and this is |
108 |
|
|
done in subroutine {\em INI\_MASKS\_ETC} which is called later by |
109 |
|
|
subroutine {\em INITIALIZE\_FIXED}. |
110 |
|
|
|
111 |
|
|
All grid descriptors are global arrays and stored in common blocks in |
112 |
|
|
{\em GRID.h} and a generally declared as {\em \_RS}. |
113 |
|
|
|
114 |
|
|
\fbox{ \begin{minipage}{4.75in} |
115 |
|
|
{\em S/R INI\_GRID} ({\em model/src/ini\_grid.F}) |
116 |
|
|
|
117 |
|
|
{\em S/R INI\_MASKS\_ETC} ({\em model/src/ini\_masks\_etc.F}) |
118 |
|
|
|
119 |
|
|
grid data: ({\em model/inc/GRID.h}) |
120 |
|
|
\end{minipage} } |
121 |
|
|
|
122 |
adcroft |
1.2 |
|
123 |
adcroft |
1.1 |
\subsection{Horizontal grid} |
124 |
cnh |
1.9 |
\label{sec:spatial_discrete_horizontal_grid} |
125 |
adcroft |
1.1 |
|
126 |
|
|
\begin{figure} |
127 |
adcroft |
1.7 |
\begin{center} |
128 |
|
|
\begin{tabular}{cc} |
129 |
jmc |
1.23 |
\raisebox{1.5in}{a)}\resizebox{!}{2in}{ \includegraphics{s_algorithm/figs/hgrid-Ac.eps}} |
130 |
|
|
& \raisebox{1.5in}{b)}\resizebox{!}{2in}{ \includegraphics{s_algorithm/figs/hgrid-Az.eps}} |
131 |
adcroft |
1.1 |
\\ |
132 |
jmc |
1.23 |
\raisebox{1.5in}{c)}\resizebox{!}{2in}{ \includegraphics{s_algorithm/figs/hgrid-Au.eps}} |
133 |
|
|
& \raisebox{1.5in}{d)}\resizebox{!}{2in}{ \includegraphics{s_algorithm/figs/hgrid-Av.eps}} |
134 |
adcroft |
1.7 |
\end{tabular} |
135 |
|
|
\end{center} |
136 |
adcroft |
1.2 |
\caption{ |
137 |
|
|
Staggering of horizontal grid descriptors (lengths and areas). The |
138 |
|
|
grid lines indicate the tracer cell boundaries and are the reference |
139 |
|
|
grid for all panels. a) The area of a tracer cell, $A_c$, is bordered |
140 |
|
|
by the lengths $\Delta x_g$ and $\Delta y_g$. b) The area of a |
141 |
|
|
vorticity cell, $A_\zeta$, is bordered by the lengths $\Delta x_c$ and |
142 |
jmc |
1.21 |
$\Delta y_c$. c) The area of a u cell, $A_w$, is bordered by the |
143 |
|
|
lengths $\Delta x_v$ and $\Delta y_f$. d) The area of a v cell, $A_s$, |
144 |
adcroft |
1.2 |
is bordered by the lengths $\Delta x_f$ and $\Delta y_u$.} |
145 |
|
|
\label{fig:hgrid} |
146 |
adcroft |
1.1 |
\end{figure} |
147 |
|
|
|
148 |
adcroft |
1.2 |
The model domain is decomposed into tiles and within each tile a |
149 |
|
|
quasi-regular grid is used. A tile is the basic unit of domain |
150 |
cnh |
1.10 |
decomposition for parallelization but may be used whether parallelized |
151 |
jmc |
1.24 |
or not; see section \ref{sec:domain_decomposition} for more details. |
152 |
jmc |
1.19 |
Although the tiles may be patched together in an unstructured manner |
153 |
adcroft |
1.2 |
(i.e. irregular or non-tessilating pattern), the interior of tiles is |
154 |
cnh |
1.10 |
a structured grid of quadrilateral cells. The horizontal coordinate |
155 |
adcroft |
1.2 |
system is orthogonal curvilinear meaning we can not necessarily treat |
156 |
cnh |
1.10 |
the two horizontal directions as separable. Instead, each cell in the |
157 |
adcroft |
1.2 |
horizontal grid is described by the length of it's sides and it's |
158 |
|
|
area. |
159 |
|
|
|
160 |
|
|
The grid information is quite general and describes any of the |
161 |
|
|
available coordinates systems, cartesian, spherical-polar or |
162 |
|
|
curvilinear. All that is necessary to distinguish between the |
163 |
cnh |
1.10 |
coordinate systems is to initialize the grid data (descriptors) |
164 |
adcroft |
1.2 |
appropriately. |
165 |
|
|
|
166 |
|
|
In the following, we refer to the orientation of quantities on the |
167 |
|
|
computational grid using geographic terminology such as points of the |
168 |
adcroft |
1.3 |
compass. |
169 |
|
|
\marginpar{Caution!} |
170 |
jmc |
1.25 |
This is purely for convenience but should not be confused |
171 |
adcroft |
1.2 |
with the actual geographic orientation of model quantities. |
172 |
|
|
|
173 |
adcroft |
1.3 |
Fig.~\ref{fig:hgrid}a shows the tracer cell (synonymous with the |
174 |
|
|
continuity cell). The length of the southern edge, $\Delta x_g$, |
175 |
|
|
western edge, $\Delta y_g$ and surface area, $A_c$, presented in the |
176 |
|
|
vertical are stored in arrays {\bf DXg}, {\bf DYg} and {\bf rAc}. |
177 |
adcroft |
1.2 |
\marginpar{$A_c$: {\bf rAc}} |
178 |
|
|
\marginpar{$\Delta x_g$: {\bf DXg}} |
179 |
|
|
\marginpar{$\Delta y_g$: {\bf DYg}} |
180 |
adcroft |
1.3 |
The ``g'' suffix indicates that the lengths are along the defining |
181 |
|
|
grid boundaries. The ``c'' suffix associates the quantity with the |
182 |
|
|
cell centers. The quantities are staggered in space and the indexing |
183 |
|
|
is such that {\bf DXg(i,j)} is positioned to the south of {\bf |
184 |
|
|
rAc(i,j)} and {\bf DYg(i,j)} positioned to the west. |
185 |
adcroft |
1.2 |
|
186 |
adcroft |
1.3 |
Fig.~\ref{fig:hgrid}b shows the vorticity cell. The length of the |
187 |
|
|
southern edge, $\Delta x_c$, western edge, $\Delta y_c$ and surface |
188 |
|
|
area, $A_\zeta$, presented in the vertical are stored in arrays {\bf |
189 |
heimbach |
1.22 |
DXc}, {\bf DYc} and {\bf rAz}. |
190 |
adcroft |
1.2 |
\marginpar{$A_\zeta$: {\bf rAz}} |
191 |
|
|
\marginpar{$\Delta x_c$: {\bf DXc}} |
192 |
|
|
\marginpar{$\Delta y_c$: {\bf DYc}} |
193 |
adcroft |
1.3 |
The ``z'' suffix indicates that the lengths are measured between the |
194 |
|
|
cell centers and the ``$\zeta$'' suffix associates points with the |
195 |
|
|
vorticity points. The quantities are staggered in space and the |
196 |
|
|
indexing is such that {\bf DXc(i,j)} is positioned to the north of |
197 |
jmc |
1.25 |
{\bf rAz(i,j)} and {\bf DYc(i,j)} positioned to the east. |
198 |
adcroft |
1.2 |
|
199 |
adcroft |
1.3 |
Fig.~\ref{fig:hgrid}c shows the ``u'' or western (w) cell. The length of |
200 |
|
|
the southern edge, $\Delta x_v$, eastern edge, $\Delta y_f$ and |
201 |
|
|
surface area, $A_w$, presented in the vertical are stored in arrays |
202 |
|
|
{\bf DXv}, {\bf DYf} and {\bf rAw}. |
203 |
adcroft |
1.2 |
\marginpar{$A_w$: {\bf rAw}} |
204 |
|
|
\marginpar{$\Delta x_v$: {\bf DXv}} |
205 |
|
|
\marginpar{$\Delta y_f$: {\bf DYf}} |
206 |
adcroft |
1.3 |
The ``v'' suffix indicates that the length is measured between the |
207 |
|
|
v-points, the ``f'' suffix indicates that the length is measured |
208 |
|
|
between the (tracer) cell faces and the ``w'' suffix associates points |
209 |
|
|
with the u-points (w stands for west). The quantities are staggered in |
210 |
|
|
space and the indexing is such that {\bf DXv(i,j)} is positioned to |
211 |
|
|
the south of {\bf rAw(i,j)} and {\bf DYf(i,j)} positioned to the east. |
212 |
adcroft |
1.2 |
|
213 |
adcroft |
1.3 |
Fig.~\ref{fig:hgrid}d shows the ``v'' or southern (s) cell. The length of |
214 |
|
|
the northern edge, $\Delta x_f$, western edge, $\Delta y_u$ and |
215 |
|
|
surface area, $A_s$, presented in the vertical are stored in arrays |
216 |
|
|
{\bf DXf}, {\bf DYu} and {\bf rAs}. |
217 |
adcroft |
1.2 |
\marginpar{$A_s$: {\bf rAs}} |
218 |
|
|
\marginpar{$\Delta x_f$: {\bf DXf}} |
219 |
|
|
\marginpar{$\Delta y_u$: {\bf DYu}} |
220 |
adcroft |
1.3 |
The ``u'' suffix indicates that the length is measured between the |
221 |
|
|
u-points, the ``f'' suffix indicates that the length is measured |
222 |
|
|
between the (tracer) cell faces and the ``s'' suffix associates points |
223 |
|
|
with the v-points (s stands for south). The quantities are staggered |
224 |
|
|
in space and the indexing is such that {\bf DXf(i,j)} is positioned to |
225 |
|
|
the north of {\bf rAs(i,j)} and {\bf DYu(i,j)} positioned to the west. |
226 |
adcroft |
1.2 |
|
227 |
|
|
\fbox{ \begin{minipage}{4.75in} |
228 |
|
|
{\em S/R INI\_CARTESIAN\_GRID} ({\em |
229 |
|
|
model/src/ini\_cartesian\_grid.F}) |
230 |
|
|
|
231 |
|
|
{\em S/R INI\_SPHERICAL\_POLAR\_GRID} ({\em |
232 |
|
|
model/src/ini\_spherical\_polar\_grid.F}) |
233 |
|
|
|
234 |
|
|
{\em S/R INI\_CURVILINEAR\_GRID} ({\em |
235 |
|
|
model/src/ini\_curvilinear\_grid.F}) |
236 |
|
|
|
237 |
|
|
$A_c$, $A_\zeta$, $A_w$, $A_s$: {\bf rAc}, {\bf rAz}, {\bf rAw}, {\bf rAs} |
238 |
|
|
({\em GRID.h}) |
239 |
|
|
|
240 |
|
|
$\Delta x_g$, $\Delta y_g$: {\bf DXg}, {\bf DYg} ({\em GRID.h}) |
241 |
|
|
|
242 |
|
|
$\Delta x_c$, $\Delta y_c$: {\bf DXc}, {\bf DYc} ({\em GRID.h}) |
243 |
|
|
|
244 |
|
|
$\Delta x_f$, $\Delta y_f$: {\bf DXf}, {\bf DYf} ({\em GRID.h}) |
245 |
|
|
|
246 |
|
|
$\Delta x_v$, $\Delta y_u$: {\bf DXv}, {\bf DYu} ({\em GRID.h}) |
247 |
|
|
|
248 |
|
|
\end{minipage} } |
249 |
|
|
|
250 |
|
|
\subsubsection{Reciprocals of horizontal grid descriptors} |
251 |
|
|
|
252 |
|
|
%\marginpar{$A_c^{-1}$: {\bf RECIP\_rAc}} |
253 |
|
|
%\marginpar{$A_\zeta^{-1}$: {\bf RECIP\_rAz}} |
254 |
|
|
%\marginpar{$A_w^{-1}$: {\bf RECIP\_rAw}} |
255 |
|
|
%\marginpar{$A_s^{-1}$: {\bf RECIP\_rAs}} |
256 |
|
|
Lengths and areas appear in the denominator of expressions as much as |
257 |
|
|
in the numerator. For efficiency and portability, we pre-calculate the |
258 |
|
|
reciprocal of the horizontal grid quantities so that in-line divisions |
259 |
|
|
can be avoided. |
260 |
|
|
|
261 |
|
|
%\marginpar{$\Delta x_g^{-1}$: {\bf RECIP\_DXg}} |
262 |
|
|
%\marginpar{$\Delta y_g^{-1}$: {\bf RECIP\_DYg}} |
263 |
|
|
%\marginpar{$\Delta x_c^{-1}$: {\bf RECIP\_DXc}} |
264 |
|
|
%\marginpar{$\Delta y_c^{-1}$: {\bf RECIP\_DYc}} |
265 |
|
|
%\marginpar{$\Delta x_f^{-1}$: {\bf RECIP\_DXf}} |
266 |
|
|
%\marginpar{$\Delta y_f^{-1}$: {\bf RECIP\_DYf}} |
267 |
|
|
%\marginpar{$\Delta x_v^{-1}$: {\bf RECIP\_DXv}} |
268 |
|
|
%\marginpar{$\Delta y_u^{-1}$: {\bf RECIP\_DYu}} |
269 |
|
|
For each grid descriptor (array) there is a reciprocal named using the |
270 |
|
|
prefix {\bf RECIP\_}. This doubles the amount of storage in {\em |
271 |
|
|
GRID.h} but they are all only 2-D descriptors. |
272 |
|
|
|
273 |
|
|
\fbox{ \begin{minipage}{4.75in} |
274 |
|
|
{\em S/R INI\_MASKS\_ETC} ({\em |
275 |
|
|
model/src/ini\_masks\_etc.F}) |
276 |
|
|
|
277 |
|
|
$A_c^{-1}$: {\bf RECIP\_Ac} ({\em GRID.h}) |
278 |
|
|
|
279 |
|
|
$A_\zeta^{-1}$: {\bf RECIP\_Az} ({\em GRID.h}) |
280 |
|
|
|
281 |
|
|
$A_w^{-1}$: {\bf RECIP\_Aw} ({\em GRID.h}) |
282 |
|
|
|
283 |
|
|
$A_s^{-1}$: {\bf RECIP\_As} ({\em GRID.h}) |
284 |
|
|
|
285 |
|
|
$\Delta x_g^{-1}$, $\Delta y_g^{-1}$: {\bf RECIP\_DXg}, {\bf RECIP\_DYg} ({\em GRID.h}) |
286 |
|
|
|
287 |
|
|
$\Delta x_c^{-1}$, $\Delta y_c^{-1}$: {\bf RECIP\_DXc}, {\bf RECIP\_DYc} ({\em GRID.h}) |
288 |
|
|
|
289 |
|
|
$\Delta x_f^{-1}$, $\Delta y_f^{-1}$: {\bf RECIP\_DXf}, {\bf RECIP\_DYf} ({\em GRID.h}) |
290 |
|
|
|
291 |
|
|
$\Delta x_v^{-1}$, $\Delta y_u^{-1}$: {\bf RECIP\_DXv}, {\bf RECIP\_DYu} ({\em GRID.h}) |
292 |
|
|
|
293 |
|
|
\end{minipage} } |
294 |
|
|
|
295 |
|
|
\subsubsection{Cartesian coordinates} |
296 |
|
|
|
297 |
|
|
Cartesian coordinates are selected when the logical flag {\bf |
298 |
|
|
using\-Cartes\-ianGrid} in namelist {\em PARM04} is set to true. The grid |
299 |
|
|
spacing can be set to uniform via scalars {\bf dXspacing} and {\bf |
300 |
|
|
dYspacing} in namelist {\em PARM04} or to variable resolution by the |
301 |
|
|
vectors {\bf DELX} and {\bf DELY}. Units are normally |
302 |
cnh |
1.10 |
meters. Non-dimensional coordinates can be used by interpreting the |
303 |
adcroft |
1.2 |
gravitational constant as the Rayleigh number. |
304 |
|
|
|
305 |
|
|
\subsubsection{Spherical-polar coordinates} |
306 |
|
|
|
307 |
|
|
Spherical coordinates are selected when the logical flag {\bf |
308 |
|
|
using\-Spherical\-PolarGrid} in namelist {\em PARM04} is set to true. The |
309 |
|
|
grid spacing can be set to uniform via scalars {\bf dXspacing} and |
310 |
|
|
{\bf dYspacing} in namelist {\em PARM04} or to variable resolution by |
311 |
|
|
the vectors {\bf DELX} and {\bf DELY}. Units of these namelist |
312 |
|
|
variables are alway degrees. The horizontal grid descriptors are |
313 |
|
|
calculated from these namelist variables have units of meters. |
314 |
|
|
|
315 |
|
|
\subsubsection{Curvilinear coordinates} |
316 |
|
|
|
317 |
|
|
Curvilinear coordinates are selected when the logical flag {\bf |
318 |
|
|
using\-Curvil\-inear\-Grid} in namelist {\em PARM04} is set to true. The |
319 |
|
|
grid spacing can not be set via the namelist. Instead, the grid |
320 |
|
|
descriptors are read from data files, one for each descriptor. As for |
321 |
|
|
other grids, the horizontal grid descriptors have units of meters. |
322 |
|
|
|
323 |
|
|
|
324 |
adcroft |
1.1 |
\subsection{Vertical grid} |
325 |
|
|
|
326 |
|
|
\begin{figure} |
327 |
adcroft |
1.7 |
\begin{center} |
328 |
|
|
\begin{tabular}{cc} |
329 |
adcroft |
1.2 |
\raisebox{4in}{a)} \resizebox{!}{4in}{ |
330 |
jmc |
1.23 |
\includegraphics{s_algorithm/figs/vgrid-cellcentered.eps}} & \raisebox{4in}{b)} |
331 |
|
|
\resizebox{!}{4in}{ \includegraphics{s_algorithm/figs/vgrid-accurate.eps}} |
332 |
adcroft |
1.7 |
\end{tabular} |
333 |
|
|
\end{center} |
334 |
adcroft |
1.1 |
\caption{Two versions of the vertical grid. a) The cell centered |
335 |
|
|
approach where the interface depths are specified and the tracer |
336 |
|
|
points centered in between the interfaces. b) The interface centered |
337 |
|
|
approach where tracer levels are specified and the w-interfaces are |
338 |
|
|
centered in between.} |
339 |
adcroft |
1.2 |
\label{fig:vgrid} |
340 |
adcroft |
1.1 |
\end{figure} |
341 |
|
|
|
342 |
adcroft |
1.2 |
As for the horizontal grid, we use the suffixes ``c'' and ``f'' to |
343 |
|
|
indicates faces and centers. Fig.~\ref{fig:vgrid}a shows the default |
344 |
|
|
vertical grid used by the model. |
345 |
|
|
\marginpar{$\Delta r_f$: {\bf DRf}} |
346 |
|
|
\marginpar{$\Delta r_c$: {\bf DRc}} |
347 |
|
|
$\Delta r_f$ is the difference in $r$ |
348 |
|
|
(vertical coordinate) between the faces (i.e. $\Delta r_f \equiv - |
349 |
|
|
\delta_k r$ where the minus sign appears due to the convention that the |
350 |
|
|
surface layer has index $k=1$.). |
351 |
|
|
|
352 |
|
|
The vertical grid is calculated in subroutine {\em |
353 |
|
|
INI\_VERTICAL\_GRID} and specified via the vector {\bf DELR} in |
354 |
|
|
namelist {\em PARM04}. The units of ``r'' are either meters or Pascals |
355 |
|
|
depending on the isomorphism being used which in turn is dependent |
356 |
cnh |
1.10 |
only on the choice of equation of state. |
357 |
adcroft |
1.2 |
|
358 |
|
|
There are alternative namelist vectors {\bf DELZ} and {\bf DELP} which |
359 |
|
|
dictate whether z- or |
360 |
|
|
\marginpar{Caution!} |
361 |
|
|
p- coordinates are to be used but we intend to |
362 |
|
|
phase this out since they are redundant. |
363 |
|
|
|
364 |
|
|
The reciprocals $\Delta r_f^{-1}$ and $\Delta r_c^{-1}$ are |
365 |
|
|
pre-calculated (also in subroutine {\em INI\_VERTICAL\_GRID}). All |
366 |
|
|
vertical grid descriptors are stored in common blocks in {\em GRID.h}. |
367 |
|
|
|
368 |
|
|
The above grid (Fig.~\ref{fig:vgrid}a) is known as the cell centered |
369 |
|
|
approach because the tracer points are at cell centers; the cell |
370 |
jmc |
1.16 |
centers are mid-way between the cell interfaces. |
371 |
jmc |
1.18 |
This discretization is selected when the thickness of the |
372 |
jmc |
1.16 |
levels are provided ({\bf delR}, parameter file {\em data}, |
373 |
|
|
namelist {\em PARM04}) |
374 |
|
|
An alternative, the vertex or interface centered approach, is shown in |
375 |
adcroft |
1.2 |
Fig.~\ref{fig:vgrid}b. Here, the interior interfaces are positioned |
376 |
|
|
mid-way between the tracer nodes (no longer cell centers). This |
377 |
|
|
approach is formally more accurate for evaluation of hydrostatic |
378 |
|
|
pressure and vertical advection but historically the cell centered |
379 |
|
|
approach has been used. An alternative form of subroutine {\em |
380 |
|
|
INI\_VERTICAL\_GRID} is used to select the interface centered approach |
381 |
jmc |
1.16 |
This form requires to specify $Nr+1$ vertical distances {\bf delRc} |
382 |
|
|
(parameter file {\em data}, namelist {\em PARM04}, e.g. |
383 |
|
|
{\em verification/ideal\_2D\_oce/input/data}) |
384 |
|
|
corresponding to surface to center, $Nr-1$ center to center, and center to |
385 |
|
|
bottom distances. |
386 |
|
|
%but no run time option is currently available. |
387 |
adcroft |
1.2 |
|
388 |
|
|
\fbox{ \begin{minipage}{4.75in} |
389 |
|
|
{\em S/R INI\_VERTICAL\_GRID} ({\em |
390 |
|
|
model/src/ini\_vertical\_grid.F}) |
391 |
|
|
|
392 |
|
|
$\Delta r_f$: {\bf DRf} ({\em GRID.h}) |
393 |
|
|
|
394 |
|
|
$\Delta r_c$: {\bf DRc} ({\em GRID.h}) |
395 |
|
|
|
396 |
|
|
$\Delta r_f^{-1}$: {\bf RECIP\_DRf} ({\em GRID.h}) |
397 |
|
|
|
398 |
|
|
$\Delta r_c^{-1}$: {\bf RECIP\_DRc} ({\em GRID.h}) |
399 |
adcroft |
1.3 |
|
400 |
|
|
\end{minipage} } |
401 |
|
|
|
402 |
|
|
|
403 |
|
|
\subsection{Topography: partially filled cells} |
404 |
edhill |
1.15 |
\begin{rawhtml} |
405 |
|
|
<!-- CMIREDIR:topo_partial_cells: --> |
406 |
|
|
\end{rawhtml} |
407 |
adcroft |
1.3 |
|
408 |
|
|
\begin{figure} |
409 |
adcroft |
1.7 |
\begin{center} |
410 |
jmc |
1.23 |
\resizebox{4.5in}{!}{\includegraphics{s_algorithm/figs/vgrid-xz.eps}} |
411 |
adcroft |
1.7 |
\end{center} |
412 |
adcroft |
1.3 |
\caption{ |
413 |
|
|
A schematic of the x-r plane showing the location of the |
414 |
|
|
non-dimensional fractions $h_c$ and $h_w$. The physical thickness of a |
415 |
|
|
tracer cell is given by $h_c(i,j,k) \Delta r_f(k)$ and the physical |
416 |
|
|
thickness of the open side is given by $h_w(i,j,k) \Delta r_f(k)$.} |
417 |
|
|
\label{fig:hfacs} |
418 |
|
|
\end{figure} |
419 |
|
|
|
420 |
adcroft |
1.11 |
\cite{adcroft:97} presented two alternatives to the step-wise finite |
421 |
adcroft |
1.3 |
difference representation of topography. The method is known to the |
422 |
|
|
engineering community as {\em intersecting boundary method}. It |
423 |
|
|
involves allowing the boundary to intersect a grid of cells thereby |
424 |
|
|
modifying the shape of those cells intersected. We suggested allowing |
425 |
cnh |
1.10 |
the topography to take on a piece-wise linear representation (shaved |
426 |
|
|
cells) or a simpler piecewise constant representation (partial step). |
427 |
adcroft |
1.3 |
Both show dramatic improvements in solution compared to the |
428 |
|
|
traditional full step representation, the piece-wise linear being the |
429 |
|
|
best. However, the storage requirements are excessive so the simpler |
430 |
|
|
piece-wise constant or partial-step method is all that is currently |
431 |
|
|
supported. |
432 |
|
|
|
433 |
|
|
Fig.~\ref{fig:hfacs} shows a schematic of the x-r plane indicating how |
434 |
|
|
the thickness of a level is determined at tracer and u points. |
435 |
|
|
\marginpar{$h_c$: {\bf hFacC}} |
436 |
|
|
\marginpar{$h_w$: {\bf hFacW}} |
437 |
|
|
\marginpar{$h_s$: {\bf hFacS}} |
438 |
|
|
The physical thickness of a tracer cell is given by $h_c(i,j,k) \Delta |
439 |
|
|
r_f(k)$ and the physical thickness of the open side is given by |
440 |
cnh |
1.10 |
$h_w(i,j,k) \Delta r_f(k)$. Three 3-D descriptors $h_c$, $h_w$ and |
441 |
adcroft |
1.3 |
$h_s$ are used to describe the geometry: {\bf hFacC}, {\bf hFacW} and |
442 |
|
|
{\bf hFacS} respectively. These are calculated in subroutine {\em |
443 |
|
|
INI\_MASKS\_ETC} along with there reciprocals {\bf RECIP\_hFacC}, {\bf |
444 |
|
|
RECIP\_hFacW} and {\bf RECIP\_hFacS}. |
445 |
|
|
|
446 |
|
|
The non-dimensional fractions (or h-facs as we call them) are |
447 |
|
|
calculated from the model depth array and then processed to avoid tiny |
448 |
|
|
volumes. The rule is that if a fraction is less than {\bf hFacMin} |
449 |
|
|
then it is rounded to the nearer of $0$ or {\bf hFacMin} or if the |
450 |
|
|
physical thickness is less than {\bf hFacMinDr} then it is similarly |
451 |
|
|
rounded. The larger of the two methods is used when there is a |
452 |
|
|
conflict. By setting {\bf hFacMinDr} equal to or larger than the |
453 |
|
|
thinnest nominal layers, $\min{(\Delta z_f)}$, but setting {\bf |
454 |
|
|
hFacMin} to some small fraction then the model will only lop thick |
455 |
|
|
layers but retain stability based on the thinnest unlopped thickness; |
456 |
|
|
$\min{(\Delta z_f,\mbox{\bf hFacMinDr})}$. |
457 |
|
|
|
458 |
|
|
\fbox{ \begin{minipage}{4.75in} |
459 |
|
|
{\em S/R INI\_MASKS\_ETC} ({\em model/src/ini\_masks\_etc.F}) |
460 |
|
|
|
461 |
|
|
$h_c$: {\bf hFacC} ({\em GRID.h}) |
462 |
|
|
|
463 |
|
|
$h_w$: {\bf hFacW} ({\em GRID.h}) |
464 |
|
|
|
465 |
|
|
$h_s$: {\bf hFacS} ({\em GRID.h}) |
466 |
|
|
|
467 |
|
|
$h_c^{-1}$: {\bf RECIP\_hFacC} ({\em GRID.h}) |
468 |
|
|
|
469 |
|
|
$h_w^{-1}$: {\bf RECIP\_hFacW} ({\em GRID.h}) |
470 |
|
|
|
471 |
|
|
$h_s^{-1}$: {\bf RECIP\_hFacS} ({\em GRID.h}) |
472 |
adcroft |
1.2 |
|
473 |
|
|
\end{minipage} } |
474 |
|
|
|
475 |
adcroft |
1.1 |
|
476 |
adcroft |
1.5 |
\section{Continuity and horizontal pressure gradient terms} |
477 |
edhill |
1.17 |
\begin{rawhtml} |
478 |
|
|
<!-- CMIREDIR:continuity_and_horizontal_pressure: --> |
479 |
|
|
\end{rawhtml} |
480 |
|
|
|
481 |
adcroft |
1.1 |
|
482 |
|
|
The core algorithm is based on the ``C grid'' discretization of the |
483 |
|
|
continuity equation which can be summarized as: |
484 |
|
|
\begin{eqnarray} |
485 |
adcroft |
1.6 |
\partial_t u + \frac{1}{\Delta x_c} \delta_i \left. \frac{ \partial \Phi}{\partial r}\right|_{s} \eta + \frac{\epsilon_{nh}}{\Delta x_c} \delta_i \Phi_{nh}' & = & G_u - \frac{1}{\Delta x_c} \delta_i \Phi_h' \label{eq:discrete-momu} \\ |
486 |
|
|
\partial_t v + \frac{1}{\Delta y_c} \delta_j \left. \frac{ \partial \Phi}{\partial r}\right|_{s} \eta + \frac{\epsilon_{nh}}{\Delta y_c} \delta_j \Phi_{nh}' & = & G_v - \frac{1}{\Delta y_c} \delta_j \Phi_h' \label{eq:discrete-momv} \\ |
487 |
|
|
\epsilon_{nh} \left( \partial_t w + \frac{1}{\Delta r_c} \delta_k \Phi_{nh}' \right) & = & \epsilon_{nh} G_w + \overline{b}^k - \frac{1}{\Delta r_c} \delta_k \Phi_{h}' \label{eq:discrete-momw} \\ |
488 |
adcroft |
1.1 |
\delta_i \Delta y_g \Delta r_f h_w u + |
489 |
|
|
\delta_j \Delta x_g \Delta r_f h_s v + |
490 |
|
|
\delta_k {\cal A}_c w & = & {\cal A}_c \delta_k (P-E)_{r=0} |
491 |
adcroft |
1.2 |
\label{eq:discrete-continuity} |
492 |
adcroft |
1.1 |
\end{eqnarray} |
493 |
|
|
where the continuity equation has been most naturally discretized by |
494 |
|
|
staggering the three components of velocity as shown in |
495 |
adcroft |
1.12 |
Fig.~\ref{fig:cgrid3d}. The grid lengths $\Delta x_c$ and $\Delta y_c$ |
496 |
adcroft |
1.1 |
are the lengths between tracer points (cell centers). The grid lengths |
497 |
|
|
$\Delta x_g$, $\Delta y_g$ are the grid lengths between cell |
498 |
|
|
corners. $\Delta r_f$ and $\Delta r_c$ are the distance (in units of |
499 |
|
|
$r$) between level interfaces (w-level) and level centers (tracer |
500 |
|
|
level). The surface area presented in the vertical is denoted ${\cal |
501 |
|
|
A}_c$. The factors $h_w$ and $h_s$ are non-dimensional fractions |
502 |
|
|
(between 0 and 1) that represent the fraction cell depth that is |
503 |
|
|
``open'' for fluid flow. |
504 |
|
|
\marginpar{$h_w$: {\bf hFacW}} |
505 |
|
|
\marginpar{$h_s$: {\bf hFacS}} |
506 |
|
|
|
507 |
|
|
The last equation, the discrete continuity equation, can be summed in |
508 |
cnh |
1.10 |
the vertical to yield the free-surface equation: |
509 |
adcroft |
1.1 |
\begin{equation} |
510 |
adcroft |
1.6 |
{\cal A}_c \partial_t \eta + \delta_i \sum_k \Delta y_g \Delta r_f h_w |
511 |
|
|
u + \delta_j \sum_k \Delta x_g \Delta r_f h_s v = {\cal |
512 |
|
|
A}_c(P-E)_{r=0} \label{eq:discrete-freesurface} |
513 |
adcroft |
1.1 |
\end{equation} |
514 |
|
|
The source term $P-E$ on the rhs of continuity accounts for the local |
515 |
|
|
addition of volume due to excess precipitation and run-off over |
516 |
|
|
evaporation and only enters the top-level of the {\em ocean} model. |
517 |
|
|
|
518 |
adcroft |
1.5 |
\section{Hydrostatic balance} |
519 |
edhill |
1.17 |
\begin{rawhtml} |
520 |
|
|
<!-- CMIREDIR:hydrostatic_balance: --> |
521 |
|
|
\end{rawhtml} |
522 |
adcroft |
1.1 |
|
523 |
|
|
The vertical momentum equation has the hydrostatic or |
524 |
|
|
quasi-hydrostatic balance on the right hand side. This discretization |
525 |
|
|
guarantees that the conversion of potential to kinetic energy as |
526 |
|
|
derived from the buoyancy equation exactly matches the form derived |
527 |
|
|
from the pressure gradient terms when forming the kinetic energy |
528 |
|
|
equation. |
529 |
|
|
|
530 |
cnh |
1.10 |
In the ocean, using z-coordinates, the hydrostatic balance terms are |
531 |
adcroft |
1.1 |
discretized: |
532 |
|
|
\begin{equation} |
533 |
|
|
\epsilon_{nh} \partial_t w |
534 |
|
|
+ g \overline{\rho'}^k + \frac{1}{\Delta z} \delta_k \Phi_h' = \ldots |
535 |
adcroft |
1.6 |
\label{eq:discrete_hydro_ocean} |
536 |
adcroft |
1.1 |
\end{equation} |
537 |
|
|
|
538 |
|
|
In the atmosphere, using p-coordinates, hydrostatic balance is |
539 |
|
|
discretized: |
540 |
|
|
\begin{equation} |
541 |
|
|
\overline{\theta'}^k + \frac{1}{\Delta \Pi} \delta_k \Phi_h' = 0 |
542 |
adcroft |
1.6 |
\label{eq:discrete_hydro_atmos} |
543 |
adcroft |
1.1 |
\end{equation} |
544 |
|
|
where $\Delta \Pi$ is the difference in Exner function between the |
545 |
|
|
pressure points. The non-hydrostatic equations are not available in |
546 |
|
|
the atmosphere. |
547 |
|
|
|
548 |
|
|
The difference in approach between ocean and atmosphere occurs because |
549 |
|
|
of the direct use of the ideal gas equation in forming the potential |
550 |
cnh |
1.10 |
energy conversion term $\alpha \omega$. The form of these conversion |
551 |
adcroft |
1.11 |
terms is discussed at length in \cite{adcroft:02}. |
552 |
adcroft |
1.1 |
|
553 |
|
|
Because of the different representation of hydrostatic balance between |
554 |
|
|
ocean and atmosphere there is no elegant way to represent both systems |
555 |
|
|
using an arbitrary coordinate. |
556 |
|
|
|
557 |
|
|
The integration for hydrostatic pressure is made in the positive $r$ |
558 |
|
|
direction (increasing k-index). For the ocean, this is from the |
559 |
|
|
free-surface down and for the atmosphere this is from the ground up. |
560 |
|
|
|
561 |
|
|
The calculations are made in the subroutine {\em |
562 |
|
|
CALC\_PHI\_HYD}. Inside this routine, one of other of the |
563 |
|
|
atmospheric/oceanic form is selected based on the string variable {\bf |
564 |
|
|
buoyancyRelation}. |
565 |
|
|
|