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/diaar5.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/diaar5.F90

    r6664 r7494  
    2424   USE phycst         ! physical constant 
    2525   USE in_out_manager  ! I/O manager 
     26   USE zdfddm 
     27   USE zdf_oce 
    2628 
    2729   IMPLICIT NONE 
     
    4244   !! * Substitutions 
    4345#  include "domzgr_substitute.h90" 
     46#  include "zdfddm_substitute.h90" 
    4447   !!---------------------------------------------------------------------- 
    4548   !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
     
    7578      INTEGER  ::   ji, jj, jk                      ! dummy loop arguments 
    7679      REAL(wp) ::   zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass 
     80      REAL(wp) ::   zaw, zbw, zrw 
    7781      ! 
    7882      REAL(wp), POINTER, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace  
     83      REAL(wp), POINTER, DIMENSION(:,:)     :: pe                         ! 2D workspace  
    7984      REAL(wp), POINTER, DIMENSION(:,:,:)   :: zrhd , zrhop               ! 3D workspace 
    8085      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace 
     
    8287      IF( nn_timing == 1 )   CALL timing_start('dia_ar5') 
    8388  
    84       CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres ) 
     89      CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres, pe ) 
    8590      CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    ) 
    8691      CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 ) 
     
    95100      CALL iom_put( 'voltot', zvol               ) 
    96101      CALL iom_put( 'sshtot', zvolssh / area_tot ) 
     102      CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) ) 
    97103 
    98104      !                      
     
    190196      CALL iom_put( 'temptot', ztemp ) 
    191197      CALL iom_put( 'saltot' , zsal  ) 
    192       ! 
    193       CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres ) 
     198 
     199      IF( iom_use( 'tnpeo' )) THEN     
     200      ! Work done against stratification by vertical mixing 
     201      ! Exclude points where rn2 is negative as convection kicks in here and 
     202      ! work is not being done against stratification 
     203          pe(:,:) = 0._wp 
     204          IF( lk_zdfddm ) THEN 
     205             DO ji=1,jpi 
     206                DO jj=1,jpj 
     207                   DO jk=1,jpk 
     208                      zrw =   ( fsdepw(ji,jj,jk  ) - fsdept(ji,jj,jk) )   & 
     209                         &  / ( fsdept(ji,jj,jk-1) - fsdept(ji,jj,jk) ) 
     210                      ! 
     211                      zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw 
     212                      zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw 
     213                      ! 
     214                      pe(ji, jj) = pe(ji, jj) - MIN(0._wp, rn2(ji,jj,jk)) * & 
     215                           &       grav * (avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  & 
     216                           &       - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) ) 
     217 
     218                   ENDDO 
     219                ENDDO 
     220             ENDDO 
     221          ELSE 
     222             DO ji=1,jpi 
     223                DO jj=1,jpj 
     224                   DO jk=1,jpk 
     225                       pe(ji,jj) = pe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk) 
     226                   ENDDO 
     227                ENDDO 
     228             ENDDO 
     229          ENDIF 
     230          CALL lbc_lnk(pe, 'T', 1._wp)          
     231          CALL iom_put( 'tnpeo', pe ) 
     232      ENDIF 
     233      ! 
     234      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres, pe ) 
    194235      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    ) 
    195236      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 ) 
Note: See TracChangeset for help on using the changeset viewer.