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

Ignore:
Timestamp:
2015-08-04T14:45:33+02:00 (9 years ago)
Author:
deazer
Message:

Upgraded branch to VERSION 3 6 STABLE

File:
1 edited

Legend:

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

    r5427 r5666  
    4848   USE iom 
    4949   USE ioipsl 
     50   USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities      
     51 
    5052#if defined key_lim2 
    5153   USE limwri_2  
     
    127129      !! 
    128130      INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
     131      INTEGER                      ::   jkbot                   ! 
    129132      REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    130133      !! 
     
    150153         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
    151154      ENDIF 
     155 
     156      CALL iom_put( "ssh" , sshn )                 ! sea surface height 
     157      if( iom_use('ssh2') )   CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height 
    152158       
    153159      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     
    156162         DO jj = 1, jpj 
    157163            DO ji = 1, jpi 
    158                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_tem) 
     164               jkbot = mbkt(ji,jj) 
     165               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
    159166            END DO 
    160167         END DO 
     
    167174         DO jj = 1, jpj 
    168175            DO ji = 1, jpi 
    169                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_sal) 
     176               jkbot = mbkt(ji,jj) 
     177               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
    170178            END DO 
    171179         END DO 
    172180         CALL iom_put( "sbs", z2d )                ! bottom salinity 
     181      ENDIF 
     182 
     183      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     184         z2d(:,:) = 0._wp 
     185         DO jj = 2, jpjm1 
     186            DO ji = fs_2, fs_jpim1   ! vector opt. 
     187               zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  & 
     188                      &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )       
     189               zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  & 
     190                      &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )  
     191               z2d(ji,jj) = rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)  
     192               ! 
     193            ENDDO 
     194         ENDDO 
     195         CALL lbc_lnk( z2d, 'T', 1. ) 
     196         CALL iom_put( "taubot", z2d )            
    173197      ENDIF 
    174198          
     
    178202         DO jj = 1, jpj 
    179203            DO ji = 1, jpi 
    180                z2d(ji,jj) = un(ji,jj,MAX(mbathy(ji,jj),1)) 
     204               jkbot = mbku(ji,jj) 
     205               z2d(ji,jj) = un(ji,jj,jkbot) 
    181206            END DO 
    182207         END DO 
    183208         CALL iom_put( "sbu", z2d )                ! bottom i-current 
    184209      ENDIF 
     210#if defined key_dynspg_ts 
     211      CALL iom_put(  "ubar", un_adv(:,:)      )    ! barotropic i-current 
     212#else 
     213      CALL iom_put(  "ubar", un_b(:,:)        )    ! barotropic i-current 
     214#endif 
    185215       
    186216      CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
     
    189219         DO jj = 1, jpj 
    190220            DO ji = 1, jpi 
    191                z2d(ji,jj) = vn(ji,jj,MAX(mbathy(ji,jj),1)) 
     221               jkbot = mbkv(ji,jj) 
     222               z2d(ji,jj) = vn(ji,jj,jkbot) 
    192223            END DO 
    193224         END DO 
    194225         CALL iom_put( "sbv", z2d )                ! bottom j-current 
     226      ENDIF 
     227#if defined key_dynspg_ts 
     228      CALL iom_put(  "vbar", vn_adv(:,:)      )    ! barotropic j-current 
     229#else 
     230      CALL iom_put(  "vbar", vn_b(:,:)        )    ! barotropic j-current 
     231#endif 
     232 
     233      CALL iom_put( "woce", wn )                   ! vertical velocity 
     234      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
     235         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
     236         z2d(:,:) = rau0 * e12t(:,:) 
     237         DO jk = 1, jpk 
     238            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) 
     239         END DO 
     240         CALL iom_put( "w_masstr" , z3d )   
     241         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) 
    195242      ENDIF 
    196243 
     
    604651         ENDIF 
    605652 
    606          IF( .NOT. lk_cpl ) THEN 
     653         IF( .NOT. ln_cpl ) THEN 
    607654            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    608655               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    613660         ENDIF 
    614661 
    615          IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     662         IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    616663            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    617664               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    636683#endif 
    637684 
    638          IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     685         IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    639686            CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
    640687               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    791838      ENDIF 
    792839 
    793       IF( .NOT. lk_cpl ) THEN 
     840      IF( .NOT. ln_cpl ) THEN 
    794841         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    795842         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    797844         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
    798845      ENDIF 
    799       IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     846      IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    800847         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    801848         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    813860#endif 
    814861 
    815       IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     862      IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    816863         CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
    817864         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.