Back to home page

MITgcm

 
 

    


File indexing completed on 2025-11-04 06:08:41 UTC

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