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 7403 for branches/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/step.F90 – NEMO

Ignore:
Timestamp:
2016-11-30T17:56:53+01:00 (8 years ago)
Author:
timgraham
Message:

Merge dev_INGV_METO_merge_2016 into branch

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/step.F90

    r6464 r7403  
    2626   !!            3.6  !  2012-07  (J. Simeon, G. Madec. C. Ethe)  Online coarsening of outputs 
    2727   !!            3.6  !  2014-04  (F. Roquet, G. Madec) New equations of state 
     28   !!            3.6  !  2014-10  (E. Clementi, P. Oddo) Add Qiao vertical mixing in case of waves 
    2829   !!            3.7  !  2014-10  (G. Madec)  LDF simplication  
    2930   !!             -   !  2014-12  (G. Madec) remove KPP scheme 
     
    7374      !!              -8- Outputs and diagnostics 
    7475      !!---------------------------------------------------------------------- 
    75       INTEGER ::   jk      ! dummy loop indice 
     76      INTEGER ::   ji,jj,jk ! dummy loop indice 
    7677      INTEGER ::   indic    ! error indicator if < 0 
    7778      INTEGER ::   kcall    ! optional integer argument (dom_vvl_sf_nxt) 
     
    128129                         CALL zdf_bfr( kstp )         ! bottom friction (if quadratic) 
    129130      !                                               ! Vertical eddy viscosity and diffusivity coefficients 
    130       IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz 
    131       IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz 
    132       IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz 
    133       IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value) 
     131      IF( lk_zdfric  )   CALL zdf_ric ( kstp )             ! Richardson number dependent Kz 
     132      IF( lk_zdftke  )   CALL zdf_tke ( kstp )             ! TKE closure scheme for Kz 
     133      IF( lk_zdfgls  )   CALL zdf_gls ( kstp )             ! GLS closure scheme for Kz 
     134      IF( ln_zdfqiao )   CALL zdf_qiao( kstp )             ! Qiao vertical mixing  
     135      ! 
     136      IF( lk_zdfcst  ) THEN                                ! Constant Kz (reset avt, avm[uv] to the background value) 
    134137         avt (:,:,:) = rn_avt0 * wmask (:,:,:) 
    135138         avmu(:,:,:) = rn_avm0 * wumask(:,:,:) 
     
    207210                         CALL dyn_adv       ( kstp )  ! advection (vector or flux form) 
    208211                         CALL dyn_vor       ( kstp )  ! vorticity term including Coriolis 
     212      IF( ln_wave .AND. ln_sdw .AND. ln_stcor)           & 
     213               &         CALL dyn_stcor     ( kstp )  ! Stokes-Coriolis forcing 
    209214                         CALL dyn_ldf       ( kstp )  ! lateral mixing 
    210215                         CALL dyn_hpg       ( kstp )  ! horizontal gradient of Hydrostatic pressure 
     
    234239      IF(.NOT.ln_cpl )   CALL dia_fwb( kstp )         ! Fresh water budget diagnostics 
    235240      IF( lk_diadct  )   CALL dia_dct( kstp )         ! Transports 
    236       IF( lk_diaar5  )   CALL dia_ar5( kstp )         ! ar5 diag 
     241                         CALL dia_ar5( kstp )         ! ar5 diag 
    237242      IF( lk_diaharm )   CALL dia_harm( kstp )        ! Tidal harmonic analysis 
    238243                         CALL dia_wri( kstp )         ! ocean model: outputs 
Note: See TracChangeset for help on using the changeset viewer.