/[MITgcm]/manual/s_examples/baroclinic_gyre/fourlayer.tex
ViewVC logotype

Contents of /manual/s_examples/baroclinic_gyre/fourlayer.tex

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.15 - (show annotations) (download) (as text)
Thu Aug 7 18:27:52 2003 UTC (21 years, 11 months ago) by edhill
Branch: MAIN
Changes since 1.14: +90 -71 lines
File MIME type: application/x-tex
Two improvements:
  - the manual build now works on a "dirty" (post-verification-runs)
    source tree
  - various small fixes to the manual (mostly Y. Hu's comments)

1 % $Header: /u/u3/gcmpack/manual/part3/case_studies/fourlayer_gyre/fourlayer.tex,v 1.14 2003/07/30 13:42:53 edhill Exp $
2 % $Name: $
3
4 \section{Four Layer Baroclinic Ocean Gyre In Spherical Coordinates}
5 \label{www:tutorials}
6 \label{sect:eg-fourlayer}
7
8 \bodytext{bgcolor="#FFFFFFFF"}
9
10 %\begin{center}
11 %{\Large \bf Using MITgcm to Simulate a Baroclinic Ocean Gyre In Spherical
12 %Polar Coordinates}
13 %
14 %\vspace*{4mm}
15 %
16 %\vspace*{3mm}
17 %{\large May 2001}
18 %\end{center}
19
20 This document describes an example experiment using MITgcm
21 to simulate a baroclinic ocean gyre in spherical
22 polar coordinates. The barotropic
23 example experiment in section \ref{sect:eg-baro}
24 illustrated how to configure the code for a single layer
25 simulation in a Cartesian grid. In this example a similar physical problem
26 is simulated, but the code is now configured
27 for four layers and in a spherical polar coordinate system.
28
29 \subsection{Overview}
30 \label{www:tutorials}
31
32 This example experiment demonstrates using the MITgcm to simulate
33 a baroclinic, wind-forced, ocean gyre circulation. The experiment
34 is a numerical rendition of the gyre circulation problem similar
35 to the problems described analytically by Stommel in 1966
36 \cite{Stommel66} and numerically in Holland et. al \cite{Holland75}.
37 \\
38
39 In this experiment the model is configured to represent a mid-latitude
40 enclosed sector of fluid on a sphere, $60^{\circ} \times 60^{\circ}$ in
41 lateral extent. The fluid is $2$~km deep and is forced
42 by a constant in time zonal wind stress, $\tau_{\lambda}$, that varies
43 sinusoidally in the north-south direction. Topologically the simulated
44 domain is a sector on a sphere and the coriolis parameter, $f$, is defined
45 according to latitude, $\varphi$
46
47 \begin{equation}
48 \label{EQ:eg-fourlayer-fcori}
49 f(\varphi) = 2 \Omega \sin( \varphi )
50 \end{equation}
51
52 \noindent with the rotation rate, $\Omega$ set to $\frac{2 \pi}{86400s}$.
53 \\
54
55 The sinusoidal wind-stress variations are defined according to
56
57 \begin{equation}
58 \label{EQ:taux}
59 \tau_{\lambda}(\varphi) = \tau_{0}\sin(\pi \frac{\varphi}{L_{\varphi}})
60 \end{equation}
61
62 \noindent where $L_{\varphi}$ is the lateral domain extent ($60^{\circ}$) and
63 $\tau_0$ is set to $0.1N m^{-2}$.
64 \\
65
66 Figure \ref{FIG:eg-fourlayer-simulation_config}
67 summarizes the configuration simulated.
68 In contrast to the example in section \ref{sect:eg-baro}, the
69 current experiment simulates a spherical polar domain. As indicated
70 by the axes in the lower left of the figure the model code works internally
71 in a locally orthogonal coordinate $(x,y,z)$. For this experiment description
72 the local orthogonal model coordinate $(x,y,z)$ is synonymous
73 with the coordinates $(\lambda,\varphi,r)$ shown in figure
74 \ref{fig:spherical-polar-coord}
75 \\
76
77 The experiment has four levels in the vertical, each of equal thickness,
78 $\Delta z = 500$~m. Initially the fluid is stratified with a reference
79 potential temperature profile,
80 $\theta_{250}=20^{\circ}$~C,
81 $\theta_{750}=10^{\circ}$~C,
82 $\theta_{1250}=8^{\circ}$~C,
83 $\theta_{1750}=6^{\circ}$~C. The equation of state used in this experiment is
84 linear
85
86 \begin{equation}
87 \label{EQ:eg-fourlayer-linear1_eos}
88 \rho = \rho_{0} ( 1 - \alpha_{\theta}\theta^{'} )
89 \end{equation}
90
91 \noindent which is implemented in the model as a density anomaly equation
92
93 \begin{equation}
94 \label{EQ:eg-fourlayer-linear1_eos_pert}
95 \rho^{'} = -\rho_{0}\alpha_{\theta}\theta^{'}
96 \end{equation}
97
98 \noindent with $\rho_{0}=999.8\,{\rm kg\,m}^{-3}$ and
99 $\alpha_{\theta}=2\times10^{-4}\,{\rm degrees}^{-1} $. Integrated forward in
100 this configuration the model state variable {\bf theta} is equivalent to
101 either in-situ temperature, $T$, or potential temperature, $\theta$. For
102 consistency with later examples, in which the equation of state is
103 non-linear, we use $\theta$ to represent temperature here. This is
104 the quantity that is carried in the model core equations.
105
106 \begin{figure}
107 \begin{center}
108 \resizebox{7.5in}{5.5in}{
109 \includegraphics*[0.2in,0.7in][10.5in,10.5in]
110 {part3/case_studies/fourlayer_gyre/simulation_config.eps} }
111 \end{center}
112 \caption{Schematic of simulation domain and wind-stress forcing function
113 for the four-layer gyre numerical experiment. The domain is enclosed by solid
114 walls at $0^{\circ}$~E, $60^{\circ}$~E, $0^{\circ}$~N and $60^{\circ}$~N.
115 An initial stratification is
116 imposed by setting the potential temperature, $\theta$, in each layer.
117 The vertical spacing, $\Delta z$, is constant and equal to $500$m.
118 }
119 \label{FIG:eg-fourlayer-simulation_config}
120 \end{figure}
121
122 \subsection{Equations solved}
123 \label{www:tutorials}
124 For this problem
125 the implicit free surface, {\bf HPE} (see section \ref{sect:hydrostatic_and_quasi-hydrostatic_forms}) form of the
126 equations described in Marshall et. al \cite{marshall:97a} are
127 employed. The flow is three-dimensional with just temperature, $\theta$, as
128 an active tracer. The equation of state is linear.
129 A horizontal Laplacian operator $\nabla_{h}^2$ provides viscous
130 dissipation and provides a diffusive sub-grid scale closure for the
131 temperature equation. A wind-stress momentum forcing is added to the momentum
132 equation for the zonal flow, $u$. Other terms in the model
133 are explicitly switched off for this experiment configuration (see section
134 \ref{SEC:eg_fourl_code_config} ). This yields an active set of equations
135 solved in this configuration, written in spherical polar coordinates as
136 follows
137
138 \begin{eqnarray}
139 \label{EQ:eg-fourlayer-model_equations}
140 \frac{Du}{Dt} - fv +
141 \frac{1}{\rho}\frac{\partial p^{\prime}}{\partial \lambda} -
142 A_{h}\nabla_{h}^2u - A_{z}\frac{\partial^{2}u}{\partial z^{2}}
143 & = &
144 \cal{F}_{\lambda}
145 \\
146 \frac{Dv}{Dt} + fu +
147 \frac{1}{\rho}\frac{\partial p^{\prime}}{\partial \varphi} -
148 A_{h}\nabla_{h}^2v - A_{z}\frac{\partial^{2}v}{\partial z^{2}}
149 & = &
150 0
151 \\
152 \frac{\partial \eta}{\partial t} + \frac{\partial H \widehat{u}}{\partial \lambda} +
153 \frac{\partial H \widehat{v}}{\partial \varphi}
154 &=&
155 0
156 \label{eq:fourl_example_continuity}
157 \\
158 \frac{D\theta}{Dt} -
159 K_{h}\nabla_{h}^2\theta - K_{z}\frac{\partial^{2}\theta}{\partial z^{2}}
160 & = &
161 0
162 \label{eq:eg_fourl_theta}
163 \\
164 p^{\prime} & = &
165 g\rho_{0} \eta + \int^{0}_{-z}\rho^{\prime} dz
166 \\
167 \rho^{\prime} & = &- \alpha_{\theta}\rho_{0}\theta^{\prime}
168 \\
169 {\cal F}_{\lambda} |_{s} & = & \frac{\tau_{\lambda}}{\rho_{0}\Delta z_{s}}
170 \\
171 {\cal F}_{\lambda} |_{i} & = & 0
172 \end{eqnarray}
173
174 \noindent where $u$ and $v$ are the components of the horizontal
175 flow vector $\vec{u}$ on the sphere ($u=\dot{\lambda},v=\dot{\varphi}$).
176 The terms $H\widehat{u}$ and $H\widehat{v}$ are the components of the vertical
177 integral term given in equation \ref{eq:free-surface} and
178 explained in more detail in section \ref{sect:pressure-method-linear-backward}.
179 However, for the problem presented here, the continuity relation (equation
180 \ref{eq:fourl_example_continuity}) differs from the general form given
181 in section \ref{sect:pressure-method-linear-backward},
182 equation \ref{eq:linear-free-surface=P-E+R}, because the source terms
183 ${\cal P}-{\cal E}+{\cal R}$
184 are all $0$.
185
186 The pressure field, $p^{\prime}$, is separated into a barotropic part
187 due to variations in sea-surface height, $\eta$, and a hydrostatic
188 part due to variations in density, $\rho^{\prime}$, integrated
189 through the water column.
190
191 The suffices ${s},{i}$ indicate surface layer and the interior of the domain.
192 The windstress forcing, ${\cal F}_{\lambda}$, is applied in the surface layer
193 by a source term in the zonal momentum equation. In the ocean interior
194 this term is zero.
195
196 In the momentum equations
197 lateral and vertical boundary conditions for the $\nabla_{h}^{2}$
198 and $\frac{\partial^{2}}{\partial z^{2}}$ operators are specified
199 when the numerical simulation is run - see section
200 \ref{SEC:eg_fourl_code_config}. For temperature
201 the boundary condition is ``zero-flux''
202 e.g. $\frac{\partial \theta}{\partial \varphi}=
203 \frac{\partial \theta}{\partial \lambda}=\frac{\partial \theta}{\partial z}=0$.
204
205
206
207 \subsection{Discrete Numerical Configuration}
208 \label{www:tutorials}
209
210 The domain is discretised with
211 a uniform grid spacing in latitude and longitude
212 $\Delta \lambda=\Delta \varphi=1^{\circ}$, so
213 that there are sixty grid cells in the zonal and meridional directions.
214 Vertically the
215 model is configured with four layers with constant depth,
216 $\Delta z$, of $500$~m. The internal, locally orthogonal, model coordinate
217 variables $x$ and $y$ are initialized from the values of
218 $\lambda$, $\varphi$, $\Delta \lambda$ and $\Delta \varphi$ in
219 radians according to
220
221 \begin{eqnarray}
222 x=r\cos(\varphi)\lambda,~\Delta x & = &r\cos(\varphi)\Delta \lambda \\
223 y=r\varphi,~\Delta y &= &r\Delta \varphi
224 \end{eqnarray}
225
226 The procedure for generating a set of internal grid variables from a
227 spherical polar grid specification is discussed in section
228 \ref{sect:spatial_discrete_horizontal_grid}.
229
230 \noindent\fbox{ \begin{minipage}{5.5in}
231 {\em S/R INI\_SPHERICAL\_POLAR\_GRID} ({\em
232 model/src/ini\_spherical\_polar\_grid.F})
233
234 $A_c$, $A_\zeta$, $A_w$, $A_s$: {\bf rAc}, {\bf rAz}, {\bf rAw}, {\bf rAs}
235 ({\em GRID.h})
236
237 $\Delta x_g$, $\Delta y_g$: {\bf DXg}, {\bf DYg} ({\em GRID.h})
238
239 $\Delta x_c$, $\Delta y_c$: {\bf DXc}, {\bf DYc} ({\em GRID.h})
240
241 $\Delta x_f$, $\Delta y_f$: {\bf DXf}, {\bf DYf} ({\em GRID.h})
242
243 $\Delta x_v$, $\Delta y_u$: {\bf DXv}, {\bf DYu} ({\em GRID.h})
244
245 \end{minipage} }\\
246
247
248
249 As described in \ref{sect:tracer_equations}, the time evolution of potential
250 temperature,
251 $\theta$, (equation \ref{eq:eg_fourl_theta})
252 is evaluated prognostically. The centered second-order scheme with
253 Adams-Bashforth time stepping described in section
254 \ref{sect:tracer_equations_abII} is used to step forward the temperature
255 equation. Prognostic terms in
256 the momentum equations are solved using flux form as
257 described in section \ref{sect:flux-form_momentum_eqautions}.
258 The pressure forces that drive the fluid motions, (
259 $\frac{\partial p^{'}}{\partial \lambda}$ and $\frac{\partial p^{'}}{\partial \varphi}$), are found by summing pressure due to surface
260 elevation $\eta$ and the hydrostatic pressure. The hydrostatic part of the
261 pressure is diagnosed explicitly by integrating density. The sea-surface
262 height, $\eta$, is diagnosed using an implicit scheme. The pressure
263 field solution method is described in sections
264 \ref{sect:pressure-method-linear-backward} and
265 \ref{sect:finding_the_pressure_field}.
266
267 \subsubsection{Numerical Stability Criteria}
268 \label{www:tutorials}
269
270 The Laplacian viscosity coefficient, $A_{h}$, is set to $400 m s^{-1}$.
271 This value is chosen to yield a Munk layer width,
272
273 \begin{eqnarray}
274 \label{EQ:eg-fourlayer-munk_layer}
275 M_{w} = \pi ( \frac { A_{h} }{ \beta } )^{\frac{1}{3}}
276 \end{eqnarray}
277
278 \noindent of $\approx 100$km. This is greater than the model
279 resolution in mid-latitudes
280 $\Delta x=r \cos(\varphi) \Delta \lambda \approx 80~{\rm km}$ at
281 $\varphi=45^{\circ}$, ensuring that the frictional
282 boundary layer is well resolved.
283 \\
284
285 \noindent The model is stepped forward with a
286 time step $\delta t=1200$secs. With this time step the stability
287 parameter to the horizontal Laplacian friction
288
289 \begin{eqnarray}
290 \label{EQ:eg-fourlayer-laplacian_stability}
291 S_{l} = 4 \frac{A_{h} \delta t}{{\Delta x}^2}
292 \end{eqnarray}
293
294 \noindent evaluates to 0.012, which is well below the 0.3 upper limit
295 for stability for this term under ABII time-stepping.
296 \\
297
298 \noindent The vertical dissipation coefficient, $A_{z}$, is set to
299 $1\times10^{-2} {\rm m}^2{\rm s}^{-1}$. The associated stability limit
300
301 \begin{eqnarray}
302 \label{EQ:eg-fourlayer-laplacian_stability_z}
303 S_{l} = 4 \frac{A_{z} \delta t}{{\Delta z}^2}
304 \end{eqnarray}
305
306 \noindent evaluates to $4.8 \times 10^{-5}$ which is again well below
307 the upper limit.
308 The values of $A_{h}$ and $A_{z}$ are also used for the horizontal ($K_{h}$)
309 and vertical ($K_{z}$) diffusion coefficients for temperature respectively.
310 \\
311
312 \noindent The numerical stability for inertial oscillations
313
314 \begin{eqnarray}
315 \label{EQ:eg-fourlayer-inertial_stability}
316 S_{i} = f^{2} {\delta t}^2
317 \end{eqnarray}
318
319 \noindent evaluates to $0.0144$, which is well below the $0.5$ upper
320 limit for stability.
321 \\
322
323 \noindent The advective CFL for a extreme maximum
324 horizontal flow
325 speed of $ | \vec{u} | = 2 ms^{-1}$
326
327 \begin{eqnarray}
328 \label{EQ:eg-fourlayer-cfl_stability}
329 C_{a} = \frac{| \vec{u} | \delta t}{ \Delta x}
330 \end{eqnarray}
331
332 \noindent evaluates to $5 \times 10^{-2}$. This is well below the stability
333 limit of 0.5.
334 \\
335
336 \noindent The stability parameter for internal gravity waves
337 propagating at $2~{\rm m}~{\rm s}^{-1}$
338
339 \begin{eqnarray}
340 \label{EQ:eg-fourlayer-igw_stability}
341 S_{c} = \frac{c_{g} \delta t}{ \Delta x}
342 \end{eqnarray}
343
344 \noindent evaluates to $\approx 5 \times 10^{-2}$. This is well below the linear
345 stability limit of 0.25.
346
347 \subsection{Code Configuration}
348 \label{www:tutorials}
349 \label{SEC:eg_fourl_code_config}
350
351 The model configuration for this experiment resides under the
352 directory {\it verification/exp2/}. The experiment files
353 \begin{itemize}
354 \item {\it input/data}
355 \item {\it input/data.pkg}
356 \item {\it input/eedata},
357 \item {\it input/windx.sin\_y},
358 \item {\it input/topog.box},
359 \item {\it code/CPP\_EEOPTIONS.h}
360 \item {\it code/CPP\_OPTIONS.h},
361 \item {\it code/SIZE.h}.
362 \end{itemize}
363 contain the code customisations and parameter settings for this
364 experiment. Below we describe the customisations to these files
365 associated with this experiment.
366
367 \subsubsection{File {\it input/data}}
368 \label{www:tutorials}
369
370 This file, reproduced completely below, specifies the main parameters
371 for the experiment. The parameters that are significant for this configuration
372 are
373
374 \begin{itemize}
375
376 \item Line 4,
377 \begin{verbatim} tRef=20.,10.,8.,6., \end{verbatim}
378 this line sets
379 the initial and reference values of potential temperature at each model
380 level in units of $^{\circ}$C.
381 The entries are ordered from surface to depth. For each
382 depth level the initial and reference profiles will be uniform in
383 $x$ and $y$. The values specified here are read into the
384 variable
385 \varlink{tRef}{tRef}
386 %{\bf
387 %\begin{rawhtml} <A href=../code_reference/vdb/names/OK.htm> \end{rawhtml}
388 %tRef
389 %\begin{rawhtml} </A>\end{rawhtml}
390 %}
391 in the model code, by procedure
392 \filelink{INI\_PARMS}{model-src-ini_parms.F}
393 %{\it
394 %\begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
395 %INI\_PARMS
396 %\begin{rawhtml} </A>\end{rawhtml}
397 %}.
398
399 %% \codelink{var:tref} tRef \endlink
400 %% \codelink{file:ini_parms} {\it INI\_PARMS } \endlink
401 %% \codelink{proc:ini_parms} {\it INI\_PARMS } \endlink
402 %% \var{tref}
403 %% \proc{ini_parms}
404 %% \file{ini_parms}
405 \newcommand{\VARtref}{
406 {\bf
407 \begin{rawhtml} <A href=../code_reference/vdb/names/OK.htm> \end{rawhtml}
408 tRef
409 \begin{rawhtml} </A>\end{rawhtml}
410 }
411 }
412
413
414
415 \fbox{
416 \begin{minipage}{5.0in}
417 {\it S/R INI\_THETA}
418 ({\it ini\_theta.F})
419 \end{minipage}
420 }
421 \filelink{ini\_theta.F}{model-src-ini_theta.F}
422 %{\bf
423 %\begin{rawhtml} <A href=../code_reference/vdb/code/98.htm> \end{rawhtml}
424 %goto code
425 %\begin{rawhtml} </A>\end{rawhtml}
426 %}
427
428
429 \item Line 6,
430 \begin{verbatim} viscAz=1.E-2, \end{verbatim}
431 this line sets the vertical Laplacian dissipation coefficient to
432 $1 \times 10^{-2} {\rm m^{2}s^{-1}}$. Boundary conditions
433 for this operator are specified later.
434 The variable
435 \varlink{viscAz}{viscAz}
436 %{\bf
437 %\begin{rawhtml} <A href=../code_reference/vdb/names/ZQ.htm> \end{rawhtml}
438 %viscAz
439 %\begin{rawhtml} </A>\end{rawhtml}
440 %}
441 is read in the routine
442 \filelink{ini\_parms.F}{model-src-ini_parms.F}
443 %{\it
444 %\begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
445 %INI\_PARMS
446 %\begin{rawhtml} </A>\end{rawhtml}
447 %}
448 and is copied into model general vertical coordinate variable
449 \varlink{viscAr}{viscAr}
450 %{\bf
451 %\begin{rawhtml} <A href=../code_reference/vdb/names/PF.htm> \end{rawhtml}
452 %viscAr
453 %\begin{rawhtml} </A>\end{rawhtml}
454 %}.
455 At each time step, the viscous term contribution to the momentum equations
456 is calculated in routine
457 %{\it S/R CALC\_DIFFUSIVITY}.
458 \varlink{CALC\_DIFFUSIVITY}{CALC_DIFFUSIVITY}
459
460 \fbox{
461 \begin{minipage}{5.0in}
462 {\it S/R CALC\_DIFFUSIVITY}({\it calc\_diffusivity.F})
463 \end{minipage}
464 }
465 %{\bf
466 %\begin{rawhtml} <A href=../code_reference/vdb/code/53.htm> \end{rawhtml}
467 %goto code
468 %\begin{rawhtml} </A>\end{rawhtml}
469 %}
470
471 \item Line 7,
472 \begin{verbatim}
473 viscAh=4.E2,
474 \end{verbatim}
475 this line sets the horizontal laplacian frictional dissipation coefficient to
476 $1 \times 10^{-2} {\rm m^{2}s^{-1}}$. Boundary conditions
477 for this operator are specified later.
478 The variable
479 \varlink{viscAh}{viscAh}
480 %{\bf
481 %\begin{rawhtml} <A href=../code_reference/vdb/names/SI.htm> \end{rawhtml}
482 %viscAh
483 %\begin{rawhtml} </A>\end{rawhtml}
484 %}
485 is read in the routine
486 \varlink{INI\_PARMS}{INI_PARMS}
487 %{\it
488 %\begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
489 %INI\_PARMS
490 %\begin{rawhtml} </A>\end{rawhtml}
491 %}
492 and applied in routines
493 %{\it CALC\_MOM\_RHS} and {\it CALC\_GW}.
494 \varlink{CALC\_MOM\_RHS}{CALC_MOM_RHS}
495 and
496 \varlink{CALC\_GW}{CALC_GW}.
497
498
499 \fbox{
500 \begin{minipage}{5.0in}
501 {\it S/R CALC\_MOM\_RHS}({\it calc\_mom\_rhs.F})
502 \end{minipage}
503 }
504 %{\bf
505 %\begin{rawhtml} <A href=../code_reference/vdb/code/60.htm> \end{rawhtml}
506 %goto code
507 %\begin{rawhtml} </A>\end{rawhtml}
508 %}
509
510 \fbox{
511 \begin{minipage}{5.0in}
512 {\it S/R CALC\_GW}({\it calc\_gw.F})
513 \end{minipage}
514 }
515 %{\bf
516 %\begin{rawhtml} <A href=../code_reference/vdb/code/58.htm> \end{rawhtml}
517 %goto code
518 %\begin{rawhtml} </A>\end{rawhtml}
519 %}
520
521 \item Lines 8,
522 \begin{verbatim}
523 no_slip_sides=.FALSE.
524 \end{verbatim}
525 this line selects a free-slip lateral boundary condition for
526 the horizontal laplacian friction operator
527 e.g. $\frac{\partial u}{\partial y}$=0 along boundaries in $y$ and
528 $\frac{\partial v}{\partial x}$=0 along boundaries in $x$.
529 The variable
530 \varlink{no\_slip\_sides}{no_slip_sides}
531 %{\bf
532 %\begin{rawhtml} <A href=../code_reference/vdb/names/UT.htm> \end{rawhtml}
533 %no\_slip\_sides
534 %\begin{rawhtml} </A>\end{rawhtml}
535 %}
536 is read in the routine
537 \varlink{INI\_PARMS}{INI_PARMS}
538 %{\it
539 %\begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
540 %INI\_PARMS
541 %\begin{rawhtml} </A>\end{rawhtml}
542 %}
543 and the boundary condition is evaluated in routine
544 %{\it S/R CALC\_MOM\_RHS}.
545
546
547 \fbox{
548 \begin{minipage}{5.0in}
549 {\it S/R CALC\_MOM\_RHS}({\it calc\_mom\_rhs.F})
550 \end{minipage}
551 }
552 {\bf
553 \begin{rawhtml} <A href=../code_reference/vdb/code/60.htm> \end{rawhtml}
554 goto code
555 \begin{rawhtml} </A>\end{rawhtml}
556 }
557
558 \item Lines 9,
559 \begin{verbatim}
560 no_slip_bottom=.TRUE.
561 \end{verbatim}
562 this line selects a no-slip boundary condition for bottom
563 boundary condition in the vertical laplacian friction operator
564 e.g. $u=v=0$ at $z=-H$, where $H$ is the local depth of the domain.
565 The variable
566 {\bf
567 \begin{rawhtml} <A href=../code_reference/vdb/names/UK.htm> \end{rawhtml}
568 no\_slip\_bottom
569 \begin{rawhtml} </A>\end{rawhtml}
570 }
571 is read in the routine
572 {\it
573 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
574 INI\_PARMS
575 \begin{rawhtml} </A>\end{rawhtml}
576 } and is applied in the routine {\it S/R CALC\_MOM\_RHS}.
577
578 \fbox{
579 \begin{minipage}{5.0in}
580 {\it S/R CALC\_MOM\_RHS}({\it calc\_mom\_rhs.F})
581 \end{minipage}
582 }
583 {\bf
584 \begin{rawhtml} <A href=../code_reference/vdb/code/60.htm> \end{rawhtml}
585 goto code
586 \begin{rawhtml} </A>\end{rawhtml}
587 }
588
589 \item Line 10,
590 \begin{verbatim}
591 diffKhT=4.E2,
592 \end{verbatim}
593 this line sets the horizontal diffusion coefficient for temperature
594 to $400\,{\rm m^{2}s^{-1}}$. The boundary condition on this
595 operator is $\frac{\partial}{\partial x}=\frac{\partial}{\partial y}=0$ at
596 all boundaries.
597 The variable
598 {\bf
599 \begin{rawhtml} <A href=../code_reference/vdb/names/RC.htm> \end{rawhtml}
600 diffKhT
601 \begin{rawhtml} </A>\end{rawhtml}
602 }
603 is read in the routine
604 {\it
605 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
606 INI\_PARMS
607 \begin{rawhtml} </A>\end{rawhtml}
608 } and used in routine {\it S/R CALC\_GT}.
609
610 \fbox{ \begin{minipage}{5.0in}
611 {\it S/R CALC\_GT}({\it calc\_gt.F})
612 \end{minipage}
613 }
614 {\bf
615 \begin{rawhtml} <A href=../code_reference/vdb/code/57.htm> \end{rawhtml}
616 goto code
617 \begin{rawhtml} </A>\end{rawhtml}
618 }
619
620 \item Line 11,
621 \begin{verbatim}
622 diffKzT=1.E-2,
623 \end{verbatim}
624 this line sets the vertical diffusion coefficient for temperature
625 to $10^{-2}\,{\rm m^{2}s^{-1}}$. The boundary condition on this
626 operator is $\frac{\partial}{\partial z}$ = 0 on all boundaries.
627 The variable
628 {\bf
629 \begin{rawhtml} <A href=../code_reference/vdb/names/ZT.htm> \end{rawhtml}
630 diffKzT
631 \begin{rawhtml} </A>\end{rawhtml}
632 }
633 is read in the routine
634 {\it
635 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
636 INI\_PARMS
637 \begin{rawhtml} </A>\end{rawhtml}
638 }.
639 It is copied into model general vertical coordinate variable
640 {\bf
641 \begin{rawhtml} <A href=../code_reference/vdb/names/PD.htm> \end{rawhtml}
642 diffKrT
643 \begin{rawhtml} </A>\end{rawhtml}
644 } which is used in routine {\it S/R CALC\_DIFFUSIVITY}.
645
646 \fbox{ \begin{minipage}{5.0in}
647 {\it S/R CALC\_DIFFUSIVITY}({\it calc\_diffusivity.F})
648 \end{minipage}
649 }
650 {\bf
651 \begin{rawhtml} <A href=../code_reference/vdb/code/53.htm> \end{rawhtml}
652 goto code
653 \begin{rawhtml} </A>\end{rawhtml}
654 }
655
656
657
658 \item Line 13,
659 \begin{verbatim}
660 tAlpha=2.E-4,
661 \end{verbatim}
662 This line sets the thermal expansion coefficient for the fluid
663 to $2 \times 10^{-4}\,{\rm degrees}^{-1}$
664 The variable
665 {\bf
666 \begin{rawhtml} <A href=../code_reference/vdb/names/ZV.htm> \end{rawhtml}
667 tAlpha
668 \begin{rawhtml} </A>\end{rawhtml}
669 }
670 is read in the routine
671 {\it
672 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
673 INI\_PARMS
674 \begin{rawhtml} </A>\end{rawhtml}
675 }. The routine {\it S/R FIND\_RHO} makes use of {\bf tAlpha}.
676
677 \fbox{
678 \begin{minipage}{5.0in}
679 {\it S/R FIND\_RHO}({\it find\_rho.F})
680 \end{minipage}
681 }
682 {\bf
683 \begin{rawhtml} <A href=../code_reference/vdb/code/79.htm> \end{rawhtml}
684 goto code
685 \begin{rawhtml} </A>\end{rawhtml}
686 }
687
688 \item Line 18,
689 \begin{verbatim}
690 eosType='LINEAR'
691 \end{verbatim}
692 This line selects the linear form of the equation of state.
693 The variable
694 {\bf
695 \begin{rawhtml} <A href=../code_reference/vdb/names/WV.htm> \end{rawhtml}
696 eosType
697 \begin{rawhtml} </A>\end{rawhtml}
698 }
699 is read in the routine
700 {\it
701 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
702 INI\_PARMS
703 \begin{rawhtml} </A>\end{rawhtml}
704 }. The values of {\bf eosType} sets which formula in routine
705 {\it FIND\_RHO} is used to calculate density.
706
707 \fbox{
708 \begin{minipage}{5.0in}
709 {\it S/R FIND\_RHO}({\it find\_rho.F})
710 \end{minipage}
711 }
712 {\bf
713 \begin{rawhtml} <A href=../code_reference/vdb/code/79.htm> \end{rawhtml}
714 goto code
715 \begin{rawhtml} </A>\end{rawhtml}
716 }
717
718
719
720 \item Line 40,
721 \begin{verbatim}
722 usingSphericalPolarGrid=.TRUE.,
723 \end{verbatim}
724 This line requests that the simulation be performed in a
725 spherical polar coordinate system. It affects the interpretation of
726 grid input parameters, for example {\bf delX} and {\bf delY} and
727 causes the grid generation routines to initialize an internal grid based
728 on spherical polar geometry.
729 The variable
730 {\bf
731 \begin{rawhtml} <A href=../code_reference/vdb/names/10T.htm> \end{rawhtml}
732 usingSphericalPolarGrid
733 \begin{rawhtml} </A>\end{rawhtml}
734 }
735 is read in the routine
736 {\it
737 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
738 INI\_PARMS
739 \begin{rawhtml} </A>\end{rawhtml}
740 }. When set to {\bf .TRUE.} the settings of {\bf delX} and {\bf delY} are
741 taken to be in degrees. These values are used in the
742 routine {\it INI\_SPEHRICAL\_POLAR\_GRID}.
743
744 \fbox{
745 \begin{minipage}{5.0in}
746 {\it S/R INI\_SPEHRICAL\_POLAR\_GRID}({\it ini\_spherical\_polar\_grid.F})
747 \end{minipage}
748 }
749 {\bf
750 \begin{rawhtml} <A href=../code_reference/vdb/code/97.htm> \end{rawhtml}
751 goto code
752 \begin{rawhtml} </A>\end{rawhtml}
753 }
754
755 \item Line 41,
756 \begin{verbatim}
757 phiMin=0.,
758 \end{verbatim}
759 This line sets the southern boundary of the modeled
760 domain to $0^{\circ}$ latitude. This value affects both the
761 generation of the locally orthogonal grid that the model
762 uses internally and affects the initialization of the coriolis force.
763 Note - it is not required to set
764 a longitude boundary, since the absolute longitude does
765 not alter the kernel equation discretisation.
766 The variable
767 {\bf
768 \begin{rawhtml} <A href=../code_reference/vdb/names/110.htm> \end{rawhtml}
769 phiMin
770 \begin{rawhtml} </A>\end{rawhtml}
771 }
772 is read in the routine
773 {\it
774 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
775 INI\_PARMS
776 \begin{rawhtml} </A>\end{rawhtml}
777 } and is used in routine {\it INI\_SPEHRICAL\_POLAR\_GRID}.
778
779 \fbox{
780 \begin{minipage}{5.0in}
781 {\it S/R INI\_SPEHRICAL\_POLAR\_GRID}({\it ini\_spherical\_polar\_grid.F})
782 \end{minipage}
783 }
784 {\bf
785 \begin{rawhtml} <A href=../code_reference/vdb/code/97.htm> \end{rawhtml}
786 goto code
787 \begin{rawhtml} </A>\end{rawhtml}
788 }
789
790 \item Line 42,
791 \begin{verbatim}
792 delX=60*1.,
793 \end{verbatim}
794 This line sets the horizontal grid spacing between each y-coordinate line
795 in the discrete grid to $1^{\circ}$ in longitude.
796 The variable
797 {\bf
798 \begin{rawhtml} <A href=../code_reference/vdb/names/10Z.htm> \end{rawhtml}
799 delX
800 \begin{rawhtml} </A>\end{rawhtml}
801 }
802 is read in the routine
803 {\it
804 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
805 INI\_PARMS
806 \begin{rawhtml} </A>\end{rawhtml}
807 } and is used in routine {\it INI\_SPEHRICAL\_POLAR\_GRID}.
808
809 \fbox{
810 \begin{minipage}{5.0in}
811 {\it S/R INI\_SPEHRICAL\_POLAR\_GRID}({\it ini\_spherical\_polar\_grid.F})
812 \end{minipage}
813 }
814 {\bf
815 \begin{rawhtml} <A href=../code_reference/vdb/code/97.htm> \end{rawhtml}
816 goto code
817 \begin{rawhtml} </A>\end{rawhtml}
818 }
819
820 \item Line 43,
821 \begin{verbatim}
822 delY=60*1.,
823 \end{verbatim}
824 This line sets the horizontal grid spacing between each y-coordinate line
825 in the discrete grid to $1^{\circ}$ in latitude.
826 The variable
827 {\bf
828 \begin{rawhtml} <A href=../code_reference/vdb/names/UB.htm> \end{rawhtml}
829 delY
830 \begin{rawhtml} </A>\end{rawhtml}
831 }
832 is read in the routine
833 {\it
834 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
835 INI\_PARMS
836 \begin{rawhtml} </A>\end{rawhtml}
837 } and is used in routine {\it INI\_SPEHRICAL\_POLAR\_GRID}.
838
839 \fbox{
840 \begin{minipage}{5.0in}
841 {\it S/R INI\_SPEHRICAL\_POLAR\_GRID}({\it ini\_spherical\_polar\_grid.F})
842 \end{minipage}
843 }
844 {\bf
845 \begin{rawhtml} <A href=../code_reference/vdb/code/97.htm> \end{rawhtml}
846 goto code
847 \begin{rawhtml} </A>\end{rawhtml}
848 }
849
850 \item Line 44,
851 \begin{verbatim}
852 delZ=500.,500.,500.,500.,
853 \end{verbatim}
854 This line sets the vertical grid spacing between each z-coordinate line
855 in the discrete grid to $500\,{\rm m}$, so that the total model depth
856 is $2\,{\rm km}$.
857 The variable
858 {\bf
859 \begin{rawhtml} <A href=../code_reference/vdb/names/10W.htm> \end{rawhtml}
860 delZ
861 \begin{rawhtml} </A>\end{rawhtml}
862 }
863 is read in the routine
864 {\it
865 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
866 INI\_PARMS
867 \begin{rawhtml} </A>\end{rawhtml}
868 }.
869 It is copied into the internal
870 model coordinate variable
871 {\bf
872 \begin{rawhtml} <A href=../code_reference/vdb/names/10Y.htm> \end{rawhtml}
873 delR
874 \begin{rawhtml} </A>\end{rawhtml}
875 } which is used in routine {\it INI\_VERTICAL\_GRID}.
876
877 \fbox{
878 \begin{minipage}{5.0in}
879 {\it S/R INI\_VERTICAL\_GRID}({\it ini\_vertical\_grid.F})
880 \end{minipage}
881 }
882 {\bf
883 \begin{rawhtml} <A href=../code_reference/vdb/code/100.htm> \end{rawhtml}
884 goto code
885 \begin{rawhtml} </A>\end{rawhtml}
886 }
887
888 \item Line 47,
889 \begin{verbatim}
890 bathyFile='topog.box'
891 \end{verbatim}
892 This line specifies the name of the file from which the domain
893 bathymetry is read. This file is a two-dimensional ($x,y$) map of
894 depths. This file is assumed to contain 64-bit binary numbers
895 giving the depth of the model at each grid cell, ordered with the x
896 coordinate varying fastest. The points are ordered from low coordinate
897 to high coordinate for both axes. The units and orientation of the
898 depths in this file are the same as used in the MITgcm code. In this
899 experiment, a depth of $0m$ indicates a solid wall and a depth
900 of $-2000m$ indicates open ocean. The matlab program
901 {\it input/gendata.m} shows an example of how to generate a
902 bathymetry file.
903 The variable
904 {\bf
905 \begin{rawhtml} <A href=../code_reference/vdb/names/179.htm> \end{rawhtml}
906 bathyFile
907 \begin{rawhtml} </A>\end{rawhtml}
908 }
909 is read in the routine
910 {\it
911 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
912 INI\_PARMS
913 \begin{rawhtml} </A>\end{rawhtml}
914 }. The bathymetry file is read in the routine {\it INI\_DEPTHS}.
915
916 \fbox{
917 \begin{minipage}{5.0in}
918 {\it S/R INI\_DEPTHS}({\it ini\_depths.F})
919 \end{minipage}
920 }
921 {\bf
922 \begin{rawhtml} <A href=../code_reference/vdb/code/88.htm> \end{rawhtml}
923 goto code
924 \begin{rawhtml} </A>\end{rawhtml}
925 }
926
927
928 \item Line 50,
929 \begin{verbatim}
930 zonalWindFile='windx.sin_y'
931 \end{verbatim}
932 This line specifies the name of the file from which the x-direction
933 (zonal) surface wind stress is read. This file is also a two-dimensional
934 ($x,y$) map and is enumerated and formatted in the same manner as the
935 bathymetry file. The matlab program {\it input/gendata.m} includes example
936 code to generate a valid
937 {\bf zonalWindFile}
938 file.
939 The variable
940 {\bf
941 \begin{rawhtml} <A href=../code_reference/vdb/names/13W.htm> \end{rawhtml}
942 zonalWindFile
943 \begin{rawhtml} </A>\end{rawhtml}
944 }
945 is read in the routine
946 {\it
947 \begin{rawhtml} <A href=../code_reference/vdb/code/94.htm> \end{rawhtml}
948 INI\_PARMS
949 \begin{rawhtml} </A>\end{rawhtml}
950 }. The wind-stress file is read in the routine
951 {\it EXTERNAL\_FIELDS\_LOAD}.
952
953 \fbox{
954 \begin{minipage}{5.0in}
955 {\it S/R EXTERNAL\_FIELDS\_LOAD}({\it external\_fields\_load.F})
956 \end{minipage}
957 }
958 {\bf
959 \begin{rawhtml} <A href=../code_reference/vdb/code/75.htm> \end{rawhtml}
960 goto code
961 \begin{rawhtml} </A>\end{rawhtml}
962 }
963
964 \end{itemize}
965
966 \noindent other lines in the file {\it input/data} are standard values.
967
968 \begin{rawhtml}<PRE>\end{rawhtml}
969 \begin{small}
970 \input{part3/case_studies/fourlayer_gyre/input/data}
971 \end{small}
972 \begin{rawhtml}</PRE>\end{rawhtml}
973
974 \subsubsection{File {\it input/data.pkg}}
975 \label{www:tutorials}
976
977 This file uses standard default values and does not contain
978 customisations for this experiment.
979
980 \subsubsection{File {\it input/eedata}}
981 \label{www:tutorials}
982
983 This file uses standard default values and does not contain
984 customisations for this experiment.
985
986 \subsubsection{File {\it input/windx.sin\_y}}
987 \label{www:tutorials}
988
989 The {\it input/windx.sin\_y} file specifies a two-dimensional ($x,y$)
990 map of wind stress ,$\tau_{x}$, values. The units used are $Nm^{-2}$ (the
991 default for MITgcm).
992 Although $\tau_{x}$ is only a function of latitude, $y$,
993 in this experiment
994 this file must still define a complete two-dimensional map in order
995 to be compatible with the standard code for loading forcing fields
996 in MITgcm (routine {\it EXTERNAL\_FIELDS\_LOAD}.
997 The included matlab program {\it input/gendata.m} gives a complete
998 code for creating the {\it input/windx.sin\_y} file.
999
1000 \subsubsection{File {\it input/topog.box}}
1001 \label{www:tutorials}
1002
1003
1004 The {\it input/topog.box} file specifies a two-dimensional ($x,y$)
1005 map of depth values. For this experiment values are either
1006 $0~{\rm m}$ or $-2000\,{\rm m}$, corresponding respectively to a wall or to deep
1007 ocean. The file contains a raw binary stream of data that is enumerated
1008 in the same way as standard MITgcm two-dimensional, horizontal arrays.
1009 The included matlab program {\it input/gendata.m} gives a complete
1010 code for creating the {\it input/topog.box} file.
1011
1012 \subsubsection{File {\it code/SIZE.h}}
1013 \label{www:tutorials}
1014
1015 Two lines are customized in this file for the current experiment
1016
1017 \begin{itemize}
1018
1019 \item Line 39,
1020 \begin{verbatim} sNx=60, \end{verbatim} this line sets
1021 the lateral domain extent in grid points for the
1022 axis aligned with the x-coordinate.
1023
1024 \item Line 40,
1025 \begin{verbatim} sNy=60, \end{verbatim} this line sets
1026 the lateral domain extent in grid points for the
1027 axis aligned with the y-coordinate.
1028
1029 \item Line 49,
1030 \begin{verbatim} Nr=4, \end{verbatim} this line sets
1031 the vertical domain extent in grid points.
1032
1033 \end{itemize}
1034
1035 \begin{small}
1036 \include{part3/case_studies/fourlayer_gyre/code/SIZE.h}
1037 \end{small}
1038
1039 \subsubsection{File {\it code/CPP\_OPTIONS.h}}
1040 \label{www:tutorials}
1041
1042 This file uses standard default values and does not contain
1043 customisations for this experiment.
1044
1045
1046 \subsubsection{File {\it code/CPP\_EEOPTIONS.h}}
1047 \label{www:tutorials}
1048
1049 This file uses standard default values and does not contain
1050 customisations for this experiment.
1051
1052 \subsubsection{Other Files }
1053 \label{www:tutorials}
1054
1055 Other files relevant to this experiment are
1056 \begin{itemize}
1057 \item {\it model/src/ini\_cori.F}. This file initializes the model
1058 coriolis variables {\bf fCorU} and {\bf fCorV}.
1059 \item {\it model/src/ini\_spherical\_polar\_grid.F} This file
1060 initializes the model grid discretisation variables {\bf
1061 dxF, dyF, dxG, dyG, dxC, dyC}.
1062 \item {\it model/src/ini\_parms.F}.
1063 \end{itemize}
1064
1065 \subsection{Running The Example}
1066 \label{www:tutorials}
1067 \label{SEC:running_the_example}
1068
1069 \subsubsection{Code Download}
1070 \label{www:tutorials}
1071
1072 In order to run the examples you must first download the code distribution.
1073 Instructions for downloading the code can be found in section
1074 \ref{sect:obtainingCode}.
1075
1076 \subsubsection{Experiment Location}
1077 \label{www:tutorials}
1078
1079 This example experiments is located under the release sub-directory
1080
1081 \vspace{5mm}
1082 {\it verification/exp2/ }
1083
1084 \subsubsection{Running the Experiment}
1085 \label{www:tutorials}
1086
1087 To run the experiment
1088
1089 \begin{enumerate}
1090 \item Set the current directory to {\it input/ }
1091
1092 \begin{verbatim}
1093 % cd input
1094 \end{verbatim}
1095
1096 \item Verify that current directory is now correct
1097
1098 \begin{verbatim}
1099 % pwd
1100 \end{verbatim}
1101
1102 You should see a response on the screen ending in
1103
1104 {\it verification/exp2/input }
1105
1106
1107 \item Run the genmake script to create the experiment {\it Makefile}
1108
1109 \begin{verbatim}
1110 % ../../../tools/genmake -mods=../code
1111 \end{verbatim}
1112
1113 \item Create a list of header file dependencies in {\it Makefile}
1114
1115 \begin{verbatim}
1116 % make depend
1117 \end{verbatim}
1118
1119 \item Build the executable file.
1120
1121 \begin{verbatim}
1122 % make
1123 \end{verbatim}
1124
1125 \item Run the {\it mitgcmuv} executable
1126
1127 \begin{verbatim}
1128 % ./mitgcmuv
1129 \end{verbatim}
1130
1131 \end{enumerate}
1132
1133

  ViewVC Help
Powered by ViewVC 1.1.22