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 14618 for NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/stprk3_stg.F90 – NEMO

Ignore:
Timestamp:
2021-03-19T15:42:32+01:00 (3 years ago)
Author:
techene
Message:

#2506 RK3 work in progess : for 2 configurag

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/stprk3_stg.F90

    r14549 r14618  
    7272      ! 
    7373      INTEGER  ::   ji, jj, jk, jtile                      ! dummy loop indices 
    74       REAL(wp) ::   ze3Ub, ze3Vb, ze3Tb, ze3Sb, z1_e3t     ! local scalars 
    75       REAL(wp) ::   ze3Ur, ze3Vr, ze3Tr, ze3Sr             !   -      - 
     74      REAL(wp) ::   ze3Tb, ze3Sb, z1_e3t     ! local scalars 
     75      REAL(wp) ::   ze3Tr, ze3Sr             !   -      - 
    7676      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zaU, zaV       ! advective horizontal velocity 
    7777      REAL(wp), DIMENSION(jpi,jpj)     ::   zub, zvb       ! advective transport  
     
    211211      CALL tra_adv( kstp, Kbb, Kmm, ts, Krhs, zaU, zaV, ww )       ! hor. + vert. advection  ==> RHS 
    212212 
    213  
    214213!===>>>>>> stg1&2:  Verify the necessity of these trends (we may need it as there are in the RHS of dynspg_ts ?) 
    215214!!gm ====>>>>   needed for heat and salt fluxes associated with mass/volume flux 
    216 !                     CALL tra_sbc( kstp,      Kmm, ts, Krhs )   ! surface boundary condition 
    217 !      IF( ln_isf )   CALL tra_isf( kstp,      Kmm, ts, Krhs )   ! ice shelf heat flux 
     215                        CALL tra_sbc( kstp,      Kmm, ts, Krhs )   ! surface boundary condition 
     216 
     217      IF( ln_isf )      CALL tra_isf( kstp,      Kmm, ts, Krhs )   ! ice shelf heat flux 
     218      IF( ln_traqsr )   CALL tra_qsr( kstp,      Kmm, ts, Krhs )   ! penetrative solar radiation qsr 
    218219!!gm 
    219220 
     
    221222!!gm ===>>>>>>  Verify the necessity of these trends  at stages 1 and 2  
    222223!           (we may need it as they are in the RHS of dynspg_ts ?) 
    223       IF(  lk_asminc .AND. ln_asmiau ) THEN               ! apply assimilation increment 
    224          IF( ln_dyninc )   CALL dyn_asm_inc( kstp, Kbb, Kmm, uu, vv, Krhs )   ! dynamics   ==> RHS 
    225          IF( ln_trainc )   CALL tra_asm_inc( kstp, Kbb, Kmm, ts    , Krhs )   ! tracers    ==> RHS 
    226       ENDIF 
     224!      IF(  lk_asminc .AND. ln_asmiau ) THEN               ! apply assimilation increment 
     225!         IF( ln_dyninc )   CALL dyn_asm_inc( kstp, Kbb, Kmm, uu, vv, Krhs )   ! dynamics   ==> RHS 
     226!         IF( ln_trainc )   CALL tra_asm_inc( kstp, Kbb, Kmm, ts    , Krhs )   ! tracers    ==> RHS 
     227!      ENDIF 
    227228!!gm  end Verif 
    228229 
     
    234235         ! 
    235236         !                                      !==  time integration  ==!   ∆t = rn_Dt/3 (stg1) or rn_Dt/2 (stg2) 
     237         IF( ln_dynadv_vec .OR. ln_linssh ) THEN   ! applied on velocity 
     238            DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     239               uu(ji,jj,jk,Kaa) = ( uu(ji,jj,jk,Kbb) + rDt * uu(ji,jj,jk,Krhs) ) * umask(ji,jj,jk) 
     240               vv(ji,jj,jk,Kaa) = ( vv(ji,jj,jk,Kbb) + rDt * vv(ji,jj,jk,Krhs) ) * vmask(ji,jj,jk) 
     241            END_3D 
     242         ELSE                                      ! applied on thickness weighted velocity 
     243            DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     244               uu(ji,jj,jk,Kaa) = (         e3u(ji,jj,jk,Kbb) * uu(ji,jj,jk,Kbb )  & 
     245                  &                 + rDt * e3u(ji,jj,jk,Kmm) * uu(ji,jj,jk,Krhs)  ) & 
     246                  &                       / e3u(ji,jj,jk,Kaa) * umask(ji,jj,jk) 
     247               vv(ji,jj,jk,Kaa) = (         e3v(ji,jj,jk,Kbb) * vv(ji,jj,jk,Kbb )  & 
     248                  &                 + rDt * e3v(ji,jj,jk,Kmm) * vv(ji,jj,jk,Krhs)  ) & 
     249                  &                       / e3v(ji,jj,jk,Kaa) * vmask(ji,jj,jk) 
     250            END_3D 
     251         ENDIF 
     252         ! 
    236253         DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
    237             ze3Ub = e3u(ji,jj,jk,Kbb) * uu(ji,jj,jk,Kbb ) 
    238             ze3Vb = e3v(ji,jj,jk,Kbb) * vv(ji,jj,jk,Kbb ) 
    239             ze3Ur = e3u(ji,jj,jk,Kmm) * uu(ji,jj,jk,Krhs) 
    240             ze3Vr = e3v(ji,jj,jk,Kmm) * vv(ji,jj,jk,Krhs) 
    241             uu(ji,jj,jk,Kaa) = ( ze3Ub + rDt * ze3Ur*umask(ji,jj,jk) ) / e3u(ji,jj,jk, Kaa) 
    242             vv(ji,jj,jk,Kaa) = ( ze3Vb + rDt * ze3Vr*vmask(ji,jj,jk) ) / e3v(ji,jj,jk, Kaa) 
    243             ! 
    244254            ze3Tb = e3t(ji,jj,jk,Kbb) * ts(ji,jj,jk,jp_tem,Kbb ) 
    245255            ze3Sb = e3t(ji,jj,jk,Kbb) * ts(ji,jj,jk,jp_sal,Kbb ) 
     
    291301         ! 
    292302                            CALL tra_ldf( kstp, Kbb, Kmm, ts, Krhs )  ! lateral mixing 
    293 !!gm ====>>>>   needed for heat and salt fluxes associated with mass/volume flux 
    294                             CALL tra_sbc( kstp,      Kmm, ts, Krhs )   ! surface boundary condition 
    295          IF( ln_isf )       CALL tra_isf( kstp,      Kmm, ts, Krhs )   ! ice shelf heat flux 
     303!!gm ====>>>>   already called in stage 1, 2 and 3 case 
     304!!st                            CALL tra_sbc( kstp,      Kmm, ts, Krhs )   ! surface boundary condition 
     305!!st         IF( ln_isf )       CALL tra_isf( kstp,      Kmm, ts, Krhs )   ! ice shelf heat flux 
    296306!!gm 
    297307         ! 
    298          IF( ln_traqsr  )   CALL tra_qsr( kstp,      Kmm, ts, Krhs )  ! penetrative solar radiation qsr 
     308!!st         IF( ln_traqsr  )   CALL tra_qsr( kstp,      Kmm, ts, Krhs )  ! penetrative solar radiation qsr 
    299309         IF( ln_trabbc  )   CALL tra_bbc( kstp,      Kmm, ts, Krhs )  ! bottom heat flux 
    300310         IF( ln_trabbl  )   CALL tra_bbl( kstp, Kbb, Kmm, ts, Krhs )  ! advective (and/or diffusive) bottom boundary layer scheme 
     
    320330         !          
    321331      END SELECT       
    322  
    323332      !                                         !==  correction of the barotropic (all stages)  ==!    at Kaa = N+1/3, N+1/2 or N+1 
    324333      !                                                           ! barotropic velocity correction 
    325       zub(:,:) = uu_b(:,:,Kaa) - SUM( e3u_0(:,:,:)*uu(:,:,:,Kaa), 3 ) * r1_hu_0(:,:) 
    326       zvb(:,:) = vv_b(:,:,Kaa) - SUM( e3v_0(:,:,:)*vv(:,:,:,Kaa), 3 ) * r1_hv_0(:,:) 
     334      zub(A2D(0)) = uu_b(A2D(0),Kaa) - SUM( e3u_0(A2D(0),:)*uu(A2D(0),:,Kaa), 3 ) * r1_hu_0(A2D(0)) 
     335      zvb(A2D(0)) = vv_b(A2D(0),Kaa) - SUM( e3v_0(A2D(0),:)*vv(A2D(0),:,Kaa), 3 ) * r1_hv_0(A2D(0)) 
     336 
     337!!st      zub(:,:) = uu_b(:,:,Kaa) - SUM( e3u_0(:,:,:)*uu(:,:,:,Kaa), 3 ) * r1_hu_0(:,:) 
     338!!st      zvb(:,:) = vv_b(:,:,Kaa) - SUM( e3v_0(:,:,:)*vv(:,:,:,Kaa), 3 ) * r1_hv_0(:,:) 
     339 
    327340      DO jk = 1, jpkm1                                            ! corrected horizontal velocity 
    328341         uu(:,:,jk,Kaa) = uu(:,:,jk,Kaa) + zub(:,:)*umask(:,:,jk) 
    329342         vv(:,:,jk,Kaa) = vv(:,:,jk,Kaa) + zvb(:,:)*vmask(:,:,jk) 
    330343      END DO 
     344 
     345!!st      DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
     346!!st         uu(ji,jj,jk,Kaa) = uu(ji,jj,jk,Kaa) + zub(ji,jj)*umask(ji,jj,jk) 
     347!!st         vv(ji,jj,jk,Kaa) = vv(ji,jj,jk,Kaa) + zvb(ji,jj)*vmask(ji,jj,jk) 
     348!!st      END_3D 
     349          
    331350 
    332351      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
Note: See TracChangeset for help on using the changeset viewer.