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 7923 for branches/UKMO/dev_r5518_GO6_package_XIOS_read/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2017-04-18T15:26:56+02:00 (7 years ago)
Author:
andmirek
Message:

merge changes up to 7573

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/dev_r5518_GO6_package_XIOS_read/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r6498 r7923  
    157157      IF( iom_use("e3tdef") )   & 
    158158         CALL iom_put( "e3tdef"  , ( ( fse3t_n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) 
     159      CALL iom_put("tpt_dep", fsdept_n(:,:,:) ) 
    159160 
    160161 
     
    323324      CALL iom_put( "hdiv", hdivn )                  ! Horizontal divergence 
    324325      ! 
    325       IF( iom_use("u_masstr") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 
     326      IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 
    326327         z3d(:,:,jpk) = 0.e0 
     328         z2d(:,:) = 0.e0 
    327329         DO jk = 1, jpkm1 
    328330            z3d(:,:,jk) = rau0 * un(:,:,jk) * e2u(:,:) * fse3u(:,:,jk) * umask(:,:,jk) 
     331            z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 
    329332         END DO 
    330333         CALL iom_put( "u_masstr", z3d )                  ! mass transport in i-direction 
     334         CALL iom_put( "u_masstr_vint", z2d )             ! mass transport in i-direction vertical sum 
    331335      ENDIF 
    332336       
     
    391395         CALL iom_put( "v_salttr", 0.5 * z2d )            !  heat transport in j-direction 
    392396      ENDIF 
     397 
     398      ! Vertical integral of temperature 
     399      IF( iom_use("tosmint") ) THEN 
     400         z2d(:,:)=0._wp 
     401         DO jk = 1, jpkm1 
     402            DO jj = 2, jpjm1 
     403               DO ji = fs_2, fs_jpim1   ! vector opt. 
     404                  z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) *  tsn(ji,jj,jk,jp_tem) 
     405               END DO 
     406            END DO 
     407         END DO 
     408         CALL lbc_lnk( z2d, 'T', -1. ) 
     409         CALL iom_put( "tosmint", z2d )  
     410      ENDIF 
     411 
     412      ! Vertical integral of salinity 
     413      IF( iom_use("somint") ) THEN 
     414         z2d(:,:)=0._wp 
     415         DO jk = 1, jpkm1 
     416            DO jj = 2, jpjm1 
     417               DO ji = fs_2, fs_jpim1   ! vector opt. 
     418                  z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 
     419               END DO 
     420            END DO 
     421         END DO 
     422         CALL lbc_lnk( z2d, 'T', -1. ) 
     423         CALL iom_put( "somint", z2d )  
     424      ENDIF 
     425 
     426      CALL iom_put( "bn2", rn2 )  !Brunt-Vaisala buoyancy frequency (N^2) 
    393427      ! 
    394428      CALL wrk_dealloc( jpi , jpj      , z2d ) 
Note: See TracChangeset for help on using the changeset viewer.