/[MITgcm]/MITgcm/pkg/diagnostics/diagnostics_main_init.F
ViewVC logotype

Contents of /MITgcm/pkg/diagnostics/diagnostics_main_init.F

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


Revision 1.43 - (show annotations) (download)
Wed Nov 30 20:58:35 2011 UTC (12 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64a, checkpoint64c, checkpoint64b, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63g, checkpoint64
Changes since 1.42: +16 -1 lines
add diags for T & S tendency which goes through Adams-Bashforth

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

  ViewVC Help
Powered by ViewVC 1.1.22