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 5980 for branches/2014/dev_r4650_UKMO7_STARTHOUR/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2015-12-02T16:20:47+01:00 (8 years ago)
Author:
timgraham
Message:

Upgraded to v3.6 revision of trunk (r5518)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2014/dev_r4650_UKMO7_STARTHOUR/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r5075 r5980  
    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  
     
    7880   !!---------------------------------------------------------------------- 
    7981   !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
    80    !! $Id $ 
     82   !! $Id$ 
    8183   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
    8284   !!---------------------------------------------------------------------- 
     
    125127      !! 
    126128      INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
     129      INTEGER                      ::   jkbot                   ! 
    127130      REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    128131      !! 
     
    142145      ENDIF 
    143146 
    144       IF( lk_vvl ) THEN 
    145          z3d(:,:,:) = tsn(:,:,:,jp_tem) * fse3t_n(:,:,:) 
    146          CALL iom_put( "toce" , z3d                        )   ! heat content 
     147      IF( .NOT.lk_vvl ) THEN 
     148         CALL iom_put( "e3t" , fse3t_n(:,:,:) ) 
     149         CALL iom_put( "e3u" , fse3u_n(:,:,:) ) 
     150         CALL iom_put( "e3v" , fse3v_n(:,:,:) ) 
     151         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
     152      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 
     156       
     157      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     158      CALL iom_put(  "sst", tsn(:,:,1,jp_tem) )    ! surface temperature 
     159      IF ( iom_use("sbt") ) THEN 
    147160         DO jj = 1, jpj 
    148161            DO ji = 1, jpi 
    149                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_tem) * fse3t_n(ji,jj,mikt(ji,jj)) 
    150             END DO 
    151          END DO   
    152          CALL iom_put( "sst"  , z2d(:,:)                 )   ! sea surface heat content       
     162               jkbot = mbkt(ji,jj) 
     163               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
     164            END DO 
     165         END DO 
     166         CALL iom_put( "sbt", z2d )                ! bottom temperature 
     167      ENDIF 
     168       
     169      CALL iom_put( "soce", tsn(:,:,:,jp_sal) )    ! 3D salinity 
     170      CALL iom_put(  "sss", tsn(:,:,1,jp_sal) )    ! surface salinity 
     171      IF ( iom_use("sbs") ) THEN 
    153172         DO jj = 1, jpj 
    154173            DO ji = 1, jpi 
    155                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_tem)**2 * fse3t_n(ji,jj,mikt(ji,jj)) 
    156             END DO 
    157          END DO   
    158          CALL iom_put( "sst2" , z2d(:,:)      )   ! sea surface content of squared temperature 
    159          z3d(:,:,:) = tsn(:,:,:,jp_sal) * fse3t_n(:,:,:)             
    160          CALL iom_put( "soce" , z3d                        )   ! salinity content 
     174               jkbot = mbkt(ji,jj) 
     175               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
     176            END DO 
     177         END DO 
     178         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 )            
     195      ENDIF 
     196          
     197      CALL iom_put( "uoce", un(:,:,:)         )    ! 3D i-current 
     198      CALL iom_put(  "ssu", un(:,:,1)         )    ! surface i-current 
     199      IF ( iom_use("sbu") ) THEN 
    161200         DO jj = 1, jpj 
    162201            DO ji = 1, jpi 
    163                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_sal) * fse3t_n(ji,jj,mikt(ji,jj)) 
    164             END DO 
    165          END DO   
    166          CALL iom_put( "sss"  , z2d(:,:)                 )   ! sea surface salinity content 
     202               jkbot = mbku(ji,jj) 
     203               z2d(ji,jj) = un(ji,jj,jkbot) 
     204            END DO 
     205         END DO 
     206         CALL iom_put( "sbu", z2d )                ! bottom i-current 
     207      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 
     213       
     214      CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
     215      CALL iom_put(  "ssv", vn(:,:,1)         )    ! surface j-current 
     216      IF ( iom_use("sbv") ) THEN 
    167217         DO jj = 1, jpj 
    168218            DO ji = 1, jpi 
    169                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_sal)**2 * fse3t_n(ji,jj,mikt(ji,jj)) 
    170             END DO 
    171          END DO   
    172          CALL iom_put( "sss2" , z2d(:,:)                 )   ! sea surface content of squared salinity 
    173       ELSE 
    174          CALL iom_put( "toce" , tsn(:,:,:,jp_tem)        )   ! temperature 
    175          IF ( iom_use("sst") ) THEN 
    176             DO jj = 1, jpj 
    177                DO ji = 1, jpi 
    178                   z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_tem) 
    179                END DO 
    180             END DO 
    181             CALL iom_put( "sst"  , z2d(:,:)            ) ! sea surface temperature 
    182          ENDIF 
    183          IF ( iom_use("sst2") )   CALL iom_put( "sst2" , z2d(:,:) * z2d(:,:) ) ! square of sea surface temperature 
    184          CALL iom_put( "soce" , tsn(:,:,:,jp_sal)          )   ! salinity 
    185          IF ( iom_use("sss") ) THEN 
    186             DO jj = 1, jpj 
    187                DO ji = 1, jpi 
    188                   z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_sal) 
    189                END DO 
    190             END DO 
    191             CALL iom_put( "sss"  , z2d(:,:)            ) ! sea surface salinity 
    192          ENDIF 
    193          CALL iom_put( "sss2" , z2d(:,:) * z2d(:,:) ) ! square of sea surface salinity 
    194       END IF 
    195       IF( lk_vvl .AND. (.NOT. ln_dynadv_vec) ) THEN 
    196          CALL iom_put( "uoce" , umask(:,:,:) * un(:,:,:) * fse3u_n(:,:,:) )    ! i-transport 
    197          CALL iom_put( "voce" , vmask(:,:,:) * vn(:,:,:) * fse3v_n(:,:,:) )    ! j-transport 
    198       ELSE 
    199          CALL iom_put( "uoce" , umask(:,:,:) * un(:,:,:)                  )    ! i-current 
    200          CALL iom_put( "voce" , vmask(:,:,:) * vn(:,:,:)                  )    ! j-current 
    201          IF ( iom_use("ssu") ) THEN 
    202             DO jj = 1, jpj 
    203                DO ji = 1, jpi 
    204                   z2d(ji,jj) = un(ji,jj,miku(ji,jj)) 
    205                END DO 
    206             END DO 
    207             CALL iom_put( "ssu"   , z2d                                    )    ! i-current 
    208          ENDIF 
    209          IF ( iom_use("ssv") ) THEN 
    210             DO jj = 1, jpj 
    211                DO ji = 1, jpi 
    212                   z2d(ji,jj) = vn(ji,jj,mikv(ji,jj)) 
    213                END DO 
    214             END DO 
    215             CALL iom_put( "ssv"   , z2d                                    )    ! j-current 
    216          ENDIF 
    217       ENDIF 
    218       CALL iom_put(    "avt"  , avt                        )    ! T vert. eddy diff. coef. 
    219       CALL iom_put(    "avm"  , avmu                       )    ! T vert. eddy visc. coef. 
    220       IF( lk_zdfddm ) THEN 
    221          CALL iom_put( "avs" , fsavs(:,:,:)                          )    ! S vert. eddy diff. coef. 
    222       ENDIF 
    223  
    224       IF ( iom_use("sstgrad2") .OR. iom_use("sstgrad2") ) THEN 
     219               jkbot = mbkv(ji,jj) 
     220               z2d(ji,jj) = vn(ji,jj,jkbot) 
     221            END DO 
     222         END DO 
     223         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(:,:,:) ) 
     240      ENDIF 
     241 
     242      CALL iom_put( "avt" , avt                        )    ! T vert. eddy diff. coef. 
     243      CALL iom_put( "avm" , avmu                       )    ! T vert. eddy visc. coef. 
     244      CALL iom_put( "avs" , fsavs(:,:,:)               )    ! S vert. eddy diff. coef. (useful only with key_zdfddm) 
     245 
     246      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 
    225247         DO jj = 2, jpjm1                                    ! sst gradient 
    226248            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    234256         CALL lbc_lnk( z2d, 'T', 1. ) 
    235257         CALL iom_put( "sstgrad2",  z2d               )    ! square of module of sst gradient 
    236          !CDIR NOVERRCHK< 
    237258         z2d(:,:) = SQRT( z2d(:,:) ) 
    238259         CALL iom_put( "sstgrad" ,  z2d               )    ! module of sst gradient 
     
    243264         z2d(:,:)  = 0._wp  
    244265         DO jk = 1, jpkm1 
    245             DO jj = 2, jpjm1 
    246                DO ji = fs_2, fs_jpim1   ! vector opt. 
     266            DO jj = 1, jpj 
     267               DO ji = 1, jpi 
    247268                  z2d(ji,jj) = z2d(ji,jj) + fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
    248269               END DO 
    249270            END DO 
    250271         END DO 
    251          CALL lbc_lnk( z2d, 'T', 1. ) 
    252272         CALL iom_put( "heatc", (rau0 * rcp) * z2d )    ! vertically integrated heat content (J/m2) 
    253273      ENDIF 
     
    256276         z2d(:,:)  = 0._wp  
    257277         DO jk = 1, jpkm1 
    258             DO jj = 2, jpjm1 
    259                DO ji = fs_2, fs_jpim1   ! vector opt. 
     278            DO jj = 1, jpj 
     279               DO ji = 1, jpi 
    260280                  z2d(ji,jj) = z2d(ji,jj) + fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
    261281               END DO 
    262282            END DO 
    263283         END DO 
    264          CALL lbc_lnk( z2d, 'T', 1. ) 
    265284         CALL iom_put( "saltc", rau0 * z2d )   ! vertically integrated salt content (PSU*kg/m2) 
    266285      ENDIF 
     
    621640         ENDIF 
    622641 
    623          IF( .NOT. lk_cpl ) THEN 
     642         IF( .NOT. ln_cpl ) THEN 
    624643            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    625644               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    630649         ENDIF 
    631650 
    632          IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     651         IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    633652            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    634653               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    653672#endif 
    654673 
    655          IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     674         IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    656675            CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
    657676               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    808827      ENDIF 
    809828 
    810       IF( .NOT. lk_cpl ) THEN 
     829      IF( .NOT. ln_cpl ) THEN 
    811830         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    812831         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    814833         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
    815834      ENDIF 
    816       IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     835      IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    817836         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    818837         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    830849#endif 
    831850 
    832       IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     851      IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    833852         CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
    834853         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.