New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 7494 for branches/2015/nemo_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2016-12-14T10:02:43+01:00 (7 years ago)
Author:
timgraham
Message:

Addition of extra diagnostic outputs in CMIP6 data request. NB. Will require update to field_def file from shaconemo

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2015/nemo_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r6348 r7494  
    156156      IF( iom_use("e3tdef") )   & 
    157157         CALL iom_put( "e3tdef"  , ( ( fse3t_n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) 
     158      CALL iom_put("tpt_dep", fsdept_n(:,:,:) ) 
     159 
    158160 
    159161 
     
    318320      CALL iom_put( "hdiv", hdivn )                  ! Horizontal divergence 
    319321      ! 
    320       IF( iom_use("u_masstr") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 
     322      IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 
    321323         z3d(:,:,jpk) = 0.e0 
     324         z2d(:,:) = 0.e0 
    322325         DO jk = 1, jpkm1 
    323326            z3d(:,:,jk) = rau0 * un(:,:,jk) * e2u(:,:) * fse3u(:,:,jk) * umask(:,:,jk) 
     327            z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 
    324328         END DO 
    325329         CALL iom_put( "u_masstr", z3d )                  ! mass transport in i-direction 
     330         CALL iom_put( "u_masstr_vint", z2d )             ! mass transport in i-direction vertical sum 
    326331      ENDIF 
    327332       
     
    386391         CALL iom_put( "v_salttr", 0.5 * z2d )            !  heat transport in j-direction 
    387392      ENDIF 
     393 
     394      ! Vertical integral of temperature 
     395      IF( iom_use("tosmint") ) THEN 
     396         z2d(:,:)=0._wp 
     397         DO jk = 1, jpkm1 
     398            DO jj = 2, jpjm1 
     399               DO ji = fs_2, fs_jpim1   ! vector opt. 
     400                  z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) *  tsn(ji,jj,jk,jp_tem) 
     401               END DO 
     402            END DO 
     403         END DO 
     404         CALL lbc_lnk( z2d, 'T', -1. ) 
     405         CALL iom_put( "tosmint", z2d )  
     406      ENDIF 
     407 
     408      ! Vertical integral of salinity 
     409      IF( iom_use("somint") ) THEN 
     410         z2d(:,:)=0._wp 
     411         DO jk = 1, jpkm1 
     412            DO jj = 2, jpjm1 
     413               DO ji = fs_2, fs_jpim1   ! vector opt. 
     414                  z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 
     415               END DO 
     416            END DO 
     417         END DO 
     418         CALL lbc_lnk( z2d, 'T', -1. ) 
     419         CALL iom_put( "somint", z2d )  
     420      ENDIF 
     421 
     422      CALL iom_put( "bn2", rn2 )  !Brunt-Vaisala buoyancy frequency (N^2) 
    388423      ! 
    389424      CALL wrk_dealloc( jpi , jpj      , z2d ) 
Note: See TracChangeset for help on using the changeset viewer.