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 6140 for trunk/NEMOGCM/NEMO/OPA_SRC/step.F90 – NEMO

Ignore:
Timestamp:
2015-12-21T12:35:23+01:00 (9 years ago)
Author:
timgraham
Message:

Merge of branches/2015/dev_merge_2015 back into trunk. Merge excludes NEMOGCM/TOOLS/OBSTOOLS/ for now due to issues with the change of file type. Will sort these manually with further commits.

Branch merged as follows:
In the working copy of branch ran:
svn merge svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk@HEAD
Small conflicts due to bug fixes applied to trunk since the dev_merge_2015 was copied. Bug fixes were applied to the branch as well so these were easy to resolve.
Branch committed at this stage

In working copy run:
svn switch svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk
to switch working copy

Run:
svn merge --reintegrate svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/2015/dev_merge_2015
to merge the branch into the trunk and then commit - no conflicts at this stage.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/NEMO/OPA_SRC/step.F90

    r5930 r6140  
    22   !!====================================================================== 
    33   !!                       ***  MODULE step  *** 
    4    !! Time-stepping    : manager of the ocean, tracer and ice time stepping 
     4   !! Time-stepping   : manager of the ocean, tracer and ice time stepping 
    55   !!====================================================================== 
    66   !! History :  OPA  !  1991-03  (G. Madec)  Original code 
     
    3535   !!---------------------------------------------------------------------- 
    3636   USE step_oce         ! time stepping definition modules 
    37    USE iom 
     37   ! 
     38   USE iom              ! xIOs server 
    3839 
    3940   IMPLICIT NONE 
     
    4243   PUBLIC   stp   ! called by nemogcm.F90 
    4344 
    44    !! * Substitutions 
    45 #  include "domzgr_substitute.h90" 
    46 !!gm   #  include "zdfddm_substitute.h90" 
    4745   !!---------------------------------------------------------------------- 
    4846   !! NEMO/OPA 3.7 , NEMO Consortium (2015) 
     
    166164 
    167165         IF( ln_zps .AND.       ln_isfcav)                               & 
    168             &            CALL zps_hde_isf( kstp, jpts, tsb, gtsu, gtsv, gtui, gtvi,   &    ! Partial steps for top cell (ISF) 
    169             &                                          rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   & 
    170             &                                   grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the first ocean level 
    171  
     166            &            CALL zps_hde_isf( kstp, jpts, tsb, gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF) 
     167            &                                          rhd, gru , grv , grui, grvi   )  ! of t, s, rd at the first ocean level 
    172168         IF( ln_traldf_triad ) THEN  
    173169                         CALL ldf_slp_triad( kstp )                       ! before slope for triad operator 
     
    183179      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    184180 
    185                          CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_hor) 
    186       IF( lk_vvl     )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors  
    187                          CALL wzv           ( kstp )  ! now cross-level velocity  
    188  
     181                            CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_hor) 
     182      IF(.NOT.ln_linssh )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors  
     183                            CALL wzv           ( kstp )  ! now cross-level velocity  
    189184!!gm : why also here ???? 
    190       IF(ln_sto_eos  )   CALL sto_pts( tsn )                             ! Random T/S fluctuations 
     185      IF( ln_sto_eos    )   CALL sto_pts( tsn )                             ! Random T/S fluctuations 
    191186!!gm 
    192                          CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) ) ! now in situ density for hpg computation 
    193  
     187                            CALL eos    ( tsn, rhd, rhop, gdept_n(:,:,:) ) ! now in situ density for hpg computation 
     188                             
    194189!!jc: fs simplification 
    195 !!jc: lines below are useless if lk_vvl=T. Keep them here (which maintains a bug if lk_vvl=F and ln_zps=T, cf ticket #1636)  
     190!!jc: lines below are useless if ln_linssh=F. Keep them here (which maintains a bug if ln_linssh=T and ln_zps=T, cf ticket #1636)  
    196191!!                                         but ensures reproductible results 
    197192!!                                         with previous versions using split-explicit free surface           
    198             IF( ln_zps .AND. .NOT. ln_isfcav)   &                           ! Partial steps: bottom before horizontal gradient 
    199                &            CALL zps_hde    ( kstp, jpts, tsn, gtsu, gtsv,  &  ! of t, s, rd at the last ocean level 
    200                &                                          rhd, gru , grv    ) 
    201             IF( ln_zps .AND.       ln_isfcav)   &                           ! Partial steps: top & bottom before horizontal gradient 
    202                &            CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv, gtui, gtvi,   &  
    203                &                                          rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   & 
    204                &                                               grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) 
     193            IF( ln_zps .AND. .NOT. ln_isfcav )                               & 
     194               &            CALL zps_hde    ( kstp, jpts, tsn, gtsu, gtsv,   &  ! Partial steps: before horizontal gradient 
     195               &                                          rhd, gru , grv     )  ! of t, s, rd at the last ocean level 
     196            IF( ln_zps .AND.       ln_isfcav )                                          & 
     197               &            CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF) 
     198               &                                          rhd, gru , grv , grui, grvi   )  ! of t, s, rd at the first ocean level 
    205199!!jc: fs simplification 
    206200                             
     
    224218      IF( ln_dynspg_ts ) THEN                         ! vertical scale factors and vertical velocity need to be updated 
    225219                            CALL div_hor    ( kstp )              ! Horizontal divergence  (2nd call in time-split case) 
    226          IF( lk_vvl )      CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors (update depth average component) 
     220         IF(.NOT.ln_linssh) CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors (update depth average component) 
    227221                            CALL wzv        ( kstp )              ! now cross-level velocity  
    228222      ENDIF 
     
    232226 
    233227      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    234       ! diagnostics and outputs              
    235       !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    236       IF( lk_floats  )   CALL flo_stp       ( kstp )  ! drifting Floats 
    237       IF( lk_diahth  )   CALL dia_hth       ( kstp )  ! Thermocline depth (20 degres isotherm depth) 
    238       IF(.NOT.ln_cpl )   CALL dia_fwb       ( kstp )  ! Fresh water budget diagnostics 
    239       IF( lk_diadct  )   CALL dia_dct       ( kstp )  ! Transports 
    240       IF( lk_diaar5  )   CALL dia_ar5       ( kstp )  ! ar5 diag 
    241       IF( lk_diaharm )   CALL dia_harm      ( kstp )  ! Tidal harmonic analysis 
    242                          CALL dia_wri       ( kstp )  ! ocean model: outputs 
     228      ! cool skin 
     229      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<       
     230      IF ( ln_diurnal )  CALL stp_diurnal( kstp ) 
     231       
     232      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     233      ! diagnostics and outputs             (ua, va, tsa used as workspace) 
     234      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     235      IF( lk_floats  )   CALL flo_stp( kstp )         ! drifting Floats 
     236      IF( nn_diacfl == 1 )   CALL dia_cfl( kstp )         ! Courant number diagnostics 
     237      IF( lk_diahth  )   CALL dia_hth( kstp )         ! Thermocline depth (20 degres isotherm depth) 
     238      IF(.NOT.ln_cpl )   CALL dia_fwb( kstp )         ! Fresh water budget diagnostics 
     239      IF( lk_diadct  )   CALL dia_dct( kstp )         ! Transports 
     240      IF( lk_diaar5  )   CALL dia_ar5( kstp )         ! ar5 diag 
     241      IF( lk_diaharm )   CALL dia_harm( kstp )        ! Tidal harmonic analysis 
     242                         CALL dia_wri( kstp )         ! ocean model: outputs 
    243243      ! 
    244244      IF( ln_crs     )   CALL crs_fld       ( kstp )  ! ocean model: online field coarsening & output 
    245  
     245       
    246246#if defined key_top 
    247247      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    293293!!    place. 
    294294!!  
    295 !!jc2: dynnxt must be the latest call. fse3t_b are indeed updated in that routine 
     295!!jc2: dynnxt must be the latest call. e3t_b are indeed updated in that routine 
    296296                         CALL tra_nxt       ( kstp )  ! finalize (bcs) tracer fields at next time step and swap 
    297297                         CALL dyn_nxt       ( kstp )  ! finalize (bcs) velocities at next time step and swap 
    298298                         CALL ssh_swp       ( kstp )  ! swap of sea surface height 
    299       IF( lk_vvl     )  CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors 
     299      IF(.NOT.ln_linssh) CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors 
    300300      ! 
    301301 
     
    313313                         CALL Agrif_Integrate_ChildGrids( stp )   
    314314 
    315       IF ( Agrif_NbStepint().EQ.0 ) THEN              ! AGRIF Update  
    316 !!jc in fact update i useless at last time step, but do it for global diagnostics 
     315      IF( Agrif_NbStepint() == 0 ) THEN               ! AGRIF Update  
     316!!jc in fact update is useless at last time step, but do it for global diagnostics 
    317317                         CALL Agrif_Update_Tra()      ! Update active tracers 
    318318                         CALL Agrif_Update_Dyn()      ! Update momentum 
    319319      ENDIF 
    320320#endif 
    321       IF( ln_diahsb  )   CALL dia_hsb       ( kstp )  ! - ML - global conservation diagnostics 
    322       IF( lk_diaobs  )   CALL dia_obs       ( kstp )  ! obs-minus-model (assimilation) diagnostics (call after dynamics update) 
     321      IF( ln_diahsb  )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics 
     322      IF( ln_diaobs  )   CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update) 
    323323 
    324324      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    326326      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    327327                         CALL stp_ctl       ( kstp, indic ) 
    328       IF( indic < 0  )   THEN 
     328      IF( indic < 0  ) THEN 
    329329                         CALL ctl_stop( 'step: indic < 0' ) 
    330330                         CALL dia_wri_state( 'output.abort', kstp ) 
    331331      ENDIF 
    332       IF( kstp == nit000 )   THEN 
     332      IF( kstp == nit000 ) THEN 
    333333                 CALL iom_close( numror )     ! close input  ocean restart file 
    334334         IF(lwm) CALL FLUSH    ( numond )     ! flush output namelist oce 
     
    353353      ! 
    354354   END SUBROUTINE stp 
    355  
    356    !!====================================================================== 
     355    
    357356END MODULE step 
Note: See TracChangeset for help on using the changeset viewer.