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

Ignore:
Timestamp:
2015-07-02T15:49:24+02:00 (9 years ago)
Author:
davestorkey
Message:

Update UKMO/dev_r5021_nn_etau_revision branch to rev 5518 of trunk
(=branching point for NEMO 3.6_stable).

File:
1 edited

Legend:

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

    r5446 r5532  
    4747   USE iom 
    4848   USE ioipsl 
     49   USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities      
     50 
    4951#if defined key_lim2 
    5052   USE limwri_2  
     
    126128      !! 
    127129      INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
     130      INTEGER                      ::   jkbot                   ! 
    128131      REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    129132      !! 
     
    149152         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
    150153      ENDIF 
     154 
     155      CALL iom_put( "ssh" , sshn )                 ! sea surface height 
     156      if( iom_use('ssh2') )   CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height 
    151157       
    152158      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     
    155161         DO jj = 1, jpj 
    156162            DO ji = 1, jpi 
    157                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_tem) 
     163               jkbot = mbkt(ji,jj) 
     164               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
    158165            END DO 
    159166         END DO 
     
    166173         DO jj = 1, jpj 
    167174            DO ji = 1, jpi 
    168                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_sal) 
     175               jkbot = mbkt(ji,jj) 
     176               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
    169177            END DO 
    170178         END DO 
    171179         CALL iom_put( "sbs", z2d )                ! bottom salinity 
     180      ENDIF 
     181 
     182      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     183         z2d(:,:) = 0._wp 
     184         DO jj = 2, jpjm1 
     185            DO ji = fs_2, fs_jpim1   ! vector opt. 
     186               zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  & 
     187                      &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )       
     188               zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  & 
     189                      &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )  
     190               z2d(ji,jj) = rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)  
     191               ! 
     192            ENDDO 
     193         ENDDO 
     194         CALL lbc_lnk( z2d, 'T', 1. ) 
     195         CALL iom_put( "taubot", z2d )            
    172196      ENDIF 
    173197          
     
    177201         DO jj = 1, jpj 
    178202            DO ji = 1, jpi 
    179                z2d(ji,jj) = un(ji,jj,MAX(mbathy(ji,jj),1)) 
     203               jkbot = mbku(ji,jj) 
     204               z2d(ji,jj) = un(ji,jj,jkbot) 
    180205            END DO 
    181206         END DO 
    182207         CALL iom_put( "sbu", z2d )                ! bottom i-current 
    183208      ENDIF 
     209#if defined key_dynspg_ts 
     210      CALL iom_put(  "ubar", un_adv(:,:)      )    ! barotropic i-current 
     211#else 
     212      CALL iom_put(  "ubar", un_b(:,:)        )    ! barotropic i-current 
     213#endif 
    184214       
    185215      CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
     
    188218         DO jj = 1, jpj 
    189219            DO ji = 1, jpi 
    190                z2d(ji,jj) = vn(ji,jj,MAX(mbathy(ji,jj),1)) 
     220               jkbot = mbkv(ji,jj) 
     221               z2d(ji,jj) = vn(ji,jj,jkbot) 
    191222            END DO 
    192223         END DO 
    193224         CALL iom_put( "sbv", z2d )                ! bottom j-current 
     225      ENDIF 
     226#if defined key_dynspg_ts 
     227      CALL iom_put(  "vbar", vn_adv(:,:)      )    ! barotropic j-current 
     228#else 
     229      CALL iom_put(  "vbar", vn_b(:,:)        )    ! barotropic j-current 
     230#endif 
     231 
     232      CALL iom_put( "woce", wn )                   ! vertical velocity 
     233      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
     234         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
     235         z2d(:,:) = rau0 * e12t(:,:) 
     236         DO jk = 1, jpk 
     237            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) 
     238         END DO 
     239         CALL iom_put( "w_masstr" , z3d )   
     240         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) 
    194241      ENDIF 
    195242 
Note: See TracChangeset for help on using the changeset viewer.