/[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.52 - (show annotations) (download)
Sat May 28 23:24:47 2016 UTC (7 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, HEAD
Changes since 1.51: +8 -1 lines
add diagnostics for surface dynamical pressure (to facilitate momemtum budget)

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

  ViewVC Help
Powered by ViewVC 1.1.22