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 4369 for branches/2013 – NEMO

Changeset 4369 for branches/2013


Ignore:
Timestamp:
2014-01-23T17:38:37+01:00 (10 years ago)
Author:
jchanut
Message:

Change time step order, ticket #1232

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/step.F90

    r4338 r4369  
    101101                         CALL sbc    ( kstp )         ! Sea Boundary Condition (including sea-ice) 
    102102                                                      ! clem: moved here for bdy ice purpose 
     103 
     104      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     105      ! Ocean physics update                (ua, va, tsa used as workspace) 
     106      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     107                         CALL bn2( tsb, rn2b )        ! before Brunt-Vaisala frequency 
     108                         CALL bn2( tsn, rn2  )        ! now    Brunt-Vaisala frequency 
     109      ! 
     110      !  VERTICAL PHYSICS 
     111                         CALL zdf_bfr( kstp )         ! bottom friction (if quadratic) 
     112      !                                               ! Vertical eddy viscosity and diffusivity coefficients 
     113      IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz 
     114      IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz 
     115      IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz 
     116      IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz 
     117      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value) 
     118         avt (:,:,:) = rn_avt0 * tmask(:,:,:) 
     119         avmu(:,:,:) = rn_avm0 * umask(:,:,:) 
     120         avmv(:,:,:) = rn_avm0 * vmask(:,:,:) 
     121      ENDIF 
     122      IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths 
     123         DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO 
     124      ENDIF 
     125      IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity 
     126 
     127      IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing 
     128 
     129      IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   & 
     130         &               CALL zdf_ddm( kstp )         ! double diffusive mixing 
     131 
     132                         CALL zdf_mxl( kstp )         ! mixed layer depth 
     133 
     134                                                      ! write TKE or GLS information in the restart file 
     135      IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' ) 
     136      IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' ) 
     137      ! 
     138      !  LATERAL  PHYSICS 
     139      ! 
     140      IF( lk_ldfslp ) THEN                            ! slope of lateral mixing 
     141                         CALL eos( tsb, rhd, gdept_0(:,:,:) )             ! before in situ density 
     142         IF( ln_zps )    CALL zps_hde( kstp, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient 
     143            &                                      rhd, gru , grv  )      ! of t, s, rd at the last ocean level 
     144         IF( ln_traldf_grif ) THEN                           ! before slope for Griffies operator 
     145                         CALL ldf_slp_grif( kstp ) 
     146         ELSE 
     147                         CALL ldf_slp( kstp, rhd, rn2b )     ! before slope for Madec operator 
     148         ENDIF 
     149      ENDIF 
     150#if defined key_traldf_c2d 
     151      IF( lk_traldf_eiv )   CALL ldf_eiv( kstp )      ! eddy induced velocity coefficient 
     152#endif 
     153#if defined key_traldf_c3d && key_traldf_smag 
     154                          CALL ldf_tra_smag( kstp )      ! eddy induced velocity coefficient 
     155#  endif 
     156#if defined key_dynldf_c3d && key_dynldf_smag 
     157                          CALL ldf_dyn_smag( kstp )      ! eddy induced velocity coefficient 
     158#  endif 
     159 
    103160      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    104161      !  Ocean dynamics : hdiv, rot, ssh, e3, wn 
    105162      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    106                          CALL zdf_bfr( kstp )         ! bottom friction (if quadratic) 
    107163                         CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_cur) 
     164 
    108165      IF( lk_dynspg_ts ) THEN 
    109166          IF( lk_vvl     )        CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors  
     
    139196      IF( lk_vvl     )   CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors 
    140197                         CALL wzv           ( kstp )  ! now cross-level velocity  
    141  
    142       !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    143       ! Ocean physics update                (ua, va, tsa used as workspace) 
    144       !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    145                          CALL bn2( tsb, rn2b )        ! before Brunt-Vaisala frequency 
    146                          CALL bn2( tsn, rn2  )        ! now    Brunt-Vaisala frequency 
    147       ! 
    148       !  VERTICAL PHYSICS 
    149       !                                               ! Vertical eddy viscosity and diffusivity coefficients 
    150       IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz 
    151       IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz 
    152       IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz 
    153       IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz 
    154       IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value) 
    155          avt (:,:,:) = rn_avt0 * tmask(:,:,:) 
    156          avmu(:,:,:) = rn_avm0 * umask(:,:,:) 
    157          avmv(:,:,:) = rn_avm0 * vmask(:,:,:) 
    158       ENDIF 
    159       IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths 
    160          DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO 
    161       ENDIF 
    162       IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity 
    163  
    164       IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing 
    165  
    166       IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   & 
    167          &               CALL zdf_ddm( kstp )         ! double diffusive mixing 
    168  
    169                          CALL zdf_mxl( kstp )         ! mixed layer depth 
    170  
    171                                                       ! write TKE or GLS information in the restart file 
    172       IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' ) 
    173       IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' ) 
    174       ! 
    175       !  LATERAL  PHYSICS 
    176       ! 
    177       IF( lk_ldfslp ) THEN                            ! slope of lateral mixing 
    178                          CALL eos( tsb, rhd, gdept_0(:,:,:) )             ! before in situ density 
    179          IF( ln_zps )    CALL zps_hde( kstp, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient 
    180             &                                      rhd, gru , grv  )      ! of t, s, rd at the last ocean level 
    181          IF( ln_traldf_grif ) THEN                           ! before slope for Griffies operator 
    182                          CALL ldf_slp_grif( kstp ) 
    183          ELSE 
    184                          CALL ldf_slp( kstp, rhd, rn2b )     ! before slope for Madec operator 
    185          ENDIF 
    186       ENDIF 
    187 #if defined key_traldf_c2d 
    188       IF( lk_traldf_eiv )   CALL ldf_eiv( kstp )      ! eddy induced velocity coefficient 
    189 #endif 
    190 #if defined key_traldf_c3d && key_traldf_smag 
    191                           CALL ldf_tra_smag( kstp )      ! eddy induced velocity coefficient 
    192 #  endif 
    193 #if defined key_dynldf_c3d && key_dynldf_smag 
    194                           CALL ldf_dyn_smag( kstp )      ! eddy induced velocity coefficient 
    195 #  endif 
    196  
    197198 
    198199      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
Note: See TracChangeset for help on using the changeset viewer.