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

Ignore:
Timestamp:
2015-07-08T17:13:59+02:00 (9 years ago)
Author:
davestorkey
Message:

Upgrade UKMO/dev_r5107_hadgem3_mct branch to trunk revision 5518
( = branching point for NEMO 3.6_stable).

File:
1 edited

Legend:

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

    r5279 r5568  
    4646   USE iom 
    4747   USE ioipsl 
     48   USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities      
     49 
    4850#if defined key_lim2 
    4951   USE limwri_2  
     
    125127      !! 
    126128      INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
     129      INTEGER                      ::   jkbot                   ! 
    127130      REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    128131      !! 
     
    148151         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
    149152      ENDIF 
     153 
     154      CALL iom_put( "ssh" , sshn )                 ! sea surface height 
     155      if( iom_use('ssh2') )   CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height 
    150156       
    151157      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     
    154160         DO jj = 1, jpj 
    155161            DO ji = 1, jpi 
    156                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_tem) 
     162               jkbot = mbkt(ji,jj) 
     163               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
    157164            END DO 
    158165         END DO 
     
    165172         DO jj = 1, jpj 
    166173            DO ji = 1, jpi 
    167                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_sal) 
     174               jkbot = mbkt(ji,jj) 
     175               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
    168176            END DO 
    169177         END DO 
    170178         CALL iom_put( "sbs", z2d )                ! bottom salinity 
     179      ENDIF 
     180 
     181      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     182         z2d(:,:) = 0._wp 
     183         DO jj = 2, jpjm1 
     184            DO ji = fs_2, fs_jpim1   ! vector opt. 
     185               zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  & 
     186                      &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )       
     187               zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  & 
     188                      &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )  
     189               z2d(ji,jj) = rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)  
     190               ! 
     191            ENDDO 
     192         ENDDO 
     193         CALL lbc_lnk( z2d, 'T', 1. ) 
     194         CALL iom_put( "taubot", z2d )            
    171195      ENDIF 
    172196          
     
    176200         DO jj = 1, jpj 
    177201            DO ji = 1, jpi 
    178                z2d(ji,jj) = un(ji,jj,MAX(mbathy(ji,jj),1)) 
     202               jkbot = mbku(ji,jj) 
     203               z2d(ji,jj) = un(ji,jj,jkbot) 
    179204            END DO 
    180205         END DO 
    181206         CALL iom_put( "sbu", z2d )                ! bottom i-current 
    182207      ENDIF 
     208#if defined key_dynspg_ts 
     209      CALL iom_put(  "ubar", un_adv(:,:)      )    ! barotropic i-current 
     210#else 
     211      CALL iom_put(  "ubar", un_b(:,:)        )    ! barotropic i-current 
     212#endif 
    183213       
    184214      CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
     
    187217         DO jj = 1, jpj 
    188218            DO ji = 1, jpi 
    189                z2d(ji,jj) = vn(ji,jj,MAX(mbathy(ji,jj),1)) 
     219               jkbot = mbkv(ji,jj) 
     220               z2d(ji,jj) = vn(ji,jj,jkbot) 
    190221            END DO 
    191222         END DO 
    192223         CALL iom_put( "sbv", z2d )                ! bottom j-current 
     224      ENDIF 
     225#if defined key_dynspg_ts 
     226      CALL iom_put(  "vbar", vn_adv(:,:)      )    ! barotropic j-current 
     227#else 
     228      CALL iom_put(  "vbar", vn_b(:,:)        )    ! barotropic j-current 
     229#endif 
     230 
     231      CALL iom_put( "woce", wn )                   ! vertical velocity 
     232      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
     233         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
     234         z2d(:,:) = rau0 * e12t(:,:) 
     235         DO jk = 1, jpk 
     236            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) 
     237         END DO 
     238         CALL iom_put( "w_masstr" , z3d )   
     239         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) 
    193240      ENDIF 
    194241 
     
    593640         ENDIF 
    594641 
    595          IF( .NOT. lk_cpl ) THEN 
     642         IF( .NOT. ln_cpl ) THEN 
    596643            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    597644               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    602649         ENDIF 
    603650 
    604          IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     651         IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    605652            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    606653               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    625672#endif 
    626673 
    627          IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     674         IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    628675            CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
    629676               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    780827      ENDIF 
    781828 
    782       IF( .NOT. lk_cpl ) THEN 
     829      IF( .NOT. ln_cpl ) THEN 
    783830         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    784831         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    786833         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
    787834      ENDIF 
    788       IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     835      IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    789836         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    790837         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    802849#endif 
    803850 
    804       IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     851      IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    805852         CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
    806853         CALL histwrite( nid_T, "soicealb", it, alb_ice(:,:,1), ndim_hT, ndex_hT )   ! ice albedo 
Note: See TracChangeset for help on using the changeset viewer.