| 1 |
zhc |
1.4 |
C $Header: /u/gcmpack/MITgcm/pkg/gmredi/gmredi_diagnostics_init.F,v 1.7 2010/01/12 21:34:09 jmc Exp $ |
| 2 |
dimitri |
1.1 |
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 |
zhc |
1.4 |
C diagMate :: diag. mate number in the (long) list of available diag. |
| 38 |
dimitri |
1.1 |
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 |
zhc |
1.4 |
INTEGER diagMate |
| 45 |
dimitri |
1.1 |
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 |
zhc |
1.4 |
diagTitle = |
| 56 |
dimitri |
1.1 |
& 'Mixing coefficient from Visbeck etal parameterization' |
| 57 |
|
|
diagUnits = 'm^2/s ' |
| 58 |
|
|
diagCode = 'SM P M1 ' |
| 59 |
zhc |
1.4 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 60 |
|
|
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
| 61 |
dimitri |
1.1 |
|
| 62 |
|
|
diagName = 'GM_hTrsL' |
| 63 |
|
|
diagTitle = 'Base depth (>0) of the Transition Layer' |
| 64 |
|
|
diagUnits = 'm ' |
| 65 |
|
|
diagCode = 'SM P M1 ' |
| 66 |
zhc |
1.4 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 67 |
|
|
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
| 68 |
dimitri |
1.1 |
|
| 69 |
|
|
diagName = 'GM_baseS' |
| 70 |
|
|
diagTitle = 'Slope at the base of the Transition Layer' |
| 71 |
|
|
diagUnits = '1 ' |
| 72 |
|
|
diagCode = 'SM P M1 ' |
| 73 |
zhc |
1.4 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 74 |
|
|
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
| 75 |
dimitri |
1.1 |
|
| 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 |
zhc |
1.4 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 82 |
|
|
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
| 83 |
dimitri |
1.1 |
|
| 84 |
|
|
diagName = 'GM_Kux ' |
| 85 |
|
|
diagTitle = 'K_11 element (U.point, X.dir) of GM-Redi tensor' |
| 86 |
|
|
diagUnits = 'm^2/s ' |
| 87 |
zhc |
1.4 |
diagCode = 'UU P MR ' |
| 88 |
|
|
diagMate = diagNum + 2 |
| 89 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 90 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 91 |
dimitri |
1.1 |
|
| 92 |
|
|
diagName = 'GM_Kvy ' |
| 93 |
|
|
diagTitle = 'K_22 element (V.point, Y.dir) of GM-Redi tensor' |
| 94 |
|
|
diagUnits = 'm^2/s ' |
| 95 |
zhc |
1.4 |
diagCode = 'VV P MR ' |
| 96 |
|
|
diagMate = diagNum |
| 97 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 98 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 99 |
dimitri |
1.1 |
|
| 100 |
|
|
diagName = 'GM_Kuz ' |
| 101 |
|
|
diagTitle = 'K_13 element (U.point, Z.dir) of GM-Redi tensor' |
| 102 |
|
|
diagUnits = 'm^2/s ' |
| 103 |
zhc |
1.4 |
diagCode = 'UU MR ' |
| 104 |
|
|
diagMate = diagNum + 2 |
| 105 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 106 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 107 |
dimitri |
1.1 |
|
| 108 |
|
|
diagName = 'GM_Kvz ' |
| 109 |
|
|
diagTitle = 'K_23 element (V.point, Z.dir) of GM-Redi tensor' |
| 110 |
|
|
diagUnits = 'm^2/s ' |
| 111 |
zhc |
1.4 |
diagCode = 'VV MR ' |
| 112 |
|
|
diagMate = diagNum |
| 113 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 114 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 115 |
dimitri |
1.1 |
|
| 116 |
|
|
diagName = 'GM_Kwx ' |
| 117 |
|
|
diagTitle = 'K_31 element (W.point, X.dir) of GM-Redi tensor' |
| 118 |
|
|
diagUnits = 'm^2/s ' |
| 119 |
zhc |
1.4 |
diagCode = 'UM LR ' |
| 120 |
|
|
diagMate = diagNum + 2 |
| 121 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 122 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 123 |
dimitri |
1.1 |
|
| 124 |
|
|
diagName = 'GM_Kwy ' |
| 125 |
|
|
diagTitle = 'K_32 element (W.point, Y.dir) of GM-Redi tensor' |
| 126 |
|
|
diagUnits = 'm^2/s ' |
| 127 |
zhc |
1.4 |
diagCode = 'VM LR ' |
| 128 |
|
|
diagMate = diagNum |
| 129 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 130 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 131 |
dimitri |
1.1 |
|
| 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 |
zhc |
1.4 |
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 137 |
|
|
I diagName, diagCode, diagUnits, diagTitle, 0, myThid ) |
| 138 |
dimitri |
1.1 |
|
| 139 |
|
|
diagName = 'GM_PsiX ' |
| 140 |
|
|
diagTitle = 'GM Bolus transport stream-function : X component' |
| 141 |
|
|
diagUnits = 'm^2/s ' |
| 142 |
zhc |
1.4 |
diagCode = 'UU LR ' |
| 143 |
|
|
diagMate = diagNum + 2 |
| 144 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 145 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 146 |
dimitri |
1.1 |
|
| 147 |
|
|
diagName = 'GM_PsiY ' |
| 148 |
|
|
diagTitle = 'GM Bolus transport stream-function : Y component' |
| 149 |
|
|
diagUnits = 'm^2/s ' |
| 150 |
zhc |
1.4 |
diagCode = 'VV LR ' |
| 151 |
|
|
diagMate = diagNum |
| 152 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 153 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 154 |
dimitri |
1.1 |
|
| 155 |
|
|
diagName = 'GM_KuzTz' |
| 156 |
|
|
diagTitle = 'Redi Off-diagonal Temperature flux: X component' |
| 157 |
|
|
diagUnits = 'degC.m^3/s ' |
| 158 |
zhc |
1.4 |
diagCode = 'UU MR ' |
| 159 |
|
|
diagMate = diagNum + 2 |
| 160 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 161 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 162 |
dimitri |
1.1 |
|
| 163 |
|
|
diagName = 'GM_KvzTz' |
| 164 |
|
|
diagTitle = 'Redi Off-diagonal Temperature flux: Y component' |
| 165 |
|
|
diagUnits = 'degC.m^3/s ' |
| 166 |
zhc |
1.4 |
diagCode = 'VV MR ' |
| 167 |
|
|
diagMate = diagNum |
| 168 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 169 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 170 |
dimitri |
1.1 |
|
| 171 |
|
|
diagName = 'GM_ubT ' |
| 172 |
|
|
diagTitle = 'Zonal Mass-Weight Bolus Transp of Pot Temp' |
| 173 |
|
|
diagUnits = 'degC.m^3/s ' |
| 174 |
zhc |
1.4 |
diagCode = 'UUr MR ' |
| 175 |
|
|
diagMate = diagNum + 2 |
| 176 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 177 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 178 |
dimitri |
1.1 |
|
| 179 |
|
|
diagName = 'GM_vbT ' |
| 180 |
|
|
diagTitle = 'Meridional Mass-Weight Bolus Transp of Pot Temp' |
| 181 |
|
|
diagUnits = 'degC.m^3/s ' |
| 182 |
zhc |
1.4 |
diagCode = 'VVr MR ' |
| 183 |
|
|
diagMate = diagNum |
| 184 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 185 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 186 |
dimitri |
1.1 |
|
| 187 |
dimitri |
1.2 |
#ifdef GM_SUBMESO |
| 188 |
|
|
CBFK These are diagnosed whether or not Bolus Advection is used... |
| 189 |
|
|
CBFK They are the quasi-streamfunction, since they are in the |
| 190 |
|
|
CBFK Visbeck -dPsix/dz=u form, not the Fox-Kemper curl(psi)=u form. |
| 191 |
|
|
CBFK They are always included in the GM part above |
| 192 |
|
|
CBFK (i.e., GM diagnostic=GM+Submeso) |
| 193 |
|
|
diagName = 'SM_PsiX ' |
| 194 |
|
|
diagTitle = 'Submeso Bolus transport quasi-streamfunction : X' |
| 195 |
|
|
diagUnits = 'm^2/s ' |
| 196 |
zhc |
1.4 |
diagCode = 'UU LR ' |
| 197 |
|
|
diagMate = diagNum + 2 |
| 198 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 199 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 200 |
dimitri |
1.2 |
|
| 201 |
|
|
diagName = 'SM_PsiY ' |
| 202 |
|
|
diagTitle = 'Submeso Bolus transport quasi-streamfunction : Y' |
| 203 |
|
|
diagUnits = 'm^2/s ' |
| 204 |
zhc |
1.4 |
diagCode = 'VV LR ' |
| 205 |
|
|
diagMate = diagNum |
| 206 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 207 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 208 |
dimitri |
1.2 |
|
| 209 |
|
|
diagName = 'SM_ubT ' |
| 210 |
|
|
diagTitle = 'Zon. Submeso Mass-Weight Bolus Transp of Pot Temp' |
| 211 |
|
|
diagUnits = 'degC.m^3/s ' |
| 212 |
zhc |
1.4 |
diagCode = 'UUr MR ' |
| 213 |
|
|
diagMate = diagNum + 2 |
| 214 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 215 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 216 |
dimitri |
1.2 |
|
| 217 |
|
|
diagName = 'SM_vbT ' |
| 218 |
|
|
diagTitle = 'Mer. Submeso Mass-Weight Bolus Transp of Pot Temp' |
| 219 |
|
|
diagUnits = 'degC.m^3/s ' |
| 220 |
zhc |
1.4 |
diagCode = 'VVr MR ' |
| 221 |
|
|
diagMate = diagNum |
| 222 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 223 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 224 |
dimitri |
1.2 |
|
| 225 |
|
|
diagName = 'SM_wbT ' |
| 226 |
|
|
diagTitle = 'Rvel Submeso Mass-Weight Bolus Transp of Pot Temp' |
| 227 |
|
|
diagUnits = 'degC.m^3/s ' |
| 228 |
zhc |
1.4 |
diagCode = 'WM LR ' |
| 229 |
|
|
diagMate = diagNum |
| 230 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 231 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 232 |
dimitri |
1.2 |
|
| 233 |
|
|
diagName = 'SM_KuzTz' |
| 234 |
|
|
diagTitle = 'Zon. Submeso Mass-Weight Kappa Transp of Pot Temp' |
| 235 |
|
|
diagUnits = 'degC.m^3/s ' |
| 236 |
zhc |
1.4 |
diagCode = 'UU MR ' |
| 237 |
|
|
diagMate = diagNum + 2 |
| 238 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 239 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 240 |
dimitri |
1.2 |
|
| 241 |
|
|
diagName = 'SM_KvzTz' |
| 242 |
|
|
diagTitle = 'Mer. Submeso Mass-Weight Kappa Transp of Pot Temp' |
| 243 |
|
|
diagUnits = 'degC.m^3/s ' |
| 244 |
zhc |
1.4 |
diagCode = 'VV MR ' |
| 245 |
|
|
diagMate = diagNum |
| 246 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 247 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 248 |
dimitri |
1.2 |
|
| 249 |
|
|
diagName = 'SM_KrddT' |
| 250 |
|
|
diagTitle = 'Rvel Submeso Mass-Weight Kappa Transp of Pot Temp' |
| 251 |
|
|
diagUnits = 'degC.m^3/s ' |
| 252 |
zhc |
1.4 |
diagCode = 'WM LR ' |
| 253 |
|
|
diagMate = diagNum |
| 254 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 255 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 256 |
dimitri |
1.2 |
#endif |
| 257 |
|
|
|
| 258 |
dimitri |
1.3 |
#ifdef ALLOW_EDDYPSI |
| 259 |
dimitri |
1.1 |
diagName = 'GMEdTauX' |
| 260 |
|
|
diagTitle = 'eddy-induced stress X-comp. estimated from Kwx' |
| 261 |
|
|
diagUnits = 'N/m^2 ' |
| 262 |
zhc |
1.4 |
diagCode = 'UM LR ' |
| 263 |
|
|
diagMate = diagNum + 2 |
| 264 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 265 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 266 |
dimitri |
1.1 |
|
| 267 |
|
|
diagName = 'GMEdTauY' |
| 268 |
|
|
diagTitle = 'eddy-induced stress Y-comp. estimated from Kwy' |
| 269 |
|
|
diagUnits = 'N/m^2 ' |
| 270 |
zhc |
1.4 |
diagCode = 'VM LR ' |
| 271 |
|
|
diagMate = diagNum |
| 272 |
|
|
CALL DIAGNOSTICS_ADDTOLIST( diagNum, |
| 273 |
|
|
I diagName, diagCode, diagUnits, diagTitle, diagMate, myThid ) |
| 274 |
dimitri |
1.1 |
#endif |
| 275 |
|
|
|
| 276 |
|
|
c ENDIF |
| 277 |
|
|
|
| 278 |
|
|
#endif /* ALLOW_DIAGNOSTICS */ |
| 279 |
|
|
|
| 280 |
|
|
RETURN |
| 281 |
|
|
END |