1 |
adcroft |
1.9 |
% $Header: /u/gcmpack/mitgcmdoc/part4/sarch.tex,v 1.8 2002/04/24 09:00:53 heimbach Exp $ |
2 |
cnh |
1.1 |
|
3 |
cnh |
1.7 |
This chapter focuses on describing the {\bf WRAPPER} environment within which |
4 |
|
|
both the core numerics and the pluggable packages operate. The description |
5 |
|
|
presented here is intended to be a detailed exposition and contains significant |
6 |
|
|
background material, as well as advanced details on working with the WRAPPER. |
7 |
adcroft |
1.9 |
The tutorial sections of this manual (see sections |
8 |
|
|
\ref{sect:tutorials} and \ref{sect:tutorialIII}) |
9 |
cnh |
1.7 |
contain more succinct, step-by-step instructions on running basic numerical |
10 |
|
|
experiments, of varous types, both sequentially and in parallel. For many |
11 |
|
|
projects simply starting from an example code and adapting it to suit a |
12 |
|
|
particular situation |
13 |
|
|
will be all that is required. |
14 |
|
|
The first part of this chapter discusses the MITgcm architecture at an |
15 |
|
|
abstract level. In the second part of the chapter we described practical |
16 |
|
|
details of the MITgcm implementation and of current tools and operating system |
17 |
|
|
features that are employed. |
18 |
cnh |
1.1 |
|
19 |
|
|
\section{Overall architectural goals} |
20 |
|
|
|
21 |
|
|
Broadly, the goals of the software architecture employed in MITgcm are |
22 |
|
|
three-fold |
23 |
|
|
|
24 |
|
|
\begin{itemize} |
25 |
|
|
\item We wish to be able to study a very broad range |
26 |
|
|
of interesting and challenging rotating fluids problems. |
27 |
|
|
\item We wish the model code to be readily targeted to |
28 |
|
|
a wide range of platforms |
29 |
|
|
\item On any given platform we would like to be |
30 |
|
|
able to achieve performance comparable to an implementation |
31 |
|
|
developed and specialized specifically for that platform. |
32 |
|
|
\end{itemize} |
33 |
|
|
|
34 |
|
|
These points are summarized in figure \ref{fig:mitgcm_architecture_goals} |
35 |
|
|
which conveys the goals of the MITgcm design. The goals lead to |
36 |
|
|
a software architecture which at the high-level can be viewed as consisting |
37 |
|
|
of |
38 |
|
|
|
39 |
|
|
\begin{enumerate} |
40 |
|
|
\item A core set of numerical and support code. This is discussed in detail in |
41 |
adcroft |
1.6 |
section \ref{sect:partII}. |
42 |
cnh |
1.1 |
\item A scheme for supporting optional "pluggable" {\bf packages} (containing |
43 |
|
|
for example mixed-layer schemes, biogeochemical schemes, atmospheric physics). |
44 |
|
|
These packages are used both to overlay alternate dynamics and to introduce |
45 |
|
|
specialized physical content onto the core numerical code. An overview of |
46 |
|
|
the {\bf package} scheme is given at the start of part \ref{part:packages}. |
47 |
|
|
\item A support framework called {\bf WRAPPER} (Wrappable Application Parallel |
48 |
|
|
Programming Environment Resource), within which the core numerics and pluggable |
49 |
|
|
packages operate. |
50 |
|
|
\end{enumerate} |
51 |
|
|
|
52 |
|
|
This chapter focuses on describing the {\bf WRAPPER} environment under which |
53 |
|
|
both the core numerics and the pluggable packages function. The description |
54 |
cnh |
1.4 |
presented here is intended to be a detailed exposition and contains significant |
55 |
cnh |
1.1 |
background material, as well as advanced details on working with the WRAPPER. |
56 |
|
|
The examples section of this manual (part \ref{part:example}) contains more |
57 |
|
|
succinct, step-by-step instructions on running basic numerical |
58 |
|
|
experiments both sequentially and in parallel. For many projects simply |
59 |
|
|
starting from an example code and adapting it to suit a particular situation |
60 |
|
|
will be all that is required. |
61 |
|
|
|
62 |
adcroft |
1.2 |
|
63 |
cnh |
1.1 |
\begin{figure} |
64 |
|
|
\begin{center} |
65 |
adcroft |
1.2 |
\resizebox{!}{2.5in}{\includegraphics{part4/mitgcm_goals.eps}} |
66 |
cnh |
1.1 |
\end{center} |
67 |
adcroft |
1.2 |
\caption{ |
68 |
|
|
The MITgcm architecture is designed to allow simulation of a wide |
69 |
cnh |
1.1 |
range of physical problems on a wide range of hardware. The computational |
70 |
|
|
resource requirements of the applications targeted range from around |
71 |
|
|
$10^7$ bytes ( $\approx 10$ megabytes ) of memory to $10^{11}$ bytes |
72 |
|
|
( $\approx 100$ gigabytes). Arithmetic operation counts for the applications of |
73 |
|
|
interest range from $10^{9}$ floating point operations to more than $10^{17}$ |
74 |
adcroft |
1.2 |
floating point operations.} |
75 |
|
|
\label{fig:mitgcm_architecture_goals} |
76 |
cnh |
1.1 |
\end{figure} |
77 |
|
|
|
78 |
|
|
\section{WRAPPER} |
79 |
|
|
|
80 |
|
|
A significant element of the software architecture utilized in |
81 |
|
|
MITgcm is a software superstructure and substructure collectively |
82 |
|
|
called the WRAPPER (Wrappable Application Parallel Programming |
83 |
|
|
Environment Resource). All numerical and support code in MITgcm is written |
84 |
|
|
to ``fit'' within the WRAPPER infrastructure. Writing code to ``fit'' within |
85 |
|
|
the WRAPPER means that coding has to follow certain, relatively |
86 |
|
|
straightforward, rules and conventions ( these are discussed further in |
87 |
adcroft |
1.6 |
section \ref{sect:specifying_a_decomposition} ). |
88 |
cnh |
1.1 |
|
89 |
|
|
The approach taken by the WRAPPER is illustrated in figure |
90 |
|
|
\ref{fig:fit_in_wrapper} which shows how the WRAPPER serves to insulate code |
91 |
|
|
that fits within it from architectural differences between hardware platforms |
92 |
|
|
and operating systems. This allows numerical code to be easily retargetted. |
93 |
adcroft |
1.2 |
|
94 |
|
|
|
95 |
cnh |
1.1 |
\begin{figure} |
96 |
|
|
\begin{center} |
97 |
adcroft |
1.2 |
\resizebox{!}{4.5in}{\includegraphics{part4/fit_in_wrapper.eps}} |
98 |
cnh |
1.1 |
\end{center} |
99 |
adcroft |
1.2 |
\caption{ |
100 |
|
|
Numerical code is written too fit within a software support |
101 |
cnh |
1.1 |
infrastructure called WRAPPER. The WRAPPER is portable and |
102 |
cnh |
1.4 |
can be specialized for a wide range of specific target hardware and |
103 |
cnh |
1.1 |
programming environments, without impacting numerical code that fits |
104 |
|
|
within the WRAPPER. Codes that fit within the WRAPPER can generally be |
105 |
|
|
made to run as fast on a particular platform as codes specially |
106 |
adcroft |
1.2 |
optimized for that platform.} |
107 |
|
|
\label{fig:fit_in_wrapper} |
108 |
cnh |
1.1 |
\end{figure} |
109 |
|
|
|
110 |
|
|
\subsection{Target hardware} |
111 |
adcroft |
1.6 |
\label{sect:target_hardware} |
112 |
cnh |
1.1 |
|
113 |
|
|
The WRAPPER is designed to target as broad as possible a range of computer |
114 |
|
|
systems. The original development of the WRAPPER took place on a |
115 |
|
|
multi-processor, CRAY Y-MP system. On that system, numerical code performance |
116 |
|
|
and scaling under the WRAPPER was in excess of that of an implementation that |
117 |
|
|
was tightly bound to the CRAY systems proprietary multi-tasking and |
118 |
|
|
micro-tasking approach. Later developments have been carried out on |
119 |
|
|
uniprocessor and multi-processor Sun systems with both uniform memory access |
120 |
|
|
(UMA) and non-uniform memory access (NUMA) designs. Significant work has also |
121 |
|
|
been undertaken on x86 cluster systems, Alpha processor based clustered SMP |
122 |
|
|
systems, and on cache-coherent NUMA (CC-NUMA) systems from Silicon Graphics. |
123 |
|
|
The MITgcm code, operating within the WRAPPER, is also used routinely used on |
124 |
|
|
large scale MPP systems (for example T3E systems and IBM SP systems). In all |
125 |
|
|
cases numerical code, operating within the WRAPPER, performs and scales very |
126 |
|
|
competitively with equivalent numerical code that has been modified to contain |
127 |
|
|
native optimizations for a particular system \ref{ref hoe and hill, ecmwf}. |
128 |
|
|
|
129 |
|
|
\subsection{Supporting hardware neutrality} |
130 |
|
|
|
131 |
adcroft |
1.6 |
The different systems listed in section \ref{sect:target_hardware} can be |
132 |
cnh |
1.1 |
categorized in many different ways. For example, one common distinction is |
133 |
|
|
between shared-memory parallel systems (SMP's, PVP's) and distributed memory |
134 |
|
|
parallel systems (for example x86 clusters and large MPP systems). This is one |
135 |
|
|
example of a difference between compute platforms that can impact an |
136 |
|
|
application. Another common distinction is between vector processing systems |
137 |
|
|
with highly specialized CPU's and memory subsystems and commodity |
138 |
|
|
microprocessor based systems. There are numerous other differences, especially |
139 |
|
|
in relation to how parallel execution is supported. To capture the essential |
140 |
|
|
differences between different platforms the WRAPPER uses a {\it machine model}. |
141 |
|
|
|
142 |
|
|
\subsection{WRAPPER machine model} |
143 |
|
|
|
144 |
|
|
Applications using the WRAPPER are not written to target just one |
145 |
|
|
particular machine (for example an IBM SP2) or just one particular family or |
146 |
|
|
class of machines (for example Parallel Vector Processor Systems). Instead the |
147 |
|
|
WRAPPER provides applications with an |
148 |
|
|
abstract {\it machine model}. The machine model is very general, however, it can |
149 |
cnh |
1.5 |
easily be specialized to fit, in a computationally efficient manner, any |
150 |
cnh |
1.1 |
computer architecture currently available to the scientific computing community. |
151 |
|
|
|
152 |
|
|
\subsection{Machine model parallelism} |
153 |
|
|
|
154 |
|
|
Codes operating under the WRAPPER target an abstract machine that is assumed to |
155 |
|
|
consist of one or more logical processors that can compute concurrently. |
156 |
cnh |
1.5 |
Computational work is divided among the logical |
157 |
cnh |
1.1 |
processors by allocating ``ownership'' to |
158 |
|
|
each processor of a certain set (or sets) of calculations. Each set of |
159 |
|
|
calculations owned by a particular processor is associated with a specific |
160 |
|
|
region of the physical space that is being simulated, only one processor will |
161 |
|
|
be associated with each such region (domain decomposition). |
162 |
|
|
|
163 |
|
|
In a strict sense the logical processors over which work is divided do not need |
164 |
|
|
to correspond to physical processors. It is perfectly possible to execute a |
165 |
|
|
configuration decomposed for multiple logical processors on a single physical |
166 |
|
|
processor. This helps ensure that numerical code that is written to fit |
167 |
|
|
within the WRAPPER will parallelize with no additional effort and is |
168 |
|
|
also useful when debugging codes. Generally, however, |
169 |
|
|
the computational domain will be subdivided over multiple logical |
170 |
|
|
processors in order to then bind those logical processors to physical |
171 |
|
|
processor resources that can compute in parallel. |
172 |
|
|
|
173 |
|
|
\subsubsection{Tiles} |
174 |
|
|
|
175 |
|
|
Computationally, associated with each region of physical |
176 |
|
|
space allocated to a particular logical processor, there will be data |
177 |
|
|
structures (arrays, scalar variables etc...) that hold the simulated state of |
178 |
|
|
that region. We refer to these data structures as being {\bf owned} by the |
179 |
cnh |
1.4 |
processor to which their |
180 |
cnh |
1.1 |
associated region of physical space has been allocated. Individual |
181 |
|
|
regions that are allocated to processors are called {\bf tiles}. A |
182 |
|
|
processor can own more |
183 |
|
|
than one tile. Figure \ref{fig:domaindecomp} shows a physical domain being |
184 |
|
|
mapped to a set of logical processors, with each processors owning a single |
185 |
|
|
region of the domain (a single tile). Except for periods of |
186 |
|
|
communication and coordination, each processor computes autonomously, working |
187 |
|
|
only with data from the tile (or tiles) that the processor owns. When multiple |
188 |
|
|
tiles are alloted to a single processor, each tile is computed on |
189 |
|
|
independently of the other tiles, in a sequential fashion. |
190 |
|
|
|
191 |
|
|
\begin{figure} |
192 |
|
|
\begin{center} |
193 |
adcroft |
1.2 |
\resizebox{5in}{!}{ |
194 |
|
|
\includegraphics{part4/domain_decomp.eps} |
195 |
cnh |
1.1 |
} |
196 |
|
|
\end{center} |
197 |
|
|
\caption{ The WRAPPER provides support for one and two dimensional |
198 |
|
|
decompositions of grid-point domains. The figure shows a hypothetical domain of |
199 |
|
|
total size $N_{x}N_{y}N_{z}$. This hypothetical domain is decomposed in |
200 |
|
|
two-dimensions along the $N_{x}$ and $N_{y}$ directions. The resulting {\bf |
201 |
|
|
tiles} are {\bf owned} by different processors. The {\bf owning} |
202 |
|
|
processors perform the |
203 |
|
|
arithmetic operations associated with a {\bf tile}. Although not illustrated |
204 |
|
|
here, a single processor can {\bf own} several {\bf tiles}. |
205 |
|
|
Whenever a processor wishes to transfer data between tiles or |
206 |
|
|
communicate with other processors it calls a WRAPPER supplied |
207 |
|
|
function. |
208 |
|
|
} \label{fig:domaindecomp} |
209 |
|
|
\end{figure} |
210 |
|
|
|
211 |
|
|
|
212 |
|
|
\subsubsection{Tile layout} |
213 |
|
|
|
214 |
|
|
Tiles consist of an interior region and an overlap region. The overlap region |
215 |
|
|
of a tile corresponds to the interior region of an adjacent tile. |
216 |
|
|
In figure \ref{fig:tiledworld} each tile would own the region |
217 |
|
|
within the black square and hold duplicate information for overlap |
218 |
|
|
regions extending into the tiles to the north, south, east and west. |
219 |
|
|
During |
220 |
|
|
computational phases a processor will reference data in an overlap region |
221 |
|
|
whenever it requires values that outside the domain it owns. Periodically |
222 |
|
|
processors will make calls to WRAPPER functions to communicate data between |
223 |
|
|
tiles, in order to keep the overlap regions up to date (see section |
224 |
adcroft |
1.6 |
\ref{sect:communication_primitives}). The WRAPPER functions can use a |
225 |
cnh |
1.1 |
variety of different mechanisms to communicate data between tiles. |
226 |
|
|
|
227 |
|
|
\begin{figure} |
228 |
|
|
\begin{center} |
229 |
adcroft |
1.2 |
\resizebox{5in}{!}{ |
230 |
|
|
\includegraphics{part4/tiled-world.eps} |
231 |
cnh |
1.1 |
} |
232 |
|
|
\end{center} |
233 |
|
|
\caption{ A global grid subdivided into tiles. |
234 |
|
|
Tiles contain a interior region and an overlap region. |
235 |
|
|
Overlap regions are periodically updated from neighboring tiles. |
236 |
|
|
} \label{fig:tiledworld} |
237 |
|
|
\end{figure} |
238 |
|
|
|
239 |
|
|
\subsection{Communication mechanisms} |
240 |
|
|
|
241 |
|
|
Logical processors are assumed to be able to exchange information |
242 |
|
|
between tiles and between each other using at least one of two possible |
243 |
|
|
mechanisms. |
244 |
|
|
|
245 |
|
|
\begin{itemize} |
246 |
|
|
\item {\bf Shared memory communication}. |
247 |
|
|
Under this mode of communication data transfers are assumed to be possible |
248 |
|
|
using direct addressing of regions of memory. In this case a CPU is able to read |
249 |
|
|
(and write) directly to regions of memory "owned" by another CPU |
250 |
|
|
using simple programming language level assignment operations of the |
251 |
|
|
the sort shown in figure \ref{fig:simple_assign}. In this way one CPU |
252 |
|
|
(CPU1 in the figure) can communicate information to another CPU (CPU2 in the |
253 |
|
|
figure) by assigning a particular value to a particular memory location. |
254 |
|
|
|
255 |
|
|
\item {\bf Distributed memory communication}. |
256 |
|
|
Under this mode of communication there is no mechanism, at the application code level, |
257 |
|
|
for directly addressing regions of memory owned and visible to another CPU. Instead |
258 |
|
|
a communication library must be used as illustrated in figure |
259 |
|
|
\ref{fig:comm_msg}. In this case CPU's must call a function in the API of the |
260 |
|
|
communication library to communicate data from a tile that it owns to a tile |
261 |
|
|
that another CPU owns. By default the WRAPPER binds to the MPI communication |
262 |
|
|
library \ref{MPI} for this style of communication. |
263 |
|
|
\end{itemize} |
264 |
|
|
|
265 |
|
|
The WRAPPER assumes that communication will use one of these two styles |
266 |
|
|
of communication. The underlying hardware and operating system support |
267 |
|
|
for the style used is not specified and can vary from system to system. |
268 |
|
|
|
269 |
|
|
\begin{figure} |
270 |
|
|
\begin{verbatim} |
271 |
|
|
|
272 |
|
|
CPU1 | CPU2 |
273 |
|
|
==== | ==== |
274 |
|
|
| |
275 |
|
|
a(3) = 8 | WHILE ( a(3) .NE. 8 ) |
276 |
|
|
| WAIT |
277 |
|
|
| END WHILE |
278 |
|
|
| |
279 |
|
|
\end{verbatim} |
280 |
|
|
\caption{ In the WRAPPER shared memory communication model, simple writes to an |
281 |
|
|
array can be made to be visible to other CPU's at the application code level. |
282 |
|
|
So that for example, if one CPU (CPU1 in the figure above) writes the value $8$ to |
283 |
|
|
element $3$ of array $a$, then other CPU's (for example CPU2 in the figure above) |
284 |
|
|
will be able to see the value $8$ when they read from $a(3)$. |
285 |
|
|
This provides a very low latency and high bandwidth communication |
286 |
|
|
mechanism. |
287 |
|
|
} \label{fig:simple_assign} |
288 |
|
|
\end{figure} |
289 |
|
|
|
290 |
|
|
\begin{figure} |
291 |
|
|
\begin{verbatim} |
292 |
|
|
|
293 |
|
|
CPU1 | CPU2 |
294 |
|
|
==== | ==== |
295 |
|
|
| |
296 |
|
|
a(3) = 8 | WHILE ( a(3) .NE. 8 ) |
297 |
|
|
CALL SEND( CPU2,a(3) ) | CALL RECV( CPU1, a(3) ) |
298 |
|
|
| END WHILE |
299 |
|
|
| |
300 |
|
|
\end{verbatim} |
301 |
|
|
\caption{ In the WRAPPER distributed memory communication model |
302 |
|
|
data can not be made directly visible to other CPU's. |
303 |
|
|
If one CPU writes the value $8$ to element $3$ of array $a$, then |
304 |
|
|
at least one of CPU1 and/or CPU2 in the figure above will need |
305 |
|
|
to call a bespoke communication library in order for the updated |
306 |
|
|
value to be communicated between CPU's. |
307 |
|
|
} \label{fig:comm_msg} |
308 |
|
|
\end{figure} |
309 |
|
|
|
310 |
|
|
\subsection{Shared memory communication} |
311 |
adcroft |
1.6 |
\label{sect:shared_memory_communication} |
312 |
cnh |
1.1 |
|
313 |
|
|
Under shared communication independent CPU's are operating |
314 |
|
|
on the exact same global address space at the application level. |
315 |
|
|
This means that CPU 1 can directly write into global |
316 |
|
|
data structures that CPU 2 ``owns'' using a simple |
317 |
|
|
assignment at the application level. |
318 |
|
|
This is the model of memory access is supported at the basic system |
319 |
|
|
design level in ``shared-memory'' systems such as PVP systems, SMP systems, |
320 |
|
|
and on distributed shared memory systems (the SGI Origin). |
321 |
|
|
On such systems the WRAPPER will generally use simple read and write statements |
322 |
|
|
to access directly application data structures when communicating between CPU's. |
323 |
|
|
|
324 |
|
|
In a system where assignments statements, like the one in figure |
325 |
|
|
\ref{fig:simple_assign} map directly to |
326 |
|
|
hardware instructions that transport data between CPU and memory banks, this |
327 |
|
|
can be a very efficient mechanism for communication. In this case two CPU's, |
328 |
|
|
CPU1 and CPU2, can communicate simply be reading and writing to an |
329 |
|
|
agreed location and following a few basic rules. The latency of this sort |
330 |
|
|
of communication is generally not that much higher than the hardware |
331 |
|
|
latency of other memory accesses on the system. The bandwidth available |
332 |
|
|
between CPU's communicating in this way can be close to the bandwidth of |
333 |
|
|
the systems main-memory interconnect. This can make this method of |
334 |
|
|
communication very efficient provided it is used appropriately. |
335 |
|
|
|
336 |
|
|
\subsubsection{Memory consistency} |
337 |
adcroft |
1.6 |
\label{sect:memory_consistency} |
338 |
cnh |
1.1 |
|
339 |
|
|
When using shared memory communication between |
340 |
|
|
multiple processors the WRAPPER level shields user applications from |
341 |
|
|
certain counter-intuitive system behaviors. In particular, one issue the |
342 |
|
|
WRAPPER layer must deal with is a systems memory model. In general the order |
343 |
|
|
of reads and writes expressed by the textual order of an application code may |
344 |
|
|
not be the ordering of instructions executed by the processor performing the |
345 |
|
|
application. The processor performing the application instructions will always |
346 |
|
|
operate so that, for the application instructions the processor is executing, |
347 |
|
|
any reordering is not apparent. However, in general machines are often |
348 |
|
|
designed so that reordering of instructions is not hidden from other second |
349 |
|
|
processors. This means that, in general, even on a shared memory system two |
350 |
|
|
processors can observe inconsistent memory values. |
351 |
|
|
|
352 |
|
|
The issue of memory consistency between multiple processors is discussed at |
353 |
|
|
length in many computer science papers, however, from a practical point of |
354 |
|
|
view, in order to deal with this issue, shared memory machines all provide |
355 |
|
|
some mechanism to enforce memory consistency when it is needed. The exact |
356 |
|
|
mechanism employed will vary between systems. For communication using shared |
357 |
|
|
memory, the WRAPPER provides a place to invoke the appropriate mechanism to |
358 |
|
|
ensure memory consistency for a particular platform. |
359 |
|
|
|
360 |
|
|
\subsubsection{Cache effects and false sharing} |
361 |
adcroft |
1.6 |
\label{sect:cache_effects_and_false_sharing} |
362 |
cnh |
1.1 |
|
363 |
|
|
Shared-memory machines often have local to processor memory caches |
364 |
|
|
which contain mirrored copies of main memory. Automatic cache-coherence |
365 |
|
|
protocols are used to maintain consistency between caches on different |
366 |
|
|
processors. These cache-coherence protocols typically enforce consistency |
367 |
|
|
between regions of memory with large granularity (typically 128 or 256 byte |
368 |
|
|
chunks). The coherency protocols employed can be expensive relative to other |
369 |
|
|
memory accesses and so care is taken in the WRAPPER (by padding synchronization |
370 |
|
|
structures appropriately) to avoid unnecessary coherence traffic. |
371 |
|
|
|
372 |
|
|
\subsubsection{Operating system support for shared memory.} |
373 |
|
|
|
374 |
|
|
Applications running under multiple threads within a single process can |
375 |
|
|
use shared memory communication. In this case {\it all} the memory locations |
376 |
|
|
in an application are potentially visible to all the compute threads. Multiple |
377 |
|
|
threads operating within a single process is the standard mechanism for |
378 |
|
|
supporting shared memory that the WRAPPER utilizes. Configuring and launching |
379 |
|
|
code to run in multi-threaded mode on specific platforms is discussed in |
380 |
adcroft |
1.6 |
section \ref{sect:running_with_threads}. However, on many systems, potentially |
381 |
cnh |
1.1 |
very efficient mechanisms for using shared memory communication between |
382 |
|
|
multiple processes (in contrast to multiple threads within a single |
383 |
|
|
process) also exist. In most cases this works by making a limited region of |
384 |
|
|
memory shared between processes. The MMAP \ref{magicgarden} and |
385 |
|
|
IPC \ref{magicgarden} facilities in UNIX systems provide this capability as do |
386 |
|
|
vendor specific tools like LAPI \ref{IBMLAPI} and IMC \ref{Memorychannel}. |
387 |
|
|
Extensions exist for the WRAPPER that allow these mechanisms |
388 |
|
|
to be used for shared memory communication. However, these mechanisms are not |
389 |
|
|
distributed with the default WRAPPER sources, because of their proprietary |
390 |
|
|
nature. |
391 |
|
|
|
392 |
|
|
\subsection{Distributed memory communication} |
393 |
adcroft |
1.6 |
\label{sect:distributed_memory_communication} |
394 |
cnh |
1.1 |
Many parallel systems are not constructed in a way where it is |
395 |
|
|
possible or practical for an application to use shared memory |
396 |
|
|
for communication. For example cluster systems consist of individual computers |
397 |
|
|
connected by a fast network. On such systems their is no notion of shared memory |
398 |
|
|
at the system level. For this sort of system the WRAPPER provides support |
399 |
|
|
for communication based on a bespoke communication library |
400 |
|
|
(see figure \ref{fig:comm_msg}). The default communication library used is MPI |
401 |
|
|
\ref{mpi}. However, it is relatively straightforward to implement bindings to |
402 |
|
|
optimized platform specific communication libraries. For example the work |
403 |
|
|
described in \ref{hoe-hill:99} substituted standard MPI communication for a |
404 |
|
|
highly optimized library. |
405 |
|
|
|
406 |
|
|
\subsection{Communication primitives} |
407 |
adcroft |
1.6 |
\label{sect:communication_primitives} |
408 |
cnh |
1.1 |
|
409 |
|
|
\begin{figure} |
410 |
|
|
\begin{center} |
411 |
adcroft |
1.2 |
\resizebox{5in}{!}{ |
412 |
|
|
\includegraphics{part4/comm-primm.eps} |
413 |
cnh |
1.1 |
} |
414 |
|
|
\end{center} |
415 |
cnh |
1.5 |
\caption{Three performance critical parallel primitives are provided |
416 |
|
|
by the WRAPPER. These primitives are always used to communicate data |
417 |
cnh |
1.1 |
between tiles. The figure shows four tiles. The curved arrows indicate |
418 |
|
|
exchange primitives which transfer data between the overlap regions at tile |
419 |
|
|
edges and interior regions for nearest-neighbor tiles. |
420 |
|
|
The straight arrows symbolize global sum operations which connect all tiles. |
421 |
|
|
The global sum operation provides both a key arithmetic primitive and can |
422 |
|
|
serve as a synchronization primitive. A third barrier primitive is also |
423 |
|
|
provided, it behaves much like the global sum primitive. |
424 |
|
|
} \label{fig:communication_primitives} |
425 |
|
|
\end{figure} |
426 |
|
|
|
427 |
|
|
|
428 |
|
|
Optimized communication support is assumed to be possibly available |
429 |
|
|
for a small number of communication operations. |
430 |
|
|
It is assumed that communication performance optimizations can |
431 |
|
|
be achieved by optimizing a small number of communication primitives. |
432 |
|
|
Three optimizable primitives are provided by the WRAPPER |
433 |
|
|
\begin{itemize} |
434 |
|
|
\item{\bf EXCHANGE} This operation is used to transfer data between interior |
435 |
|
|
and overlap regions of neighboring tiles. A number of different forms of this |
436 |
|
|
operation are supported. These different forms handle |
437 |
|
|
\begin{itemize} |
438 |
|
|
\item Data type differences. Sixty-four bit and thirty-two bit fields may be handled |
439 |
|
|
separately. |
440 |
|
|
\item Bindings to different communication methods. |
441 |
|
|
Exchange primitives select between using shared memory or distributed |
442 |
|
|
memory communication. |
443 |
|
|
\item Transformation operations required when transporting |
444 |
|
|
data between different grid regions. Transferring data |
445 |
|
|
between faces of a cube-sphere grid, for example, involves a rotation |
446 |
|
|
of vector components. |
447 |
|
|
\item Forward and reverse mode computations. Derivative calculations require |
448 |
|
|
tangent linear and adjoint forms of the exchange primitives. |
449 |
|
|
|
450 |
|
|
\end{itemize} |
451 |
|
|
|
452 |
|
|
\item{\bf GLOBAL SUM} The global sum operation is a central arithmetic |
453 |
|
|
operation for the pressure inversion phase of the MITgcm algorithm. |
454 |
|
|
For certain configurations scaling can be highly sensitive to |
455 |
|
|
the performance of the global sum primitive. This operation is a collective |
456 |
|
|
operation involving all tiles of the simulated domain. Different forms |
457 |
|
|
of the global sum primitive exist for handling |
458 |
|
|
\begin{itemize} |
459 |
|
|
\item Data type differences. Sixty-four bit and thirty-two bit fields may be handled |
460 |
|
|
separately. |
461 |
|
|
\item Bindings to different communication methods. |
462 |
|
|
Exchange primitives select between using shared memory or distributed |
463 |
|
|
memory communication. |
464 |
|
|
\item Forward and reverse mode computations. Derivative calculations require |
465 |
|
|
tangent linear and adjoint forms of the exchange primitives. |
466 |
|
|
\end{itemize} |
467 |
|
|
|
468 |
|
|
\item{\bf BARRIER} The WRAPPER provides a global synchronization function |
469 |
|
|
called barrier. This is used to synchronize computations over all tiles. |
470 |
|
|
The {\bf BARRIER} and {\bf GLOBAL SUM} primitives have much in common and in |
471 |
|
|
some cases use the same underlying code. |
472 |
|
|
\end{itemize} |
473 |
|
|
|
474 |
|
|
|
475 |
|
|
\subsection{Memory architecture} |
476 |
|
|
|
477 |
|
|
The WRAPPER machine model is aimed to target efficiently systems with |
478 |
|
|
highly pipelined memory architectures and systems with deep memory |
479 |
|
|
hierarchies that favor memory reuse. This is achieved by supporting a |
480 |
|
|
flexible tiling strategy as shown in figure \ref{fig:tiling-strategy}. |
481 |
|
|
Within a CPU computations are carried out sequentially on each tile |
482 |
|
|
in turn. By reshaping tiles according to the target platform it is |
483 |
|
|
possible to automatically tune code to improve memory performance. |
484 |
|
|
On a vector machine a given domain might be sub-divided into a few |
485 |
|
|
long, thin regions. On a commodity microprocessor based system, however, |
486 |
|
|
the same region could be simulated use many more smaller |
487 |
|
|
sub-domains. |
488 |
|
|
|
489 |
|
|
|
490 |
|
|
\begin{figure} |
491 |
|
|
\begin{center} |
492 |
adcroft |
1.2 |
\resizebox{5in}{!}{ |
493 |
|
|
\includegraphics{part4/tiling_detail.eps} |
494 |
cnh |
1.1 |
} |
495 |
|
|
\end{center} |
496 |
|
|
\caption{The tiling strategy that the WRAPPER supports allows tiles |
497 |
|
|
to be shaped to suit the underlying system memory architecture. |
498 |
|
|
Compact tiles that lead to greater memory reuse can be used on cache |
499 |
|
|
based systems (upper half of figure) with deep memory hierarchies, long tiles |
500 |
|
|
with large inner loops can be used to exploit vector systems having |
501 |
|
|
highly pipelined memory systems. |
502 |
|
|
} \label{fig:tiling-strategy} |
503 |
|
|
\end{figure} |
504 |
|
|
|
505 |
|
|
\newpage |
506 |
|
|
\subsection{Summary} |
507 |
|
|
Following the discussion above, the machine model that the WRAPPER |
508 |
|
|
presents to an application has the following characteristics |
509 |
|
|
|
510 |
|
|
\begin{itemize} |
511 |
|
|
\item The machine consists of one or more logical processors. \vspace{-3mm} |
512 |
|
|
\item Each processor operates on tiles that it owns.\vspace{-3mm} |
513 |
|
|
\item A processor may own more than one tile.\vspace{-3mm} |
514 |
|
|
\item Processors may compute concurrently.\vspace{-3mm} |
515 |
|
|
\item Exchange of information between tiles is handled by the |
516 |
|
|
machine (WRAPPER) not by the application. |
517 |
|
|
\end{itemize} |
518 |
|
|
Behind the scenes this allows the WRAPPER to adapt the machine model |
519 |
|
|
functions to exploit hardware on which |
520 |
|
|
\begin{itemize} |
521 |
|
|
\item Processors may be able to communicate very efficiently with each other |
522 |
|
|
using shared memory. \vspace{-3mm} |
523 |
|
|
\item An alternative communication mechanism based on a relatively |
524 |
|
|
simple inter-process communication API may be required.\vspace{-3mm} |
525 |
|
|
\item Shared memory may not necessarily obey sequential consistency, |
526 |
|
|
however some mechanism will exist for enforcing memory consistency. |
527 |
|
|
\vspace{-3mm} |
528 |
|
|
\item Memory consistency that is enforced at the hardware level |
529 |
|
|
may be expensive. Unnecessary triggering of consistency protocols |
530 |
|
|
should be avoided. \vspace{-3mm} |
531 |
|
|
\item Memory access patterns may need to either repetitive or highly |
532 |
|
|
pipelined for optimum hardware performance. \vspace{-3mm} |
533 |
|
|
\end{itemize} |
534 |
|
|
|
535 |
|
|
This generic model captures the essential hardware ingredients |
536 |
|
|
of almost all successful scientific computer systems designed in the |
537 |
|
|
last 50 years. |
538 |
|
|
|
539 |
|
|
\section{Using the WRAPPER} |
540 |
|
|
|
541 |
|
|
In order to support maximum portability the WRAPPER is implemented primarily |
542 |
|
|
in sequential Fortran 77. At a practical level the key steps provided by the |
543 |
|
|
WRAPPER are |
544 |
|
|
\begin{enumerate} |
545 |
|
|
\item specifying how a domain will be decomposed |
546 |
|
|
\item starting a code in either sequential or parallel modes of operations |
547 |
|
|
\item controlling communication between tiles and between concurrently |
548 |
|
|
computing CPU's. |
549 |
|
|
\end{enumerate} |
550 |
|
|
This section describes the details of each of these operations. |
551 |
adcroft |
1.6 |
Section \ref{sect:specifying_a_decomposition} explains how the way in which |
552 |
cnh |
1.1 |
a domain is decomposed (or composed) is expressed. Section |
553 |
adcroft |
1.6 |
\ref{sect:starting_a_code} describes practical details of running codes |
554 |
cnh |
1.1 |
in various different parallel modes on contemporary computer systems. |
555 |
adcroft |
1.6 |
Section \ref{sect:controlling_communication} explains the internal information |
556 |
cnh |
1.1 |
that the WRAPPER uses to control how information is communicated between |
557 |
|
|
tiles. |
558 |
|
|
|
559 |
|
|
\subsection{Specifying a domain decomposition} |
560 |
adcroft |
1.6 |
\label{sect:specifying_a_decomposition} |
561 |
cnh |
1.1 |
|
562 |
|
|
At its heart much of the WRAPPER works only in terms of a collection of tiles |
563 |
|
|
which are interconnected to each other. This is also true of application |
564 |
|
|
code operating within the WRAPPER. Application code is written as a series |
565 |
|
|
of compute operations, each of which operates on a single tile. If |
566 |
|
|
application code needs to perform operations involving data |
567 |
|
|
associated with another tile, it uses a WRAPPER function to obtain |
568 |
|
|
that data. |
569 |
|
|
The specification of how a global domain is constructed from tiles or alternatively |
570 |
|
|
how a global domain is decomposed into tiles is made in the file {\em SIZE.h}. |
571 |
|
|
This file defines the following parameters \\ |
572 |
|
|
|
573 |
|
|
\fbox{ |
574 |
|
|
\begin{minipage}{4.75in} |
575 |
|
|
Parameters: {\em sNx, sNy, OLx, OLy, nSx, nSy, nPx, nPy} \\ |
576 |
|
|
File: {\em model/inc/SIZE.h} |
577 |
|
|
\end{minipage} |
578 |
|
|
} \\ |
579 |
|
|
|
580 |
|
|
Together these parameters define a tiling decomposition of the style shown in |
581 |
|
|
figure \ref{fig:labelled_tile}. The parameters {\em sNx} and {\em sNy} define |
582 |
|
|
the size of an individual tile. The parameters {\em OLx} and {\em OLy} define the |
583 |
|
|
maximum size of the overlap extent. This must be set to the maximum width |
584 |
|
|
of the computation stencil that the numerical code finite-difference operations |
585 |
|
|
require between overlap region updates. The maximum overlap required |
586 |
|
|
by any of the operations in the MITgcm code distributed with this release is three grid |
587 |
|
|
points. This is set by the requirements of the $\nabla^4$ dissipation and |
588 |
|
|
diffusion operator. Code modifications and enhancements that involve adding wide |
589 |
|
|
finite-difference stencils may require increasing {\em OLx} and {\em OLy}. |
590 |
|
|
Setting {\em OLx} and {\em OLy} to a too large value will decrease code |
591 |
|
|
performance (because redundant computations will be performed), however it will |
592 |
|
|
not cause any other problems. |
593 |
|
|
|
594 |
|
|
\begin{figure} |
595 |
|
|
\begin{center} |
596 |
adcroft |
1.2 |
\resizebox{5in}{!}{ |
597 |
|
|
\includegraphics{part4/size_h.eps} |
598 |
cnh |
1.1 |
} |
599 |
|
|
\end{center} |
600 |
|
|
\caption{ The three level domain decomposition hierarchy employed by the |
601 |
|
|
WRAPPER. A domain is composed of tiles. Multiple tiles can be allocated |
602 |
|
|
to a single process. Multiple processes can exist, each with multiple tiles. |
603 |
|
|
Tiles within a process can be spread over multiple compute threads. |
604 |
|
|
} \label{fig:labelled_tile} |
605 |
|
|
\end{figure} |
606 |
|
|
|
607 |
|
|
The parameters {\em nSx} and {\em nSy} specify the number of tiles that will |
608 |
|
|
be created within a single process. Each of these tiles will have internal |
609 |
|
|
dimensions of {\em sNx} and {\em sNy}. If, when the code is executed, these tiles are |
610 |
|
|
allocated to different threads of a process that are then bound to |
611 |
|
|
different physical processors ( see the multi-threaded |
612 |
adcroft |
1.6 |
execution discussion in section \ref{sect:starting_the_code} ) then |
613 |
cnh |
1.1 |
computation will be performed concurrently on each tile. However, it is also |
614 |
|
|
possible to run the same decomposition within a process running a single thread on |
615 |
|
|
a single processor. In this case the tiles will be computed over sequentially. |
616 |
|
|
If the decomposition is run in a single process running multiple threads |
617 |
|
|
but attached to a single physical processor, then, in general, the computation |
618 |
|
|
for different tiles will be interleaved by system level software. |
619 |
|
|
This too is a valid mode of operation. |
620 |
|
|
|
621 |
|
|
The parameters {\em sNx, sNy, OLx, OLy, nSx} and{\em nSy} are used extensively by |
622 |
|
|
numerical code. The settings of {\em sNx, sNy, OLx} and {\em OLy} |
623 |
|
|
are used to form the loop ranges for many numerical calculations and |
624 |
|
|
to provide dimensions for arrays holding numerical state. |
625 |
|
|
The {\em nSx} and{\em nSy} are used in conjunction with the thread number |
626 |
|
|
parameter {\em myThid}. Much of the numerical code operating within the |
627 |
|
|
WRAPPER takes the form |
628 |
|
|
\begin{verbatim} |
629 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
630 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
631 |
|
|
: |
632 |
|
|
a block of computations ranging |
633 |
|
|
over 1,sNx +/- OLx and 1,sNy +/- OLy grid points |
634 |
|
|
: |
635 |
|
|
ENDDO |
636 |
|
|
ENDDO |
637 |
|
|
|
638 |
|
|
communication code to sum a number or maybe update |
639 |
|
|
tile overlap regions |
640 |
|
|
|
641 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
642 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
643 |
|
|
: |
644 |
|
|
another block of computations ranging |
645 |
|
|
over 1,sNx +/- OLx and 1,sNy +/- OLy grid points |
646 |
|
|
: |
647 |
|
|
ENDDO |
648 |
|
|
ENDDO |
649 |
|
|
\end{verbatim} |
650 |
|
|
The variables {\em myBxLo(myThid), myBxHi(myThid), myByLo(myThid)} and {\em |
651 |
|
|
myByHi(myThid)} set the bounds of the loops in {\em bi} and {\em bj} in this |
652 |
|
|
schematic. These variables specify the subset of the tiles in |
653 |
|
|
the range 1,{\em nSx} and 1,{\em nSy} that the logical processor bound to |
654 |
|
|
thread number {\em myThid} owns. The thread number variable {\em myThid} |
655 |
|
|
ranges from 1 to the total number of threads requested at execution time. |
656 |
|
|
For each value of {\em myThid} the loop scheme above will step sequentially |
657 |
|
|
through the tiles owned by that thread. However, different threads will |
658 |
|
|
have different ranges of tiles assigned to them, so that separate threads can |
659 |
|
|
compute iterations of the {\em bi}, {\em bj} loop concurrently. |
660 |
|
|
Within a {\em bi}, {\em bj} loop |
661 |
|
|
computation is performed concurrently over as many processes and threads |
662 |
|
|
as there are physical processors available to compute. |
663 |
|
|
|
664 |
|
|
The amount of computation that can be embedded |
665 |
|
|
a single loop over {\em bi} and {\em bj} varies for different parts of the |
666 |
|
|
MITgcm algorithm. Figure \ref{fig:bibj_extract} shows a code extract |
667 |
|
|
from the two-dimensional implicit elliptic solver. This portion of the |
668 |
|
|
code computes the l2Norm of a vector whose elements are held in |
669 |
|
|
the array {\em cg2d\_r} writing the final result to scalar variable |
670 |
|
|
{\em err}. In this case, because the l2norm requires a global reduction, |
671 |
|
|
the {\em bi},{\em bj} loop only contains one statement. This computation |
672 |
|
|
phase is then followed by a communication phase in which all threads and |
673 |
|
|
processes must participate. However, |
674 |
|
|
in other areas of the MITgcm code entries subsections of code are within |
675 |
|
|
a single {\em bi},{\em bj} loop. For example the evaluation of all |
676 |
|
|
the momentum equation prognostic terms ( see {\em S/R DYNAMICS()}) |
677 |
|
|
is within a single {\em bi},{\em bj} loop. |
678 |
|
|
|
679 |
|
|
\begin{figure} |
680 |
|
|
\begin{verbatim} |
681 |
|
|
REAL*8 cg2d_r(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
682 |
|
|
REAL*8 err |
683 |
|
|
: |
684 |
|
|
: |
685 |
|
|
other computations |
686 |
|
|
: |
687 |
|
|
: |
688 |
|
|
err = 0. |
689 |
|
|
DO bj=myByLo(myThid),myByHi(myThid) |
690 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
691 |
|
|
DO J=1,sNy |
692 |
|
|
DO I=1,sNx |
693 |
|
|
err = err + |
694 |
|
|
& cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
695 |
|
|
ENDDO |
696 |
|
|
ENDDO |
697 |
|
|
ENDDO |
698 |
|
|
ENDDO |
699 |
|
|
|
700 |
|
|
CALL GLOBAL_SUM_R8( err , myThid ) |
701 |
|
|
err = SQRT(err) |
702 |
|
|
|
703 |
|
|
\end{verbatim} |
704 |
|
|
\caption{Example of numerical code for calculating |
705 |
|
|
the l2-norm of a vector within the WRAPPER. Notice that |
706 |
|
|
under the WRAPPER arrays such as {\em cg2d\_r} have two extra trailing |
707 |
|
|
dimensions. These right most indices are tile indexes. Different |
708 |
|
|
threads with a single process operate on different ranges of tile |
709 |
|
|
index, as controlled by the settings of |
710 |
|
|
{\em myByLo, myByHi, myBxLo} and {\em myBxHi}. |
711 |
|
|
} \label{fig:bibj_extract} |
712 |
|
|
\end{figure} |
713 |
|
|
|
714 |
|
|
The final decomposition parameters are {\em nPx} and {\em nPy}. These parameters |
715 |
|
|
are used to indicate to the WRAPPER level how many processes (each with |
716 |
|
|
{\em nSx}$\times${\em nSy} tiles) will be used for this simulation. |
717 |
|
|
This information is needed during initialization and during I/O phases. |
718 |
|
|
However, unlike the variables {\em sNx, sNy, OLx, OLy, nSx} and {\em nSy} |
719 |
|
|
the values of {\em nPx} and {\em nPy} are absent |
720 |
|
|
from the core numerical and support code. |
721 |
|
|
|
722 |
|
|
\subsubsection{Examples of {\em SIZE.h} specifications} |
723 |
|
|
|
724 |
|
|
The following different {\em SIZE.h} parameter setting illustrate how to |
725 |
|
|
interpret the values of {\em sNx, sNy, OLx, OLy, nSx, nSy, nPx} |
726 |
|
|
and {\em nPy}. |
727 |
|
|
\begin{enumerate} |
728 |
|
|
\item |
729 |
|
|
\begin{verbatim} |
730 |
|
|
PARAMETER ( |
731 |
|
|
& sNx = 90, |
732 |
|
|
& sNy = 40, |
733 |
|
|
& OLx = 3, |
734 |
|
|
& OLy = 3, |
735 |
|
|
& nSx = 1, |
736 |
|
|
& nSy = 1, |
737 |
|
|
& nPx = 1, |
738 |
|
|
& nPy = 1) |
739 |
|
|
\end{verbatim} |
740 |
|
|
This sets up a single tile with x-dimension of ninety grid points, y-dimension of |
741 |
|
|
forty grid points, and x and y overlaps of three grid points each. |
742 |
|
|
\item |
743 |
|
|
\begin{verbatim} |
744 |
|
|
PARAMETER ( |
745 |
|
|
& sNx = 45, |
746 |
|
|
& sNy = 20, |
747 |
|
|
& OLx = 3, |
748 |
|
|
& OLy = 3, |
749 |
|
|
& nSx = 1, |
750 |
|
|
& nSy = 1, |
751 |
|
|
& nPx = 2, |
752 |
|
|
& nPy = 2) |
753 |
|
|
\end{verbatim} |
754 |
|
|
This sets up tiles with x-dimension of forty-five grid points, y-dimension of |
755 |
|
|
twenty grid points, and x and y overlaps of three grid points each. There are |
756 |
|
|
four tiles allocated to four separate processes ({\em nPx=2,nPy=2}) and |
757 |
|
|
arranged so that the global domain size is again ninety grid points in x and |
758 |
|
|
forty grid points in y. In general the formula for global grid size (held in |
759 |
|
|
model variables {\em Nx} and {\em Ny}) is |
760 |
|
|
\begin{verbatim} |
761 |
|
|
Nx = sNx*nSx*nPx |
762 |
|
|
Ny = sNy*nSy*nPy |
763 |
|
|
\end{verbatim} |
764 |
|
|
\item |
765 |
|
|
\begin{verbatim} |
766 |
|
|
PARAMETER ( |
767 |
|
|
& sNx = 90, |
768 |
|
|
& sNy = 10, |
769 |
|
|
& OLx = 3, |
770 |
|
|
& OLy = 3, |
771 |
|
|
& nSx = 1, |
772 |
|
|
& nSy = 2, |
773 |
|
|
& nPx = 1, |
774 |
|
|
& nPy = 2) |
775 |
|
|
\end{verbatim} |
776 |
|
|
This sets up tiles with x-dimension of ninety grid points, y-dimension of |
777 |
|
|
ten grid points, and x and y overlaps of three grid points each. There are |
778 |
|
|
four tiles allocated to two separate processes ({\em nPy=2}) each of which |
779 |
|
|
has two separate sub-domains {\em nSy=2}, |
780 |
|
|
The global domain size is again ninety grid points in x and |
781 |
|
|
forty grid points in y. The two sub-domains in each process will be computed |
782 |
|
|
sequentially if they are given to a single thread within a single process. |
783 |
|
|
Alternatively if the code is invoked with multiple threads per process |
784 |
|
|
the two domains in y may be computed on concurrently. |
785 |
|
|
\item |
786 |
|
|
\begin{verbatim} |
787 |
|
|
PARAMETER ( |
788 |
|
|
& sNx = 32, |
789 |
|
|
& sNy = 32, |
790 |
|
|
& OLx = 3, |
791 |
|
|
& OLy = 3, |
792 |
|
|
& nSx = 6, |
793 |
|
|
& nSy = 1, |
794 |
|
|
& nPx = 1, |
795 |
|
|
& nPy = 1) |
796 |
|
|
\end{verbatim} |
797 |
|
|
This sets up tiles with x-dimension of thirty-two grid points, y-dimension of |
798 |
|
|
thirty-two grid points, and x and y overlaps of three grid points each. |
799 |
|
|
There are six tiles allocated to six separate logical processors ({\em nSx=6}). |
800 |
|
|
This set of values can be used for a cube sphere calculation. |
801 |
|
|
Each tile of size $32 \times 32$ represents a face of the |
802 |
cnh |
1.4 |
cube. Initializing the tile connectivity correctly ( see section |
803 |
adcroft |
1.6 |
\ref{sect:cube_sphere_communication}. allows the rotations associated with |
804 |
cnh |
1.1 |
moving between the six cube faces to be embedded within the |
805 |
|
|
tile-tile communication code. |
806 |
|
|
\end{enumerate} |
807 |
|
|
|
808 |
|
|
|
809 |
|
|
\subsection{Starting the code} |
810 |
adcroft |
1.6 |
\label{sect:starting_the_code} |
811 |
cnh |
1.1 |
When code is started under the WRAPPER, execution begins in a main routine {\em |
812 |
|
|
eesupp/src/main.F} that is owned by the WRAPPER. Control is transferred |
813 |
|
|
to the application through a routine called {\em THE\_MODEL\_MAIN()} |
814 |
|
|
once the WRAPPER has initialized correctly and has created the necessary variables |
815 |
|
|
to support subsequent calls to communication routines |
816 |
|
|
by the application code. The startup calling sequence followed by the |
817 |
|
|
WRAPPER is shown in figure \ref{fig:wrapper_startup}. |
818 |
|
|
|
819 |
|
|
\begin{figure} |
820 |
heimbach |
1.8 |
{\footnotesize |
821 |
cnh |
1.1 |
\begin{verbatim} |
822 |
|
|
|
823 |
|
|
MAIN |
824 |
|
|
| |
825 |
|
|
|--EEBOOT :: WRAPPER initialization |
826 |
|
|
| | |
827 |
|
|
| |-- EEBOOT_MINMAL :: Minimal startup. Just enough to |
828 |
|
|
| | allow basic I/O. |
829 |
|
|
| |-- EEINTRO_MSG :: Write startup greeting. |
830 |
|
|
| | |
831 |
|
|
| |-- EESET_PARMS :: Set WRAPPER parameters |
832 |
|
|
| | |
833 |
|
|
| |-- EEWRITE_EEENV :: Print WRAPPER parameter settings |
834 |
|
|
| | |
835 |
|
|
| |-- INI_PROCS :: Associate processes with grid regions. |
836 |
|
|
| | |
837 |
|
|
| |-- INI_THREADING_ENVIRONMENT :: Associate threads with grid regions. |
838 |
|
|
| | |
839 |
|
|
| |--INI_COMMUNICATION_PATTERNS :: Initialize between tile |
840 |
|
|
| :: communication data structures |
841 |
|
|
| |
842 |
|
|
| |
843 |
|
|
|--CHECK_THREADS :: Validate multiple thread start up. |
844 |
|
|
| |
845 |
|
|
|--THE_MODEL_MAIN :: Numerical code top-level driver routine |
846 |
|
|
|
847 |
|
|
|
848 |
|
|
\end{verbatim} |
849 |
heimbach |
1.8 |
} |
850 |
cnh |
1.1 |
\caption{Main stages of the WRAPPER startup procedure. |
851 |
|
|
This process proceeds transfer of control to application code, which |
852 |
|
|
occurs through the procedure {\em THE\_MODEL\_MAIN()}. |
853 |
|
|
} \label{fig:wrapper_startup} |
854 |
|
|
\end{figure} |
855 |
|
|
|
856 |
|
|
\subsubsection{Multi-threaded execution} |
857 |
adcroft |
1.6 |
\label{sect:multi-threaded-execution} |
858 |
cnh |
1.1 |
Prior to transferring control to the procedure {\em THE\_MODEL\_MAIN()} the |
859 |
|
|
WRAPPER may cause several coarse grain threads to be initialized. The routine |
860 |
|
|
{\em THE\_MODEL\_MAIN()} is called once for each thread and is passed a single |
861 |
|
|
stack argument which is the thread number, stored in the |
862 |
|
|
variable {\em myThid}. In addition to specifying a decomposition with |
863 |
adcroft |
1.6 |
multiple tiles per process ( see section \ref{sect:specifying_a_decomposition}) |
864 |
cnh |
1.1 |
configuring and starting a code to run using multiple threads requires the following |
865 |
|
|
steps.\\ |
866 |
|
|
|
867 |
|
|
\paragraph{Compilation} |
868 |
|
|
First the code must be compiled with appropriate multi-threading directives |
869 |
|
|
active in the file {\em main.F} and with appropriate compiler flags |
870 |
|
|
to request multi-threading support. The header files |
871 |
|
|
{\em MAIN\_PDIRECTIVES1.h} and {\em MAIN\_PDIRECTIVES2.h} |
872 |
|
|
contain directives compatible with compilers for Sun, Compaq, SGI, |
873 |
|
|
Hewlett-Packard SMP systems and CRAY PVP systems. These directives can be |
874 |
|
|
activated by using compile time |
875 |
|
|
directives {\em -DTARGET\_SUN}, |
876 |
|
|
{\em -DTARGET\_DEC}, {\em -DTARGET\_SGI}, {\em -DTARGET\_HP} |
877 |
|
|
or {\em -DTARGET\_CRAY\_VECTOR} respectively. Compiler options |
878 |
|
|
for invoking multi-threaded compilation vary from system to system |
879 |
|
|
and from compiler to compiler. The options will be described |
880 |
|
|
in the individual compiler documentation. For the Fortran compiler |
881 |
|
|
from Sun the following options are needed to correctly compile |
882 |
|
|
multi-threaded code |
883 |
|
|
\begin{verbatim} |
884 |
|
|
-stackvar -explicitpar -vpara -noautopar |
885 |
|
|
\end{verbatim} |
886 |
|
|
These options are specific to the Sun compiler. Other compilers |
887 |
|
|
will use different syntax that will be described in their |
888 |
|
|
documentation. The effect of these options is as follows |
889 |
|
|
\begin{enumerate} |
890 |
|
|
\item {\bf -stackvar} Causes all local variables to be allocated in stack |
891 |
|
|
storage. This is necessary for local variables to ensure that they are private |
892 |
|
|
to their thread. Note, when using this option it may be necessary to override |
893 |
|
|
the default limit on stack-size that the operating system assigns to a process. |
894 |
|
|
This can normally be done by changing the settings of the command shells |
895 |
|
|
{\em stack-size} limit variable. However, on some systems changing this limit |
896 |
|
|
will require privileged administrator access to modify system parameters. |
897 |
|
|
|
898 |
|
|
\item {\bf -explicitpar} Requests that multiple threads be spawned |
899 |
|
|
in response to explicit directives in the application code. These |
900 |
|
|
directives are inserted with syntax appropriate to the particular target |
901 |
|
|
platform when, for example, the {\em -DTARGET\_SUN} flag is selected. |
902 |
|
|
|
903 |
|
|
\item {\bf -vpara} This causes the compiler to describe the multi-threaded |
904 |
|
|
configuration it is creating. This is not required |
905 |
|
|
but it can be useful when trouble shooting. |
906 |
|
|
|
907 |
|
|
\item {\bf -noautopar} This inhibits any automatic multi-threaded |
908 |
|
|
parallelization the compiler may otherwise generate. |
909 |
|
|
|
910 |
|
|
\end{enumerate} |
911 |
|
|
|
912 |
|
|
|
913 |
|
|
An example of valid settings for the {\em eedata} file for a |
914 |
|
|
domain with two subdomains in y and running with two threads is shown |
915 |
|
|
below |
916 |
|
|
\begin{verbatim} |
917 |
|
|
nTx=1,nTy=2 |
918 |
|
|
\end{verbatim} |
919 |
|
|
This set of values will cause computations to stay within a single |
920 |
|
|
thread when moving across the {\em nSx} sub-domains. In the y-direction, |
921 |
|
|
however, sub-domains will be split equally between two threads. |
922 |
|
|
|
923 |
|
|
\paragraph{Multi-threading files and parameters} The following |
924 |
|
|
files and variables are used in setting up multi-threaded execution.\\ |
925 |
|
|
|
926 |
|
|
\fbox{ |
927 |
|
|
\begin{minipage}{4.75in} |
928 |
|
|
File: {\em eesupp/inc/MAIN\_PDIRECTIVES1.h}\\ |
929 |
|
|
File: {\em eesupp/inc/MAIN\_PDIRECTIVES2.h}\\ |
930 |
|
|
File: {\em model/src/THE\_MODEL\_MAIN.F}\\ |
931 |
|
|
File: {\em eesupp/src/MAIN.F}\\ |
932 |
|
|
File: {\em tools/genmake}\\ |
933 |
|
|
File: {\em eedata}\\ |
934 |
|
|
CPP: {\em TARGET\_SUN}\\ |
935 |
|
|
CPP: {\em TARGET\_DEC}\\ |
936 |
|
|
CPP: {\em TARGET\_HP }\\ |
937 |
|
|
CPP: {\em TARGET\_SGI}\\ |
938 |
|
|
CPP: {\em TARGET\_CRAY\_VECTOR}\\ |
939 |
|
|
Parameter: {\em nTx}\\ |
940 |
|
|
Parameter: {\em nTy} |
941 |
|
|
\end{minipage} |
942 |
|
|
} \\ |
943 |
|
|
|
944 |
|
|
\subsubsection{Multi-process execution} |
945 |
adcroft |
1.6 |
\label{sect:multi-process-execution} |
946 |
cnh |
1.1 |
|
947 |
|
|
Despite its appealing programming model, multi-threaded execution remains |
948 |
|
|
less common then multi-process execution. One major reason for this |
949 |
|
|
is that many system libraries are still not ``thread-safe''. This means that for |
950 |
|
|
example on some systems it is not safe to call system routines to |
951 |
|
|
do I/O when running in multi-threaded mode, except for in a limited set of |
952 |
|
|
circumstances. Another reason is that support for multi-threaded programming |
953 |
|
|
models varies between systems. |
954 |
|
|
|
955 |
|
|
Multi-process execution is more ubiquitous. |
956 |
|
|
In order to run code in a multi-process configuration a decomposition |
957 |
adcroft |
1.6 |
specification ( see section \ref{sect:specifying_a_decomposition}) |
958 |
cnh |
1.3 |
is given ( in which the at least one of the |
959 |
cnh |
1.1 |
parameters {\em nPx} or {\em nPy} will be greater than one) |
960 |
|
|
and then, as for multi-threaded operation, |
961 |
|
|
appropriate compile time and run time steps must be taken. |
962 |
|
|
|
963 |
|
|
\paragraph{Compilation} Multi-process execution under the WRAPPER |
964 |
|
|
assumes that the portable, MPI libraries are available |
965 |
|
|
for controlling the start-up of multiple processes. The MPI libraries |
966 |
|
|
are not required, although they are usually used, for performance |
967 |
|
|
critical communication. However, in order to simplify the task |
968 |
|
|
of controlling and coordinating the start up of a large number |
969 |
|
|
(hundreds and possibly even thousands) of copies of the same |
970 |
|
|
program, MPI is used. The calls to the MPI multi-process startup |
971 |
|
|
routines must be activated at compile time. This is done |
972 |
|
|
by setting the {\em ALLOW\_USE\_MPI} and {\em ALWAYS\_USE\_MPI} |
973 |
|
|
flags in the {\em CPP\_EEOPTIONS.h} file.\\ |
974 |
|
|
|
975 |
|
|
\fbox{ |
976 |
|
|
\begin{minipage}{4.75in} |
977 |
|
|
File: {\em eesupp/inc/CPP\_EEOPTIONS.h}\\ |
978 |
|
|
CPP: {\em ALLOW\_USE\_MPI}\\ |
979 |
|
|
CPP: {\em ALWAYS\_USE\_MPI}\\ |
980 |
|
|
Parameter: {\em nPx}\\ |
981 |
|
|
Parameter: {\em nPy} |
982 |
|
|
\end{minipage} |
983 |
|
|
} \\ |
984 |
|
|
|
985 |
|
|
Additionally, compile time options are required to link in the |
986 |
|
|
MPI libraries and header files. Examples of these options |
987 |
|
|
can be found in the {\em genmake} script that creates makefiles |
988 |
|
|
for compilation. When this script is executed with the {bf -mpi} |
989 |
|
|
flag it will generate a makefile that includes |
990 |
|
|
paths for search for MPI head files and for linking in |
991 |
|
|
MPI libraries. For example the {\bf -mpi} flag on a |
992 |
|
|
Silicon Graphics IRIX system causes a |
993 |
|
|
Makefile with the compilation command |
994 |
|
|
Graphics IRIX system \begin{verbatim} |
995 |
|
|
mpif77 -I/usr/local/mpi/include -DALLOW_USE_MPI -DALWAYS_USE_MPI |
996 |
|
|
\end{verbatim} |
997 |
|
|
to be generated. |
998 |
|
|
This is the correct set of options for using the MPICH open-source |
999 |
|
|
version of MPI, when it has been installed under the subdirectory |
1000 |
|
|
/usr/local/mpi. |
1001 |
|
|
However, on many systems there may be several |
1002 |
|
|
versions of MPI installed. For example many systems have both |
1003 |
|
|
the open source MPICH set of libraries and a vendor specific native form |
1004 |
|
|
of the MPI libraries. The correct setup to use will depend on the |
1005 |
|
|
local configuration of your system.\\ |
1006 |
|
|
|
1007 |
|
|
\fbox{ |
1008 |
|
|
\begin{minipage}{4.75in} |
1009 |
|
|
File: {\em tools/genmake} |
1010 |
|
|
\end{minipage} |
1011 |
|
|
} \\ |
1012 |
|
|
\paragraph{\bf Execution} The mechanics of starting a program in |
1013 |
|
|
multi-process mode under MPI is not standardized. Documentation |
1014 |
|
|
associated with the distribution of MPI installed on a system will |
1015 |
|
|
describe how to start a program using that distribution. |
1016 |
|
|
For the free, open-source MPICH system the MITgcm program is started |
1017 |
|
|
using a command such as |
1018 |
|
|
\begin{verbatim} |
1019 |
|
|
mpirun -np 64 -machinefile mf ./mitgcmuv |
1020 |
|
|
\end{verbatim} |
1021 |
cnh |
1.5 |
In this example the text {\em -np 64} specifies the number of processes |
1022 |
cnh |
1.1 |
that will be created. The numeric value {\em 64} must be equal to the |
1023 |
|
|
product of the processor grid settings of {\em nPx} and {\em nPy} |
1024 |
|
|
in the file {\em SIZE.h}. The parameter {\em mf} specifies that a text file |
1025 |
|
|
called ``mf'' will be read to get a list of processor names on |
1026 |
|
|
which the sixty-four processes will execute. The syntax of this file |
1027 |
|
|
is specified by the MPI distribution |
1028 |
|
|
\\ |
1029 |
|
|
|
1030 |
|
|
\fbox{ |
1031 |
|
|
\begin{minipage}{4.75in} |
1032 |
|
|
File: {\em SIZE.h}\\ |
1033 |
|
|
Parameter: {\em nPx}\\ |
1034 |
|
|
Parameter: {\em nPy} |
1035 |
|
|
\end{minipage} |
1036 |
|
|
} \\ |
1037 |
|
|
|
1038 |
adcroft |
1.2 |
|
1039 |
|
|
\paragraph{Environment variables} |
1040 |
|
|
On most systems multi-threaded execution also requires the setting |
1041 |
|
|
of a special environment variable. On many machines this variable |
1042 |
|
|
is called PARALLEL and its values should be set to the number |
1043 |
|
|
of parallel threads required. Generally the help pages associated |
1044 |
|
|
with the multi-threaded compiler on a machine will explain |
1045 |
|
|
how to set the required environment variables for that machines. |
1046 |
|
|
|
1047 |
|
|
\paragraph{Runtime input parameters} |
1048 |
|
|
Finally the file {\em eedata} needs to be configured to indicate |
1049 |
|
|
the number of threads to be used in the x and y directions. |
1050 |
|
|
The variables {\em nTx} and {\em nTy} in this file are used to |
1051 |
|
|
specify the information required. The product of {\em nTx} and |
1052 |
|
|
{\em nTy} must be equal to the number of threads spawned i.e. |
1053 |
|
|
the setting of the environment variable PARALLEL. |
1054 |
|
|
The value of {\em nTx} must subdivide the number of sub-domains |
1055 |
|
|
in x ({\em nSx}) exactly. The value of {\em nTy} must subdivide the |
1056 |
|
|
number of sub-domains in y ({\em nSy}) exactly. |
1057 |
cnh |
1.1 |
The multiprocess startup of the MITgcm executable {\em mitgcmuv} |
1058 |
|
|
is controlled by the routines {\em EEBOOT\_MINIMAL()} and |
1059 |
|
|
{\em INI\_PROCS()}. The first routine performs basic steps required |
1060 |
|
|
to make sure each process is started and has a textual output |
1061 |
|
|
stream associated with it. By default two output files are opened |
1062 |
|
|
for each process with names {\bf STDOUT.NNNN} and {\bf STDERR.NNNN}. |
1063 |
|
|
The {\bf NNNNN} part of the name is filled in with the process |
1064 |
|
|
number so that process number 0 will create output files |
1065 |
|
|
{\bf STDOUT.0000} and {\bf STDERR.0000}, process number 1 will create |
1066 |
|
|
output files {\bf STDOUT.0001} and {\bf STDERR.0001} etc... These files |
1067 |
|
|
are used for reporting status and configuration information and |
1068 |
|
|
for reporting error conditions on a process by process basis. |
1069 |
adcroft |
1.2 |
The {\em EEBOOT\_MINIMAL()} procedure also sets the variables |
1070 |
cnh |
1.1 |
{\em myProcId} and {\em MPI\_COMM\_MODEL}. |
1071 |
|
|
These variables are related |
1072 |
|
|
to processor identification are are used later in the routine |
1073 |
|
|
{\em INI\_PROCS()} to allocate tiles to processes. |
1074 |
|
|
|
1075 |
|
|
Allocation of processes to tiles in controlled by the routine |
1076 |
|
|
{\em INI\_PROCS()}. For each process this routine sets |
1077 |
|
|
the variables {\em myXGlobalLo} and {\em myYGlobalLo}. |
1078 |
|
|
These variables specify (in index space) the coordinate |
1079 |
|
|
of the southern most and western most corner of the |
1080 |
|
|
southern most and western most tile owned by this process. |
1081 |
|
|
The variables {\em pidW}, {\em pidE}, {\em pidS} and {\em pidN} |
1082 |
|
|
are also set in this routine. These are used to identify |
1083 |
|
|
processes holding tiles to the west, east, south and north |
1084 |
|
|
of this process. These values are stored in global storage |
1085 |
|
|
in the header file {\em EESUPPORT.h} for use by |
1086 |
|
|
communication routines. |
1087 |
|
|
\\ |
1088 |
|
|
|
1089 |
|
|
\fbox{ |
1090 |
|
|
\begin{minipage}{4.75in} |
1091 |
|
|
File: {\em eesupp/src/eeboot\_minimal.F}\\ |
1092 |
|
|
File: {\em eesupp/src/ini\_procs.F} \\ |
1093 |
|
|
File: {\em eesupp/inc/EESUPPORT.h} \\ |
1094 |
|
|
Parameter: {\em myProcId} \\ |
1095 |
|
|
Parameter: {\em MPI\_COMM\_MODEL} \\ |
1096 |
|
|
Parameter: {\em myXGlobalLo} \\ |
1097 |
|
|
Parameter: {\em myYGlobalLo} \\ |
1098 |
|
|
Parameter: {\em pidW } \\ |
1099 |
|
|
Parameter: {\em pidE } \\ |
1100 |
|
|
Parameter: {\em pidS } \\ |
1101 |
|
|
Parameter: {\em pidN } |
1102 |
|
|
\end{minipage} |
1103 |
|
|
} \\ |
1104 |
|
|
|
1105 |
|
|
|
1106 |
|
|
\subsection{Controlling communication} |
1107 |
|
|
The WRAPPER maintains internal information that is used for communication |
1108 |
|
|
operations and that can be customized for different platforms. This section |
1109 |
|
|
describes the information that is held and used. |
1110 |
adcroft |
1.2 |
|
1111 |
cnh |
1.1 |
\begin{enumerate} |
1112 |
|
|
\item {\bf Tile-tile connectivity information} For each tile the WRAPPER |
1113 |
|
|
sets a flag that sets the tile number to the north, south, east and |
1114 |
|
|
west of that tile. This number is unique over all tiles in a |
1115 |
|
|
configuration. The number is held in the variables {\em tileNo} |
1116 |
|
|
( this holds the tiles own number), {\em tileNoN}, {\em tileNoS}, |
1117 |
|
|
{\em tileNoE} and {\em tileNoW}. A parameter is also stored with each tile |
1118 |
|
|
that specifies the type of communication that is used between tiles. |
1119 |
|
|
This information is held in the variables {\em tileCommModeN}, |
1120 |
|
|
{\em tileCommModeS}, {\em tileCommModeE} and {\em tileCommModeW}. |
1121 |
|
|
This latter set of variables can take one of the following values |
1122 |
|
|
{\em COMM\_NONE}, {\em COMM\_MSG}, {\em COMM\_PUT} and {\em COMM\_GET}. |
1123 |
|
|
A value of {\em COMM\_NONE} is used to indicate that a tile has no |
1124 |
cnh |
1.4 |
neighbor to communicate with on a particular face. A value |
1125 |
cnh |
1.1 |
of {\em COMM\_MSG} is used to indicated that some form of distributed |
1126 |
|
|
memory communication is required to communicate between |
1127 |
adcroft |
1.6 |
these tile faces ( see section \ref{sect:distributed_memory_communication}). |
1128 |
cnh |
1.1 |
A value of {\em COMM\_PUT} or {\em COMM\_GET} is used to indicate |
1129 |
|
|
forms of shared memory communication ( see section |
1130 |
adcroft |
1.6 |
\ref{sect:shared_memory_communication}). The {\em COMM\_PUT} value indicates |
1131 |
cnh |
1.1 |
that a CPU should communicate by writing to data structures owned by another |
1132 |
|
|
CPU. A {\em COMM\_GET} value indicates that a CPU should communicate by reading |
1133 |
|
|
from data structures owned by another CPU. These flags affect the behavior |
1134 |
|
|
of the WRAPPER exchange primitive |
1135 |
|
|
(see figure \ref{fig:communication_primitives}). The routine |
1136 |
|
|
{\em ini\_communication\_patterns()} is responsible for setting the |
1137 |
|
|
communication mode values for each tile. |
1138 |
|
|
\\ |
1139 |
|
|
|
1140 |
|
|
\fbox{ |
1141 |
|
|
\begin{minipage}{4.75in} |
1142 |
|
|
File: {\em eesupp/src/ini\_communication\_patterns.F}\\ |
1143 |
|
|
File: {\em eesupp/inc/EESUPPORT.h} \\ |
1144 |
|
|
Parameter: {\em tileNo} \\ |
1145 |
|
|
Parameter: {\em tileNoE} \\ |
1146 |
|
|
Parameter: {\em tileNoW} \\ |
1147 |
|
|
Parameter: {\em tileNoN} \\ |
1148 |
|
|
Parameter: {\em tileNoS} \\ |
1149 |
|
|
Parameter: {\em tileCommModeE} \\ |
1150 |
|
|
Parameter: {\em tileCommModeW} \\ |
1151 |
|
|
Parameter: {\em tileCommModeN} \\ |
1152 |
|
|
Parameter: {\em tileCommModeS} \\ |
1153 |
|
|
\end{minipage} |
1154 |
|
|
} \\ |
1155 |
|
|
|
1156 |
|
|
\item {\bf MP directives} |
1157 |
|
|
The WRAPPER transfers control to numerical application code through |
1158 |
|
|
the routine {\em THE\_MODEL\_MAIN}. This routine is called in a way |
1159 |
|
|
that allows for it to be invoked by several threads. Support for this |
1160 |
|
|
is based on using multi-processing (MP) compiler directives. |
1161 |
|
|
Most commercially available Fortran compilers support the generation |
1162 |
|
|
of code to spawn multiple threads through some form of compiler directives. |
1163 |
|
|
As this is generally much more convenient than writing code to interface |
1164 |
|
|
to operating system libraries to explicitly spawn threads, and on some systems |
1165 |
|
|
this may be the only method available the WRAPPER is distributed with |
1166 |
|
|
template MP directives for a number of systems. |
1167 |
|
|
|
1168 |
|
|
These directives are inserted into the code just before and after the |
1169 |
|
|
transfer of control to numerical algorithm code through the routine |
1170 |
|
|
{\em THE\_MODEL\_MAIN}. Figure \ref{fig:mp_directives} shows an example of |
1171 |
|
|
the code that performs this process for a Silicon Graphics system. |
1172 |
|
|
This code is extracted from the files {\em main.F} and |
1173 |
|
|
{\em MAIN\_PDIRECTIVES1.h}. The variable {\em nThreads} specifies |
1174 |
|
|
how many instances of the routine {\em THE\_MODEL\_MAIN} will |
1175 |
|
|
be created. The value of {\em nThreads} is set in the routine |
1176 |
|
|
{\em INI\_THREADING\_ENVIRONMENT}. The value is set equal to the |
1177 |
|
|
the product of the parameters {\em nTx} and {\em nTy} that |
1178 |
|
|
are read from the file {\em eedata}. If the value of {\em nThreads} |
1179 |
|
|
is inconsistent with the number of threads requested from the |
1180 |
|
|
operating system (for example by using an environment |
1181 |
adcroft |
1.6 |
variable as described in section \ref{sect:multi_threaded_execution}) |
1182 |
cnh |
1.1 |
then usually an error will be reported by the routine |
1183 |
|
|
{\em CHECK\_THREADS}.\\ |
1184 |
|
|
|
1185 |
|
|
\fbox{ |
1186 |
|
|
\begin{minipage}{4.75in} |
1187 |
|
|
File: {\em eesupp/src/ini\_threading\_environment.F}\\ |
1188 |
|
|
File: {\em eesupp/src/check\_threads.F} \\ |
1189 |
|
|
File: {\em eesupp/src/main.F} \\ |
1190 |
|
|
File: {\em eesupp/inc/MAIN\_PDIRECTIVES1.h} \\ |
1191 |
|
|
File: {\em eedata } \\ |
1192 |
|
|
Parameter: {\em nThreads} \\ |
1193 |
|
|
Parameter: {\em nTx} \\ |
1194 |
|
|
Parameter: {\em nTy} \\ |
1195 |
|
|
\end{minipage} |
1196 |
|
|
} |
1197 |
|
|
|
1198 |
|
|
\item {\bf memsync flags} |
1199 |
adcroft |
1.6 |
As discussed in section \ref{sect:memory_consistency}, when using shared memory, |
1200 |
cnh |
1.1 |
a low-level system function may be need to force memory consistency. |
1201 |
|
|
The routine {\em MEMSYNC()} is used for this purpose. This routine should |
1202 |
|
|
not need modifying and the information below is only provided for |
1203 |
|
|
completeness. A logical parameter {\em exchNeedsMemSync} set |
1204 |
|
|
in the routine {\em INI\_COMMUNICATION\_PATTERNS()} controls whether |
1205 |
|
|
the {\em MEMSYNC()} primitive is called. In general this |
1206 |
|
|
routine is only used for multi-threaded execution. |
1207 |
|
|
The code that goes into the {\em MEMSYNC()} |
1208 |
|
|
routine is specific to the compiler and |
1209 |
|
|
processor being used for multi-threaded execution and in general |
1210 |
|
|
must be written using a short code snippet of assembly language. |
1211 |
|
|
For an Ultra Sparc system the following code snippet is used |
1212 |
|
|
\begin{verbatim} |
1213 |
|
|
asm("membar #LoadStore|#StoreStore"); |
1214 |
|
|
\end{verbatim} |
1215 |
cnh |
1.4 |
for an Alpha based system the equivalent code reads |
1216 |
cnh |
1.1 |
\begin{verbatim} |
1217 |
|
|
asm("mb"); |
1218 |
|
|
\end{verbatim} |
1219 |
|
|
while on an x86 system the following code is required |
1220 |
|
|
\begin{verbatim} |
1221 |
|
|
asm("lock; addl $0,0(%%esp)": : :"memory") |
1222 |
|
|
\end{verbatim} |
1223 |
|
|
|
1224 |
|
|
\item {\bf Cache line size} |
1225 |
adcroft |
1.6 |
As discussed in section \ref{sect:cache_effects_and_false_sharing}, |
1226 |
cnh |
1.1 |
milti-threaded codes explicitly avoid penalties associated with excessive |
1227 |
cnh |
1.5 |
coherence traffic on an SMP system. To do this the shared memory data structures |
1228 |
cnh |
1.1 |
used by the {\em GLOBAL\_SUM}, {\em GLOBAL\_MAX} and {\em BARRIER} routines |
1229 |
|
|
are padded. The variables that control the padding are set in the |
1230 |
|
|
header file {\em EEPARAMS.h}. These variables are called |
1231 |
|
|
{\em cacheLineSize}, {\em lShare1}, {\em lShare4} and |
1232 |
|
|
{\em lShare8}. The default values should not normally need changing. |
1233 |
|
|
\item {\bf \_BARRIER} |
1234 |
|
|
This is a CPP macro that is expanded to a call to a routine |
1235 |
cnh |
1.5 |
which synchronizes all the logical processors running under the |
1236 |
cnh |
1.1 |
WRAPPER. Using a macro here preserves flexibility to insert |
1237 |
|
|
a specialized call in-line into application code. By default this |
1238 |
|
|
resolves to calling the procedure {\em BARRIER()}. The default |
1239 |
|
|
setting for the \_BARRIER macro is given in the file {\em CPP\_EEMACROS.h}. |
1240 |
|
|
|
1241 |
|
|
\item {\bf \_GSUM} |
1242 |
|
|
This is a CPP macro that is expanded to a call to a routine |
1243 |
cnh |
1.5 |
which sums up a floating point number |
1244 |
cnh |
1.1 |
over all the logical processors running under the |
1245 |
|
|
WRAPPER. Using a macro here provides extra flexibility to insert |
1246 |
|
|
a specialized call in-line into application code. By default this |
1247 |
cnh |
1.5 |
resolves to calling the procedure {\em GLOBAL\_SUM\_R8()} ( for |
1248 |
|
|
64-bit floating point operands) |
1249 |
|
|
or {\em GLOBAL\_SUM\_R4()} (for 32-bit floating point operands). The default |
1250 |
cnh |
1.1 |
setting for the \_GSUM macro is given in the file {\em CPP\_EEMACROS.h}. |
1251 |
|
|
The \_GSUM macro is a performance critical operation, especially for |
1252 |
|
|
large processor count, small tile size configurations. |
1253 |
adcroft |
1.6 |
The custom communication example discussed in section \ref{sect:jam_example} |
1254 |
cnh |
1.1 |
shows how the macro is used to invoke a custom global sum routine |
1255 |
|
|
for a specific set of hardware. |
1256 |
|
|
|
1257 |
|
|
\item {\bf \_EXCH} |
1258 |
|
|
The \_EXCH CPP macro is used to update tile overlap regions. |
1259 |
|
|
It is qualified by a suffix indicating whether overlap updates are for |
1260 |
|
|
two-dimensional ( \_EXCH\_XY ) or three dimensional ( \_EXCH\_XYZ ) |
1261 |
|
|
physical fields and whether fields are 32-bit floating point |
1262 |
|
|
( \_EXCH\_XY\_R4, \_EXCH\_XYZ\_R4 ) or 64-bit floating point |
1263 |
|
|
( \_EXCH\_XY\_R8, \_EXCH\_XYZ\_R8 ). The macro mappings are defined |
1264 |
|
|
in the header file {\em CPP\_EEMACROS.h}. As with \_GSUM, the |
1265 |
|
|
\_EXCH operation plays a crucial role in scaling to small tile, |
1266 |
|
|
large logical and physical processor count configurations. |
1267 |
adcroft |
1.6 |
The example in section \ref{sect:jam_example} discusses defining an |
1268 |
cnh |
1.5 |
optimized and specialized form on the \_EXCH operation. |
1269 |
cnh |
1.1 |
|
1270 |
|
|
The \_EXCH operation is also central to supporting grids such as |
1271 |
|
|
the cube-sphere grid. In this class of grid a rotation may be required |
1272 |
|
|
between tiles. Aligning the coordinate requiring rotation with the |
1273 |
cnh |
1.5 |
tile decomposition, allows the coordinate transformation to |
1274 |
cnh |
1.1 |
be embedded within a custom form of the \_EXCH primitive. |
1275 |
|
|
|
1276 |
|
|
\item {\bf Reverse Mode} |
1277 |
|
|
The communication primitives \_EXCH and \_GSUM both employ |
1278 |
|
|
hand-written adjoint forms (or reverse mode) forms. |
1279 |
|
|
These reverse mode forms can be found in the |
1280 |
cnh |
1.5 |
source code directory {\em pkg/autodiff}. |
1281 |
cnh |
1.1 |
For the global sum primitive the reverse mode form |
1282 |
|
|
calls are to {\em GLOBAL\_ADSUM\_R4} and |
1283 |
|
|
{\em GLOBAL\_ADSUM\_R8}. The reverse mode form of the |
1284 |
cnh |
1.5 |
exchange primitives are found in routines |
1285 |
cnh |
1.1 |
prefixed {\em ADEXCH}. The exchange routines make calls to |
1286 |
|
|
the same low-level communication primitives as the forward mode |
1287 |
|
|
operations. However, the routine argument {\em simulationMode} |
1288 |
|
|
is set to the value {\em REVERSE\_SIMULATION}. This signifies |
1289 |
|
|
ti the low-level routines that the adjoint forms of the |
1290 |
|
|
appropriate communication operation should be performed. |
1291 |
|
|
\item {\bf MAX\_NO\_THREADS} |
1292 |
|
|
The variable {\em MAX\_NO\_THREADS} is used to indicate the |
1293 |
|
|
maximum number of OS threads that a code will use. This |
1294 |
|
|
value defaults to thirty-two and is set in the file {\em EEPARAMS.h}. |
1295 |
|
|
For single threaded execution it can be reduced to one if required. |
1296 |
cnh |
1.5 |
The value; is largely private to the WRAPPER and application code |
1297 |
cnh |
1.1 |
will nor normally reference the value, except in the following scenario. |
1298 |
|
|
|
1299 |
|
|
For certain physical parametrization schemes it is necessary to have |
1300 |
|
|
a substantial number of work arrays. Where these arrays are allocated |
1301 |
|
|
in heap storage ( for example COMMON blocks ) multi-threaded |
1302 |
|
|
execution will require multiple instances of the COMMON block data. |
1303 |
|
|
This can be achieved using a Fortran 90 module construct, however, |
1304 |
|
|
if this might be unavailable then the work arrays can be extended |
1305 |
|
|
with dimensions use the tile dimensioning scheme of {\em nSx} |
1306 |
|
|
and {\em nSy} ( as described in section |
1307 |
adcroft |
1.6 |
\ref{sect:specifying_a_decomposition}). However, if the configuration |
1308 |
cnh |
1.1 |
being specified involves many more tiles than OS threads then |
1309 |
|
|
it can save memory resources to reduce the variable |
1310 |
|
|
{\em MAX\_NO\_THREADS} to be equal to the actual number of threads that |
1311 |
cnh |
1.5 |
will be used and to declare the physical parameterization |
1312 |
|
|
work arrays with a single {\em MAX\_NO\_THREADS} extra dimension. |
1313 |
cnh |
1.1 |
An example of this is given in the verification experiment |
1314 |
|
|
{\em aim.5l\_cs}. Here the default setting of |
1315 |
|
|
{\em MAX\_NO\_THREADS} is altered to |
1316 |
|
|
\begin{verbatim} |
1317 |
|
|
INTEGER MAX_NO_THREADS |
1318 |
|
|
PARAMETER ( MAX_NO_THREADS = 6 ) |
1319 |
|
|
\end{verbatim} |
1320 |
|
|
and several work arrays for storing intermediate calculations are |
1321 |
|
|
created with declarations of the form. |
1322 |
|
|
\begin{verbatim} |
1323 |
|
|
common /FORCIN/ sst1(ngp,MAX_NO_THREADS) |
1324 |
|
|
\end{verbatim} |
1325 |
cnh |
1.5 |
This declaration scheme is not used widely, because most global data |
1326 |
cnh |
1.1 |
is used for permanent not temporary storage of state information. |
1327 |
|
|
In the case of permanent state information this approach cannot be used |
1328 |
|
|
because there has to be enough storage allocated for all tiles. |
1329 |
|
|
However, the technique can sometimes be a useful scheme for reducing memory |
1330 |
cnh |
1.5 |
requirements in complex physical parameterizations. |
1331 |
adcroft |
1.2 |
\end{enumerate} |
1332 |
|
|
|
1333 |
|
|
\begin{figure} |
1334 |
|
|
\begin{verbatim} |
1335 |
|
|
C-- |
1336 |
|
|
C-- Parallel directives for MIPS Pro Fortran compiler |
1337 |
|
|
C-- |
1338 |
|
|
C Parallel compiler directives for SGI with IRIX |
1339 |
|
|
C$PAR PARALLEL DO |
1340 |
|
|
C$PAR& CHUNK=1,MP_SCHEDTYPE=INTERLEAVE, |
1341 |
|
|
C$PAR& SHARE(nThreads),LOCAL(myThid,I) |
1342 |
|
|
C |
1343 |
|
|
DO I=1,nThreads |
1344 |
|
|
myThid = I |
1345 |
|
|
|
1346 |
|
|
C-- Invoke nThreads instances of the numerical model |
1347 |
|
|
CALL THE_MODEL_MAIN(myThid) |
1348 |
|
|
|
1349 |
|
|
ENDDO |
1350 |
|
|
\end{verbatim} |
1351 |
|
|
\caption{Prior to transferring control to |
1352 |
|
|
the procedure {\em THE\_MODEL\_MAIN()} the WRAPPER may use |
1353 |
|
|
MP directives to spawn multiple threads. |
1354 |
|
|
} \label{fig:mp_directives} |
1355 |
|
|
\end{figure} |
1356 |
cnh |
1.1 |
|
1357 |
|
|
|
1358 |
|
|
\subsubsection{Specializing the Communication Code} |
1359 |
|
|
|
1360 |
|
|
The isolation of performance critical communication primitives and the |
1361 |
|
|
sub-division of the simulation domain into tiles is a powerful tool. |
1362 |
|
|
Here we show how it can be used to improve application performance and |
1363 |
cnh |
1.5 |
how it can be used to adapt to new griding approaches. |
1364 |
cnh |
1.1 |
|
1365 |
|
|
\subsubsection{JAM example} |
1366 |
adcroft |
1.6 |
\label{sect:jam_example} |
1367 |
cnh |
1.1 |
On some platforms a big performance boost can be obtained by |
1368 |
|
|
binding the communication routines {\em \_EXCH} and |
1369 |
|
|
{\em \_GSUM} to specialized native libraries ) fro example the |
1370 |
|
|
shmem library on CRAY T3E systems). The {\em LETS\_MAKE\_JAM} CPP flag |
1371 |
|
|
is used as an illustration of a specialized communication configuration |
1372 |
|
|
that substitutes for standard, portable forms of {\em \_EXCH} and |
1373 |
|
|
{\em \_GSUM}. It affects three source files {\em eeboot.F}, |
1374 |
|
|
{\em CPP\_EEMACROS.h} and {\em cg2d.F}. When the flag is defined |
1375 |
|
|
is has the following effects. |
1376 |
|
|
\begin{itemize} |
1377 |
|
|
\item An extra phase is included at boot time to initialize the custom |
1378 |
|
|
communications library ( see {\em ini\_jam.F}). |
1379 |
|
|
\item The {\em \_GSUM} and {\em \_EXCH} macro definitions are replaced |
1380 |
|
|
with calls to custom routines ( see {\em gsum\_jam.F} and {\em exch\_jam.F}) |
1381 |
|
|
\item a highly specialized form of the exchange operator (optimized |
1382 |
cnh |
1.4 |
for overlap regions of width one) is substituted into the elliptic |
1383 |
cnh |
1.1 |
solver routine {\em cg2d.F}. |
1384 |
|
|
\end{itemize} |
1385 |
|
|
Developing specialized code for other libraries follows a similar |
1386 |
|
|
pattern. |
1387 |
|
|
|
1388 |
|
|
\subsubsection{Cube sphere communication} |
1389 |
adcroft |
1.6 |
\label{sect:cube_sphere_communication} |
1390 |
cnh |
1.1 |
Actual {\em \_EXCH} routine code is generated automatically from |
1391 |
|
|
a series of template files, for example {\em exch\_rx.template}. |
1392 |
|
|
This is done to allow a large number of variations on the exchange |
1393 |
|
|
process to be maintained. One set of variations supports the |
1394 |
cnh |
1.4 |
cube sphere grid. Support for a cube sphere grid in MITgcm is based |
1395 |
cnh |
1.1 |
on having each face of the cube as a separate tile (or tiles). |
1396 |
cnh |
1.4 |
The exchange routines are then able to absorb much of the |
1397 |
cnh |
1.1 |
detailed rotation and reorientation required when moving around the |
1398 |
|
|
cube grid. The set of {\em \_EXCH} routines that contain the |
1399 |
|
|
word cube in their name perform these transformations. |
1400 |
|
|
They are invoked when the run-time logical parameter |
1401 |
|
|
{\em useCubedSphereExchange} is set true. To facilitate the |
1402 |
|
|
transformations on a staggered C-grid, exchange operations are defined |
1403 |
cnh |
1.4 |
separately for both vector and scalar quantities and for |
1404 |
cnh |
1.1 |
grid-centered and for grid-face and corner quantities. |
1405 |
|
|
Three sets of exchange routines are defined. Routines |
1406 |
|
|
with names of the form {\em exch\_rx} are used to exchange |
1407 |
|
|
cell centered scalar quantities. Routines with names of the form |
1408 |
|
|
{\em exch\_uv\_rx} are used to exchange vector quantities located at |
1409 |
|
|
the C-grid velocity points. The vector quantities exchanged by the |
1410 |
|
|
{\em exch\_uv\_rx} routines can either be signed (for example velocity |
1411 |
|
|
components) or un-signed (for example grid-cell separations). |
1412 |
|
|
Routines with names of the form {\em exch\_z\_rx} are used to exchange |
1413 |
|
|
quantities at the C-grid vorticity point locations. |
1414 |
|
|
|
1415 |
|
|
|
1416 |
|
|
|
1417 |
|
|
|
1418 |
|
|
\section{MITgcm execution under WRAPPER} |
1419 |
|
|
|
1420 |
|
|
Fitting together the WRAPPER elements, package elements and |
1421 |
|
|
MITgcm core equation elements of the source code produces calling |
1422 |
adcroft |
1.6 |
sequence shown in section \ref{sect:calling_sequence} |
1423 |
cnh |
1.1 |
|
1424 |
|
|
\subsection{Annotated call tree for MITgcm and WRAPPER} |
1425 |
adcroft |
1.6 |
\label{sect:calling_sequence} |
1426 |
cnh |
1.1 |
|
1427 |
|
|
WRAPPER layer. |
1428 |
|
|
|
1429 |
heimbach |
1.8 |
{\footnotesize |
1430 |
cnh |
1.1 |
\begin{verbatim} |
1431 |
|
|
|
1432 |
|
|
MAIN |
1433 |
|
|
| |
1434 |
|
|
|--EEBOOT :: WRAPPER initialization |
1435 |
|
|
| | |
1436 |
|
|
| |-- EEBOOT_MINMAL :: Minimal startup. Just enough to |
1437 |
|
|
| | allow basic I/O. |
1438 |
|
|
| |-- EEINTRO_MSG :: Write startup greeting. |
1439 |
|
|
| | |
1440 |
|
|
| |-- EESET_PARMS :: Set WRAPPER parameters |
1441 |
|
|
| | |
1442 |
|
|
| |-- EEWRITE_EEENV :: Print WRAPPER parameter settings |
1443 |
|
|
| | |
1444 |
|
|
| |-- INI_PROCS :: Associate processes with grid regions. |
1445 |
|
|
| | |
1446 |
|
|
| |-- INI_THREADING_ENVIRONMENT :: Associate threads with grid regions. |
1447 |
|
|
| | |
1448 |
|
|
| |--INI_COMMUNICATION_PATTERNS :: Initialize between tile |
1449 |
|
|
| :: communication data structures |
1450 |
|
|
| |
1451 |
|
|
| |
1452 |
|
|
|--CHECK_THREADS :: Validate multiple thread start up. |
1453 |
|
|
| |
1454 |
|
|
|--THE_MODEL_MAIN :: Numerical code top-level driver routine |
1455 |
|
|
|
1456 |
|
|
\end{verbatim} |
1457 |
heimbach |
1.8 |
} |
1458 |
cnh |
1.1 |
|
1459 |
|
|
Core equations plus packages. |
1460 |
|
|
|
1461 |
heimbach |
1.8 |
{\footnotesize |
1462 |
cnh |
1.1 |
\begin{verbatim} |
1463 |
|
|
C |
1464 |
|
|
C |
1465 |
|
|
C Invocation from WRAPPER level... |
1466 |
|
|
C : |
1467 |
|
|
C : |
1468 |
|
|
C | |
1469 |
|
|
C |-THE_MODEL_MAIN :: Primary driver for the MITgcm algorithm |
1470 |
|
|
C | :: Called from WRAPPER level numerical |
1471 |
cnh |
1.4 |
C | :: code invocation routine. On entry |
1472 |
cnh |
1.1 |
C | :: to THE_MODEL_MAIN separate thread and |
1473 |
|
|
C | :: separate processes will have been established. |
1474 |
|
|
C | :: Each thread and process will have a unique ID |
1475 |
|
|
C | :: but as yet it will not be associated with a |
1476 |
|
|
C | :: specific region in decomposed discrete space. |
1477 |
|
|
C | |
1478 |
|
|
C |-INITIALISE_FIXED :: Set fixed model arrays such as topography, |
1479 |
|
|
C | | :: grid, solver matrices etc.. |
1480 |
|
|
C | | |
1481 |
|
|
C | |-INI_PARMS :: Routine to set kernel model parameters. |
1482 |
|
|
C | | :: By default kernel parameters are read from file |
1483 |
|
|
C | | :: "data" in directory in which code executes. |
1484 |
|
|
C | | |
1485 |
cnh |
1.4 |
C | |-MON_INIT :: Initializes monitor package ( see pkg/monitor ) |
1486 |
cnh |
1.1 |
C | | |
1487 |
cnh |
1.4 |
C | |-INI_GRID :: Control grid array (vert. and hori.) initialization. |
1488 |
cnh |
1.1 |
C | | | :: Grid arrays are held and described in GRID.h. |
1489 |
|
|
C | | | |
1490 |
cnh |
1.4 |
C | | |-INI_VERTICAL_GRID :: Initialize vertical grid arrays. |
1491 |
cnh |
1.1 |
C | | | |
1492 |
cnh |
1.4 |
C | | |-INI_CARTESIAN_GRID :: Cartesian horiz. grid initialization |
1493 |
cnh |
1.1 |
C | | | :: (calculate grid from kernel parameters). |
1494 |
|
|
C | | | |
1495 |
|
|
C | | |-INI_SPHERICAL_POLAR_GRID :: Spherical polar horiz. grid |
1496 |
cnh |
1.4 |
C | | | :: initialization (calculate grid from |
1497 |
cnh |
1.1 |
C | | | :: kernel parameters). |
1498 |
|
|
C | | | |
1499 |
|
|
C | | |-INI_CURVILINEAR_GRID :: General orthogonal, structured horiz. |
1500 |
cnh |
1.4 |
C | | :: grid initializations. ( input from raw |
1501 |
cnh |
1.1 |
C | | :: grid files, LONC.bin, DXF.bin etc... ) |
1502 |
|
|
C | | |
1503 |
|
|
C | |-INI_DEPTHS :: Read (from "bathyFile") or set bathymetry/orgography. |
1504 |
|
|
C | | |
1505 |
|
|
C | |-INI_MASKS_ETC :: Derive horizontal and vertical cell fractions and |
1506 |
|
|
C | | :: land masking for solid-fluid boundaries. |
1507 |
|
|
C | | |
1508 |
|
|
C | |-INI_LINEAR_PHSURF :: Set ref. surface Bo_surf |
1509 |
|
|
C | | |
1510 |
|
|
C | |-INI_CORI :: Set coriolis term. zero, f-plane, beta-plane, |
1511 |
cnh |
1.4 |
C | | :: sphere options are coded. |
1512 |
cnh |
1.1 |
C | | |
1513 |
|
|
C | |-PACAKGES_BOOT :: Start up the optional package environment. |
1514 |
|
|
C | | :: Runtime selection of active packages. |
1515 |
|
|
C | | |
1516 |
|
|
C | |-PACKAGES_READPARMS :: Call active package internal parameter load. |
1517 |
|
|
C | | | |
1518 |
|
|
C | | |-GMREDI_READPARMS :: GM Package. see pkg/gmredi |
1519 |
|
|
C | | |-KPP_READPARMS :: KPP Package. see pkg/kpp |
1520 |
|
|
C | | |-SHAP_FILT_READPARMS :: Shapiro filter package. see pkg/shap_filt |
1521 |
|
|
C | | |-OBCS_READPARMS :: Open bndy package. see pkg/obcs |
1522 |
|
|
C | | |-AIM_READPARMS :: Intermediate Atmos. pacakage. see pkg/aim |
1523 |
|
|
C | | |-COST_READPARMS :: Cost function package. see pkg/cost |
1524 |
|
|
C | | |-CTRL_INIT :: Control vector support package. see pkg/ctrl |
1525 |
|
|
C | | |-OPTIM_READPARMS :: Optimisation support package. see pkg/ctrl |
1526 |
|
|
C | | |-GRDCHK_READPARMS :: Gradient check package. see pkg/grdchk |
1527 |
|
|
C | | |-ECCO_READPARMS :: ECCO Support Package. see pkg/ecco |
1528 |
|
|
C | | |
1529 |
|
|
C | |-PACKAGES_CHECK |
1530 |
|
|
C | | | |
1531 |
|
|
C | | |-KPP_CHECK :: KPP Package. pkg/kpp |
1532 |
cnh |
1.4 |
C | | |-OBCS_CHECK :: Open bndy Package. pkg/obcs |
1533 |
cnh |
1.1 |
C | | |-GMREDI_CHECK :: GM Package. pkg/gmredi |
1534 |
|
|
C | | |
1535 |
|
|
C | |-PACKAGES_INIT_FIXED |
1536 |
|
|
C | | |-OBCS_INIT_FIXED :: Open bndy Package. see pkg/obcs |
1537 |
|
|
C | | |-FLT_INIT :: Floats Package. see pkg/flt |
1538 |
|
|
C | | |
1539 |
|
|
C | |-ZONAL_FILT_INIT :: FFT filter Package. see pkg/zonal_filt |
1540 |
|
|
C | | |
1541 |
|
|
C | |-INI_CG2D :: 2d con. grad solver initialisation. |
1542 |
|
|
C | | |
1543 |
|
|
C | |-INI_CG3D :: 3d con. grad solver initialisation. |
1544 |
|
|
C | | |
1545 |
|
|
C | |-CONFIG_SUMMARY :: Provide synopsis of kernel setup. |
1546 |
|
|
C | :: Includes annotated table of kernel |
1547 |
|
|
C | :: parameter settings. |
1548 |
|
|
C | |
1549 |
|
|
C |-CTRL_UNPACK :: Control vector support package. see pkg/ctrl |
1550 |
|
|
C | |
1551 |
|
|
C |-ADTHE_MAIN_LOOP :: Derivative evaluating form of main time stepping loop |
1552 |
cnh |
1.4 |
C ! :: Auotmatically generated by TAMC/TAF. |
1553 |
cnh |
1.1 |
C | |
1554 |
|
|
C |-CTRL_PACK :: Control vector support package. see pkg/ctrl |
1555 |
|
|
C | |
1556 |
|
|
C |-GRDCHK_MAIN :: Gradient check package. see pkg/grdchk |
1557 |
|
|
C | |
1558 |
|
|
C |-THE_MAIN_LOOP :: Main timestepping loop routine. |
1559 |
|
|
C | | |
1560 |
|
|
C | |-INITIALISE_VARIA :: Set the initial conditions for time evolving |
1561 |
|
|
C | | | :: variables |
1562 |
|
|
C | | | |
1563 |
|
|
C | | |-INI_LINEAR_PHISURF :: Set ref. surface Bo_surf |
1564 |
|
|
C | | | |
1565 |
|
|
C | | |-INI_CORI :: Set coriolis term. zero, f-plane, beta-plane, |
1566 |
cnh |
1.4 |
C | | | :: sphere options are coded. |
1567 |
cnh |
1.1 |
C | | | |
1568 |
|
|
C | | |-INI_CG2D :: 2d con. grad solver initialisation. |
1569 |
|
|
C | | |-INI_CG3D :: 3d con. grad solver initialisation. |
1570 |
|
|
C | | |-INI_MIXING :: Initialise diapycnal diffusivity. |
1571 |
|
|
C | | |-INI_DYNVARS :: Initialise to zero all DYNVARS.h arrays (dynamical |
1572 |
|
|
C | | | :: fields). |
1573 |
|
|
C | | | |
1574 |
cnh |
1.4 |
C | | |-INI_FIELDS :: Control initializing model fields to non-zero |
1575 |
cnh |
1.1 |
C | | | |-INI_VEL :: Initialize 3D flow field. |
1576 |
|
|
C | | | |-INI_THETA :: Set model initial temperature field. |
1577 |
|
|
C | | | |-INI_SALT :: Set model initial salinity field. |
1578 |
|
|
C | | | |-INI_PSURF :: Set model initial free-surface height/pressure. |
1579 |
|
|
C | | | |
1580 |
|
|
C | | |-INI_TR1 :: Set initial tracer 1 distribution. |
1581 |
|
|
C | | | |
1582 |
|
|
C | | |-THE_CORRECTION_STEP :: Step forward to next time step. |
1583 |
|
|
C | | | | :: Here applied to move restart conditions |
1584 |
|
|
C | | | | :: (saved in mid timestep) to correct level in |
1585 |
|
|
C | | | | :: time (only used for pre-c35). |
1586 |
|
|
C | | | | |
1587 |
|
|
C | | | |-CALC_GRAD_PHI_SURF :: Return DDx and DDy of surface pressure |
1588 |
|
|
C | | | |-CORRECTION_STEP :: Pressure correction to momentum |
1589 |
|
|
C | | | |-CYCLE_TRACER :: Move tracers forward in time. |
1590 |
|
|
C | | | |-OBCS_APPLY :: Open bndy package. see pkg/obcs |
1591 |
|
|
C | | | |-SHAP_FILT_APPLY :: Shapiro filter package. see pkg/shap_filt |
1592 |
|
|
C | | | |-ZONAL_FILT_APPLY :: FFT filter package. see pkg/zonal_filt |
1593 |
|
|
C | | | |-CONVECTIVE_ADJUSTMENT :: Control static instability mixing. |
1594 |
|
|
C | | | | |-FIND_RHO :: Find adjacent densities. |
1595 |
|
|
C | | | | |-CONVECT :: Mix static instability. |
1596 |
|
|
C | | | | |-TIMEAVE_CUMULATE :: Update convection statistics. |
1597 |
|
|
C | | | | |
1598 |
|
|
C | | | |-CALC_EXACT_ETA :: Change SSH to flow divergence. |
1599 |
|
|
C | | | |
1600 |
|
|
C | | |-CONVECTIVE_ADJUSTMENT_INI :: Control static instability mixing |
1601 |
|
|
C | | | | :: Extra time history interactions. |
1602 |
|
|
C | | | | |
1603 |
|
|
C | | | |-FIND_RHO :: Find adjacent densities. |
1604 |
|
|
C | | | |-CONVECT :: Mix static instability. |
1605 |
|
|
C | | | |-TIMEAVE_CUMULATE :: Update convection statistics. |
1606 |
|
|
C | | | |
1607 |
|
|
C | | |-PACKAGES_INIT_VARIABLES :: Does initialisation of time evolving |
1608 |
|
|
C | | | | :: package data. |
1609 |
|
|
C | | | | |
1610 |
|
|
C | | | |-GMREDI_INIT :: GM package. ( see pkg/gmredi ) |
1611 |
|
|
C | | | |-KPP_INIT :: KPP package. ( see pkg/kpp ) |
1612 |
|
|
C | | | |-KPP_OPEN_DIAGS |
1613 |
|
|
C | | | |-OBCS_INIT_VARIABLES :: Open bndy. package. ( see pkg/obcs ) |
1614 |
|
|
C | | | |-AIM_INIT :: Interm. atmos package. ( see pkg/aim ) |
1615 |
|
|
C | | | |-CTRL_MAP_INI :: Control vector package.( see pkg/ctrl ) |
1616 |
|
|
C | | | |-COST_INIT :: Cost function package. ( see pkg/cost ) |
1617 |
|
|
C | | | |-ECCO_INIT :: ECCO support package. ( see pkg/ecco ) |
1618 |
|
|
C | | | |-INI_FORCING :: Set model initial forcing fields. |
1619 |
|
|
C | | | | :: Either set in-line or from file as shown. |
1620 |
|
|
C | | | |-READ_FLD_XY_RS(zonalWindFile) |
1621 |
|
|
C | | | |-READ_FLD_XY_RS(meridWindFile) |
1622 |
|
|
C | | | |-READ_FLD_XY_RS(surfQFile) |
1623 |
|
|
C | | | |-READ_FLD_XY_RS(EmPmRfile) |
1624 |
|
|
C | | | |-READ_FLD_XY_RS(thetaClimFile) |
1625 |
|
|
C | | | |-READ_FLD_XY_RS(saltClimFile) |
1626 |
|
|
C | | | |-READ_FLD_XY_RS(surfQswFile) |
1627 |
|
|
C | | | |
1628 |
|
|
C | | |-CALC_SURF_DR :: Calculate the new surface level thickness. |
1629 |
|
|
C | | |-UPDATE_SURF_DR :: Update the surface-level thickness fraction. |
1630 |
|
|
C | | |-UPDATE_CG2D :: Update 2d conjugate grad. for Free-Surf. |
1631 |
|
|
C | | |-STATE_SUMMARY :: Summarize model prognostic variables. |
1632 |
|
|
C | | |-TIMEAVE_STATVARS :: Time averaging package ( see pkg/timeave ). |
1633 |
|
|
C | | |
1634 |
|
|
C | |-WRITE_STATE :: Controlling routine for IO to dump model state. |
1635 |
|
|
C | | |-WRITE_REC_XYZ_RL :: Single file I/O |
1636 |
|
|
C | | |-WRITE_FLD_XYZ_RL :: Multi-file I/O |
1637 |
|
|
C | | |
1638 |
|
|
C | |-MONITOR :: Monitor state ( see pkg/monitor ) |
1639 |
|
|
C | |-CTRL_MAP_FORCING :: Control vector support package. ( see pkg/ctrl ) |
1640 |
|
|
C====|>| |
1641 |
|
|
C====|>| **************************** |
1642 |
|
|
C====|>| BEGIN MAIN TIMESTEPPING LOOP |
1643 |
|
|
C====|>| **************************** |
1644 |
|
|
C====|>| |
1645 |
|
|
C/\ | |-FORWARD_STEP :: Step forward a time-step ( AT LAST !!! ) |
1646 |
|
|
C/\ | | | |
1647 |
|
|
C/\ | | |-DUMMY_IN_STEPPING :: autodiff package ( pkg/autoduff ). |
1648 |
|
|
C/\ | | |-CALC_EXACT_ETA :: Change SSH to flow divergence. |
1649 |
|
|
C/\ | | |-CALC_SURF_DR :: Calculate the new surface level thickness. |
1650 |
|
|
C/\ | | |-EXF_GETFORCING :: External forcing package. ( pkg/exf ) |
1651 |
|
|
C/\ | | |-EXTERNAL_FIELDS_LOAD :: Control loading time dep. external data. |
1652 |
cnh |
1.4 |
C/\ | | | | :: Simple interpolation between end-points |
1653 |
cnh |
1.1 |
C/\ | | | | :: for forcing datasets. |
1654 |
|
|
C/\ | | | | |
1655 |
|
|
C/\ | | | |-EXCH :: Sync forcing. in overlap regions. |
1656 |
|
|
C/\ | | | |
1657 |
|
|
C/\ | | |-THERMODYNAMICS :: theta, salt + tracer equations driver. |
1658 |
|
|
C/\ | | | | |
1659 |
|
|
C/\ | | | |-INTEGRATE_FOR_W :: Integrate for vertical velocity. |
1660 |
|
|
C/\ | | | |-OBCS_APPLY_W :: Open bndy. package ( see pkg/obcs ). |
1661 |
|
|
C/\ | | | |-FIND_RHO :: Calculates [rho(S,T,z)-Rhonil] of a slice |
1662 |
|
|
C/\ | | | |-GRAD_SIGMA :: Calculate isoneutral gradients |
1663 |
|
|
C/\ | | | |-CALC_IVDC :: Set Implicit Vertical Diffusivity for Convection |
1664 |
|
|
C/\ | | | | |
1665 |
|
|
C/\ | | | |-OBCS_CALC :: Open bndy. package ( see pkg/obcs ). |
1666 |
|
|
C/\ | | | |-EXTERNAL_FORCING_SURF:: Accumulates appropriately dimensioned |
1667 |
|
|
C/\ | | | | :: forcing terms. |
1668 |
|
|
C/\ | | | | |
1669 |
|
|
C/\ | | | |-GMREDI_CALC_TENSOR :: GM package ( see pkg/gmredi ). |
1670 |
|
|
C/\ | | | |-GMREDI_CALC_TENSOR_DUMMY :: GM package ( see pkg/gmredi ). |
1671 |
|
|
C/\ | | | |-KPP_CALC :: KPP package ( see pkg/kpp ). |
1672 |
|
|
C/\ | | | |-KPP_CALC_DUMMY :: KPP package ( see pkg/kpp ). |
1673 |
|
|
C/\ | | | |-AIM_DO_ATMOS_PHYSICS :: Intermed. atmos package ( see pkg/aim ). |
1674 |
|
|
C/\ | | | |-GAD_ADVECTION :: Generalised advection driver (multi-dim |
1675 |
|
|
C/\ | | | | advection case) (see pkg/gad). |
1676 |
|
|
C/\ | | | |-CALC_COMMON_FACTORS :: Calculate common data (such as volume flux) |
1677 |
|
|
C/\ | | | |-CALC_DIFFUSIVITY :: Calculate net vertical diffusivity |
1678 |
|
|
C/\ | | | | | |
1679 |
|
|
C/\ | | | | |-GMREDI_CALC_DIFF :: GM package ( see pkg/gmredi ). |
1680 |
|
|
C/\ | | | | |-KPP_CALC_DIFF :: KPP package ( see pkg/kpp ). |
1681 |
|
|
C/\ | | | | |
1682 |
|
|
C/\ | | | |-CALC_GT :: Calculate the temperature tendency terms |
1683 |
|
|
C/\ | | | | | |
1684 |
|
|
C/\ | | | | |-GAD_CALC_RHS :: Generalised advection package |
1685 |
|
|
C/\ | | | | | :: ( see pkg/gad ) |
1686 |
|
|
C/\ | | | | |-EXTERNAL_FORCING_T :: Problem specific forcing for temperature. |
1687 |
|
|
C/\ | | | | |-ADAMS_BASHFORTH2 :: Extrapolate tendencies forward in time. |
1688 |
|
|
C/\ | | | | |-FREESURF_RESCALE_G :: Re-scale Gt for free-surface height. |
1689 |
|
|
C/\ | | | | |
1690 |
|
|
C/\ | | | |-TIMESTEP_TRACER :: Step tracer field forward in time |
1691 |
|
|
C/\ | | | | |
1692 |
|
|
C/\ | | | |-CALC_GS :: Calculate the salinity tendency terms |
1693 |
|
|
C/\ | | | | | |
1694 |
|
|
C/\ | | | | |-GAD_CALC_RHS :: Generalised advection package |
1695 |
|
|
C/\ | | | | | :: ( see pkg/gad ) |
1696 |
|
|
C/\ | | | | |-EXTERNAL_FORCING_S :: Problem specific forcing for salt. |
1697 |
|
|
C/\ | | | | |-ADAMS_BASHFORTH2 :: Extrapolate tendencies forward in time. |
1698 |
|
|
C/\ | | | | |-FREESURF_RESCALE_G :: Re-scale Gs for free-surface height. |
1699 |
|
|
C/\ | | | | |
1700 |
|
|
C/\ | | | |-TIMESTEP_TRACER :: Step tracer field forward in time |
1701 |
|
|
C/\ | | | | |
1702 |
|
|
C/\ | | | |-CALC_GTR1 :: Calculate other tracer(s) tendency terms |
1703 |
|
|
C/\ | | | | | |
1704 |
|
|
C/\ | | | | |-GAD_CALC_RHS :: Generalised advection package |
1705 |
|
|
C/\ | | | | | :: ( see pkg/gad ) |
1706 |
|
|
C/\ | | | | |-EXTERNAL_FORCING_TR:: Problem specific forcing for tracer. |
1707 |
|
|
C/\ | | | | |-ADAMS_BASHFORTH2 :: Extrapolate tendencies forward in time. |
1708 |
|
|
C/\ | | | | |-FREESURF_RESCALE_G :: Re-scale Gs for free-surface height. |
1709 |
|
|
C/\ | | | | |
1710 |
|
|
C/\ | | | |-TIMESTEP_TRACER :: Step tracer field forward in time |
1711 |
|
|
C/\ | | | |-OBCS_APPLY_TS :: Open bndy. package (see pkg/obcs ). |
1712 |
|
|
C/\ | | | |-FREEZE :: Limit range of temperature. |
1713 |
|
|
C/\ | | | | |
1714 |
|
|
C/\ | | | |-IMPLDIFF :: Solve vertical implicit diffusion equation. |
1715 |
|
|
C/\ | | | |-OBCS_APPLY_TS :: Open bndy. package (see pkg/obcs ). |
1716 |
|
|
C/\ | | | | |
1717 |
|
|
C/\ | | | |-AIM_AIM2DYN_EXCHANGES :: Inetermed. atmos (see pkg/aim). |
1718 |
|
|
C/\ | | | |-EXCH :: Update overlaps |
1719 |
|
|
C/\ | | | |
1720 |
|
|
C/\ | | |-DYNAMICS :: Momentum equations driver. |
1721 |
|
|
C/\ | | | | |
1722 |
|
|
C/\ | | | |-CALC_GRAD_PHI_SURF :: Calculate the gradient of the surface |
1723 |
|
|
C/\ | | | | Potential anomaly. |
1724 |
|
|
C/\ | | | |-CALC_VISCOSITY :: Calculate net vertical viscosity |
1725 |
|
|
C/\ | | | | |-KPP_CALC_VISC :: KPP package ( see pkg/kpp ). |
1726 |
|
|
C/\ | | | | |
1727 |
|
|
C/\ | | | |-CALC_PHI_HYD :: Integrate the hydrostatic relation. |
1728 |
|
|
C/\ | | | |-MOM_FLUXFORM :: Flux form mom eqn. package ( see |
1729 |
|
|
C/\ | | | | pkg/mom_fluxform ). |
1730 |
|
|
C/\ | | | |-MOM_VECINV :: Vector invariant form mom eqn. package ( see |
1731 |
|
|
C/\ | | | | pkg/mom_vecinv ). |
1732 |
|
|
C/\ | | | |-TIMESTEP :: Step momentum fields forward in time |
1733 |
|
|
C/\ | | | |-OBCS_APPLY_UV :: Open bndy. package (see pkg/obcs ). |
1734 |
|
|
C/\ | | | | |
1735 |
|
|
C/\ | | | |-IMPLDIFF :: Solve vertical implicit diffusion equation. |
1736 |
|
|
C/\ | | | |-OBCS_APPLY_UV :: Open bndy. package (see pkg/obcs ). |
1737 |
|
|
C/\ | | | | |
1738 |
|
|
C/\ | | | |-TIMEAVE_CUMUL_1T :: Time averaging package ( see pkg/timeave ). |
1739 |
|
|
C/\ | | | |-TIMEAVE_CUMUATE :: Time averaging package ( see pkg/timeave ). |
1740 |
|
|
C/\ | | | |-DEBUG_STATS_RL :: Quick debug package ( see pkg/debug ). |
1741 |
|
|
C/\ | | | |
1742 |
|
|
C/\ | | |-CALC_GW :: vert. momentum tendency terms ( NH, QH only ). |
1743 |
|
|
C/\ | | | |
1744 |
|
|
C/\ | | |-UPDATE_SURF_DR :: Update the surface-level thickness fraction. |
1745 |
|
|
C/\ | | | |
1746 |
|
|
C/\ | | |-UPDATE_CG2D :: Update 2d conjugate grad. for Free-Surf. |
1747 |
|
|
C/\ | | | |
1748 |
|
|
C/\ | | |-SOLVE_FOR_PRESSURE :: Find surface pressure. |
1749 |
|
|
C/\ | | | |-CALC_DIV_GHAT :: Form the RHS of the surface pressure eqn. |
1750 |
|
|
C/\ | | | |-CG2D :: Two-dim pre-con. conjugate-gradient. |
1751 |
|
|
C/\ | | | |-CG3D :: Three-dim pre-con. conjugate-gradient solver. |
1752 |
|
|
C/\ | | | |
1753 |
|
|
C/\ | | |-THE_CORRECTION_STEP :: Step forward to next time step. |
1754 |
|
|
C/\ | | | | |
1755 |
|
|
C/\ | | | |-CALC_GRAD_PHI_SURF :: Return DDx and DDy of surface pressure |
1756 |
|
|
C/\ | | | |-CORRECTION_STEP :: Pressure correction to momentum |
1757 |
|
|
C/\ | | | |-CYCLE_TRACER :: Move tracers forward in time. |
1758 |
|
|
C/\ | | | |-OBCS_APPLY :: Open bndy package. see pkg/obcs |
1759 |
|
|
C/\ | | | |-SHAP_FILT_APPLY :: Shapiro filter package. see pkg/shap_filt |
1760 |
|
|
C/\ | | | |-ZONAL_FILT_APPLY :: FFT filter package. see pkg/zonal_filt |
1761 |
|
|
C/\ | | | |-CONVECTIVE_ADJUSTMENT :: Control static instability mixing. |
1762 |
|
|
C/\ | | | | |-FIND_RHO :: Find adjacent densities. |
1763 |
|
|
C/\ | | | | |-CONVECT :: Mix static instability. |
1764 |
|
|
C/\ | | | | |-TIMEAVE_CUMULATE :: Update convection statistics. |
1765 |
|
|
C/\ | | | | |
1766 |
|
|
C/\ | | | |-CALC_EXACT_ETA :: Change SSH to flow divergence. |
1767 |
|
|
C/\ | | | |
1768 |
|
|
C/\ | | |-DO_FIELDS_BLOCKING_EXCHANGES :: Sync up overlap regions. |
1769 |
|
|
C/\ | | | |-EXCH |
1770 |
|
|
C/\ | | | |
1771 |
|
|
C/\ | | |-FLT_MAIN :: Float package ( pkg/flt ). |
1772 |
|
|
C/\ | | | |
1773 |
|
|
C/\ | | |-MONITOR :: Monitor package ( pkg/monitor ). |
1774 |
|
|
C/\ | | | |
1775 |
|
|
C/\ | | |-DO_THE_MODEL_IO :: Standard diagnostic I/O. |
1776 |
|
|
C/\ | | | |-WRITE_STATE :: Core state I/O |
1777 |
|
|
C/\ | | | |-TIMEAVE_STATV_WRITE :: Time averages. see pkg/timeave |
1778 |
|
|
C/\ | | | |-AIM_WRITE_DIAGS :: Intermed. atmos diags. see pkg/aim |
1779 |
|
|
C/\ | | | |-GMREDI_DIAGS :: GM diags. see pkg/gmredi |
1780 |
|
|
C/\ | | | |-KPP_DO_DIAGS :: KPP diags. see pkg/kpp |
1781 |
|
|
C/\ | | | |
1782 |
|
|
C/\ | | |-WRITE_CHECKPOINT :: Do I/O for restart files. |
1783 |
|
|
C/\ | | |
1784 |
|
|
C/\ | |-COST_TILE :: Cost function package. ( see pkg/cost ) |
1785 |
|
|
C<===|=| |
1786 |
|
|
C<===|=| ************************** |
1787 |
|
|
C<===|=| END MAIN TIMESTEPPING LOOP |
1788 |
|
|
C<===|=| ************************** |
1789 |
|
|
C<===|=| |
1790 |
|
|
C | |-COST_FINAL :: Cost function package. ( see pkg/cost ) |
1791 |
|
|
C | |
1792 |
|
|
C |-WRITE_CHECKPOINT :: Final state storage, for restart. |
1793 |
|
|
C | |
1794 |
|
|
C |-TIMER_PRINTALL :: Computational timing summary |
1795 |
|
|
C | |
1796 |
|
|
C |-COMM_STATS :: Summarise inter-proc and inter-thread communication |
1797 |
|
|
C :: events. |
1798 |
|
|
C |
1799 |
|
|
\end{verbatim} |
1800 |
heimbach |
1.8 |
} |
1801 |
cnh |
1.1 |
|
1802 |
|
|
\subsection{Measuring and Characterizing Performance} |
1803 |
|
|
|
1804 |
|
|
TO BE DONE (CNH) |
1805 |
|
|
|
1806 |
|
|
\subsection{Estimating Resource Requirements} |
1807 |
|
|
|
1808 |
|
|
TO BE DONE (CNH) |
1809 |
|
|
|
1810 |
|
|
\subsubsection{Atlantic 1/6 degree example} |
1811 |
|
|
\subsubsection{Dry Run testing} |
1812 |
|
|
\subsubsection{Adjoint Resource Requirements} |
1813 |
|
|
\subsubsection{State Estimation Environment Resources} |
1814 |
|
|
|