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

Ignore:
Timestamp:
2015-04-29T12:17:12+02:00 (9 years ago)
Author:
davestorkey
Message:

Update UKMO nn_etau_revision branch with trunk changes to rev 5107.

File:
1 edited

Legend:

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

    r5239 r5240  
    143143      ENDIF 
    144144 
    145       IF( lk_vvl ) THEN 
    146          z3d(:,:,:) = tsn(:,:,:,jp_tem) * fse3t_n(:,:,:) 
    147          CALL iom_put( "toce" , z3d                        )   ! heat content 
     145      IF( .NOT.lk_vvl ) THEN 
     146         CALL iom_put( "e3t" , fse3t_n(:,:,:) ) 
     147         CALL iom_put( "e3u" , fse3u_n(:,:,:) ) 
     148         CALL iom_put( "e3v" , fse3v_n(:,:,:) ) 
     149         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
     150      ENDIF 
     151       
     152      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     153      CALL iom_put(  "sst", tsn(:,:,1,jp_tem) )    ! surface temperature 
     154      IF ( iom_use("sbt") ) THEN 
    148155         DO jj = 1, jpj 
    149156            DO ji = 1, jpi 
    150                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_tem) * fse3t_n(ji,jj,mikt(ji,jj)) 
    151             END DO 
    152          END DO   
    153          CALL iom_put( "sst"  , z2d(:,:)                 )   ! sea surface heat content       
     157               z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_tem) 
     158            END DO 
     159         END DO 
     160         CALL iom_put( "sbt", z2d )                ! bottom temperature 
     161      ENDIF 
     162       
     163      CALL iom_put( "soce", tsn(:,:,:,jp_sal) )    ! 3D salinity 
     164      CALL iom_put(  "sss", tsn(:,:,1,jp_sal) )    ! surface salinity 
     165      IF ( iom_use("sbs") ) THEN 
    154166         DO jj = 1, jpj 
    155167            DO ji = 1, jpi 
    156                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_tem)**2 * fse3t_n(ji,jj,mikt(ji,jj)) 
    157             END DO 
    158          END DO   
    159          CALL iom_put( "sst2" , z2d(:,:)      )   ! sea surface content of squared temperature 
    160          z3d(:,:,:) = tsn(:,:,:,jp_sal) * fse3t_n(:,:,:)             
    161          CALL iom_put( "soce" , z3d                        )   ! salinity content 
     168               z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_sal) 
     169            END DO 
     170         END DO 
     171         CALL iom_put( "sbs", z2d )                ! bottom salinity 
     172      ENDIF 
     173          
     174      CALL iom_put( "uoce", un(:,:,:)         )    ! 3D i-current 
     175      CALL iom_put(  "ssu", un(:,:,1)         )    ! surface i-current 
     176      IF ( iom_use("sbu") ) THEN 
    162177         DO jj = 1, jpj 
    163178            DO ji = 1, jpi 
    164                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_sal) * fse3t_n(ji,jj,mikt(ji,jj)) 
    165             END DO 
    166          END DO   
    167          CALL iom_put( "sss"  , z2d(:,:)                 )   ! sea surface salinity content 
     179               z2d(ji,jj) = un(ji,jj,MAX(mbathy(ji,jj),1)) 
     180            END DO 
     181         END DO 
     182         CALL iom_put( "sbu", z2d )                ! bottom i-current 
     183      ENDIF 
     184       
     185      CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
     186      CALL iom_put(  "ssv", vn(:,:,1)         )    ! surface j-current 
     187      IF ( iom_use("sbv") ) THEN 
    168188         DO jj = 1, jpj 
    169189            DO ji = 1, jpi 
    170                z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_sal)**2 * fse3t_n(ji,jj,mikt(ji,jj)) 
    171             END DO 
    172          END DO   
    173          CALL iom_put( "sss2" , z2d(:,:)                 )   ! sea surface content of squared salinity 
    174       ELSE 
    175          CALL iom_put( "toce" , tsn(:,:,:,jp_tem)        )   ! temperature 
    176          IF ( iom_use("sst") ) THEN 
    177             DO jj = 1, jpj 
    178                DO ji = 1, jpi 
    179                   z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_tem) 
    180                END DO 
    181             END DO 
    182             CALL iom_put( "sst"  , z2d(:,:)            ) ! sea surface temperature 
    183          ENDIF 
    184          IF ( iom_use("sst2") )   CALL iom_put( "sst2" , z2d(:,:) * z2d(:,:) ) ! square of sea surface temperature 
    185          CALL iom_put( "soce" , tsn(:,:,:,jp_sal)          )   ! salinity 
    186          IF ( iom_use("sss") ) THEN 
    187             DO jj = 1, jpj 
    188                DO ji = 1, jpi 
    189                   z2d(ji,jj) = tsn(ji,jj,mikt(ji,jj),jp_sal) 
    190                END DO 
    191             END DO 
    192             CALL iom_put( "sss"  , z2d(:,:)            ) ! sea surface salinity 
    193          ENDIF 
    194          CALL iom_put( "sss2" , z2d(:,:) * z2d(:,:) ) ! square of sea surface salinity 
    195       END IF 
    196       IF( lk_vvl .AND. (.NOT. ln_dynadv_vec) ) THEN 
    197          CALL iom_put( "uoce" , umask(:,:,:) * un(:,:,:) * fse3u_n(:,:,:) )    ! i-transport 
    198          CALL iom_put( "voce" , vmask(:,:,:) * vn(:,:,:) * fse3v_n(:,:,:) )    ! j-transport 
    199       ELSE 
    200          CALL iom_put( "uoce" , umask(:,:,:) * un(:,:,:)                  )    ! i-current 
    201          CALL iom_put( "voce" , vmask(:,:,:) * vn(:,:,:)                  )    ! j-current 
    202          IF ( iom_use("ssu") ) THEN 
    203             DO jj = 1, jpj 
    204                DO ji = 1, jpi 
    205                   z2d(ji,jj) = un(ji,jj,miku(ji,jj)) 
    206                END DO 
    207             END DO 
    208             CALL iom_put( "ssu"   , z2d                                    )    ! i-current 
    209          ENDIF 
    210          IF ( iom_use("ssv") ) THEN 
    211             DO jj = 1, jpj 
    212                DO ji = 1, jpi 
    213                   z2d(ji,jj) = vn(ji,jj,mikv(ji,jj)) 
    214                END DO 
    215             END DO 
    216             CALL iom_put( "ssv"   , z2d                                    )    ! j-current 
    217          ENDIF 
    218       ENDIF 
    219       CALL iom_put(    "avt"  , avt                        )    ! T vert. eddy diff. coef. 
    220       CALL iom_put(    "avm"  , avmu                       )    ! T vert. eddy visc. coef. 
    221       IF( lk_zdftke ) THEN   
    222          CALL iom_put( "tke"      , en                               )    ! TKE budget: Turbulent Kinetic Energy   
    223          CALL iom_put( "tke_niw"  , e_niw                            )    ! TKE budget: Near-inertial waves   
    224       ENDIF 
    225       IF( lk_zdfddm ) THEN 
    226          CALL iom_put( "avs" , fsavs(:,:,:)                          )    ! S vert. eddy diff. coef. 
    227       ENDIF 
    228  
    229       IF ( iom_use("sstgrad2") .OR. iom_use("sstgrad2") ) THEN 
     190               z2d(ji,jj) = vn(ji,jj,MAX(mbathy(ji,jj),1)) 
     191            END DO 
     192         END DO 
     193         CALL iom_put( "sbv", z2d )                ! bottom j-current 
     194      ENDIF 
     195 
     196      CALL iom_put( "avt" , avt                        )    ! T vert. eddy diff. coef. 
     197      CALL iom_put( "avm" , avmu                       )    ! T vert. eddy visc. coef. 
     198      IF( lk_zdftke ) THEN    
     199         CALL iom_put( "tke"      , en                               )    ! TKE budget: Turbulent Kinetic Energy    
     200         CALL iom_put( "tke_niw"  , e_niw                            )    ! TKE budget: Near-inertial waves    
     201      ENDIF  
     202      CALL iom_put( "avs" , fsavs(:,:,:)               )    ! S vert. eddy diff. coef. (useful only with key_zdfddm) 
     203 
     204      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 
    230205         DO jj = 2, jpjm1                                    ! sst gradient 
    231206            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    239214         CALL lbc_lnk( z2d, 'T', 1. ) 
    240215         CALL iom_put( "sstgrad2",  z2d               )    ! square of module of sst gradient 
    241          !CDIR NOVERRCHK< 
    242216         z2d(:,:) = SQRT( z2d(:,:) ) 
    243217         CALL iom_put( "sstgrad" ,  z2d               )    ! module of sst gradient 
     
    248222         z2d(:,:)  = 0._wp  
    249223         DO jk = 1, jpkm1 
    250             DO jj = 2, jpjm1 
    251                DO ji = fs_2, fs_jpim1   ! vector opt. 
     224            DO jj = 1, jpj 
     225               DO ji = 1, jpi 
    252226                  z2d(ji,jj) = z2d(ji,jj) + fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
    253227               END DO 
    254228            END DO 
    255229         END DO 
    256          CALL lbc_lnk( z2d, 'T', 1. ) 
    257230         CALL iom_put( "heatc", (rau0 * rcp) * z2d )    ! vertically integrated heat content (J/m2) 
    258231      ENDIF 
     
    261234         z2d(:,:)  = 0._wp  
    262235         DO jk = 1, jpkm1 
    263             DO jj = 2, jpjm1 
    264                DO ji = fs_2, fs_jpim1   ! vector opt. 
     236            DO jj = 1, jpj 
     237               DO ji = 1, jpi 
    265238                  z2d(ji,jj) = z2d(ji,jj) + fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
    266239               END DO 
    267240            END DO 
    268241         END DO 
    269          CALL lbc_lnk( z2d, 'T', 1. ) 
    270242         CALL iom_put( "saltc", rau0 * z2d )   ! vertically integrated salt content (PSU*kg/m2) 
    271243      ENDIF 
Note: See TracChangeset for help on using the changeset viewer.