1 |
C $Header: /u/gcmpack/MITgcm/pkg/diagnostics/diagnostics_main_init.F,v 1.36 2010/01/11 19:45:11 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "DIAG_OPTIONS.h" |
5 |
|
6 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
7 |
CBOP 0 |
8 |
C !ROUTINE: DIAGNOSTICS_MAIN_INIT |
9 |
|
10 |
C !INTERFACE: |
11 |
SUBROUTINE DIAGNOSTICS_MAIN_INIT( myThid ) |
12 |
|
13 |
C !DESCRIPTION: |
14 |
C Initialize available diagnostics list for variables of the main code |
15 |
C (not part of a package): set the following attributes: |
16 |
C name (=cdiag), parsing code (=gdiag), units (=udiag), and title (=tdiag) |
17 |
C Notes: 1) diagnostics defined here do not require any EQUIVALENCE |
18 |
C since they get filled-in with S/R FILL_DIAGNOSTICS |
19 |
C 2) GDIAG is defined as character*16 and can be to character*1 |
20 |
C parse(16) with the following codes currently defined: |
21 |
|
22 |
C \begin{center} |
23 |
C \begin{tabular}[h]{|c|c|}\hline |
24 |
C \textbf{Positions} & \textbf{Characters} |
25 |
C & \textbf{Meanings} \\\hline |
26 |
C parse(1) & S & scalar \\ |
27 |
C & U & vector component in X direction \\ |
28 |
C & V & vector component in Y direction \\ |
29 |
C & W & vector component in vertical direction \\ |
30 |
C parse(2) & U & C-grid U-Point \\ |
31 |
C & V & C-grid V-Point \\ |
32 |
C & M & C-grid Mass Point \\ |
33 |
C & Z & C-grid Corner Point \\ |
34 |
C parse(3) & & Used for Level Integrated output: cumulate levels \\ |
35 |
C & r & same but cumulate product by model level thickness \\ |
36 |
C & R & same but cumulate product by hFac & level thickness \\ |
37 |
C parse(4) & P & positive definite \\ |
38 |
C parse(5 ) & C & with counter array \\ |
39 |
C & D & disable an array for output \\ |
40 |
C parse(6--8) & '123' & retired, formerly: 3-digit mate number \\ |
41 |
C parse(9) & U & model-level plus 1/2 \\ |
42 |
C & M & model-level middle \\ |
43 |
C & L & model-level minus 1/2 \\ |
44 |
C parse(10) & 0 & levels = 0 \\ |
45 |
C & 1 & levels = 1 \\ |
46 |
C & R & levels = Nr \\ |
47 |
C & L & levels = MAX(Nr,NrPhys) \\ |
48 |
C & M & levels = MAX(Nr,NrPhys) - 1 \\ |
49 |
C & G & levels = Ground_level Number \\ |
50 |
C & I & levels = sea-Ice_level Number \\ |
51 |
C & X & free levels option (need to be set explicitly) \\ |
52 |
C \end{tabular} |
53 |
C \end{center} |
54 |
|
55 |
C !USES: |
56 |
IMPLICIT NONE |
57 |
#include "SIZE.h" |
58 |
#include "EEPARAMS.h" |
59 |
#include "PARAMS.h" |
60 |
|
61 |
C !INPUT PARAMETERS: |
62 |
INTEGER myThid |
63 |
CEOP |
64 |
|
65 |
C !LOCAL VARIABLES: |
66 |
C rTitle :: r-coordinate title |
67 |
C eTitle :: free-surface title |
68 |
C fTitle :: fixed boundary title |
69 |
C pTitle :: "Phi" title |
70 |
C sTitle :: "salt" title |
71 |
INTEGER diagNum |
72 |
INTEGER diagMate |
73 |
CHARACTER*8 diagName |
74 |
CHARACTER*16 diagCode |
75 |
CHARACTER*16 diagUnits |
76 |
CHARACTER*(80) diagTitle |
77 |
CHARACTER*2 rUnit2c |
78 |
CHARACTER*4 tUnit4c, sUnit4c |
79 |
CHARACTER*(10) rTitle, eTitle, fTitle |
80 |
CHARACTER*(20) pTitle, sTitle |
81 |
|
82 |
CHARACTER*(16) DIAGS_MK_UNITS |
83 |
EXTERNAL DIAGS_MK_UNITS |
84 |
CHARACTER*(80) DIAGS_MK_TITLE |
85 |
EXTERNAL DIAGS_MK_TITLE |
86 |
|
87 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
88 |
C For each output variable, |
89 |
C specify Name (cdiag, 8c), Descriptions (tdiag, *c), Units (udiag, 16c) |
90 |
C and Type/Parms (location on C grid, 2D/3D, ...) (gdiag, 16c) |
91 |
C---------------------------------------------------------------------- |
92 |
|
93 |
IF ( usingPCoords ) THEN |
94 |
rUnit2c= 'Pa' |
95 |
rTitle = ' Pressure ' |
96 |
pTitle = ' Geopotential ' |
97 |
ELSE |
98 |
rUnit2c= 'm ' |
99 |
rTitle = ' Height ' |
100 |
pTitle = 'Pressure Pot.(p/rho)' |
101 |
ENDIF |
102 |
IF ( fluidIsAir ) THEN |
103 |
tUnit4c= 'K ' |
104 |
sUnit4c= 'g/kg' |
105 |
sTitle = ' Specific Humidity ' |
106 |
#ifdef ALLOW_FIZHI |
107 |
IF (useFIZHI) sUnit4c= 'kg/kg' |
108 |
#endif /* ALLOW_FIZHI */ |
109 |
ELSE |
110 |
tUnit4c= 'degC' |
111 |
sUnit4c= 'psu ' |
112 |
sTitle = ' Salinity ' |
113 |
ENDIF |
114 |
C- free-surface (eTitle) and fixed-boundary (fTitle) position: |
115 |
IF ( fluidIsAir ) THEN |
116 |
eTitle = ' Surface ' |
117 |
fTitle = ' Top ' |
118 |
ELSEIF ( usingPCoords ) THEN |
119 |
eTitle = ' Bottom ' |
120 |
fTitle = ' Surface ' |
121 |
ELSE |
122 |
eTitle = ' Surface ' |
123 |
fTitle = ' Bottom ' |
124 |
ENDIF |
125 |
|
126 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
127 |
C- state variables of the main code (and related quadratic var): |
128 |
|
129 |
|
130 |
diagName = 'ETAN ' |
131 |
diagTitle = DIAGS_MK_TITLE( eTitle//rTitle//' Anomaly', myThid ) |
132 |
c IF ( fluidIsWater .AND. usingZCoords ) |
133 |
c &diagTitle = 'Sea Surface Elevation' |
134 |
diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid ) |
135 |
diagCode = 'SM M1 ' |
136 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
137 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
138 |
|
139 |
diagName = 'ETANSQ ' |
140 |
diagTitle = DIAGS_MK_TITLE( 'Square of '//eTitle//rTitle |
141 |
I //' Anomaly', myThid ) |
142 |
diagUnits = DIAGS_MK_UNITS( rUnit2c//'^2', myThid ) |
143 |
diagCode = 'SM P M1 ' |
144 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
145 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
146 |
|
147 |
diagName = 'DETADT2 ' |
148 |
diagTitle = DIAGS_MK_TITLE( 'Square of '//eTitle//rTitle |
149 |
I //' Anomaly Tendency', myThid ) |
150 |
diagUnits = DIAGS_MK_UNITS( rUnit2c//'^2/s^2', myThid ) |
151 |
diagCode = 'SM M1 ' |
152 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
153 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
154 |
|
155 |
diagName = 'THETA ' |
156 |
diagTitle = 'Potential Temperature' |
157 |
diagUnits = DIAGS_MK_UNITS( tUnit4c, myThid ) |
158 |
diagCode = 'SMR MR ' |
159 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
160 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
161 |
|
162 |
c diagName = 'SST ' |
163 |
c diagTitle = 'Sea Surface Temperature (degC,K)' |
164 |
c diagUnits = DIAGS_MK_UNITS( tUnit4c, myThid ) |
165 |
c diagCode = 'SM M1 ' |
166 |
c CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
167 |
c I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
168 |
|
169 |
diagName = 'SALT ' |
170 |
diagTitle = DIAGS_MK_TITLE( sTitle, myThid ) |
171 |
diagUnits = DIAGS_MK_UNITS( sUnit4c, myThid ) |
172 |
diagCode = 'SMR MR ' |
173 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
174 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
175 |
|
176 |
diagName = 'RELHUM ' |
177 |
diagTitle = 'Relative Humidity' |
178 |
diagUnits = 'percent ' |
179 |
diagCode = 'SMR MR ' |
180 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
181 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
182 |
|
183 |
c diagName = 'SSS ' |
184 |
c diagTitle = 'Sea Surface Salinity ' |
185 |
c diagUnits = DIAGS_MK_UNITS( sUnit4c, myThid ) |
186 |
c diagCode = 'SM M1 ' |
187 |
c CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
188 |
c I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
189 |
|
190 |
IF ( fluidIsWater ) THEN |
191 |
diagName = 'SALTanom' |
192 |
diagTitle = 'Salt anomaly (=SALT-35; g/kg)' |
193 |
diagUnits = DIAGS_MK_UNITS( sUnit4c, myThid ) |
194 |
diagCode = 'SMR MR ' |
195 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
196 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
197 |
ENDIF |
198 |
|
199 |
diagName = 'UVEL ' |
200 |
diagTitle = 'Zonal Component of Velocity (m/s)' |
201 |
diagUnits = 'm/s ' |
202 |
diagCode = 'UUR MR ' |
203 |
diagMate = diagNum + 2 |
204 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
205 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
206 |
|
207 |
diagName = 'VVEL ' |
208 |
diagTitle = 'Meridional Component of Velocity (m/s)' |
209 |
diagUnits = 'm/s ' |
210 |
diagCode = 'VVR MR ' |
211 |
diagMate = diagNum |
212 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
213 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
214 |
|
215 |
diagName = 'WVEL ' |
216 |
diagTitle = 'Vertical Component of Velocity (r_units/s)' |
217 |
diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s', myThid ) |
218 |
diagCode = 'WM LR ' |
219 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
220 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
221 |
|
222 |
diagName = 'THETASQ ' |
223 |
diagTitle = 'Square of Potential Temperature' |
224 |
diagUnits = DIAGS_MK_UNITS( tUnit4c//'^2', myThid ) |
225 |
diagCode = 'SMRP MR ' |
226 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
227 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
228 |
|
229 |
diagName = 'SALTSQ ' |
230 |
diagTitle = DIAGS_MK_TITLE( 'Square of '//sTitle, myThid ) |
231 |
diagUnits = DIAGS_MK_UNITS( '('//sUnit4c//')^2', myThid ) |
232 |
diagCode = 'SMRP MR ' |
233 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
234 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
235 |
|
236 |
IF ( fluidIsWater ) THEN |
237 |
diagName = 'SALTSQan' |
238 |
diagTitle = 'Square of Salt anomaly (=(SALT-35)^2 (g^2/kg^2)' |
239 |
diagUnits = DIAGS_MK_UNITS( '('//sUnit4c//')^2', myThid ) |
240 |
diagCode = 'SMRP MR ' |
241 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
242 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
243 |
ENDIF |
244 |
|
245 |
diagName = 'UVELSQ ' |
246 |
diagTitle = 'Square of Zonal Comp of Velocity (m^2/s^2)' |
247 |
diagUnits = 'm^2/s^2 ' |
248 |
diagCode = 'UURP MR ' |
249 |
diagMate = diagNum + 2 |
250 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
251 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
252 |
|
253 |
diagName = 'VVELSQ ' |
254 |
diagTitle = 'Square of Meridional Comp of Velocity (m^2/s^2)' |
255 |
diagUnits = 'm^2/s^2 ' |
256 |
diagCode = 'VVRP MR ' |
257 |
diagMate = diagNum |
258 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
259 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
260 |
|
261 |
diagName = 'WVELSQ ' |
262 |
diagTitle = 'Square of Vertical Comp of Velocity' |
263 |
diagUnits = DIAGS_MK_UNITS( rUnit2c//'^2/s^2', myThid ) |
264 |
diagCode = 'WM P LR ' |
265 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
266 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
267 |
|
268 |
diagName = 'UE_VEL_C' |
269 |
diagTitle = 'Eastward Velocity (m/s) (cell center)' |
270 |
diagUnits = 'm/s ' |
271 |
diagCode = 'UMR MR ' |
272 |
diagMate = diagNum + 2 |
273 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
274 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
275 |
|
276 |
diagName = 'VN_VEL_C' |
277 |
diagTitle = 'Northward Velocity (m/s) (cell center)' |
278 |
diagUnits = 'm/s ' |
279 |
diagCode = 'VMR MR ' |
280 |
diagMate = diagNum |
281 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
282 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
283 |
|
284 |
diagName = 'UV_VEL_C' |
285 |
diagTitle ='Product of horizontal Comp of velocity (cell center)' |
286 |
diagUnits = 'm^2/s^2 ' |
287 |
diagCode = 'UMR MR ' |
288 |
diagMate = diagNum + 1 |
289 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
290 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
291 |
|
292 |
diagName = 'UV_VEL_Z' |
293 |
diagTitle = 'Meridional Transport of Zonal Momentum (m^2/s^2)' |
294 |
diagUnits = 'm^2/s^2 ' |
295 |
diagCode = 'UZR MR ' |
296 |
diagMate = diagNum + 1 |
297 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
298 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
299 |
|
300 |
diagName = 'WU_VEL ' |
301 |
diagTitle = 'Vertical Transport of Zonal Momentum' |
302 |
diagUnits = DIAGS_MK_UNITS( 'm.'//rUnit2c//'/s^2', myThid ) |
303 |
diagCode = 'WU LR ' |
304 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
305 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
306 |
|
307 |
diagName = 'WV_VEL ' |
308 |
diagTitle ='Vertical Transport of Meridional Momentum' |
309 |
diagUnits = DIAGS_MK_UNITS( 'm.'//rUnit2c//'/s^2', myThid ) |
310 |
diagCode = 'WV LR ' |
311 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
312 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
313 |
|
314 |
diagName = 'UVELMASS' |
315 |
diagTitle = 'Zonal Mass-Weighted Comp of Velocity (m/s)' |
316 |
diagUnits = 'm/s ' |
317 |
diagCode = 'UUr MR ' |
318 |
diagMate = diagNum + 2 |
319 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
320 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
321 |
|
322 |
diagName = 'VVELMASS' |
323 |
diagTitle = 'Meridional Mass-Weighted Comp of Velocity (m/s)' |
324 |
diagUnits = 'm/s ' |
325 |
diagCode = 'VVr MR ' |
326 |
diagMate = diagNum |
327 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
328 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
329 |
|
330 |
diagName = 'WVELMASS' |
331 |
diagTitle = 'Vertical Mass-Weighted Comp of Velocity' |
332 |
diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s', myThid ) |
333 |
diagCode = 'WM LR ' |
334 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
335 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
336 |
|
337 |
diagName = 'UTHMASS ' |
338 |
diagTitle = 'Zonal Mass-Weight Transp of Pot Temp' |
339 |
diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid ) |
340 |
diagCode = 'UUr MR ' |
341 |
diagMate = diagNum + 2 |
342 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
343 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
344 |
|
345 |
diagName = 'VTHMASS ' |
346 |
diagTitle = 'Meridional Mass-Weight Transp of Pot Temp' |
347 |
diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid ) |
348 |
diagCode = 'VVr MR ' |
349 |
diagMate = diagNum |
350 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
351 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
352 |
|
353 |
diagName = 'WTHMASS ' |
354 |
diagTitle = 'Vertical Mass-Weight Transp of Pot Temp (K.m/s)' |
355 |
diagUnits = DIAGS_MK_UNITS(tUnit4c//'.'//rUnit2c//'/s', myThid ) |
356 |
diagCode = 'WM LR ' |
357 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
358 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
359 |
|
360 |
diagName = 'USLTMASS' |
361 |
diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of ' |
362 |
I //sTitle, myThid ) |
363 |
diagUnits = DIAGS_MK_UNITS(sUnit4c//'.m/s', myThid ) |
364 |
diagCode = 'UUr MR ' |
365 |
diagMate = diagNum + 2 |
366 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
367 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
368 |
|
369 |
diagName = 'VSLTMASS' |
370 |
diagTitle = DIAGS_MK_TITLE( 'Meridional Mass-Weight Transp of ' |
371 |
I //sTitle, myThid ) |
372 |
diagUnits = DIAGS_MK_UNITS(sUnit4c//'.m/s', myThid ) |
373 |
diagCode = 'VVr MR ' |
374 |
diagMate = diagNum |
375 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
376 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
377 |
|
378 |
diagName = 'WSLTMASS' |
379 |
diagTitle = DIAGS_MK_TITLE( 'Vertical Mass-Weight Transp of ' |
380 |
I //sTitle, myThid ) |
381 |
diagUnits = DIAGS_MK_UNITS(sUnit4c//'.'//rUnit2c//'/s', myThid ) |
382 |
diagCode = 'WM LR ' |
383 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
384 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
385 |
|
386 |
diagName = 'UVELTH ' |
387 |
diagTitle = 'Zonal Transport of Pot Temp' |
388 |
diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid ) |
389 |
diagCode = 'UUR MR ' |
390 |
diagMate = diagNum + 2 |
391 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
392 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
393 |
|
394 |
diagName = 'VVELTH ' |
395 |
diagTitle = 'Meridional Transport of Pot Temp' |
396 |
diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid ) |
397 |
diagCode = 'VVR MR ' |
398 |
diagMate = diagNum |
399 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
400 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
401 |
|
402 |
diagName = 'WVELTH ' |
403 |
diagTitle = 'Vertical Transport of Pot Temp' |
404 |
diagUnits = DIAGS_MK_UNITS(tUnit4c//'.'//rUnit2c//'/s', myThid ) |
405 |
diagCode = 'WM LR ' |
406 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
407 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
408 |
|
409 |
diagName = 'UVELSLT ' |
410 |
diagTitle = DIAGS_MK_TITLE( 'Zonal Transport of ' |
411 |
I //sTitle, myThid ) |
412 |
diagUnits = DIAGS_MK_UNITS( sUnit4c//'.m/s', myThid ) |
413 |
diagCode = 'UUR MR ' |
414 |
diagMate = diagNum + 2 |
415 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
416 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
417 |
|
418 |
diagName = 'VVELSLT ' |
419 |
diagTitle = DIAGS_MK_TITLE( 'Meridional Transport of ' |
420 |
I //sTitle, myThid ) |
421 |
diagUnits = DIAGS_MK_UNITS( sUnit4c//'.m/s', myThid ) |
422 |
diagCode = 'VVR MR ' |
423 |
diagMate = diagNum |
424 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
425 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
426 |
|
427 |
diagName = 'WVELSLT ' |
428 |
diagTitle = DIAGS_MK_TITLE( 'Vertical Transport of ' |
429 |
I //sTitle, myThid ) |
430 |
diagUnits = DIAGS_MK_UNITS(sUnit4c//'.'//rUnit2c//'/s', myThid ) |
431 |
diagCode = 'WM LR ' |
432 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
433 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
434 |
|
435 |
diagName = 'UVELPHI ' |
436 |
diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of ' |
437 |
I //pTitle//' Anomaly', myThid ) |
438 |
diagUnits = 'm^3/s^3 ' |
439 |
diagCode = 'UUr MR ' |
440 |
diagMate = diagNum + 2 |
441 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
442 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
443 |
|
444 |
diagName = 'VVELPHI ' |
445 |
diagTitle = DIAGS_MK_TITLE( 'Meridional Mass-Weight Transp of ' |
446 |
I //pTitle//' Anomaly', myThid ) |
447 |
diagUnits = 'm^3/s^3 ' |
448 |
diagCode = 'VVr MR ' |
449 |
diagMate = diagNum |
450 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
451 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
452 |
|
453 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
454 |
|
455 |
diagName = 'RHOAnoma' |
456 |
diagTitle = 'Density Anomaly (=Rho-rhoConst)' |
457 |
diagUnits = 'kg/m^3 ' |
458 |
diagCode = 'SMR MR ' |
459 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
460 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
461 |
|
462 |
diagName = 'RHOANOSQ' |
463 |
diagTitle = 'Square of Density Anomaly (=(Rho-rhoConst)^2)' |
464 |
diagUnits = 'kg^2/m^6 ' |
465 |
diagCode = 'SMRP MR ' |
466 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
467 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
468 |
|
469 |
diagName = 'URHOMASS' |
470 |
diagTitle = 'Zonal Transport of Density' |
471 |
diagUnits = 'kg/m^2/s ' |
472 |
diagCode = 'UUr MR ' |
473 |
diagMate = diagNum + 2 |
474 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
475 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
476 |
|
477 |
diagName = 'VRHOMASS' |
478 |
diagTitle = 'Meridional Transport of Density' |
479 |
diagUnits = 'kg/m^2/s ' |
480 |
diagCode = 'VVr MR ' |
481 |
diagMate = diagNum |
482 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
483 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
484 |
|
485 |
diagName = 'WRHOMASS' |
486 |
diagTitle = 'Vertical Transport of Potential Density' |
487 |
diagUnits = 'kg/m^2/s ' |
488 |
diagCode = 'WM LR ' |
489 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
490 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
491 |
|
492 |
diagName = 'PHIHYD ' |
493 |
diagTitle = DIAGS_MK_TITLE( 'Hydrostatic ' |
494 |
I //pTitle//' Anomaly', myThid ) |
495 |
diagUnits = 'm^2/s^2 ' |
496 |
diagCode = 'SMR MR ' |
497 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
498 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
499 |
|
500 |
diagName = 'PHIHYDSQ' |
501 |
diagTitle = DIAGS_MK_TITLE( 'Square of Hyd. ' |
502 |
I //pTitle//' Anomaly', myThid ) |
503 |
diagUnits = 'm^4/s^4 ' |
504 |
diagCode = 'SMRP MR ' |
505 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
506 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
507 |
|
508 |
diagName = 'PHIBOT ' |
509 |
c diagTitle = 'ocean bottom pressure / top. atmos geo-Potential' |
510 |
diagTitle = DIAGS_MK_TITLE( fTitle |
511 |
I //pTitle//' Anomaly', myThid ) |
512 |
diagUnits = 'm^2/s^2 ' |
513 |
diagCode = 'SM M1 ' |
514 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
515 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
516 |
|
517 |
diagName = 'PHIBOTSQ' |
518 |
c diagTitle = 'Square of ocean bottom pressure / top. geo-Potential' |
519 |
diagTitle = DIAGS_MK_TITLE( 'Square of '//fTitle |
520 |
I //pTitle//' Anomaly', myThid ) |
521 |
diagUnits = 'm^4/s^4 ' |
522 |
diagCode = 'SM P M1 ' |
523 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
524 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
525 |
|
526 |
#ifdef ALLOW_NONHYDROSTATIC |
527 |
diagName = 'PHI_NH ' |
528 |
diagTitle = DIAGS_MK_TITLE( 'Non-Hydrostatic '//pTitle, myThid ) |
529 |
diagUnits = 'm^2/s^2 ' |
530 |
diagCode = 'SMR MR ' |
531 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
532 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
533 |
#endif /* ALLOW_NONHYDROSTATIC */ |
534 |
|
535 |
diagName = 'MXLDEPTH' |
536 |
diagTitle = 'Mixed-Layer Depth (>0)' |
537 |
diagUnits = 'm ' |
538 |
diagCode = 'SM M1 ' |
539 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
540 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
541 |
|
542 |
diagName = 'DRHODR ' |
543 |
diagTitle = 'Stratification: d.Sigma/dr (kg/m3/r_unit)' |
544 |
diagUnits = 'kg/m^4 ' |
545 |
IF ( usingPCoords ) diagUnits = 's^2/m^2 ' |
546 |
diagCode = 'SM LR ' |
547 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
548 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
549 |
|
550 |
diagName = 'CONVADJ ' |
551 |
diagTitle = 'Convective Adjustment Index [0-1] ' |
552 |
diagUnits = 'fraction ' |
553 |
diagCode = 'SMR LR ' |
554 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
555 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
556 |
|
557 |
C-- surface fluxes: |
558 |
diagName = 'oceTAUX ' |
559 |
diagTitle = 'zonal surface wind stress, >0 increases uVel' |
560 |
diagUnits = 'N/m^2 ' |
561 |
diagCode = 'UU U1 ' |
562 |
diagMate = diagNum + 2 |
563 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
564 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
565 |
|
566 |
diagName = 'oceTAUY ' |
567 |
diagTitle = 'meridional surf. wind stress, >0 increases vVel' |
568 |
diagUnits = 'N/m^2 ' |
569 |
diagCode = 'VV U1 ' |
570 |
diagMate = diagNum |
571 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
572 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
573 |
|
574 |
diagName = 'atmPload' |
575 |
diagTitle = 'Atmospheric pressure loading' |
576 |
diagUnits = 'Pa ' |
577 |
diagCode = 'SM U1 ' |
578 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
579 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
580 |
|
581 |
diagName = 'sIceLoad' |
582 |
diagTitle = 'sea-ice loading (in Mass of ice+snow / area unit)' |
583 |
diagUnits = 'kg/m^2 ' |
584 |
diagCode = 'SM U1 ' |
585 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
586 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
587 |
|
588 |
diagName = 'oceFWflx' |
589 |
diagTitle = 'net surface Fresh-Water flux into the ocean' |
590 |
& //' (+=down), >0 decreases salinity' |
591 |
diagUnits = 'kg/m^2/s ' |
592 |
diagCode = 'SM U1 ' |
593 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
594 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
595 |
|
596 |
diagName = 'oceSflux' |
597 |
diagTitle = 'net surface Salt flux into the ocean (+=down),' |
598 |
& //' >0 increases salinity' |
599 |
diagUnits = 'g/m^2/s ' |
600 |
diagCode = 'SM U1 ' |
601 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
602 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
603 |
|
604 |
diagName = 'oceQnet ' |
605 |
diagTitle = 'net surface heat flux into the ocean (+=down),' |
606 |
& //' >0 increases theta' |
607 |
diagUnits = 'W/m^2 ' |
608 |
diagCode = 'SM U1 ' |
609 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
610 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
611 |
|
612 |
diagName = 'oceQsw ' |
613 |
diagTitle = 'net Short-Wave radiation (+=down),' |
614 |
& //' >0 increases theta' |
615 |
diagUnits = 'W/m^2 ' |
616 |
diagCode = 'SM U1 ' |
617 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
618 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
619 |
|
620 |
diagName = 'oceFreez' |
621 |
diagTitle = 'heating from freezing of sea-water (allowFreezing=T)' |
622 |
diagUnits = 'W/m^2 ' |
623 |
diagCode = 'SM U1 ' |
624 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
625 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
626 |
|
627 |
diagName = 'TRELAX ' |
628 |
diagTitle = 'surface temperature relaxation, >0 increases theta' |
629 |
diagUnits = 'W/m^2 ' |
630 |
diagCode = 'SM U1 ' |
631 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
632 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
633 |
|
634 |
diagName = 'SRELAX ' |
635 |
diagTitle = 'surface salinity relaxation, >0 increases salt' |
636 |
diagUnits = 'g/m^2/s ' |
637 |
diagCode = 'SM U1 ' |
638 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
639 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
640 |
|
641 |
diagName = 'surForcT' |
642 |
diagTitle = 'model surface forcing for Temperature,' |
643 |
& //' >0 increases theta' |
644 |
diagUnits = 'W/m^2 ' |
645 |
diagCode = 'SM U1 ' |
646 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
647 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
648 |
|
649 |
diagName = 'surForcS' |
650 |
diagTitle = 'model surface forcing for Salinity,' |
651 |
& //' >0 increases salinity' |
652 |
diagUnits = 'g/m^2/s ' |
653 |
diagCode = 'SM U1 ' |
654 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
655 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
656 |
|
657 |
diagName = 'TFLUX ' |
658 |
diagTitle = 'total heat flux (match heat-content variations),' |
659 |
& //' >0 increases theta' |
660 |
diagUnits = 'W/m^2 ' |
661 |
diagCode = 'SM U1 ' |
662 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
663 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
664 |
|
665 |
diagName = 'SFLUX ' |
666 |
diagTitle = 'total salt flux (match salt-content variations),' |
667 |
& //' >0 increases salt' |
668 |
diagUnits = 'g/m^2/s ' |
669 |
diagCode = 'SM U1 ' |
670 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
671 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
672 |
|
673 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
674 |
|
675 |
diagName = 'RCENTER ' |
676 |
c diagTitle = 'Cell-Center r-Position (Pressure, Height) (Pa,m)' |
677 |
diagTitle = DIAGS_MK_TITLE( 'Cell-Center '//rTitle, myThid ) |
678 |
diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid ) |
679 |
diagCode = 'SM MR ' |
680 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
681 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
682 |
|
683 |
diagName = 'RSURF ' |
684 |
c diagTitle = 'Free-Surface r-Position (Pressure, Height) (Pa,m)' |
685 |
diagTitle = DIAGS_MK_TITLE( eTitle//rTitle, myThid ) |
686 |
diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid ) |
687 |
diagCode = 'SM M1 ' |
688 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
689 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
690 |
|
691 |
diagName = 'TOTUTEND' |
692 |
diagTitle = 'Tendency of Zonal Component of Velocity' |
693 |
diagUnits = 'm/s/day ' |
694 |
diagCode = 'UUR MR ' |
695 |
diagMate = diagNum + 2 |
696 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
697 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
698 |
|
699 |
diagName = 'TOTVTEND' |
700 |
diagTitle = 'Tendency of Meridional Component of Velocity' |
701 |
diagUnits = 'm/s/day ' |
702 |
diagCode = 'VVR MR ' |
703 |
diagMate = diagNum |
704 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
705 |
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
706 |
|
707 |
diagName = 'TOTTTEND' |
708 |
diagTitle = 'Tendency of Potential Temperature' |
709 |
diagUnits = DIAGS_MK_UNITS( tUnit4c//'/day', myThid ) |
710 |
diagCode = 'SMR MR ' |
711 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
712 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
713 |
|
714 |
diagName = 'TOTSTEND' |
715 |
diagTitle = DIAGS_MK_TITLE('Tendency of '//sTitle, myThid ) |
716 |
diagUnits = DIAGS_MK_UNITS( sUnit4c//'/day', myThid ) |
717 |
diagCode = 'SMR MR ' |
718 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
719 |
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
720 |
|
721 |
|
722 |
RETURN |
723 |
END |