Back to home page

MITgcm

 
 

    


File indexing completed on 2024-07-17 05:10:40 UTC

view on githubraw file Latest commit acacc28f on 2024-07-17 03:59:01 UTC
09ceb40cd6 Jean*0001 #include "DIAG_OPTIONS.h"
                0002 
                0003 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
acacc28f7f Jean*0004 CBOP
09ceb40cd6 Jean*0005 C     !ROUTINE: DIAGNOSTICS_MAIN_INIT
                0006 
                0007 C     !INTERFACE:
                0008       SUBROUTINE DIAGNOSTICS_MAIN_INIT( myThid )
                0009 
                0010 C     !DESCRIPTION:
                0011 C     Initialize available diagnostics list for variables of the main code
                0012 C     (not part of a package): set the following attributes:
                0013 C     name (=cdiag), parsing code (=gdiag), units (=udiag), and title (=tdiag)
                0014 C     Notes: 1) diagnostics defined here do not require any EQUIVALENCE
                0015 C            since they get filled-in with S/R FILL_DIAGNOSTICS
                0016 C            2) GDIAG is defined as character*16 and can be to character*1
                0017 C            parse(16) with the following codes currently defined:
                0018 
acacc28f7f Jean*0019 C \begin{center}
                0020 C   \begin{tabular}[h]{|c|c|}\hline
                0021 C     \textbf{Positions} & \textbf{Characters} & \textbf{Meanings} \\\hline
                0022 C     parse(1)  &  S  &  scalar \\
                0023 C               &  U  &  vector component in X direction \\
                0024 C               &  V  &  vector component in Y direction \\
                0025 C               &  W  &  vector component in vertical direction \\
                0026 C     parse(2)  &  U  &  C-grid U-Point  \\
                0027 C               &  V  &  C-grid V-Point  \\
                0028 C               &  M  &  C-grid Mass Point  \\
                0029 C               &  Z  &  C-grid Corner Point  \\
                0030 C     parse(3)  &     &  Used for Level Integrated output: cumulate levels \\
                0031 C               &  r  &  same but cumulate product by model level thickness \\
                0032 C               &  R  &  same but cumulate product by hFac & level thickness \\
                0033 C     parse(4)  &  P  &  positive definite  \\
                0034 C               &  A  &  Adjoint variable diagnostics \\
                0035 C     parse(5 ) &  C  &  with counter array  \\
                0036 C               &  P  &  post-processed (not filled up) from other diags  \\
                0037 C               &  D  &  disable an array for output  \\
                0038 C     parse(6--8) & '123'  &  retired, formerly: 3-digit mate number \\
                0039 C     parse(9)  &  U  &  model-level plus 1/2  \\
                0040 C               &  M  &  model-level middle  \\
                0041 C               &  L  &  model-level minus 1/2  \\
                0042 C     parse(10) &  0  &  levels = 0  \\
                0043 C               &  1  &  levels = 1  \\
                0044 C               &  R  &  levels = Nr  \\
                0045 C               &  L  &  levels = MAX(Nr,NrPhys)  \\
                0046 C               &  M  &  levels = MAX(Nr,NrPhys) - 1  \\
                0047 C               &  G  &  levels = Ground_level Number \\
                0048 C               &  I  &  levels = sea-Ice_level Number \\
                0049 C               &  X  &  free levels option (need to be set explicitly) \\
                0050 C   \end{tabular}
                0051 C \end{center}
09ceb40cd6 Jean*0052 
                0053 C     !USES:
                0054       IMPLICIT NONE
                0055 #include "SIZE.h"
                0056 #include "EEPARAMS.h"
                0057 #include "PARAMS.h"
                0058 
                0059 C     !INPUT PARAMETERS:
                0060       INTEGER myThid
                0061 CEOP
                0062 
                0063 C     !LOCAL VARIABLES:
9ba8ef70b1 Jean*0064 C     rTitle     :: r-coordinate title
                0065 C     eTitle     :: free-surface title
                0066 C     fTitle     :: fixed boundary title
                0067 C     pTitle     :: "Phi"  title
                0068 C     sTitle     :: "salt" title
                0069       INTEGER        diagNum
27569449c9 Jean*0070       INTEGER        diagMate
9ba8ef70b1 Jean*0071       CHARACTER*8    diagName
                0072       CHARACTER*16   diagCode
                0073       CHARACTER*16   diagUnits
09ceb40cd6 Jean*0074       CHARACTER*(80) diagTitle
9ba8ef70b1 Jean*0075       CHARACTER*2    rUnit2c
3ac63c2262 Jean*0076       CHARACTER*4    tUnit4c
                0077       CHARACTER*5    sUnit5c
9ba8ef70b1 Jean*0078       CHARACTER*(10) rTitle, eTitle, fTitle
                0079       CHARACTER*(20) pTitle, sTitle
db26b4dd29 Jean*0080 
                0081       CHARACTER*(16) DIAGS_MK_UNITS
                0082       EXTERNAL DIAGS_MK_UNITS
9ba8ef70b1 Jean*0083       CHARACTER*(80) DIAGS_MK_TITLE
                0084       EXTERNAL DIAGS_MK_TITLE
09ceb40cd6 Jean*0085 
                0086 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0087 C     For each output variable,
9ba8ef70b1 Jean*0088 C     specify Name (cdiag, 8c), Descriptions (tdiag, *c), Units (udiag, 16c)
09ceb40cd6 Jean*0089 C         and Type/Parms (location on C grid, 2D/3D, ...) (gdiag, 16c)
                0090 C----------------------------------------------------------------------
                0091 
db26b4dd29 Jean*0092       IF ( usingPCoords ) THEN
                0093         rUnit2c= 'Pa'
9ba8ef70b1 Jean*0094         rTitle = ' Pressure '
                0095         pTitle = ' Geopotential       '
db26b4dd29 Jean*0096       ELSE
                0097         rUnit2c= 'm '
9ba8ef70b1 Jean*0098         rTitle = ' Height   '
                0099         pTitle = 'Pressure Pot.(p/rho)'
db26b4dd29 Jean*0100       ENDIF
                0101       IF ( fluidIsAir ) THEN
                0102         tUnit4c= 'K   '
3ac63c2262 Jean*0103         sUnit5c= 'kg/kg'
9ba8ef70b1 Jean*0104         sTitle = ' Specific Humidity  '
3ac63c2262 Jean*0105         IF (useAIM) sUnit5c= 'g/kg '
                0106       ELSEIF ( eosType.EQ.'TEOS10' ) THEN
                0107         tUnit4c= 'degC'
                0108         sUnit5c= 'g/kg '
                0109 c       tTitle = 'Conservative Temp.   '
                0110         sTitle = ' Absolute Salinity  '
db26b4dd29 Jean*0111       ELSE
                0112         tUnit4c= 'degC'
ba0b047096 Mart*0113         sUnit5c= 'g/kg '
3ac63c2262 Jean*0114 c       tTitle = 'Potential Temperature'
9ba8ef70b1 Jean*0115         sTitle = ' Salinity           '
                0116       ENDIF
                0117 C-    free-surface (eTitle) and fixed-boundary (fTitle) position:
                0118       IF ( fluidIsAir ) THEN
                0119         eTitle = ' Surface  '
                0120         fTitle = ' Top      '
                0121       ELSEIF ( usingPCoords ) THEN
                0122         eTitle = ' Bottom   '
                0123         fTitle = ' Surface  '
                0124       ELSE
                0125         eTitle = ' Surface  '
                0126         fTitle = ' Bottom   '
db26b4dd29 Jean*0127       ENDIF
                0128 
8c664c8b58 Jean*0129 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0130 C-    state variables of the main code (and related quadratic var):
                0131 
                0132       diagName  = 'ETAN    '
9ba8ef70b1 Jean*0133       diagTitle = DIAGS_MK_TITLE( eTitle//rTitle//' Anomaly', myThid )
                0134 c     IF ( fluidIsWater .AND. usingZCoords )
                0135 c    &diagTitle = 'Sea Surface Elevation'
db26b4dd29 Jean*0136       diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid )
8c664c8b58 Jean*0137       diagCode  = 'SM      M1      '
27569449c9 Jean*0138       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0139      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0140 
                0141       diagName  = 'ETANSQ  '
9ba8ef70b1 Jean*0142       diagTitle = DIAGS_MK_TITLE( 'Square of '//eTitle//rTitle
                0143      I                          //' Anomaly', myThid )
db26b4dd29 Jean*0144       diagUnits = DIAGS_MK_UNITS( rUnit2c//'^2', myThid )
27569449c9 Jean*0145       diagCode  = 'SM P    M1      '
                0146       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0147      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0148 
                0149       diagName  = 'DETADT2 '
9ba8ef70b1 Jean*0150       diagTitle = DIAGS_MK_TITLE( 'Square of '//eTitle//rTitle
                0151      I                          //' Anomaly Tendency', myThid )
db26b4dd29 Jean*0152       diagUnits = DIAGS_MK_UNITS( rUnit2c//'^2/s^2', myThid )
8c664c8b58 Jean*0153       diagCode  = 'SM      M1      '
27569449c9 Jean*0154       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0155      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0156 
                0157       diagName  = 'THETA   '
9ba8ef70b1 Jean*0158       diagTitle = 'Potential Temperature'
db26b4dd29 Jean*0159       diagUnits = DIAGS_MK_UNITS( tUnit4c, myThid )
5644fea420 Jean*0160       diagCode  = 'SMR     MR      '
27569449c9 Jean*0161       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0162      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0163 
b21f2d376f Jean*0164 c     diagName  = 'SST     '
                0165 c     diagTitle = 'Sea Surface Temperature (degC,K)'
                0166 c     diagUnits = DIAGS_MK_UNITS( tUnit4c, myThid )
                0167 c     diagCode  = 'SM      M1      '
27569449c9 Jean*0168 c     CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0169 c    I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
f4f8850bd9 Dimi*0170 
8c664c8b58 Jean*0171       diagName  = 'SALT    '
9ba8ef70b1 Jean*0172       diagTitle = DIAGS_MK_TITLE( sTitle,  myThid )
3ac63c2262 Jean*0173       diagUnits = DIAGS_MK_UNITS( sUnit5c, myThid )
5644fea420 Jean*0174       diagCode  = 'SMR     MR      '
27569449c9 Jean*0175       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0176      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
d316948822 Andr*0177 
                0178       diagName  = 'RELHUM  '
9ba8ef70b1 Jean*0179       diagTitle = 'Relative Humidity'
d316948822 Andr*0180       diagUnits = 'percent         '
5644fea420 Jean*0181       diagCode  = 'SMR     MR      '
27569449c9 Jean*0182       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0183      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0184 
b21f2d376f Jean*0185 c     diagName  = 'SSS     '
                0186 c     diagTitle = 'Sea Surface Salinity '
3ac63c2262 Jean*0187 c     diagUnits = DIAGS_MK_UNITS( sUnit5c, myThid )
b21f2d376f Jean*0188 c     diagCode  = 'SM      M1      '
27569449c9 Jean*0189 c     CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0190 c    I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
f4f8850bd9 Dimi*0191 
9ba8ef70b1 Jean*0192       IF ( fluidIsWater ) THEN
23753a76a9 Dimi*0193       diagName  = 'SALTanom'
b21f2d376f Jean*0194       diagTitle = 'Salt anomaly (=SALT-35; g/kg)'
3ac63c2262 Jean*0195       diagUnits = DIAGS_MK_UNITS( sUnit5c, myThid )
5644fea420 Jean*0196       diagCode  = 'SMR     MR      '
27569449c9 Jean*0197       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0198      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
9ba8ef70b1 Jean*0199       ENDIF
23753a76a9 Dimi*0200 
8c664c8b58 Jean*0201       diagName  = 'UVEL    '
                0202       diagTitle = 'Zonal Component of Velocity (m/s)'
db26b4dd29 Jean*0203       diagUnits = 'm/s             '
5644fea420 Jean*0204       diagCode  = 'UUR     MR      '
27569449c9 Jean*0205       diagMate  = diagNum + 2
                0206       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0207      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0208 
                0209       diagName  = 'VVEL    '
                0210       diagTitle = 'Meridional Component of Velocity (m/s)'
db26b4dd29 Jean*0211       diagUnits = 'm/s             '
5644fea420 Jean*0212       diagCode  = 'VVR     MR      '
27569449c9 Jean*0213       diagMate  = diagNum
                0214       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0215      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0216 
                0217       diagName  = 'WVEL    '
db26b4dd29 Jean*0218       diagTitle = 'Vertical Component of Velocity (r_units/s)'
                0219       diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s', myThid )
8c664c8b58 Jean*0220       diagCode  = 'WM      LR      '
27569449c9 Jean*0221       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0222      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0223 
                0224       diagName  = 'THETASQ '
9ba8ef70b1 Jean*0225       diagTitle = 'Square of Potential Temperature'
db26b4dd29 Jean*0226       diagUnits = DIAGS_MK_UNITS( tUnit4c//'^2', myThid )
5644fea420 Jean*0227       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0228       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0229      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0230 
                0231       diagName  = 'SALTSQ  '
9ba8ef70b1 Jean*0232       diagTitle = DIAGS_MK_TITLE( 'Square of '//sTitle, myThid )
3ac63c2262 Jean*0233       diagUnits = DIAGS_MK_UNITS( '('//sUnit5c//')^2', myThid )
5644fea420 Jean*0234       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0235       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0236      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
23753a76a9 Dimi*0237 
9ba8ef70b1 Jean*0238       IF ( fluidIsWater ) THEN
23753a76a9 Dimi*0239       diagName  = 'SALTSQan'
6d1c87ad40 Dimi*0240       diagTitle = 'Square of Salt anomaly (=(SALT-35)^2 (g^2/kg^2)'
3ac63c2262 Jean*0241       diagUnits = DIAGS_MK_UNITS( '('//sUnit5c//')^2', myThid )
5644fea420 Jean*0242       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0243       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0244      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
9ba8ef70b1 Jean*0245       ENDIF
8c664c8b58 Jean*0246 
                0247       diagName  = 'UVELSQ  '
                0248       diagTitle = 'Square of Zonal Comp of Velocity (m^2/s^2)'
db26b4dd29 Jean*0249       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0250       diagCode  = 'UURP    MR      '
27569449c9 Jean*0251       diagMate  = diagNum + 2
                0252       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0253      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0254 
                0255       diagName  = 'VVELSQ  '
                0256       diagTitle = 'Square of Meridional Comp of Velocity (m^2/s^2)'
db26b4dd29 Jean*0257       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0258       diagCode  = 'VVRP    MR      '
27569449c9 Jean*0259       diagMate  = diagNum
                0260       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0261      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0262 
                0263       diagName  = 'WVELSQ  '
9ba8ef70b1 Jean*0264       diagTitle = 'Square of Vertical Comp of Velocity'
db26b4dd29 Jean*0265       diagUnits = DIAGS_MK_UNITS( rUnit2c//'^2/s^2', myThid )
27569449c9 Jean*0266       diagCode  = 'WM P    LR      '
                0267       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0268      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0269 
06349ecfd7 Jean*0270       diagName  = 'UE_VEL_C'
                0271       diagTitle = 'Eastward Velocity (m/s) (cell center)'
                0272       diagUnits = 'm/s             '
5644fea420 Jean*0273       diagCode  = 'UMR     MR      '
06349ecfd7 Jean*0274       diagMate  = diagNum + 2
                0275       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0276      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0277 
                0278       diagName  = 'VN_VEL_C'
                0279       diagTitle = 'Northward Velocity (m/s) (cell center)'
                0280       diagUnits = 'm/s             '
5644fea420 Jean*0281       diagCode  = 'VMR     MR      '
06349ecfd7 Jean*0282       diagMate  = diagNum
                0283       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0284      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0285 
db26b4dd29 Jean*0286       diagName  = 'UV_VEL_C'
                0287       diagTitle ='Product of horizontal Comp of velocity (cell center)'
                0288       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0289       diagCode  = 'UMR     MR      '
27569449c9 Jean*0290       diagMate  = diagNum + 1
                0291       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0292      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
db26b4dd29 Jean*0293 
                0294       diagName  = 'UV_VEL_Z'
8c664c8b58 Jean*0295       diagTitle = 'Meridional Transport of Zonal Momentum (m^2/s^2)'
db26b4dd29 Jean*0296       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0297       diagCode  = 'UZR     MR      '
27569449c9 Jean*0298       diagMate  = diagNum + 1
                0299       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0300      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0301 
6df3568a20 Jean*0302       diagName  = 'WU_VEL  '
9ba8ef70b1 Jean*0303       diagTitle = 'Vertical Transport of Zonal Momentum'
6df3568a20 Jean*0304       diagUnits = DIAGS_MK_UNITS( 'm.'//rUnit2c//'/s^2', myThid )
                0305       diagCode  = 'WU      LR      '
27569449c9 Jean*0306       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0307      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
6df3568a20 Jean*0308 
                0309       diagName  = 'WV_VEL  '
9ba8ef70b1 Jean*0310       diagTitle ='Vertical Transport of Meridional Momentum'
6df3568a20 Jean*0311       diagUnits = DIAGS_MK_UNITS( 'm.'//rUnit2c//'/s^2', myThid )
                0312       diagCode  = 'WV      LR      '
27569449c9 Jean*0313       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0314      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
db26b4dd29 Jean*0315 
8c664c8b58 Jean*0316       diagName  = 'UVELMASS'
                0317       diagTitle = 'Zonal Mass-Weighted Comp of Velocity (m/s)'
db26b4dd29 Jean*0318       diagUnits = 'm/s             '
5644fea420 Jean*0319       diagCode  = 'UUr     MR      '
27569449c9 Jean*0320       diagMate  = diagNum + 2
                0321       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0322      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0323 
                0324       diagName  = 'VVELMASS'
                0325       diagTitle = 'Meridional Mass-Weighted Comp of Velocity (m/s)'
db26b4dd29 Jean*0326       diagUnits = 'm/s             '
5644fea420 Jean*0327       diagCode  = 'VVr     MR      '
27569449c9 Jean*0328       diagMate  = diagNum
                0329       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0330      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0331 
                0332       diagName  = 'WVELMASS'
9ba8ef70b1 Jean*0333       diagTitle = 'Vertical Mass-Weighted Comp of Velocity'
db26b4dd29 Jean*0334       diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s', myThid )
fe67d630cb Jean*0335       diagCode  = 'WM      LR      '
27569449c9 Jean*0336       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0337      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0338 
d5ada4102b Jean*0339       diagName  = 'PhiVEL  '
                0340       diagTitle = 'Horizontal Velocity Potential (m^2/s)'
                0341       diagUnits = 'm^2/s           '
                0342       diagCode  = 'SMR P   MR      '
                0343 C-    use 'UVELMASS' as mate.
                0344       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0345      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0346 
e823c5094b Jean*0347       diagName  = 'PsiVEL  '
                0348       diagTitle = 'Horizontal Velocity Stream-Function'
                0349       diagUnits = DIAGS_MK_UNITS( rUnit2c//'.m^2/s', myThid )
                0350       diagCode  = 'SZ  P   MR      '
d5ada4102b Jean*0351 C-    use 'PhiVEL' as mate.
e823c5094b Jean*0352       diagMate  = diagNum
                0353       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0354      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
d5ada4102b Jean*0355 
8c664c8b58 Jean*0356       diagName  = 'UTHMASS '
9ba8ef70b1 Jean*0357       diagTitle = 'Zonal Mass-Weight Transp of Pot Temp'
db26b4dd29 Jean*0358       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0359       diagCode  = 'UUr     MR      '
27569449c9 Jean*0360       diagMate  = diagNum + 2
                0361       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0362      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0363 
                0364       diagName  = 'VTHMASS '
9ba8ef70b1 Jean*0365       diagTitle = 'Meridional Mass-Weight Transp of Pot Temp'
db26b4dd29 Jean*0366       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0367       diagCode  = 'VVr     MR      '
27569449c9 Jean*0368       diagMate  = diagNum
                0369       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0370      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0371 
                0372       diagName  = 'WTHMASS '
db26b4dd29 Jean*0373       diagTitle = 'Vertical Mass-Weight Transp of Pot Temp (K.m/s)'
                0374       diagUnits = DIAGS_MK_UNITS(tUnit4c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0375       diagCode  = 'WM      LR      '
27569449c9 Jean*0376       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0377      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0378 
                0379       diagName  = 'USLTMASS'
9ba8ef70b1 Jean*0380       diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of '
                0381      I                           //sTitle, myThid )
3ac63c2262 Jean*0382       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.m/s', myThid )
5644fea420 Jean*0383       diagCode  = 'UUr     MR      '
27569449c9 Jean*0384       diagMate  = diagNum + 2
                0385       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0386      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0387 
                0388       diagName  = 'VSLTMASS'
9ba8ef70b1 Jean*0389       diagTitle = DIAGS_MK_TITLE( 'Meridional Mass-Weight Transp of '
                0390      I                           //sTitle, myThid )
3ac63c2262 Jean*0391       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.m/s', myThid )
5644fea420 Jean*0392       diagCode  = 'VVr     MR      '
27569449c9 Jean*0393       diagMate  = diagNum
                0394       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0395      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0396 
                0397       diagName  = 'WSLTMASS'
9ba8ef70b1 Jean*0398       diagTitle = DIAGS_MK_TITLE( 'Vertical Mass-Weight Transp of '
                0399      I                           //sTitle, myThid )
3ac63c2262 Jean*0400       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0401       diagCode  = 'WM      LR      '
27569449c9 Jean*0402       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0403      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0404 
                0405       diagName  = 'UVELTH  '
9ba8ef70b1 Jean*0406       diagTitle = 'Zonal Transport of Pot Temp'
db26b4dd29 Jean*0407       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0408       diagCode  = 'UUR     MR      '
27569449c9 Jean*0409       diagMate  = diagNum + 2
                0410       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0411      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0412 
                0413       diagName  = 'VVELTH  '
9ba8ef70b1 Jean*0414       diagTitle = 'Meridional Transport of Pot Temp'
db26b4dd29 Jean*0415       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0416       diagCode  = 'VVR     MR      '
27569449c9 Jean*0417       diagMate  = diagNum
                0418       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0419      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0420 
                0421       diagName  = 'WVELTH  '
9ba8ef70b1 Jean*0422       diagTitle = 'Vertical Transport of Pot Temp'
db26b4dd29 Jean*0423       diagUnits = DIAGS_MK_UNITS(tUnit4c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0424       diagCode  = 'WM      LR      '
27569449c9 Jean*0425       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0426      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0427 
                0428       diagName  = 'UVELSLT '
9ba8ef70b1 Jean*0429       diagTitle = DIAGS_MK_TITLE( 'Zonal Transport of '
                0430      I                          //sTitle, myThid )
3ac63c2262 Jean*0431       diagUnits = DIAGS_MK_UNITS( sUnit5c//'.m/s', myThid )
5644fea420 Jean*0432       diagCode  = 'UUR     MR      '
27569449c9 Jean*0433       diagMate  = diagNum + 2
                0434       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0435      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0436 
                0437       diagName  = 'VVELSLT '
9ba8ef70b1 Jean*0438       diagTitle = DIAGS_MK_TITLE( 'Meridional Transport of '
                0439      I                          //sTitle, myThid )
3ac63c2262 Jean*0440       diagUnits = DIAGS_MK_UNITS( sUnit5c//'.m/s', myThid )
5644fea420 Jean*0441       diagCode  = 'VVR     MR      '
27569449c9 Jean*0442       diagMate  = diagNum
                0443       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0444      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0445 
                0446       diagName  = 'WVELSLT '
9ba8ef70b1 Jean*0447       diagTitle = DIAGS_MK_TITLE( 'Vertical Transport of '
                0448      I                          //sTitle, myThid )
3ac63c2262 Jean*0449       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0450       diagCode  = 'WM      LR      '
27569449c9 Jean*0451       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0452      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0453 
fefa9e1c60 Andr*0454       diagName  = 'UVELPHI '
f038101e44 Davi*0455       diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of '
9ba8ef70b1 Jean*0456      I                 //pTitle//' Anomaly', myThid )
fefa9e1c60 Andr*0457       diagUnits = 'm^3/s^3         '
5644fea420 Jean*0458       diagCode  = 'UUr     MR      '
27569449c9 Jean*0459       diagMate  = diagNum + 2
                0460       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0461      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fefa9e1c60 Andr*0462 
                0463       diagName  = 'VVELPHI '
8507f4494e Jean*0464       diagTitle = DIAGS_MK_TITLE( 'Merid. Mass-Weight Transp of '
9ba8ef70b1 Jean*0465      I                 //pTitle//' Anomaly', myThid )
fefa9e1c60 Andr*0466       diagUnits = 'm^3/s^3         '
5644fea420 Jean*0467       diagCode  = 'VVr     MR      '
27569449c9 Jean*0468       diagMate  = diagNum
                0469       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0470      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fefa9e1c60 Andr*0471 
8c664c8b58 Jean*0472 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0473 
c5c831ccb0 Jean*0474       diagName  = 'RHOAnoma'
                0475       diagTitle = 'Density Anomaly (=Rho-rhoConst)'
                0476       diagUnits = 'kg/m^3          '
5644fea420 Jean*0477       diagCode  = 'SMR     MR      '
27569449c9 Jean*0478       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0479      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
c5c831ccb0 Jean*0480 
fe67d630cb Jean*0481       diagName  = 'RHOANOSQ'
                0482       diagTitle = 'Square of Density Anomaly (=(Rho-rhoConst)^2)'
1d220297dd Jean*0483       diagUnits = 'kg^2/m^6        '
5644fea420 Jean*0484       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0485       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0486      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
fe67d630cb Jean*0487 
                0488       diagName  = 'URHOMASS'
                0489       diagTitle = 'Zonal Transport of Density'
                0490       diagUnits = 'kg/m^2/s        '
5644fea420 Jean*0491       diagCode  = 'UUr     MR      '
27569449c9 Jean*0492       diagMate  = diagNum + 2
                0493       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0494      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fe67d630cb Jean*0495 
                0496       diagName  = 'VRHOMASS'
                0497       diagTitle = 'Meridional Transport of Density'
                0498       diagUnits = 'kg/m^2/s        '
5644fea420 Jean*0499       diagCode  = 'VVr     MR      '
27569449c9 Jean*0500       diagMate  = diagNum
                0501       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0502      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fe67d630cb Jean*0503 
                0504       diagName  = 'WRHOMASS'
91331dd8db Jean*0505       diagTitle = 'Vertical Transport of Density'
                0506       diagUnits = 'kg/m^2/s        '
                0507       diagCode  = 'WM      LR      '
                0508       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0509      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0510 
                0511       diagName  = 'WdRHO_P '
                0512       diagTitle = 'Vertical velocity times delta^k(Rho)_at-const-P'
                0513       diagUnits = 'kg/m^2/s        '
                0514       diagCode  = 'WM      LR      '
                0515       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0516      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0517 
                0518       diagName  = 'WdRHOdP '
                0519       diagTitle = 'Vertical velocity times delta^k(Rho)_at-const-T,S'
fe67d630cb Jean*0520       diagUnits = 'kg/m^2/s        '
                0521       diagCode  = 'WM      LR      '
27569449c9 Jean*0522       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0523      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
fe67d630cb Jean*0524 
8ab13aa3ca Jean*0525       diagName  = 'PHIHYD  '
9ba8ef70b1 Jean*0526       diagTitle = DIAGS_MK_TITLE( 'Hydrostatic '
                0527      I                           //pTitle//' Anomaly', myThid )
db26b4dd29 Jean*0528       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0529       diagCode  = 'SMR     MR      '
27569449c9 Jean*0530       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0531      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8ab13aa3ca Jean*0532 
598b42a4a0 Andr*0533       diagName  = 'PHIHYDSQ'
9ba8ef70b1 Jean*0534       diagTitle = DIAGS_MK_TITLE( 'Square of Hyd. '
                0535      I                           //pTitle//' Anomaly', myThid )
598b42a4a0 Andr*0536       diagUnits = 'm^4/s^4         '
5644fea420 Jean*0537       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0538       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0539      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
598b42a4a0 Andr*0540 
8ab13aa3ca Jean*0541       diagName  = 'PHIBOT  '
9ba8ef70b1 Jean*0542 c     diagTitle = 'ocean bottom pressure / top. atmos geo-Potential'
                0543       diagTitle = DIAGS_MK_TITLE( fTitle
                0544      I                           //pTitle//' Anomaly', myThid )
db26b4dd29 Jean*0545       diagUnits = 'm^2/s^2         '
cdca88753e Dimi*0546       diagCode  = 'SM      M1      '
27569449c9 Jean*0547       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0548      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8ab13aa3ca Jean*0549 
                0550       diagName  = 'PHIBOTSQ'
9ba8ef70b1 Jean*0551 c     diagTitle = 'Square of ocean bottom pressure / top. geo-Potential'
                0552       diagTitle = DIAGS_MK_TITLE( 'Square of '//fTitle
                0553      I                           //pTitle//' Anomaly', myThid )
db26b4dd29 Jean*0554       diagUnits = 'm^4/s^4         '
27569449c9 Jean*0555       diagCode  = 'SM P    M1      '
                0556       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0557      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
ec20054932 Jean*0558 
6dc894f7a3 Jean*0559       diagName  = 'PHI_SURF'
                0560       diagTitle = DIAGS_MK_TITLE('Surface Dynamical '//pTitle,myThid)
                0561       diagUnits = 'm^2/s^2         '
                0562       diagCode  = 'SM      M1      '
                0563       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0564      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0565 
30cd0dda73 Jean*0566 #ifdef NONLIN_FRSURF
                0567       diagName  = 'PHIHYDcR'
                0568       diagTitle = DIAGS_MK_TITLE( 'Hydrostatic '
                0569      I                       //pTitle//' Anomaly @ const r', myThid )
                0570       diagUnits = 'm^2/s^2         '
                0571       diagCode  = 'SMR     MR      '
                0572       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0573      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0574 #endif
                0575 
5b0d80d77c Jean*0576 #ifdef ALLOW_NONHYDROSTATIC
                0577       diagName  = 'PHI_NH  '
                0578       diagTitle = DIAGS_MK_TITLE( 'Non-Hydrostatic '//pTitle, myThid )
                0579       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0580       diagCode  = 'SMR     MR      '
5b0d80d77c Jean*0581       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0582      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0583 #endif /* ALLOW_NONHYDROSTATIC */
                0584 
ec20054932 Jean*0585       diagName  = 'MXLDEPTH'
                0586       diagTitle = 'Mixed-Layer Depth (>0)'
                0587       diagUnits = 'm               '
                0588       diagCode  = 'SM      M1      '
27569449c9 Jean*0589       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0590      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8ab13aa3ca Jean*0591 
09ceb40cd6 Jean*0592       diagName  = 'DRHODR  '
db26b4dd29 Jean*0593       diagTitle = 'Stratification: d.Sigma/dr (kg/m3/r_unit)'
                0594       diagUnits = 'kg/m^4          '
                0595       IF ( usingPCoords ) diagUnits = 's^2/m^2         '
09ceb40cd6 Jean*0596       diagCode  = 'SM      LR      '
27569449c9 Jean*0597       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0598      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
09ceb40cd6 Jean*0599 
45a2b2d8df Jean*0600       diagName  = 'CONVADJ '
                0601       diagTitle = 'Convective Adjustment Index [0-1] '
                0602       diagUnits = 'fraction        '
5644fea420 Jean*0603       diagCode  = 'SMR     LR      '
27569449c9 Jean*0604       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0605      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
09ceb40cd6 Jean*0606 
8c664c8b58 Jean*0607 C--   surface fluxes:
ee8b184348 Jean*0608       diagName  = 'oceTAUX '
9ba8ef70b1 Jean*0609       diagTitle = 'zonal surface wind stress, >0 increases uVel'
db26b4dd29 Jean*0610       diagUnits = 'N/m^2           '
27569449c9 Jean*0611       diagCode  = 'UU      U1      '
                0612       diagMate  = diagNum + 2
                0613       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0614      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
e07cbe4110 Dimi*0615 
ee8b184348 Jean*0616       diagName  = 'oceTAUY '
9ba8ef70b1 Jean*0617       diagTitle = 'meridional surf. wind stress, >0 increases vVel'
db26b4dd29 Jean*0618       diagUnits = 'N/m^2           '
27569449c9 Jean*0619       diagCode  = 'VV      U1      '
                0620       diagMate  = diagNum
                0621       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0622      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
e07cbe4110 Dimi*0623 
b0bbe4b4f5 Jean*0624       diagName  = 'atmPload'
118f5617eb Jean*0625       diagTitle = 'Atmospheric pressure loading anomaly (vs surf_pRef)'
b0bbe4b4f5 Jean*0626       diagUnits = 'Pa              '
                0627       diagCode  = 'SM      U1      '
27569449c9 Jean*0628       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0629      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0630 
                0631       diagName  = 'sIceLoad'
                0632       diagTitle = 'sea-ice loading (in Mass of ice+snow / area unit)'
                0633       diagUnits = 'kg/m^2          '
                0634       diagCode  = 'SM      U1      '
27569449c9 Jean*0635       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0636      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0637 
                0638       diagName  = 'oceFWflx'
                0639       diagTitle = 'net surface Fresh-Water flux into the ocean'
                0640      &          //' (+=down), >0 decreases salinity'
                0641       diagUnits = 'kg/m^2/s        '
                0642       diagCode  = 'SM      U1      '
27569449c9 Jean*0643       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0644      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0645 
                0646       diagName  = 'oceSflux'
                0647       diagTitle = 'net surface Salt flux into the ocean (+=down),'
                0648      &          //' >0 increases salinity'
                0649       diagUnits = 'g/m^2/s         '
                0650       diagCode  = 'SM      U1      '
27569449c9 Jean*0651       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0652      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0653 
                0654       diagName  = 'oceQnet '
                0655       diagTitle = 'net surface heat flux into the ocean (+=down),'
                0656      &          //' >0 increases theta'
db26b4dd29 Jean*0657       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0658       diagCode  = 'SM      U1      '
27569449c9 Jean*0659       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0660      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0661 
b0bbe4b4f5 Jean*0662       diagName  = 'oceQsw  '
                0663       diagTitle = 'net Short-Wave radiation (+=down),'
                0664      &          //' >0 increases theta'
db26b4dd29 Jean*0665       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0666       diagCode  = 'SM      U1      '
27569449c9 Jean*0667       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0668      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0669 
ee8b184348 Jean*0670       diagName  = 'oceFreez'
b0bbe4b4f5 Jean*0671       diagTitle = 'heating from freezing of sea-water (allowFreezing=T)'
db26b4dd29 Jean*0672       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0673       diagCode  = 'SM      U1      '
27569449c9 Jean*0674       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0675      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0676 
b0bbe4b4f5 Jean*0677       diagName  = 'TRELAX  '
                0678       diagTitle = 'surface temperature relaxation, >0 increases theta'
                0679       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0680       diagCode  = 'SM      U1      '
27569449c9 Jean*0681       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0682      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0683 
                0684       diagName  = 'SRELAX  '
9ba8ef70b1 Jean*0685       diagTitle = 'surface salinity relaxation, >0 increases salt'
db26b4dd29 Jean*0686       diagUnits = 'g/m^2/s         '
e07cbe4110 Dimi*0687       diagCode  = 'SM      U1      '
27569449c9 Jean*0688       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0689      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0690 
b0bbe4b4f5 Jean*0691       diagName  = 'surForcT'
                0692       diagTitle = 'model surface forcing for Temperature,'
                0693      &          //' >0 increases theta'
256b55ab75 Patr*0694       diagUnits = 'W/m^2           '
                0695       diagCode  = 'SM      U1      '
27569449c9 Jean*0696       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0697      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
256b55ab75 Patr*0698 
b0bbe4b4f5 Jean*0699       diagName  = 'surForcS'
                0700       diagTitle = 'model surface forcing for Salinity,'
                0701      &          //' >0 increases salinity'
                0702       diagUnits = 'g/m^2/s         '
                0703       diagCode  = 'SM      U1      '
27569449c9 Jean*0704       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0705      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0706 
                0707       diagName  = 'TFLUX   '
                0708       diagTitle = 'total heat flux (match heat-content variations),'
                0709      &          //' >0 increases theta'
                0710       diagUnits = 'W/m^2           '
                0711       diagCode  = 'SM      U1      '
27569449c9 Jean*0712       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0713      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0714 
                0715       diagName  = 'SFLUX   '
                0716       diagTitle = 'total salt flux (match salt-content variations),'
                0717      &          //' >0 increases salt'
                0718       diagUnits = 'g/m^2/s         '
                0719       diagCode  = 'SM      U1      '
27569449c9 Jean*0720       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0721      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0722 
09ceb40cd6 Jean*0723 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b0bbe4b4f5 Jean*0724 
ee8b184348 Jean*0725       diagName  = 'RCENTER '
9ba8ef70b1 Jean*0726 c     diagTitle = 'Cell-Center r-Position (Pressure, Height) (Pa,m)'
                0727       diagTitle = DIAGS_MK_TITLE( 'Cell-Center '//rTitle, myThid )
b3f1143d77 Andr*0728       diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid )
                0729       diagCode  = 'SM      MR      '
27569449c9 Jean*0730       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0731      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
09ceb40cd6 Jean*0732 
b5d5d65bd6 Andr*0733       diagName  = 'RSURF   '
9ba8ef70b1 Jean*0734 c     diagTitle = 'Free-Surface r-Position (Pressure, Height) (Pa,m)'
                0735       diagTitle = DIAGS_MK_TITLE( eTitle//rTitle, myThid )
b5d5d65bd6 Andr*0736       diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid )
                0737       diagCode  = 'SM      M1      '
27569449c9 Jean*0738       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0739      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b5d5d65bd6 Andr*0740 
ff78c25591 Jean*0741       IF ( nonlinFreeSurf.GT.0 ) THEN
                0742        diagName  = 'hFactorC'
                0743        diagTitle = 'Center cell-thickness fraction [-]'
                0744        diagUnits = '1               '
                0745        diagCode  = 'SMr     MR      '
                0746        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0747      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0748 
                0749        diagName  = 'hFactorW'
                0750        diagTitle = 'Western-Edge cell-thickness fraction [-]'
                0751        diagUnits = '1               '
                0752        diagCode  = 'SUr     MR      '
                0753        diagMate  = diagNum + 2
                0754        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0755      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0756 
                0757        diagName  = 'hFactorS'
                0758        diagTitle = 'Southern-Edge cell-thickness fraction [-]'
                0759        diagUnits = '1               '
                0760        diagCode  = 'SVr     MR      '
                0761        diagMate  = diagNum
                0762        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0763      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0764       ENDIF
                0765 
cd2de3b19d Andr*0766       diagName  = 'TOTUTEND'
32bf9eae1e Jean*0767       diagTitle = 'Tendency of Zonal Component of Velocity'
                0768       diagUnits = 'm/s/day         '
5644fea420 Jean*0769       diagCode  = 'UUR     MR      '
27569449c9 Jean*0770       diagMate  = diagNum + 2
                0771       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0772      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
cd2de3b19d Andr*0773 
                0774       diagName  = 'TOTVTEND'
32bf9eae1e Jean*0775       diagTitle = 'Tendency of Meridional Component of Velocity'
                0776       diagUnits = 'm/s/day         '
5644fea420 Jean*0777       diagCode  = 'VVR     MR      '
27569449c9 Jean*0778       diagMate  = diagNum
                0779       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0780      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
cd2de3b19d Andr*0781 
f87dbbf625 Andr*0782       diagName  = 'TOTTTEND'
9ba8ef70b1 Jean*0783       diagTitle = 'Tendency of Potential Temperature'
32bf9eae1e Jean*0784       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/day', myThid )
5644fea420 Jean*0785       diagCode  = 'SMR     MR      '
27569449c9 Jean*0786       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0787      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
cd2de3b19d Andr*0788 
                0789       diagName  = 'TOTSTEND'
9ba8ef70b1 Jean*0790       diagTitle = DIAGS_MK_TITLE('Tendency of '//sTitle, myThid )
3ac63c2262 Jean*0791       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/day', myThid )
5644fea420 Jean*0792       diagCode  = 'SMR     MR      '
27569449c9 Jean*0793       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0794      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
cd2de3b19d Andr*0795 
bec73443e7 Jean*0796       diagName  = 'MoistCor'
                0797       diagTitle = 'Heating correction due to moist thermodynamics'
588c33df54 Jean*0798       diagUnits = 'W/m^2           '
bec73443e7 Jean*0799       diagCode  = 'SM      MR      '
                0800       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0801      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0802 
1a019c043f Jean*0803 #ifdef ALLOW_FRICTION_HEATING
                0804       diagName  = 'HeatDiss'
                0805       diagTitle = 'Heating from frictional dissipation'
                0806       diagUnits = 'W/m^2           '
                0807       diagCode  = 'SM      MR      '
                0808       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0809      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0810 #endif /* ALLOW_FRICTION_HEATING */
                0811 
8507f4494e Jean*0812 #ifdef ALLOW_GENERIC_ADVDIFF
                0813       diagName  = 'gT_Forc '
                0814       diagTitle = 'Potential Temp. forcing tendency'
                0815       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/s', myThid )
                0816       diagCode  = 'SMR     MR      '
9d9b466d92 Jean*0817       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
8507f4494e Jean*0818      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
9d9b466d92 Jean*0819 
8507f4494e Jean*0820       diagName  = 'gS_Forc '
                0821       diagTitle = DIAGS_MK_TITLE(
                0822      &            sTitle//'forcing tendency', myThid )
                0823       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/s', myThid )
                0824       diagCode  = 'SMR     MR      '
9d9b466d92 Jean*0825       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0826      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0827 
                0828       diagName  = 'AB_gT   '
                0829       diagTitle = 'Potential Temp. tendency from Adams-Bashforth'
                0830       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/s', myThid )
                0831       diagCode  = 'SMR     MR      '
                0832       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0833      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0834 
                0835       diagName  = 'AB_gS   '
                0836       diagTitle = DIAGS_MK_TITLE(
                0837      &            sTitle//'tendency from Adams-Bashforth', myThid )
3ac63c2262 Jean*0838       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/s', myThid )
9d9b466d92 Jean*0839       diagCode  = 'SMR     MR      '
                0840       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0841      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
004ff77f8b Jean*0842 
                0843       diagName  = 'gTinAB  '
                0844       diagTitle = 'Potential Temp. tendency going in Adams-Bashforth'
                0845       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/s', myThid )
                0846       diagCode  = 'SMR     MR      '
                0847       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0848      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0849 
                0850       diagName  = 'gSinAB  '
                0851       diagTitle = DIAGS_MK_TITLE(
                0852      &            sTitle//'tendency going in Adams-Bashforth', myThid )
3ac63c2262 Jean*0853       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/s', myThid )
004ff77f8b Jean*0854       diagCode  = 'SMR     MR      '
                0855       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0856      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
9d9b466d92 Jean*0857 #endif /* ALLOW_GENERIC_ADVDIFF */
b5d5d65bd6 Andr*0858 
8507f4494e Jean*0859       diagName  = 'AB_gU   '
                0860       diagTitle = 'U momentum tendency from Adams-Bashforth'
                0861       diagUnits = 'm/s^2           '
                0862       diagCode  = 'UUR     MR      '
                0863       diagMate  = diagNum + 2
                0864       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0865      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0866 
                0867       diagName  = 'AB_gV   '
                0868       diagTitle = 'V momentum tendency from Adams-Bashforth'
                0869       diagUnits = 'm/s^2           '
                0870       diagCode  = 'VVR     MR      '
                0871       diagMate  = diagNum
                0872       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0873      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0874 
                0875 #ifdef ALLOW_NONHYDROSTATIC
                0876       diagName  = 'AB_gW   '
                0877       diagTitle = 'W momentum tendency from Adams-Bashforth'
                0878       diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s^2', myThid )
                0879       diagCode  = 'WM      LR      '
                0880       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0881      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0882 #endif /* ALLOW_NONHYDROSTATIC */
                0883 
3f7db86ed7 Mich*0884 #ifdef ALLOW_EDDYPSI
                0885       diagName  = 'TAUXEDDY'
                0886       diagTitle = 'Zonal Eddy Stress'
1a019c043f Jean*0887       diagUnits = 'N/m^2           '
                0888       diagCode  = 'UU      LR      '
                0889       diagMate  = diagNum + 2
3f7db86ed7 Mich*0890       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0891      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8507f4494e Jean*0892 
3f7db86ed7 Mich*0893       diagName  = 'TAUYEDDY'
                0894       diagTitle = 'Meridional Eddy Stress'
1a019c043f Jean*0895       diagUnits = 'N/m^2           '
                0896       diagCode  = 'VV      LR      '
3f7db86ed7 Mich*0897       diagMate  = diagNum
                0898       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0899      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0900 
8507f4494e Jean*0901 # ifdef ALLOW_GMREDI
1a019c043f Jean*0902       diagName  = 'U_EulerM'
                0903       diagTitle = 'Zonal Eulerian-Mean Velocity (m/s)'
3f7db86ed7 Mich*0904       diagUnits = 'm/s             '
                0905       diagCode  = 'UUR     MR      '
                0906       diagMate  = diagNum + 2
                0907       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0908      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8507f4494e Jean*0909 
1a019c043f Jean*0910       diagName  = 'V_EulerM'
                0911       diagTitle = 'Meridional Eulerian-Mean Velocity (m/s)'
3f7db86ed7 Mich*0912       diagUnits = 'm/s             '
                0913       diagCode  = 'VVR     MR      '
                0914       diagMate  = diagNum
                0915       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0916      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8507f4494e Jean*0917 # endif /* ALLOW_GMREDI */
3f7db86ed7 Mich*0918 #endif /* ALLOW_EDDYPSI */
0b5a2cab3a Mich*0919 
fb2c52cfd4 aver*0920 #ifdef INCLUDE_SOUNDSPEED_CALC_CODE
                0921       diagName  = 'CSound  '
                0922       diagTitle = 'Speed of sound in seawater'
                0923       diagUnits = 'm/s             '
                0924       diagCode  = 'SMR     MR      '
                0925       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0926      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0927 #endif
                0928 
a5ec81ed49 Timo*0929 #ifdef ALLOW_AUTODIFF
                0930 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0931 C     Adjoint state variables
                0932 
                0933       diagName  = 'ADJetan '
                0934       diagTitle = 'dJ/dEtaN: Sensitivity to sea surface height anomaly'
                0935       diagUnits = 'dJ/m            '
                0936       diagCode  = 'SM A    M1      '
                0937       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0938      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0939 
                0940       diagName  = 'ADJuvel '
                0941       diagTitle = 'dJ/dU: Sensitivity to zonal velocity'
                0942       diagUnits = 'dJ/(m/s)        '
                0943       diagCode  = 'UURA    MR      '
                0944       diagMate  = diagNum + 2
                0945       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0946      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0947 
                0948       diagName  = 'ADJvvel '
                0949       diagTitle = 'dJ/dV: Sensitivity to meridional velocity'
                0950       diagUnits = 'dJ/(m/s)        '
                0951       diagCode  = 'VVRA    MR      '
                0952       diagMate  = diagNum
                0953       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0954      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0955 
                0956       diagName  = 'ADJwvel '
                0957       diagTitle = 'dJ/dW: Sensitivity to vertical velocity'
                0958       diagUnits = 'dJ/(m/s)        '
                0959       diagCode  = 'WM A    LR      '
                0960       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0961      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0962 
                0963       diagName  = 'ADJtheta'
                0964       diagTitle = 'dJ/dTheta: Sensitivity to potential temperature'
                0965       diagUnits = 'dJ/degC         '
                0966       diagCode  = 'SMRA    MR      '
                0967       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0968      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0969 
                0970       diagName  = 'ADJsalt '
                0971       diagTitle = 'dJ/dSalt: Sensitivity to salinity'
a10c595eb6 Timo*0972       diagUnits = DIAGS_MK_UNITS( 'dJ/('//sUnit5c//')', myThid )
a5ec81ed49 Timo*0973       diagCode  = 'SMRA    MR      '
                0974       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0975      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0976 
a10c595eb6 Timo*0977 C--   surface fluxes:
                0978       diagName  = 'ADJtaux '
                0979       diagTitle = 'dJ/dTaux: Senstivity to zonal surface wind stress'
                0980       diagUnits = 'dJ/(N/m^2)      '
                0981       diagCode  = 'UU A    U1      '
                0982       diagMate  = diagNum + 2
                0983       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0984      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0985 
                0986       diagName  = 'ADJtauy '
                0987       diagTitle = 'dJ/dTauy: Sensitivity to merid. surface wind stress'
                0988       diagUnits = 'dJ/(N/m^2)      '
                0989       diagCode  = 'VV A    U1      '
                0990       diagMate  = diagNum
                0991       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0992      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0993 
                0994       diagName  = 'ADJempmr'
1ac7194fda plla*0995       diagTitle = 'dJ/dEmPmR: Sensitivity to net surface freshwater'
                0996      &          //' flux'
a10c595eb6 Timo*0997       diagUnits = 'dJ/(kg/m^2/s)   '
                0998       diagCode  = 'SM A    U1      '
                0999       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1000      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1001 
                1002       diagName  = 'ADJqnet'
                1003       diagTitle = 'dJ/dQnet: Sensitivity to net surface heat flux'
                1004       diagUnits = 'dJ/(W/m^2)      '
                1005       diagCode  = 'SM A    U1      '
                1006       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1007      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1008 
                1009       diagName  = 'ADJqsw'
                1010       diagTitle = 'dJ/dQsw: Sensitivitiy to net Short-Wave radiation'
                1011       diagUnits = 'dJ/(W/m^2)      '
                1012       diagCode  = 'SM A    U1      '
                1013       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1014      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1015 
                1016       diagName  = 'ADJsst'
                1017       diagTitle = 'dJ/dSST: Sensitivity to Sea Surface Temperature'
                1018       diagUnits = 'dJ/K'
                1019       diagCode  = 'SM A    M1      '
                1020       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1021      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1022 
                1023       diagName  = 'ADJsss'
                1024       diagTitle = 'dJ/dSSS: Sensitivity to Sea Surface Salinity'
                1025       diagUnits = DIAGS_MK_UNITS( 'dJ/('//sUnit5c//')', myThid )
                1026       diagCode  = 'SM A    M1      '
                1027       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1028      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1029 
                1030 C Miscellanious
                1031       diagName  = 'ADJbtdrg'
                1032       diagTitle = 'dJ/dCd: Sensitivity to bottom drag coefficient'
                1033       diagUnits = 'dJ/d()'
                1034       diagCode  = 'SM A    M1      '
                1035       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1036      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1037 
                1038       diagName  = 'ADJdifkr'
                1039       diagTitle = 'dJ/dKr: Sensitivity to vertical diffusivity'
                1040       diagUnits = 'dJ/d(m^2/s))'
                1041       diagCode  = 'SMRA    MR      '
                1042       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1043      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1044 
                1045       diagName  = 'ADJepsix'
                1046       diagTitle = 'dJ/dEddyPsiX: Sensitivity to zonal eddy'
                1047      &           //'streamfunction'
                1048       diagUnits = 'dJ/(m^2/s)      '
                1049       diagCode  = 'UURA    UR      '
                1050       diagMate  = diagNum + 2
                1051       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1052      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                1053 
                1054       diagName  = 'ADJepsiy'
                1055       diagTitle = 'dJ/dEddyPsiY: Sensitivity to meridional eddy'
                1056      &           //'streamfunction'
                1057       diagUnits = 'dJ/(m^2/s)      '
                1058       diagCode  = 'VVRA    UR      '
                1059       diagMate  = diagNum
                1060       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1061      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
a5ec81ed49 Timo*1062 #endif /* ALLOW_AUTODIFF */
                1063 
09ceb40cd6 Jean*1064       RETURN
                1065       END