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 7237 – NEMO

Changeset 7237


Ignore:
Timestamp:
2016-11-16T13:48:20+01:00 (7 years ago)
Author:
timgraham
Message:

Changes in diaar5.F90

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2016/dev_r7233_CMIP6_diags_trunk_version/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90

    r6665 r7237  
    2424   USE phycst         ! physical constant 
    2525   USE in_out_manager  ! I/O manager 
     26   USE zdfddm 
     27   USE zdf_oce 
    2628 
    2729   IMPLICIT NONE 
     
    4042   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   sn0          ! initial salinity 
    4143       
     44   !! * Substitutions 
     45#  include "zdfddm_substitute.h90" 
    4246   !!---------------------------------------------------------------------- 
    4347   !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
     
    8084      IF( nn_timing == 1 )   CALL timing_start('dia_ar5') 
    8185  
    82       CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres ) 
     86      CALL wrk_alloc( jpi , jpj              , zarea_ssh , zbotpres, pe ) 
    8387      CALL wrk_alloc( jpi , jpj , jpk        , zrhd      , zrhop    ) 
    8488      CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn                 ) 
     
    191195      CALL iom_put( 'temptot', ztemp ) 
    192196      CALL iom_put( 'saltot' , zsal  ) 
    193       ! 
    194       CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres ) 
     197 
     198      IF( iom_use( 'tnpeo' )) THEN     
     199      ! Work done against stratification by vertical mixing 
     200      ! Exclude points where rn2 is negative as convection kicks in here and 
     201      ! work is not being done against stratification 
     202          pe(:,:) = 0._wp 
     203          IF( lk_zdfddm ) THEN 
     204             DO ji=1,jpi 
     205                DO jj=1,jpj 
     206                   DO jk=1,jpk 
     207                      zrw =   ( gdepw_n(ji,jj,jk  ) - gdept_n(ji,jj,jk) )   & 
     208                         &  / ( gdept_n(ji,jj,jk-1) - gdept_n(ji,jj,jk) ) 
     209                      ! 
     210                      zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw 
     211                      zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw 
     212                      ! 
     213                      pe(ji, jj) = pe(ji, jj) - MIN(0._wp, rn2(ji,jj,jk)) * & 
     214                           &       grav * (avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) )  & 
     215                           &       - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) ) 
     216 
     217                   ENDDO 
     218                ENDDO 
     219             ENDDO 
     220          ELSE 
     221             DO ji=1,jpi 
     222                DO jj=1,jpj 
     223                   DO jk=1,jpk 
     224                       pe(ji,jj) = pe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk) 
     225                   ENDDO 
     226                ENDDO 
     227             ENDDO 
     228          ENDIF 
     229          CALL lbc_lnk(pe, 'T', 1._wp)          
     230          CALL iom_put( 'tnpeo', pe ) 
     231      ENDIF 
     232      ! 
     233      CALL wrk_dealloc( jpi , jpj              , zarea_ssh , zbotpres, pe ) 
    195234      CALL wrk_dealloc( jpi , jpj , jpk        , zrhd      , zrhop    ) 
    196235      CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn                 ) 
Note: See TracChangeset for help on using the changeset viewer.