/[MITgcm]/MITgcm/pkg/gmredi/gmredi_diagnostics_init.F
ViewVC logotype

Contents of /MITgcm/pkg/gmredi/gmredi_diagnostics_init.F

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


Revision 1.7 - (show annotations) (download)
Tue Jan 12 21:34:09 2010 UTC (14 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62r, checkpoint62q, checkpoint62p
Changes since 1.6: +70 -54 lines
replace remaining calls to old S/R DIAGNOSTICS_ADD2LIST with newer S/R
 DIAGNOSTICS_ADDTOLIST ;
 fix add few parser code 3 ; fix few vertical position flag (U -> L).

1 C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_diagnostics_init.F,v 1.6 2008/05/30 02:50:16 gforget Exp $
2 C $Name: $
3
4 #include "GMREDI_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GMREDI_DIAGNOSTICS_INIT
8 C !INTERFACE:
9 SUBROUTINE GMREDI_DIAGNOSTICS_INIT( myThid )
10
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | SUBROUTINE GMREDI_DIAGNOSTICS_INIT
14 C | o Routine to initialize list of all available diagnostics
15 C | for GM/Redi package
16 C *==========================================================*
17 C \ev
18 C !USES:
19 IMPLICIT NONE
20
21 C === Global variables ===
22 #include "EEPARAMS.h"
23 c #include "SIZE.h"
24 c #include "PARAMS.h"
25 c #include "GMREDI.h"
26
27 C !INPUT/OUTPUT PARAMETERS:
28 C === Routine arguments ===
29 C myThid :: my Thread Id number
30 INTEGER myThid
31 CEOP
32
33 #ifdef ALLOW_DIAGNOSTICS
34 C !LOCAL VARIABLES:
35 C === Local variables ===
36 C diagNum :: diagnostics number in the (long) list of available diag.
37 C diagMate :: diag. mate number in the (long) list of available diag.
38 C diagName :: local short name (8c) of a diagnostics
39 C diagCode :: local parser field with characteristics of the diagnostics
40 C cf head of S/R DIAGNOSTICS_INIT_EARLY or DIAGNOSTICS_MAIN_INIT
41 C diagUnits :: local string (16c): physical units of a diagnostic field
42 C diagTitle :: local string (80c): description of field in diagnostic
43 INTEGER diagNum
44 INTEGER diagMate
45 CHARACTER*8 diagName
46 CHARACTER*16 diagCode
47 CHARACTER*16 diagUnits
48 CHARACTER*(80) diagTitle
49
50 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
51
52 c IF ( useDiagnotics ) THEN
53
54 diagName = 'GM_VisbK'
55 diagTitle =
56 & 'Mixing coefficient from Visbeck etal parameterization'
57 diagUnits = 'm^2/s '
58 diagCode = 'SM P M1 '
59 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
60 I diagName, diagCode, diagUnits, diagTitle, 0, myThid )
61
62 diagName = 'GM_hTrsL'
63 diagTitle = 'Base depth (>0) of the Transition Layer'
64 diagUnits = 'm '
65 diagCode = 'SM P M1 '
66 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
67 I diagName, diagCode, diagUnits, diagTitle, 0, myThid )
68
69 diagName = 'GM_baseS'
70 diagTitle = 'Slope at the base of the Transition Layer'
71 diagUnits = '1 '
72 diagCode = 'SM P M1 '
73 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
74 I diagName, diagCode, diagUnits, diagTitle, 0, myThid )
75
76 diagName = 'GM_rLamb'
77 diagTitle =
78 & 'Slope vertical gradient at Trans. Layer Base (=recip.Lambda)'
79 diagUnits = '1/m '
80 diagCode = 'SM P M1 '
81 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
82 I diagName, diagCode, diagUnits, diagTitle, 0, myThid )
83
84 diagName = 'GM_Kux '
85 diagTitle = 'K_11 element (U.point, X.dir) of GM-Redi tensor'
86 diagUnits = 'm^2/s '
87 diagCode = 'UU P MR '
88 diagMate = diagNum + 2
89 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
90 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
91
92 diagName = 'GM_Kvy '
93 diagTitle = 'K_22 element (V.point, Y.dir) of GM-Redi tensor'
94 diagUnits = 'm^2/s '
95 diagCode = 'VV P MR '
96 diagMate = diagNum
97 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
98 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
99
100 diagName = 'GM_Kuz '
101 diagTitle = 'K_13 element (U.point, Z.dir) of GM-Redi tensor'
102 diagUnits = 'm^2/s '
103 diagCode = 'UU MR '
104 diagMate = diagNum + 2
105 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
106 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
107
108 diagName = 'GM_Kvz '
109 diagTitle = 'K_23 element (V.point, Z.dir) of GM-Redi tensor'
110 diagUnits = 'm^2/s '
111 diagCode = 'VV MR '
112 diagMate = diagNum
113 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
114 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
115
116 diagName = 'GM_Kwx '
117 diagTitle = 'K_31 element (W.point, X.dir) of GM-Redi tensor'
118 diagUnits = 'm^2/s '
119 diagCode = 'UM LR '
120 diagMate = diagNum + 2
121 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
122 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
123
124 diagName = 'GM_Kwy '
125 diagTitle = 'K_32 element (W.point, Y.dir) of GM-Redi tensor'
126 diagUnits = 'm^2/s '
127 diagCode = 'VM LR '
128 diagMate = diagNum
129 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
130 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
131
132 diagName = 'GM_Kwz '
133 diagTitle = 'K_33 element (W.point, Z.dir) of GM-Redi tensor'
134 diagUnits = 'm^2/s '
135 diagCode = 'WM P LR '
136 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
137 I diagName, diagCode, diagUnits, diagTitle, 0, myThid )
138
139 diagName = 'GM_PsiX '
140 diagTitle = 'GM Bolus transport stream-function : X component'
141 diagUnits = 'm^2/s '
142 diagCode = 'UU LR '
143 diagMate = diagNum + 2
144 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
145 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
146
147 diagName = 'GM_PsiY '
148 diagTitle = 'GM Bolus transport stream-function : Y component'
149 diagUnits = 'm^2/s '
150 diagCode = 'VV LR '
151 diagMate = diagNum
152 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
153 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
154
155 diagName = 'GM_KuzTz'
156 diagTitle = 'Redi Off-diagonal Temperature flux: X component'
157 diagUnits = 'degC.m^3/s '
158 diagCode = 'UU MR '
159 diagMate = diagNum + 2
160 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
161 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
162
163 diagName = 'GM_KvzTz'
164 diagTitle = 'Redi Off-diagonal Temperature flux: Y component'
165 diagUnits = 'degC.m^3/s '
166 diagCode = 'VV MR '
167 diagMate = diagNum
168 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
169 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
170
171 diagName = 'GM_ubT '
172 diagTitle = 'Zonal Mass-Weight Bolus Transp of Pot Temp'
173 diagUnits = 'degC.m^3/s '
174 diagCode = 'UUr MR '
175 diagMate = diagNum + 2
176 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
177 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
178
179 diagName = 'GM_vbT '
180 diagTitle = 'Meridional Mass-Weight Bolus Transp of Pot Temp'
181 diagUnits = 'degC.m^3/s '
182 diagCode = 'VVr MR '
183 diagMate = diagNum
184 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
185 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
186
187 #ifdef ALLOW_EDDYPSI
188 diagName = 'GMEdTauX'
189 diagTitle = 'eddy-induced stress X-comp. estimated from Kwx'
190 diagUnits = 'N/m^2 '
191 diagCode = 'UM LR '
192 diagMate = diagNum + 2
193 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
194 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
195
196 diagName = 'GMEdTauY'
197 diagTitle = 'eddy-induced stress Y-comp. estimated from Kwy'
198 diagUnits = 'N/m^2 '
199 diagCode = 'VM LR '
200 diagMate = diagNum
201 CALL DIAGNOSTICS_ADDTOLIST( diagNum,
202 I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
203 #endif
204
205 c ENDIF
206
207 #endif /* ALLOW_DIAGNOSTICS */
208
209 RETURN
210 END

  ViewVC Help
Powered by ViewVC 1.1.22