Back to home page

MITgcm

 
 

    


File indexing completed on 2025-02-02 06:10:59 UTC

view on githubraw file Latest commit 701e10a9 on 2025-02-01 19:15:20 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 
06349ecfd7 Jean*0280       diagName  = 'UE_VEL_C'
                0281       diagTitle = 'Eastward Velocity (m/s) (cell center)'
                0282       diagUnits = 'm/s             '
5644fea420 Jean*0283       diagCode  = 'UMR     MR      '
06349ecfd7 Jean*0284       diagMate  = diagNum + 2
                0285       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0286      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0287 
                0288       diagName  = 'VN_VEL_C'
                0289       diagTitle = 'Northward Velocity (m/s) (cell center)'
                0290       diagUnits = 'm/s             '
5644fea420 Jean*0291       diagCode  = 'VMR     MR      '
06349ecfd7 Jean*0292       diagMate  = diagNum
                0293       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0294      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0295 
db26b4dd29 Jean*0296       diagName  = 'UV_VEL_C'
                0297       diagTitle ='Product of horizontal Comp of velocity (cell center)'
                0298       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0299       diagCode  = 'UMR     MR      '
27569449c9 Jean*0300       diagMate  = diagNum + 1
                0301       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0302      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
db26b4dd29 Jean*0303 
                0304       diagName  = 'UV_VEL_Z'
8c664c8b58 Jean*0305       diagTitle = 'Meridional Transport of Zonal Momentum (m^2/s^2)'
db26b4dd29 Jean*0306       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0307       diagCode  = 'UZR     MR      '
27569449c9 Jean*0308       diagMate  = diagNum + 1
                0309       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0310      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0311 
6df3568a20 Jean*0312       diagName  = 'WU_VEL  '
9ba8ef70b1 Jean*0313       diagTitle = 'Vertical Transport of Zonal Momentum'
6df3568a20 Jean*0314       diagUnits = DIAGS_MK_UNITS( 'm.'//rUnit2c//'/s^2', myThid )
                0315       diagCode  = 'WU      LR      '
27569449c9 Jean*0316       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0317      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
6df3568a20 Jean*0318 
                0319       diagName  = 'WV_VEL  '
9ba8ef70b1 Jean*0320       diagTitle ='Vertical Transport of Meridional Momentum'
6df3568a20 Jean*0321       diagUnits = DIAGS_MK_UNITS( 'm.'//rUnit2c//'/s^2', myThid )
                0322       diagCode  = 'WV      LR      '
27569449c9 Jean*0323       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0324      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
db26b4dd29 Jean*0325 
8c664c8b58 Jean*0326       diagName  = 'UVELMASS'
                0327       diagTitle = 'Zonal Mass-Weighted Comp of Velocity (m/s)'
db26b4dd29 Jean*0328       diagUnits = 'm/s             '
5644fea420 Jean*0329       diagCode  = 'UUr     MR      '
27569449c9 Jean*0330       diagMate  = diagNum + 2
                0331       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0332      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0333 
                0334       diagName  = 'VVELMASS'
                0335       diagTitle = 'Meridional Mass-Weighted Comp of Velocity (m/s)'
db26b4dd29 Jean*0336       diagUnits = 'm/s             '
5644fea420 Jean*0337       diagCode  = 'VVr     MR      '
27569449c9 Jean*0338       diagMate  = diagNum
                0339       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0340      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0341 
                0342       diagName  = 'WVELMASS'
9ba8ef70b1 Jean*0343       diagTitle = 'Vertical Mass-Weighted Comp of Velocity'
db26b4dd29 Jean*0344       diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s', myThid )
fe67d630cb Jean*0345       diagCode  = 'WM      LR      '
27569449c9 Jean*0346       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0347      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0348 
d5ada4102b Jean*0349       diagName  = 'PhiVEL  '
                0350       diagTitle = 'Horizontal Velocity Potential (m^2/s)'
                0351       diagUnits = 'm^2/s           '
                0352       diagCode  = 'SMR P   MR      '
                0353 C-    use 'UVELMASS' as mate.
                0354       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0355      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0356 
e823c5094b Jean*0357       diagName  = 'PsiVEL  '
                0358       diagTitle = 'Horizontal Velocity Stream-Function'
                0359       diagUnits = DIAGS_MK_UNITS( rUnit2c//'.m^2/s', myThid )
                0360       diagCode  = 'SZ  P   MR      '
d5ada4102b Jean*0361 C-    use 'PhiVEL' as mate.
e823c5094b Jean*0362       diagMate  = diagNum
                0363       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0364      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
d5ada4102b Jean*0365 
8c664c8b58 Jean*0366       diagName  = 'UTHMASS '
701e10a905 Mart*0367       diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of '
                0368      I                           //tTitl2, myThid )
db26b4dd29 Jean*0369       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0370       diagCode  = 'UUr     MR      '
27569449c9 Jean*0371       diagMate  = diagNum + 2
                0372       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0373      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0374 
                0375       diagName  = 'VTHMASS '
701e10a905 Mart*0376       diagTitle = DIAGS_MK_TITLE( 'Meridional Mass-Weight Transp of '
                0377      I                           //tTitl2, myThid )
db26b4dd29 Jean*0378       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0379       diagCode  = 'VVr     MR      '
27569449c9 Jean*0380       diagMate  = diagNum
                0381       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0382      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0383 
                0384       diagName  = 'WTHMASS '
701e10a905 Mart*0385       diagTitle = DIAGS_MK_TITLE( 'Vertical Mass-Weight Transp of '
                0386      I                           //tTitl2, myThid )
db26b4dd29 Jean*0387       diagUnits = DIAGS_MK_UNITS(tUnit4c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0388       diagCode  = 'WM      LR      '
27569449c9 Jean*0389       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0390      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0391 
                0392       diagName  = 'USLTMASS'
9ba8ef70b1 Jean*0393       diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of '
                0394      I                           //sTitle, myThid )
3ac63c2262 Jean*0395       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.m/s', myThid )
5644fea420 Jean*0396       diagCode  = 'UUr     MR      '
27569449c9 Jean*0397       diagMate  = diagNum + 2
                0398       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0399      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0400 
                0401       diagName  = 'VSLTMASS'
9ba8ef70b1 Jean*0402       diagTitle = DIAGS_MK_TITLE( 'Meridional Mass-Weight Transp of '
                0403      I                           //sTitle, myThid )
3ac63c2262 Jean*0404       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.m/s', myThid )
5644fea420 Jean*0405       diagCode  = 'VVr     MR      '
27569449c9 Jean*0406       diagMate  = diagNum
                0407       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0408      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0409 
                0410       diagName  = 'WSLTMASS'
9ba8ef70b1 Jean*0411       diagTitle = DIAGS_MK_TITLE( 'Vertical Mass-Weight Transp of '
                0412      I                           //sTitle, myThid )
3ac63c2262 Jean*0413       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0414       diagCode  = 'WM      LR      '
27569449c9 Jean*0415       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0416      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0417 
                0418       diagName  = 'UVELTH  '
701e10a905 Mart*0419       diagTitle = DIAGS_MK_TITLE( 'Zonal Transport of '
                0420      I                          //tTitl2, myThid )
db26b4dd29 Jean*0421       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0422       diagCode  = 'UUR     MR      '
27569449c9 Jean*0423       diagMate  = diagNum + 2
                0424       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0425      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0426 
                0427       diagName  = 'VVELTH  '
701e10a905 Mart*0428       diagTitle = DIAGS_MK_TITLE( 'Meridional Transport of '
                0429      I                          //tTitl2, myThid )
db26b4dd29 Jean*0430       diagUnits = DIAGS_MK_UNITS( tUnit4c//'.m/s', myThid )
5644fea420 Jean*0431       diagCode  = 'VVR     MR      '
27569449c9 Jean*0432       diagMate  = diagNum
                0433       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0434      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0435 
                0436       diagName  = 'WVELTH  '
701e10a905 Mart*0437       diagTitle = DIAGS_MK_TITLE( 'Vertical Transport of '
                0438      I                          //tTitl2, myThid )
db26b4dd29 Jean*0439       diagUnits = DIAGS_MK_UNITS(tUnit4c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0440       diagCode  = 'WM      LR      '
27569449c9 Jean*0441       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0442      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0443 
                0444       diagName  = 'UVELSLT '
9ba8ef70b1 Jean*0445       diagTitle = DIAGS_MK_TITLE( 'Zonal Transport of '
                0446      I                          //sTitle, myThid )
3ac63c2262 Jean*0447       diagUnits = DIAGS_MK_UNITS( sUnit5c//'.m/s', myThid )
5644fea420 Jean*0448       diagCode  = 'UUR     MR      '
27569449c9 Jean*0449       diagMate  = diagNum + 2
                0450       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0451      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0452 
                0453       diagName  = 'VVELSLT '
9ba8ef70b1 Jean*0454       diagTitle = DIAGS_MK_TITLE( 'Meridional Transport of '
                0455      I                          //sTitle, myThid )
3ac63c2262 Jean*0456       diagUnits = DIAGS_MK_UNITS( sUnit5c//'.m/s', myThid )
5644fea420 Jean*0457       diagCode  = 'VVR     MR      '
27569449c9 Jean*0458       diagMate  = diagNum
                0459       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0460      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8c664c8b58 Jean*0461 
                0462       diagName  = 'WVELSLT '
9ba8ef70b1 Jean*0463       diagTitle = DIAGS_MK_TITLE( 'Vertical Transport of '
                0464      I                          //sTitle, myThid )
3ac63c2262 Jean*0465       diagUnits = DIAGS_MK_UNITS(sUnit5c//'.'//rUnit2c//'/s', myThid )
fe67d630cb Jean*0466       diagCode  = 'WM      LR      '
27569449c9 Jean*0467       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0468      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8c664c8b58 Jean*0469 
fefa9e1c60 Andr*0470       diagName  = 'UVELPHI '
f038101e44 Davi*0471       diagTitle = DIAGS_MK_TITLE( 'Zonal Mass-Weight Transp of '
9ba8ef70b1 Jean*0472      I                 //pTitle//' Anomaly', myThid )
fefa9e1c60 Andr*0473       diagUnits = 'm^3/s^3         '
5644fea420 Jean*0474       diagCode  = 'UUr     MR      '
27569449c9 Jean*0475       diagMate  = diagNum + 2
                0476       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0477      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fefa9e1c60 Andr*0478 
                0479       diagName  = 'VVELPHI '
8507f4494e Jean*0480       diagTitle = DIAGS_MK_TITLE( 'Merid. Mass-Weight Transp of '
9ba8ef70b1 Jean*0481      I                 //pTitle//' Anomaly', myThid )
fefa9e1c60 Andr*0482       diagUnits = 'm^3/s^3         '
5644fea420 Jean*0483       diagCode  = 'VVr     MR      '
27569449c9 Jean*0484       diagMate  = diagNum
                0485       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0486      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fefa9e1c60 Andr*0487 
8c664c8b58 Jean*0488 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0489 
c5c831ccb0 Jean*0490       diagName  = 'RHOAnoma'
                0491       diagTitle = 'Density Anomaly (=Rho-rhoConst)'
                0492       diagUnits = 'kg/m^3          '
5644fea420 Jean*0493       diagCode  = 'SMR     MR      '
27569449c9 Jean*0494       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0495      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
c5c831ccb0 Jean*0496 
fe67d630cb Jean*0497       diagName  = 'RHOANOSQ'
                0498       diagTitle = 'Square of Density Anomaly (=(Rho-rhoConst)^2)'
1d220297dd Jean*0499       diagUnits = 'kg^2/m^6        '
5644fea420 Jean*0500       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0501       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0502      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
fe67d630cb Jean*0503 
                0504       diagName  = 'URHOMASS'
                0505       diagTitle = 'Zonal Transport of Density'
                0506       diagUnits = 'kg/m^2/s        '
5644fea420 Jean*0507       diagCode  = 'UUr     MR      '
27569449c9 Jean*0508       diagMate  = diagNum + 2
                0509       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0510      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fe67d630cb Jean*0511 
                0512       diagName  = 'VRHOMASS'
                0513       diagTitle = 'Meridional Transport of Density'
                0514       diagUnits = 'kg/m^2/s        '
5644fea420 Jean*0515       diagCode  = 'VVr     MR      '
27569449c9 Jean*0516       diagMate  = diagNum
                0517       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0518      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
fe67d630cb Jean*0519 
                0520       diagName  = 'WRHOMASS'
91331dd8db Jean*0521       diagTitle = 'Vertical Transport of Density'
                0522       diagUnits = 'kg/m^2/s        '
                0523       diagCode  = 'WM      LR      '
                0524       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0525      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0526 
                0527       diagName  = 'WdRHO_P '
                0528       diagTitle = 'Vertical velocity times delta^k(Rho)_at-const-P'
                0529       diagUnits = 'kg/m^2/s        '
                0530       diagCode  = 'WM      LR      '
                0531       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0532      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0533 
                0534       diagName  = 'WdRHOdP '
                0535       diagTitle = 'Vertical velocity times delta^k(Rho)_at-const-T,S'
fe67d630cb Jean*0536       diagUnits = 'kg/m^2/s        '
                0537       diagCode  = 'WM      LR      '
27569449c9 Jean*0538       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0539      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
fe67d630cb Jean*0540 
8ab13aa3ca Jean*0541       diagName  = 'PHIHYD  '
9ba8ef70b1 Jean*0542       diagTitle = DIAGS_MK_TITLE( 'Hydrostatic '
                0543      I                           //pTitle//' Anomaly', myThid )
db26b4dd29 Jean*0544       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0545       diagCode  = 'SMR     MR      '
27569449c9 Jean*0546       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0547      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8ab13aa3ca Jean*0548 
598b42a4a0 Andr*0549       diagName  = 'PHIHYDSQ'
9ba8ef70b1 Jean*0550       diagTitle = DIAGS_MK_TITLE( 'Square of Hyd. '
                0551      I                           //pTitle//' Anomaly', myThid )
598b42a4a0 Andr*0552       diagUnits = 'm^4/s^4         '
5644fea420 Jean*0553       diagCode  = 'SMRP    MR      '
27569449c9 Jean*0554       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0555      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
598b42a4a0 Andr*0556 
8ab13aa3ca Jean*0557       diagName  = 'PHIBOT  '
9ba8ef70b1 Jean*0558 c     diagTitle = 'ocean bottom pressure / top. atmos geo-Potential'
                0559       diagTitle = DIAGS_MK_TITLE( fTitle
                0560      I                           //pTitle//' Anomaly', myThid )
db26b4dd29 Jean*0561       diagUnits = 'm^2/s^2         '
cdca88753e Dimi*0562       diagCode  = 'SM      M1      '
27569449c9 Jean*0563       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0564      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8ab13aa3ca Jean*0565 
                0566       diagName  = 'PHIBOTSQ'
9ba8ef70b1 Jean*0567 c     diagTitle = 'Square of ocean bottom pressure / top. geo-Potential'
                0568       diagTitle = DIAGS_MK_TITLE( 'Square of '//fTitle
                0569      I                           //pTitle//' Anomaly', myThid )
db26b4dd29 Jean*0570       diagUnits = 'm^4/s^4         '
27569449c9 Jean*0571       diagCode  = 'SM P    M1      '
                0572       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0573      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
ec20054932 Jean*0574 
6dc894f7a3 Jean*0575       diagName  = 'PHI_SURF'
                0576       diagTitle = DIAGS_MK_TITLE('Surface Dynamical '//pTitle,myThid)
                0577       diagUnits = 'm^2/s^2         '
                0578       diagCode  = 'SM      M1      '
                0579       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0580      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0581 
30cd0dda73 Jean*0582 #ifdef NONLIN_FRSURF
                0583       diagName  = 'PHIHYDcR'
                0584       diagTitle = DIAGS_MK_TITLE( 'Hydrostatic '
                0585      I                       //pTitle//' Anomaly @ const r', myThid )
                0586       diagUnits = 'm^2/s^2         '
                0587       diagCode  = 'SMR     MR      '
                0588       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0589      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0590 #endif
                0591 
5b0d80d77c Jean*0592 #ifdef ALLOW_NONHYDROSTATIC
                0593       diagName  = 'PHI_NH  '
                0594       diagTitle = DIAGS_MK_TITLE( 'Non-Hydrostatic '//pTitle, myThid )
                0595       diagUnits = 'm^2/s^2         '
5644fea420 Jean*0596       diagCode  = 'SMR     MR      '
5b0d80d77c Jean*0597       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0598      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0599 #endif /* ALLOW_NONHYDROSTATIC */
                0600 
ec20054932 Jean*0601       diagName  = 'MXLDEPTH'
                0602       diagTitle = 'Mixed-Layer Depth (>0)'
                0603       diagUnits = 'm               '
                0604       diagCode  = 'SM      M1      '
27569449c9 Jean*0605       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0606      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
8ab13aa3ca Jean*0607 
09ceb40cd6 Jean*0608       diagName  = 'DRHODR  '
db26b4dd29 Jean*0609       diagTitle = 'Stratification: d.Sigma/dr (kg/m3/r_unit)'
                0610       diagUnits = 'kg/m^4          '
                0611       IF ( usingPCoords ) diagUnits = 's^2/m^2         '
09ceb40cd6 Jean*0612       diagCode  = 'SM      LR      '
27569449c9 Jean*0613       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0614      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
09ceb40cd6 Jean*0615 
45a2b2d8df Jean*0616       diagName  = 'CONVADJ '
                0617       diagTitle = 'Convective Adjustment Index [0-1] '
                0618       diagUnits = 'fraction        '
5644fea420 Jean*0619       diagCode  = 'SMR     LR      '
27569449c9 Jean*0620       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0621      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
09ceb40cd6 Jean*0622 
701e10a905 Mart*0623       IF ( fluidIsWater ) THEN
                0624        diagName  = 'GCM_SST '
                0625        diagTitle = 'Ocean model Sea Surface Temperature (degC)'
                0626        diagUnits = 'degC            '
                0627        diagCode  = 'SM      M1      '
                0628        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0629      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0630       ENDIF
                0631 
8c664c8b58 Jean*0632 C--   surface fluxes:
ee8b184348 Jean*0633       diagName  = 'oceTAUX '
9ba8ef70b1 Jean*0634       diagTitle = 'zonal surface wind stress, >0 increases uVel'
db26b4dd29 Jean*0635       diagUnits = 'N/m^2           '
27569449c9 Jean*0636       diagCode  = 'UU      U1      '
                0637       diagMate  = diagNum + 2
                0638       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0639      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
e07cbe4110 Dimi*0640 
ee8b184348 Jean*0641       diagName  = 'oceTAUY '
9ba8ef70b1 Jean*0642       diagTitle = 'meridional surf. wind stress, >0 increases vVel'
db26b4dd29 Jean*0643       diagUnits = 'N/m^2           '
27569449c9 Jean*0644       diagCode  = 'VV      U1      '
                0645       diagMate  = diagNum
                0646       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0647      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
e07cbe4110 Dimi*0648 
b0bbe4b4f5 Jean*0649       diagName  = 'atmPload'
118f5617eb Jean*0650       diagTitle = 'Atmospheric pressure loading anomaly (vs surf_pRef)'
b0bbe4b4f5 Jean*0651       diagUnits = 'Pa              '
                0652       diagCode  = 'SM      U1      '
27569449c9 Jean*0653       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0654      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0655 
                0656       diagName  = 'sIceLoad'
                0657       diagTitle = 'sea-ice loading (in Mass of ice+snow / area unit)'
                0658       diagUnits = 'kg/m^2          '
                0659       diagCode  = 'SM      U1      '
27569449c9 Jean*0660       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0661      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0662 
                0663       diagName  = 'oceFWflx'
                0664       diagTitle = 'net surface Fresh-Water flux into the ocean'
                0665      &          //' (+=down), >0 decreases salinity'
                0666       diagUnits = 'kg/m^2/s        '
                0667       diagCode  = 'SM      U1      '
27569449c9 Jean*0668       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0669      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0670 
                0671       diagName  = 'oceSflux'
                0672       diagTitle = 'net surface Salt flux into the ocean (+=down),'
                0673      &          //' >0 increases salinity'
                0674       diagUnits = 'g/m^2/s         '
                0675       diagCode  = 'SM      U1      '
27569449c9 Jean*0676       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0677      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0678 
                0679       diagName  = 'oceQnet '
                0680       diagTitle = 'net surface heat flux into the ocean (+=down),'
                0681      &          //' >0 increases theta'
db26b4dd29 Jean*0682       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0683       diagCode  = 'SM      U1      '
27569449c9 Jean*0684       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0685      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0686 
b0bbe4b4f5 Jean*0687       diagName  = 'oceQsw  '
                0688       diagTitle = 'net Short-Wave radiation (+=down),'
                0689      &          //' >0 increases theta'
db26b4dd29 Jean*0690       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0691       diagCode  = 'SM      U1      '
27569449c9 Jean*0692       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0693      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0694 
ee8b184348 Jean*0695       diagName  = 'oceFreez'
b0bbe4b4f5 Jean*0696       diagTitle = 'heating from freezing of sea-water (allowFreezing=T)'
db26b4dd29 Jean*0697       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0698       diagCode  = 'SM      U1      '
27569449c9 Jean*0699       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0700      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0701 
b0bbe4b4f5 Jean*0702       diagName  = 'TRELAX  '
                0703       diagTitle = 'surface temperature relaxation, >0 increases theta'
                0704       diagUnits = 'W/m^2           '
e07cbe4110 Dimi*0705       diagCode  = 'SM      U1      '
27569449c9 Jean*0706       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0707      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0708 
                0709       diagName  = 'SRELAX  '
9ba8ef70b1 Jean*0710       diagTitle = 'surface salinity relaxation, >0 increases salt'
db26b4dd29 Jean*0711       diagUnits = 'g/m^2/s         '
e07cbe4110 Dimi*0712       diagCode  = 'SM      U1      '
27569449c9 Jean*0713       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0714      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
e07cbe4110 Dimi*0715 
b0bbe4b4f5 Jean*0716       diagName  = 'surForcT'
                0717       diagTitle = 'model surface forcing for Temperature,'
                0718      &          //' >0 increases theta'
256b55ab75 Patr*0719       diagUnits = 'W/m^2           '
                0720       diagCode  = 'SM      U1      '
27569449c9 Jean*0721       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0722      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
256b55ab75 Patr*0723 
b0bbe4b4f5 Jean*0724       diagName  = 'surForcS'
                0725       diagTitle = 'model surface forcing for Salinity,'
                0726      &          //' >0 increases salinity'
                0727       diagUnits = 'g/m^2/s         '
                0728       diagCode  = 'SM      U1      '
27569449c9 Jean*0729       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0730      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0731 
                0732       diagName  = 'TFLUX   '
                0733       diagTitle = 'total heat flux (match heat-content variations),'
                0734      &          //' >0 increases theta'
                0735       diagUnits = 'W/m^2           '
                0736       diagCode  = 'SM      U1      '
27569449c9 Jean*0737       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0738      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0739 
                0740       diagName  = 'SFLUX   '
                0741       diagTitle = 'total salt flux (match salt-content variations),'
                0742      &          //' >0 increases salt'
                0743       diagUnits = 'g/m^2/s         '
                0744       diagCode  = 'SM      U1      '
27569449c9 Jean*0745       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0746      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b0bbe4b4f5 Jean*0747 
09ceb40cd6 Jean*0748 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
b0bbe4b4f5 Jean*0749 
ee8b184348 Jean*0750       diagName  = 'RCENTER '
9ba8ef70b1 Jean*0751 c     diagTitle = 'Cell-Center r-Position (Pressure, Height) (Pa,m)'
                0752       diagTitle = DIAGS_MK_TITLE( 'Cell-Center '//rTitle, myThid )
b3f1143d77 Andr*0753       diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid )
                0754       diagCode  = 'SM      MR      '
27569449c9 Jean*0755       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0756      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
09ceb40cd6 Jean*0757 
b5d5d65bd6 Andr*0758       diagName  = 'RSURF   '
9ba8ef70b1 Jean*0759 c     diagTitle = 'Free-Surface r-Position (Pressure, Height) (Pa,m)'
                0760       diagTitle = DIAGS_MK_TITLE( eTitle//rTitle, myThid )
b5d5d65bd6 Andr*0761       diagUnits = DIAGS_MK_UNITS( rUnit2c, myThid )
                0762       diagCode  = 'SM      M1      '
27569449c9 Jean*0763       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0764      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
b5d5d65bd6 Andr*0765 
ff78c25591 Jean*0766       IF ( nonlinFreeSurf.GT.0 ) THEN
                0767        diagName  = 'hFactorC'
                0768        diagTitle = 'Center cell-thickness fraction [-]'
                0769        diagUnits = '1               '
                0770        diagCode  = 'SMr     MR      '
                0771        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0772      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0773 
                0774        diagName  = 'hFactorW'
                0775        diagTitle = 'Western-Edge cell-thickness fraction [-]'
                0776        diagUnits = '1               '
                0777        diagCode  = 'SUr     MR      '
                0778        diagMate  = diagNum + 2
                0779        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0780      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0781 
                0782        diagName  = 'hFactorS'
                0783        diagTitle = 'Southern-Edge cell-thickness fraction [-]'
                0784        diagUnits = '1               '
                0785        diagCode  = 'SVr     MR      '
                0786        diagMate  = diagNum
                0787        CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0788      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0789       ENDIF
                0790 
cd2de3b19d Andr*0791       diagName  = 'TOTUTEND'
32bf9eae1e Jean*0792       diagTitle = 'Tendency of Zonal Component of Velocity'
                0793       diagUnits = 'm/s/day         '
5644fea420 Jean*0794       diagCode  = 'UUR     MR      '
27569449c9 Jean*0795       diagMate  = diagNum + 2
                0796       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0797      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
cd2de3b19d Andr*0798 
                0799       diagName  = 'TOTVTEND'
32bf9eae1e Jean*0800       diagTitle = 'Tendency of Meridional Component of Velocity'
                0801       diagUnits = 'm/s/day         '
5644fea420 Jean*0802       diagCode  = 'VVR     MR      '
27569449c9 Jean*0803       diagMate  = diagNum
                0804       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0805      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
cd2de3b19d Andr*0806 
f87dbbf625 Andr*0807       diagName  = 'TOTTTEND'
701e10a905 Mart*0808       diagTitle = DIAGS_MK_TITLE( 'Tendency of '
                0809      I                  //tTitl1//' Temperature', myThid )
32bf9eae1e Jean*0810       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/day', myThid )
5644fea420 Jean*0811       diagCode  = 'SMR     MR      '
27569449c9 Jean*0812       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0813      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
cd2de3b19d Andr*0814 
                0815       diagName  = 'TOTSTEND'
9ba8ef70b1 Jean*0816       diagTitle = DIAGS_MK_TITLE('Tendency of '//sTitle, myThid )
3ac63c2262 Jean*0817       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/day', myThid )
5644fea420 Jean*0818       diagCode  = 'SMR     MR      '
27569449c9 Jean*0819       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0820      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
cd2de3b19d Andr*0821 
bec73443e7 Jean*0822       diagName  = 'MoistCor'
                0823       diagTitle = 'Heating correction due to moist thermodynamics'
588c33df54 Jean*0824       diagUnits = 'W/m^2           '
bec73443e7 Jean*0825       diagCode  = 'SM      MR      '
                0826       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0827      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0828 
1a019c043f Jean*0829 #ifdef ALLOW_FRICTION_HEATING
                0830       diagName  = 'HeatDiss'
                0831       diagTitle = 'Heating from frictional dissipation'
                0832       diagUnits = 'W/m^2           '
                0833       diagCode  = 'SM      MR      '
                0834       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0835      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0836 #endif /* ALLOW_FRICTION_HEATING */
                0837 
8507f4494e Jean*0838 #ifdef ALLOW_GENERIC_ADVDIFF
                0839       diagName  = 'gT_Forc '
701e10a905 Mart*0840       diagTitle = DIAGS_MK_TITLE( tTitl1//' Temp. '
                0841      I                  //'forcing tendency', myThid )
8507f4494e Jean*0842       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/s', myThid )
                0843       diagCode  = 'SMR     MR      '
9d9b466d92 Jean*0844       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
8507f4494e Jean*0845      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
9d9b466d92 Jean*0846 
8507f4494e Jean*0847       diagName  = 'gS_Forc '
701e10a905 Mart*0848       diagTitle = DIAGS_MK_TITLE( sTitle
                0849      &                  //'forcing tendency', myThid )
8507f4494e Jean*0850       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/s', myThid )
                0851       diagCode  = 'SMR     MR      '
9d9b466d92 Jean*0852       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0853      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0854 
                0855       diagName  = 'AB_gT   '
701e10a905 Mart*0856       diagTitle = DIAGS_MK_TITLE( tTitl1//' Temp. '
                0857      I                  //'tendency from Adams-Bashforth', myThid )
9d9b466d92 Jean*0858       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/s', myThid )
                0859       diagCode  = 'SMR     MR      '
                0860       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0861      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0862 
                0863       diagName  = 'AB_gS   '
701e10a905 Mart*0864       diagTitle = DIAGS_MK_TITLE( sTitle
                0865      &                  //'tendency from Adams-Bashforth', myThid )
3ac63c2262 Jean*0866       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/s', myThid )
9d9b466d92 Jean*0867       diagCode  = 'SMR     MR      '
                0868       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0869      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
004ff77f8b Jean*0870 
                0871       diagName  = 'gTinAB  '
701e10a905 Mart*0872       diagTitle = DIAGS_MK_TITLE( tTitl1//' Temp. '
                0873      I                  //'tendency going in Adams-Bashforth', myThid )
004ff77f8b Jean*0874       diagUnits = DIAGS_MK_UNITS( tUnit4c//'/s', myThid )
                0875       diagCode  = 'SMR     MR      '
                0876       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0877      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0878 
                0879       diagName  = 'gSinAB  '
701e10a905 Mart*0880       diagTitle = DIAGS_MK_TITLE( sTitle
                0881      &                  //'tendency going in Adams-Bashforth', myThid )
3ac63c2262 Jean*0882       diagUnits = DIAGS_MK_UNITS( sUnit5c//'/s', myThid )
004ff77f8b Jean*0883       diagCode  = 'SMR     MR      '
                0884       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0885      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
9d9b466d92 Jean*0886 #endif /* ALLOW_GENERIC_ADVDIFF */
b5d5d65bd6 Andr*0887 
8507f4494e Jean*0888       diagName  = 'AB_gU   '
                0889       diagTitle = 'U momentum tendency from Adams-Bashforth'
                0890       diagUnits = 'm/s^2           '
                0891       diagCode  = 'UUR     MR      '
                0892       diagMate  = diagNum + 2
                0893       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0894      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0895 
                0896       diagName  = 'AB_gV   '
                0897       diagTitle = 'V momentum tendency from Adams-Bashforth'
                0898       diagUnits = 'm/s^2           '
                0899       diagCode  = 'VVR     MR      '
                0900       diagMate  = diagNum
                0901       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0902      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0903 
                0904 #ifdef ALLOW_NONHYDROSTATIC
                0905       diagName  = 'AB_gW   '
                0906       diagTitle = 'W momentum tendency from Adams-Bashforth'
                0907       diagUnits = DIAGS_MK_UNITS( rUnit2c//'/s^2', myThid )
                0908       diagCode  = 'WM      LR      '
                0909       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0910      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0911 #endif /* ALLOW_NONHYDROSTATIC */
                0912 
3f7db86ed7 Mich*0913 #ifdef ALLOW_EDDYPSI
                0914       diagName  = 'TAUXEDDY'
                0915       diagTitle = 'Zonal Eddy Stress'
1a019c043f Jean*0916       diagUnits = 'N/m^2           '
                0917       diagCode  = 'UU      LR      '
                0918       diagMate  = diagNum + 2
3f7db86ed7 Mich*0919       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0920      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8507f4494e Jean*0921 
3f7db86ed7 Mich*0922       diagName  = 'TAUYEDDY'
                0923       diagTitle = 'Meridional Eddy Stress'
1a019c043f Jean*0924       diagUnits = 'N/m^2           '
                0925       diagCode  = 'VV      LR      '
3f7db86ed7 Mich*0926       diagMate  = diagNum
                0927       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0928      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0929 
8507f4494e Jean*0930 # ifdef ALLOW_GMREDI
1a019c043f Jean*0931       diagName  = 'U_EulerM'
                0932       diagTitle = 'Zonal Eulerian-Mean Velocity (m/s)'
3f7db86ed7 Mich*0933       diagUnits = 'm/s             '
                0934       diagCode  = 'UUR     MR      '
                0935       diagMate  = diagNum + 2
                0936       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0937      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8507f4494e Jean*0938 
1a019c043f Jean*0939       diagName  = 'V_EulerM'
                0940       diagTitle = 'Meridional Eulerian-Mean Velocity (m/s)'
3f7db86ed7 Mich*0941       diagUnits = 'm/s             '
                0942       diagCode  = 'VVR     MR      '
                0943       diagMate  = diagNum
                0944       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0945      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
8507f4494e Jean*0946 # endif /* ALLOW_GMREDI */
3f7db86ed7 Mich*0947 #endif /* ALLOW_EDDYPSI */
0b5a2cab3a Mich*0948 
fb2c52cfd4 aver*0949 #ifdef INCLUDE_SOUNDSPEED_CALC_CODE
                0950       diagName  = 'CSound  '
                0951       diagTitle = 'Speed of sound in seawater'
                0952       diagUnits = 'm/s             '
                0953       diagCode  = 'SMR     MR      '
                0954       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0955      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0956 #endif
                0957 
a5ec81ed49 Timo*0958 #ifdef ALLOW_AUTODIFF
                0959 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
                0960 C     Adjoint state variables
                0961 
                0962       diagName  = 'ADJetan '
                0963       diagTitle = 'dJ/dEtaN: Sensitivity to sea surface height anomaly'
                0964       diagUnits = 'dJ/m            '
                0965       diagCode  = 'SM A    M1      '
                0966       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0967      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0968 
                0969       diagName  = 'ADJuvel '
                0970       diagTitle = 'dJ/dU: Sensitivity to zonal velocity'
                0971       diagUnits = 'dJ/(m/s)        '
                0972       diagCode  = 'UURA    MR      '
                0973       diagMate  = diagNum + 2
                0974       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0975      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0976 
                0977       diagName  = 'ADJvvel '
                0978       diagTitle = 'dJ/dV: Sensitivity to meridional velocity'
                0979       diagUnits = 'dJ/(m/s)        '
                0980       diagCode  = 'VVRA    MR      '
                0981       diagMate  = diagNum
                0982       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0983      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                0984 
                0985       diagName  = 'ADJwvel '
                0986       diagTitle = 'dJ/dW: Sensitivity to vertical velocity'
                0987       diagUnits = 'dJ/(m/s)        '
                0988       diagCode  = 'WM A    LR      '
                0989       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0990      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0991 
                0992       diagName  = 'ADJtheta'
                0993       diagTitle = 'dJ/dTheta: Sensitivity to potential temperature'
                0994       diagUnits = 'dJ/degC         '
                0995       diagCode  = 'SMRA    MR      '
                0996       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                0997      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                0998 
                0999       diagName  = 'ADJsalt '
                1000       diagTitle = 'dJ/dSalt: Sensitivity to salinity'
a10c595eb6 Timo*1001       diagUnits = DIAGS_MK_UNITS( 'dJ/('//sUnit5c//')', myThid )
a5ec81ed49 Timo*1002       diagCode  = 'SMRA    MR      '
                1003       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1004      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1005 
a10c595eb6 Timo*1006 C--   surface fluxes:
                1007       diagName  = 'ADJtaux '
                1008       diagTitle = 'dJ/dTaux: Senstivity to zonal surface wind stress'
                1009       diagUnits = 'dJ/(N/m^2)      '
                1010       diagCode  = 'UU A    U1      '
                1011       diagMate  = diagNum + 2
                1012       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1013      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                1014 
                1015       diagName  = 'ADJtauy '
                1016       diagTitle = 'dJ/dTauy: Sensitivity to merid. surface wind stress'
                1017       diagUnits = 'dJ/(N/m^2)      '
                1018       diagCode  = 'VV A    U1      '
                1019       diagMate  = diagNum
                1020       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1021      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                1022 
                1023       diagName  = 'ADJempmr'
1ac7194fda plla*1024       diagTitle = 'dJ/dEmPmR: Sensitivity to net surface freshwater'
                1025      &          //' flux'
a10c595eb6 Timo*1026       diagUnits = 'dJ/(kg/m^2/s)   '
                1027       diagCode  = 'SM A    U1      '
                1028       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1029      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1030 
                1031       diagName  = 'ADJqnet'
                1032       diagTitle = 'dJ/dQnet: Sensitivity to net surface heat flux'
                1033       diagUnits = 'dJ/(W/m^2)      '
                1034       diagCode  = 'SM A    U1      '
                1035       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1036      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1037 
                1038       diagName  = 'ADJqsw'
                1039       diagTitle = 'dJ/dQsw: Sensitivitiy to net Short-Wave radiation'
                1040       diagUnits = 'dJ/(W/m^2)      '
                1041       diagCode  = 'SM A    U1      '
                1042       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1043      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1044 
                1045       diagName  = 'ADJsst'
                1046       diagTitle = 'dJ/dSST: Sensitivity to Sea Surface Temperature'
                1047       diagUnits = 'dJ/K'
                1048       diagCode  = 'SM A    M1      '
                1049       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1050      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1051 
                1052       diagName  = 'ADJsss'
                1053       diagTitle = 'dJ/dSSS: Sensitivity to Sea Surface Salinity'
                1054       diagUnits = DIAGS_MK_UNITS( 'dJ/('//sUnit5c//')', myThid )
                1055       diagCode  = 'SM A    M1      '
                1056       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1057      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1058 
                1059 C Miscellanious
                1060       diagName  = 'ADJbtdrg'
                1061       diagTitle = 'dJ/dCd: Sensitivity to bottom drag coefficient'
                1062       diagUnits = 'dJ/d()'
                1063       diagCode  = 'SM A    M1      '
                1064       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1065      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1066 
                1067       diagName  = 'ADJdifkr'
                1068       diagTitle = 'dJ/dKr: Sensitivity to vertical diffusivity'
                1069       diagUnits = 'dJ/d(m^2/s))'
                1070       diagCode  = 'SMRA    MR      '
                1071       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1072      I          diagName, diagCode, diagUnits, diagTitle, 0, myThid )
                1073 
                1074       diagName  = 'ADJepsix'
                1075       diagTitle = 'dJ/dEddyPsiX: Sensitivity to zonal eddy'
                1076      &           //'streamfunction'
                1077       diagUnits = 'dJ/(m^2/s)      '
                1078       diagCode  = 'UURA    UR      '
                1079       diagMate  = diagNum + 2
                1080       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1081      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
                1082 
                1083       diagName  = 'ADJepsiy'
                1084       diagTitle = 'dJ/dEddyPsiY: Sensitivity to meridional eddy'
                1085      &           //'streamfunction'
                1086       diagUnits = 'dJ/(m^2/s)      '
                1087       diagCode  = 'VVRA    UR      '
                1088       diagMate  = diagNum
                1089       CALL DIAGNOSTICS_ADDTOLIST( diagNum,
                1090      I   diagName, diagCode, diagUnits, diagTitle, diagMate, myThid )
a5ec81ed49 Timo*1091 #endif /* ALLOW_AUTODIFF */
                1092 
09ceb40cd6 Jean*1093       RETURN
                1094       END