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 9019 for branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/C1D/step_c1d.F90 – NEMO

Ignore:
Timestamp:
2017-12-13T15:58:53+01:00 (6 years ago)
Author:
timgraham
Message:

Merge of dev_CNRS_2017 into branch

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/C1D/step_c1d.F90

    r6140 r9019  
    2727   PUBLIC stp_c1d      ! called by opa.F90 
    2828 
    29    !! * Substitutions 
    30 #  include "zdfddm_substitute.h90" 
    3129   !!---------------------------------------------------------------------- 
    3230   !! NEMO/C1D 3.7 , NEMO Consortium (2015) 
     
    7674                         CALL bn2( tsb, rab_b, rn2b ) ! before Brunt-Vaisala frequency 
    7775                         CALL bn2( tsn, rab_n, rn2  ) ! now    Brunt-Vaisala frequency 
    78       !  VERTICAL PHYSICS    
    79                          CALL zdf_bfr( kstp )         ! bottom friction 
    80       !                                               ! Vertical eddy viscosity and diffusivity coefficients 
    81       IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz 
    82       IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz 
    83       IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz 
    84       IF( lk_zdfcst  )   THEN                            ! Constant Kz (reset avt, avm[uv] to the background value) 
    85          avt (:,:,:) = rn_avt0 * tmask(:,:,:) 
    86          avmu(:,:,:) = rn_avm0 * umask(:,:,:) 
    87          avmv(:,:,:) = rn_avm0 * vmask(:,:,:) 
    88       ENDIF 
     76       
     77      !  VERTICAL PHYSICS 
     78                         CALL zdf_phy( kstp )         ! vertical physics update (bfr, avt, avs, avm + MLD) 
    8979 
    90       IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths 
    91          DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:)   ;   END DO 
    92       ENDIF 
    93       IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity 
    94       IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing 
    95       IF( lk_zdfddm  )   CALL zdf_ddm( kstp )         ! double diffusive mixing 
    96                          CALL zdf_mxl( kstp )         ! mixed layer depth 
     80      IF(.NOT.ln_linssh )   CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_hor) 
     81      IF(.NOT.ln_linssh )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors  
    9782 
    98                                                       ! write tke information in the restart file 
    99       IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' ) 
    100                                                       ! write gls information in the restart file 
    101       IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' ) 
    102  
     83      IF(.NOT.ln_linssh )   CALL wzv           ( kstp )  ! now cross-level velocity  
    10384      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    10485      ! diagnostics and outputs             (ua, va, ta, sa used as workspace) 
     
    123104      IF( ln_traqsr )   CALL tra_qsr( kstp )       ! penetrative solar radiation qsr 
    124105      IF( ln_tradmp )   CALL tra_dmp( kstp )       ! internal damping trends- tracers 
     106      IF(.NOT.ln_linssh)CALL tra_adv( kstp )       ! horizontal & vertical advection 
     107      IF( ln_zdfosm  )  CALL tra_osm( kstp )       ! OSMOSIS non-local tracer fluxes 
    125108                        CALL tra_zdf( kstp )       ! vertical mixing 
    126109                        CALL eos( tsn, rhd, rhop, gdept_0(:,:,:) )   ! now potential density for zdfmxl 
     
    138121      IF( ln_dyndmp )   CALL dyn_dmp    ( kstp )   ! internal damping trends- momentum 
    139122                        CALL dyn_cor_c1d( kstp )   ! vorticity term including Coriolis 
     123      IF( ln_zdfosm  )  CALL dyn_osm    ( kstp )   ! OSMOSIS non-local velocity fluxes 
    140124                        CALL dyn_zdf    ( kstp )   ! vertical diffusion 
    141125                        CALL dyn_nxt    ( kstp )   ! lateral velocity at next time step 
     126      IF(.NOT.ln_linssh)CALL ssh_swp    ( kstp )   ! swap of sea surface height 
     127 
     128      IF(.NOT.ln_linssh)CALL dom_vvl_sf_swp( kstp )! swap of vertical scale factors 
    142129 
    143130      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
Note: See TracChangeset for help on using the changeset viewer.