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 12732 for NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/stpMLF.F90 – NEMO

Ignore:
Timestamp:
2020-04-09T21:06:01+02:00 (4 years ago)
Author:
techene
Message:

some cleaning and proper module/routine name, mini bug introduced and corrected in sbcice_cice

File:
1 moved

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/stpMLF.F90

    r12731 r12732  
    1 MODULE steplf 
     1MODULE stepMLF 
    22   !!====================================================================== 
    33   !!                       ***  MODULE step  *** 
     
    3535 
    3636   !!---------------------------------------------------------------------- 
    37    !!   stplf             : OPA system time-stepping 
    38    !!---------------------------------------------------------------------- 
    39    USE step_oce         ! time stepping definition modules 
     37   !!   stp_MLF       : OPA system time-stepping 
     38   !!---------------------------------------------------------------------- 
     39   USE step_oce       ! time stepping definition modules 
    4040   ! 
    41    USE iom              ! xIOs server 
    42    USE domqe 
    43    USE traatfqco         ! time filtering                   (tra_atf_qco routine) 
    44    USE dynatfqco         ! time filtering                   (dyn_atf_qco routine) 
    45    USE dynspg_ts         ! surface pressure gradient: split-explicit scheme (define un_adv) 
    46    USE bdydyn            ! ocean open boundary conditions (define bdy_dyn) 
     41   USE iom            ! xIOs server 
     42   USE domqco 
     43   USE traatfqco      ! time filtering                   (tra_atf_qco routine) 
     44   USE dynatfqco      ! time filtering                   (dyn_atf_qco routine) 
     45   USE dynspg_ts      ! surface pressure gradient: split-explicit scheme (define un_adv) 
     46   USE bdydyn         ! ocean open boundary conditions (define bdy_dyn) 
    4747 
    4848   IMPLICIT NONE 
    4949   PRIVATE 
    5050 
    51    PUBLIC   stplf   ! called by nemogcm.F90 
     51   PUBLIC   stp_MLF   ! called by nemogcm.F90 
    5252 
    5353   !!---------------------------------------------------------------------- 
     
    6464 
    6565#if defined key_agrif 
    66    RECURSIVE SUBROUTINE stplf( ) 
     66   RECURSIVE SUBROUTINE stp_MLF( ) 
    6767      INTEGER             ::   kstp   ! ocean time-step index 
    6868#else 
    69    SUBROUTINE stplf( kstp ) 
     69   SUBROUTINE stp_MLF( kstp ) 
    7070      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index 
    7171#endif 
    7272      !!---------------------------------------------------------------------- 
    73       !!                     ***  ROUTINE stplf  *** 
     73      !!                     ***  ROUTINE stp_MLF  *** 
    7474      !! 
    7575      !! ** Purpose : - Time stepping of OPA  (momentum and active tracer eqs.) 
     
    9191      INTEGER ::   kcall        ! optional integer argument (dom_vvl_sf_nxt) 
    9292!!st patch 
    93       REAL(wp), DIMENSION(jpi,jpj,jpk) :: zgdept  ! st patch 
     93      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zgdept  
    9494      !! --------------------------------------------------------------------- 
    9595#if defined key_agrif 
     
    106106#endif 
    107107      ! 
    108       IF( ln_timing )   CALL timing_start('stplf') 
     108      IF( ln_timing )   CALL timing_start('stp_MLF') 
    109109      ! 
    110110      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    142142      IF( ln_isf     )   CALL isf_stp ( kstp, Nnn ) 
    143143                         CALL sbc     ( kstp, Nbb, Nnn )              ! Sea Boundary Condition (including sea-ice) 
    144  
    145       !             !!st      !!!!!!!!!!!!!!!!!!!!!!! 
    146       ! 
    147       ! emp = 0._wp 
    148       ! emp_b = 0._wp 
    149       ! qns = 0._wp 
    150       ! qsr = 0._wp 
    151       ! qns_b = 0._wp 
    152       ! 
    153       !             !!st 
    154144 
    155145      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    196186      !  Ocean dynamics : hdiv, ssh, e3, u, v, w 
    197187      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    198 !!st patch 
    199 !!st patch to be able to use substitution 
    200188      DO jk = 1, jpk 
    201189         zgdept(:,:,jk) = gdept(:,:,jk,Nnn) 
    202190      END DO 
    203 !!st end 
    204                             CALL ssh_nxt       ( kstp, Nbb, Nnn, ssh,  Naa )    ! after ssh (includes call to div_hor) 
    205       IF( .NOT.ln_linssh )  CALL dom_qe_r3c    ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa) ) 
    206 !!st      IF( .NOT.ln_linssh )  CALL dom_h_nxt     ( kstp, Nbb, Nnn,       Naa )    ! after vertical scale factors 
    207       !IF( .NOT.ln_linssh )  CALL dom_qe_sf_nxt ( kstp, Nbb, Nnn,      Naa )    ! after vertical scale factors 
     191                            CALL ssh_nxt       ( kstp, Nbb, Nnn, ssh,  Naa )   ! after ssh (includes call to div_hor) 
     192      IF( .NOT.ln_linssh )  CALL dom_qco_r3c    ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa) )   ! "after" ssh./h._0 ratio 
    208193                            CALL wzv           ( kstp, Nbb, Nnn, Naa, ww  )    ! Nnn cross-level velocity 
    209       IF( ln_zad_Aimp )     CALL wAimp         ( kstp,      Nnn           )  ! Adaptive-implicit vertical advection partitioning 
    210                             CALL eos    ( ts(:,:,:,:,Nnn), rhd, rhop, zgdept )  ! now in situ density for hpg computation 
     194      IF( ln_zad_Aimp )     CALL wAimp         ( kstp,      Nnn           )    ! Adaptive-implicit vertical advection partitioning 
     195                            CALL eos    ( ts(:,:,:,:,Nnn), rhd, rhop, zgdept ) ! now in situ density for hpg computation 
    211196 
    212197 
     
    231216      IF( ln_dynspg_ts ) THEN                         ! vertical scale factors and vertical velocity need to be updated 
    232217                            CALL div_hor    ( kstp, Nbb, Nnn )                ! Horizontal divergence  (2nd call in time-split case) 
    233          IF(.NOT.ln_linssh) CALL dom_qe_r3c ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) ) 
    234          !IF(.NOT.ln_linssh) CALL dom_qe_sf_nxt ( kstp, Nbb, Nnn, Naa, kcall=2 )  ! after vertical scale factors (update depth average component) 
    235 !!st         IF(.NOT.ln_linssh) CALL dom_h_nxt  ( kstp, Nbb, Nnn, Naa, kcall=2 )  ! after vertical scale factors (update depth average component) 
     218         IF(.NOT.ln_linssh) CALL dom_qco_r3c ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )   ! "after" ssh./h._0 ratio 
    236219      ENDIF 
    237220                            CALL dyn_zdf    ( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa  )  ! vertical diffusion 
     
    311294!!    place. 
    312295!! 
    313                          CALL zdyn_ts       ( Nnn, Naa, uu, vv )                   ! barotrope ajustment 
    314                          CALL finalize_sbc  ( kstp, Nbb, Naa, uu, vv, ts )                   ! boundary condifions 
    315                          CALL ssh_atf       ( kstp, Nbb, Nnn, Naa, ssh )                     ! time filtering of "now" sea surface height 
    316                          CALL dom_qe_r3c    ( ssh(:,:,Nnn), r3t_f, r3u_f, r3v_f )            ! "now" ssh/h_0 ratio from filtrered ssh 
    317                          CALL tra_atf_qco   ( kstp, Nbb, Nnn, Naa, ts )                      ! time filtering of "now" tracer arrays 
    318                          CALL dyn_atf_qco   ( kstp, Nbb, Nnn, Naa, uu, vv  )  ! time filtering of "now" velocities and scale factors 
     296                         CALL zdyn_ts       ( Nnn, Naa, uu, vv )                    ! barotrope ajustment 
     297                         CALL finalize_sbc  ( kstp, Nbb, Naa, uu, vv, ts )          ! boundary condifions 
     298                         CALL ssh_atf       ( kstp, Nbb, Nnn, Naa, ssh )            ! time filtering of "now" sea surface height 
     299                         CALL dom_qco_r3c    ( ssh(:,:,Nnn), r3t_f, r3u_f, r3v_f )   ! "now" ssh/h_0 ratio from filtrered ssh 
     300                         CALL tra_atf_qco   ( kstp, Nbb, Nnn, Naa, ts )             ! time filtering of "now" tracer arrays 
     301                         CALL dyn_atf_qco   ( kstp, Nbb, Nnn, Naa, uu, vv  )        ! time filtering of "now" velocities and scale factors 
    319302                         r3t(:,:,Nnn) = r3t_f(:,:) 
    320303                         r3u(:,:,Nnn) = r3u_f(:,:) 
     
    328311      Naa = Nrhs 
    329312      ! 
    330       !IF(.NOT.ln_linssh) CALL dom_qe_sf_update( kstp, Nbb, Nnn, Naa )  ! recompute vertical scale factors 
    331 !!st      IF(.NOT.ln_linssh) CALL dom_h_update  ( kstp, Nbb, Nnn, Naa )  ! recompute vertical scale factors 
    332313      ! 
    333314      IF( ln_diahsb  )   CALL dia_hsb       ( kstp, Nbb, Nnn )  ! - ML - global conservation diagnostics 
     
    346327      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    347328                         Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs      ! agrif_oce module copies of time level indices 
    348                          CALL Agrif_Integrate_ChildGrids( stplf )       ! allows to finish all the Child Grids before updating 
     329                         CALL Agrif_Integrate_ChildGrids( stp_MLF )       ! allows to finish all the Child Grids before updating 
    349330 
    350331                         IF( Agrif_NbStepint() == 0 ) THEN 
     
    385366      ENDIF 
    386367      ! 
    387       IF( ln_timing )   CALL timing_stop('stplf') 
    388       ! 
    389    END SUBROUTINE stplf 
     368      IF( ln_timing )   CALL timing_stop('stp_MLF') 
     369      ! 
     370   END SUBROUTINE stp_MLF 
    390371 
    391372 
     
    481462   ! 
    482463   !!====================================================================== 
    483 END MODULE steplf 
     464END MODULE stepMLF 
Note: See TracChangeset for help on using the changeset viewer.