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 14205 – NEMO

Changeset 14205


Ignore:
Timestamp:
2020-12-17T19:23:48+01:00 (3 years ago)
Author:
techene
Message:

#2385 cosmetics and cleaning

Location:
NEMO/trunk/src/OCE
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/DYN/dynspg_ts.F90

    r14064 r14205  
    162162      REAL(wp), DIMENSION(jpi,jpj) :: zCdU_u, zCdU_v   ! top/bottom stress at u- & v-points 
    163163      REAL(wp), DIMENSION(jpi,jpj) :: zhU, zhV         ! fluxes 
    164 #if defined key_qco  
    165       REAL(wp), DIMENSION(jpi, jpj, jpk) :: ze3u, ze3v 
    166 #endif 
     164!!st#if defined key_qco  
     165!!st      REAL(wp), DIMENSION(jpi, jpj, jpk) :: ze3u, ze3v 
     166!!st#endif 
    167167      ! 
    168168      REAL(wp) ::   zwdramp                     ! local scalar - only used if ln_wd_dl = .True.  
     
    232232      !                                   !  ---------------------------  ! 
    233233#if defined key_qco  
    234       DO jk = 1 , jpk 
    235          ze3u(:,:,jk) = e3u(:,:,jk,Kmm) 
    236          ze3v(:,:,jk) = e3v(:,:,jk,Kmm) 
    237       END DO 
    238       ! 
    239       zu_frc(:,:) = SUM( ze3u(:,:,:) * uu(:,:,:,Krhs) * umask(:,:,:) , DIM=3 ) * r1_hu(:,:,Kmm) 
    240       zv_frc(:,:) = SUM( ze3v(:,:,:) * vv(:,:,:,Krhs) * vmask(:,:,:) , DIM=3 ) * r1_hv(:,:,Kmm) 
     234!!st      DO jk = 1 , jpk 
     235!!st         ze3u(:,:,jk) = e3u(:,:,jk,Kmm) 
     236!!st         ze3v(:,:,jk) = e3v(:,:,jk,Kmm) 
     237!!st      END DO 
     238!!st      ! 
     239!!st      zu_frc(:,:) = SUM( ze3u(:,:,:) * uu(:,:,:,Krhs) * umask(:,:,:) , DIM=3 ) * r1_hu(:,:,Kmm) 
     240!!st      zv_frc(:,:) = SUM( ze3v(:,:,:) * vv(:,:,:,Krhs) * vmask(:,:,:) , DIM=3 ) * r1_hv(:,:,Kmm) 
     241 
     242      zu_frc(:,:) = SUM( e3u_0(:,:,:) * uu(:,:,:,Krhs) * umask(:,:,:) , DIM=3 ) * r1_hu_0(:,:) 
     243      zv_frc(:,:) = SUM( e3v_0(:,:,:) * vv(:,:,:,Krhs) * vmask(:,:,:) , DIM=3 ) * r1_hv_0(:,:) 
     244       
    241245#else 
    242246      zu_frc(:,:) = SUM( e3u(:,:,:,Kmm) * uu(:,:,:,Krhs) * umask(:,:,:) , DIM=3 ) * r1_hu(:,:,Kmm) 
  • NEMO/trunk/src/OCE/DYN/sshwzv.F90

    r14139 r14205  
    261261 
    262262 
    263    SUBROUTINE ssh_atf( kt, Kbb, Kmm, Kaa, pssh, pssh_f ) 
     263   SUBROUTINE ssh_atf( kt, Kbb, Kmm, Kaa, pssh ) 
    264264      !!---------------------------------------------------------------------- 
    265265      !!                    ***  ROUTINE ssh_atf  *** 
     
    278278      INTEGER                         , INTENT(in   ) ::   kt             ! ocean time-step index 
    279279      INTEGER                         , INTENT(in   ) ::   Kbb, Kmm, Kaa  ! ocean time level indices 
    280       REAL(wp), DIMENSION(jpi,jpj,jpt)          , TARGET, INTENT(inout) ::   pssh           ! SSH field 
    281       REAL(wp), DIMENSION(jpi,jpj    ), OPTIONAL, TARGET, INTENT(  out) ::   pssh_f         ! filtered SSH field 
     280      REAL(wp), DIMENSION(jpi,jpj,jpt), INTENT(inout) ::   pssh           ! SSH field 
    282281      ! 
    283282      REAL(wp) ::   zcoef   ! local scalar 
    284       REAL(wp), POINTER, DIMENSION(:,:) ::   zssh   ! pointer for filtered SSH  
    285283      !!---------------------------------------------------------------------- 
    286284      ! 
     
    292290         IF(lwp) WRITE(numout,*) '~~~~~~~ ' 
    293291      ENDIF 
    294       !              !==  Euler time-stepping: no filter, just swap  ==! 
    295       IF ( .NOT.( l_1st_euler ) ) THEN   ! Only do time filtering for leapfrog timesteps 
    296          IF( PRESENT( pssh_f ) ) THEN   ;   zssh => pssh_f 
    297          ELSE                           ;   zssh => pssh(:,:,Kmm) 
    298          ENDIF 
    299          !                                                  ! filtered "now" field 
    300          pssh(:,:,Kmm) = pssh(:,:,Kmm) + rn_atfp * ( pssh(:,:,Kbb) - 2 * pssh(:,:,Kmm) + pssh(:,:,Kaa) ) 
     292      ! 
     293      IF( .NOT.l_1st_euler ) THEN   ! Apply Asselin time filter on Kmm field (not on euler 1st) 
    301294         ! 
    302          IF( .NOT.ln_linssh ) THEN                          ! "now" <-- with forcing removed 
     295         IF( ln_linssh ) THEN                ! filtered "now" field 
     296            pssh(:,:,Kmm) = pssh(:,:,Kmm) + rn_atfp * ( pssh(:,:,Kbb) - 2 * pssh(:,:,Kmm) + pssh(:,:,Kaa) ) 
     297            ! 
     298         ELSE                                ! filtered "now" field with forcing removed 
    303299            zcoef = rn_atfp * rn_Dt * r1_rho0 
    304             pssh(:,:,Kmm) = pssh(:,:,Kmm) - zcoef * (     emp_b(:,:) - emp   (:,:)   & 
    305                &                             - rnf_b(:,:)        + rnf   (:,:)       & 
    306                &                             + fwfisf_cav_b(:,:) - fwfisf_cav(:,:)   & 
    307                &                             + fwfisf_par_b(:,:) - fwfisf_par(:,:)   ) * ssmask(:,:) 
     300            pssh(:,:,Kmm) = pssh(:,:,Kmm) + rn_atfp * ( pssh(:,:,Kbb) - 2 * pssh(:,:,Kmm) + pssh(:,:,Kaa) )   & 
     301               &                          - zcoef   * (         emp_b(:,:) -        emp(:,:)   & 
     302               &                                              - rnf_b(:,:) +        rnf(:,:)   & 
     303               &                                       + fwfisf_cav_b(:,:) - fwfisf_cav(:,:)   & 
     304               &                                       + fwfisf_par_b(:,:) - fwfisf_par(:,:)   ) * ssmask(:,:) 
    308305 
    309306            ! ice sheet coupling 
    310307            IF( ln_isf .AND. ln_isfcpl .AND. kt == nit000+1 )   & 
    311                &   pssh(:,:,Kbb) = pssh(:,:,Kbb) - rn_atfp * rn_Dt * ( risfcpl_ssh(:,:) - 0.0 ) * ssmask(:,:) 
     308               &   pssh(:,:,Kbb) = pssh(:,:,Kbb) - rn_atfp * rn_Dt * ( risfcpl_ssh(:,:) - 0._wp ) * ssmask(:,:) 
    312309 
    313310         ENDIF 
  • NEMO/trunk/src/OCE/stpmlf.F90

    r14201 r14205  
    5050   USE traatf_qco     ! time filtering                 (tra_atf_qco routine) 
    5151   USE dynatf_qco     ! time filtering                 (dyn_atf_qco routine) 
    52     
    53    USE bdydyn         ! ocean open boundary conditions (define bdy_dyn) 
    54  
    55 #if defined key_agrif 
    56    USE agrif_oce_interp 
    57 #endif 
    5852 
    5953   IMPLICIT NONE 
     
    6862#  include "do_loop_substitute.h90" 
    6963#  include "domzgr_substitute.h90" 
     64#  include "do_loop_substitute.h90" 
    7065   !!---------------------------------------------------------------------- 
    7166   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    9994      !!---------------------------------------------------------------------- 
    10095      INTEGER ::   ji, jj, jk, jtile   ! dummy loop indice 
    101       REAL(wp),              DIMENSION(jpi,jpj,jpk) ::   zgdept 
    102       REAL(wp), ALLOCATABLE, DIMENSION(:,:)         ::   zssh_f 
     96      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   zgdept 
    10397      !! --------------------------------------------------------------------- 
    10498#if defined key_agrif 
     
    210204      !  Ocean dynamics : hdiv, ssh, e3, u, v, w 
    211205      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    212       DO jk = 1, jpk 
    213          zgdept(:,:,jk) = gdept(:,:,jk,Nnn) 
    214       END DO 
     206       
    215207                         CALL ssh_nxt    ( kstp, Nbb, Nnn, ssh,  Naa )   ! after ssh (includes call to div_hor) 
    216208      IF( .NOT.lk_linssh ) THEN 
     
    221213                         CALL wzv        ( kstp, Nbb, Nnn, Naa, ww  )    ! Nnn cross-level velocity 
    222214      IF( ln_zad_Aimp )  CALL wAimp      ( kstp,      Nnn           )    ! Adaptive-implicit vertical advection partitioning 
     215                         ALLOCATE( zgdept(jpi,jpj,jpk) ) 
     216                         DO jk = 1, jpk 
     217                            zgdept(:,:,jk) = gdept(:,:,jk,Nnn) 
     218                         END DO 
    223219                         CALL eos        ( ts(:,:,:,:,Nnn), rhd, rhop, zgdept ) ! now in situ density for hpg computation 
    224  
     220                         DEALLOCATE( zgdept ) 
    225221 
    226222                         uu(:,:,:,Nrhs) = 0._wp            ! set dynamics trends to zero 
     
    240236                         CALL dyn_hpg( kstp,      Nnn      , uu, vv, Nrhs )  ! horizontal gradient of Hydrostatic pressure 
    241237                         CALL dyn_spg( kstp, Nbb, Nnn, Nrhs, uu, vv, ssh, uu_b, vv_b, Naa )  ! surface pressure gradient 
    242  
    243                                                       ! With split-explicit free surface, since now transports have been updated and ssh(:,:,Nrhs) as well 
    244       IF( ln_dynspg_ts ) THEN                         ! vertical scale factors and vertical velocity need to be updated 
     238                          
     239      IF( ln_dynspg_ts ) THEN      ! With split-explicit free surface, since now transports have been updated and ssh(:,:,Nrhs) 
     240                                   ! as well as vertical scale factors and vertical velocity need to be updated 
    245241                            CALL div_hor    ( kstp, Nbb, Nnn )                ! Horizontal divergence  (2nd call in time-split case) 
    246          IF(.NOT.lk_linssh) CALL dom_qco_r3c ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )   ! update ssh/h_0 ratio at t,u,v,f pts  
     242         IF(.NOT.lk_linssh) CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )   ! update ssh/h_0 ratio at t,u,v,f pts  
    247243      ENDIF 
    248244                            CALL dyn_zdf    ( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa  )  ! vertical diffusion 
     
    307303 
    308304#if defined key_agrif 
    309       IF(.NOT. Agrif_Root()) THEN 
     305      IF(.NOT. Agrif_Root() ) THEN 
    310306         IF( ln_tile    )   CALL dom_tile( ntsi, ntsj, ntei, ntej, ktile = 0 ) 
    311307                            CALL Agrif_Sponge_tra        ! tracers sponge 
     
    319315                            CALL tra_adv    ( kstp, Nbb, Nnn, ts, Nrhs )  ! hor. + vert. advection ==> RHS 
    320316         IF( ln_zdfmfc  )   CALL tra_mfc    ( kstp, Nbb,      ts, Nrhs )  ! Mass Flux Convection 
    321          IF( ln_zdfosm  )   CALL tra_osm    ( kstp,      Nnn, ts, Nrhs )  ! OSMOSIS non-local tracer fluxes ==> RHS 
    322          IF( lrst_oce .AND. ln_zdfosm )   & 
    323             &               CALL osm_rst    ( kstp,      Nnn, 'WRITE'  )  ! write OSMOSIS outputs + ww (so must do here) to restarts 
     317         IF( ln_zdfosm  ) THEN 
     318                            CALL tra_osm    ( kstp,      Nnn, ts, Nrhs )  ! OSMOSIS non-local tracer fluxes ==> RHS 
     319            IF( lrst_oce )  CALL osm_rst    ( kstp,      Nnn, 'WRITE'  )  ! write OSMOSIS outputs + ww (so must do here) to restarts 
     320         ENDIF 
    324321                            CALL tra_ldf    ( kstp, Nbb, Nnn, ts, Nrhs )  ! lateral mixing 
    325322 
     
    391388      IF( Agrif_NbStepint() == 0 .AND. nstop == 0 )   & 
    392389         &               CALL Agrif_update_all( )                  ! Update all components 
    393       ENDIF 
    394390 
    395391#endif 
     
    494490      !! ** Action :   puu(Kaa),pvv(Kaa)   after horizontal velocity and tracers 
    495491      !!---------------------------------------------------------------------- 
     492#if defined key_agrif 
     493      USE agrif_oce_interp 
     494#endif 
     495      USE bdydyn         ! ocean open boundary conditions (define bdy_dyn) 
     496      !! 
    496497      INTEGER                                  , INTENT(in   ) ::   kt         ! ocean time-step index 
    497498      INTEGER                                  , INTENT(in   ) ::   Kbb, Kaa   ! before and after time level indices 
Note: See TracChangeset for help on using the changeset viewer.