/[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.46 - (show annotations) (download)
Wed Sep 18 01:41:05 2013 UTC (10 years, 7 months ago) by m_bates
Branch: MAIN
Changes since 1.45: +19 -1 lines
Initial diagnostics for the eddy stress (available for ALLOW_EDDYPSI)

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

  ViewVC Help
Powered by ViewVC 1.1.22