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

Changeset 14618


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

#2506 RK3 work in progess : for 2 configurag

Location:
NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE
Files:
8 edited

Legend:

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

    r14072 r14618  
    9494      ! ------------------------- ! 
    9595      z_frc_trd_v = r1_rho0 * glob_sum( 'diahsb', - ( emp(:,:) - rnf(:,:) + fwfisf_cav(:,:) + fwfisf_par(:,:) ) * surf(:,:) )   ! volume fluxes 
     96#if defined key_RK3 
     97      CALL ctl_stop( 'dia_hsb: not yet instrumented for RK3' ) 
     98#else 
    9699      z_frc_trd_t =           glob_sum( 'diahsb', sbc_tsc(:,:,jp_tem) * surf(:,:) )                       ! heat fluxes 
    97100      z_frc_trd_s =           glob_sum( 'diahsb', sbc_tsc(:,:,jp_sal) * surf(:,:) )                       ! salt fluxes 
     101#endif 
    98102      !                    !  Add runoff    heat & salt input 
    99103      IF( ln_rnf    )   z_frc_trd_t = z_frc_trd_t + glob_sum( 'diahsb', rnf_tsc(:,:,jp_tem) * surf(:,:) ) 
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/DYN/sshwzv.F90

    r14547 r14618  
    216216!               &                                       - e3t(:,:,jk,Kbb)  )   ) * tmask(:,:,jk) 
    217217!!gm slightly faster : 
     218#if defined key_qco 
    218219            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk)                            & 
    219220               &                            + r1_Dt * e3t_0(:,:,jk) * ( r3t(:,:,Kaa) - r3t(:,:,Kbb) )  ) * tmask(:,:,jk) 
    220 !!gm end 
     221#else 
     222            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk)   & 
     223               &                            + r1_Dt * (  e3t(:,:,jk,Kaa)      & 
     224               &                                       - e3t(:,:,jk,Kbb)  )   ) * tmask(:,:,jk) 
     225#endif 
     226!!gm end !!st need to be there when not using qco  
    221227 
    222228 
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/TRA/traatf_qco.F90

    r14072 r14618  
    22   !!====================================================================== 
    33   !!                       ***  MODULE  traatf_qco  *** 
    4    !! Ocean active tracers:  Asselin time filtering for temperature and salinity 
     4   !! Ocean active tracers:  MLF, Asselin time filtering for temperature and salinity 
    55   !!====================================================================== 
    66   !! History :  OPA  !  1991-11  (G. Madec)  Original code 
     
    1717   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  semi-implicit hpg with asselin filter + modified LF-RA 
    1818   !!             -   !  2010-05  (C. Ethe, G. Madec)  merge TRC-TRA 
    19    !!            4.1  !  2019-08  (A. Coward, D. Storkey) rename tranxt.F90 -> traatfLF.F90. Now only does time filtering. 
    20    !!---------------------------------------------------------------------- 
    21  
     19   !!            4.1  !  2019-08  (A. Coward, D. Storkey) rename tranxt.F90 -> traatf.F90. Now only does time filtering. 
     20   !!            4.2  !  2020-06  (S. Techene, G. Madec) qco version of traatf.F90 
     21   !!---------------------------------------------------------------------- 
     22#if defined key_RK3 
     23   !!---------------------------------------------------------------------- 
     24   !!   'key_RK3'             EMPTY MODULE            3rd order Runge-Kutta  
     25   !!---------------------------------------------------------------------- 
     26#else 
    2227   !!---------------------------------------------------------------------- 
    2328   !!   tra_atf       : time filtering on tracers 
     
    2530   !!   tra_atf_vvl   : time filtering on tracers : variable volume case 
    2631   !!---------------------------------------------------------------------- 
    27    USE oce             ! ocean dynamics and tracers variables 
    28    USE dom_oce         ! ocean space and time domain variables 
    29    USE sbc_oce         ! surface boundary condition: ocean 
    30    USE sbcrnf          ! river runoffs 
    31    USE isf_oce         ! ice shelf melting 
    32    USE zdf_oce         ! ocean vertical mixing 
    33    USE domvvl          ! variable volume 
    34    USE trd_oce         ! trends: ocean variables 
    35    USE trdtra          ! trends manager: tracers 
    36    USE traqsr          ! penetrative solar radiation (needed for nksr) 
    37    USE phycst          ! physical constant 
    38    USE ldftra          ! lateral physics : tracers 
    39    USE ldfslp          ! lateral physics : slopes 
    40    USE bdy_oce  , ONLY : ln_bdy 
    41    USE bdytra          ! open boundary condition (bdy_tra routine) 
     32   USE oce            ! ocean dynamics and tracers variables 
     33   USE dom_oce        ! ocean space and time domain variables 
     34   USE sbc_oce        ! surface boundary condition: ocean 
     35   USE sbcrnf         ! river runoffs 
     36   USE isf_oce        ! ice shelf melting 
     37   USE zdf_oce        ! ocean vertical mixing 
     38   USE domvvl         ! variable volume 
     39   USE trd_oce        ! trends: ocean variables 
     40   USE trdtra         ! trends manager: tracers 
     41   USE traqsr         ! penetrative solar radiation (needed for nksr) 
     42   USE phycst         ! physical constant 
     43   USE ldftra         ! lateral physics : tracers 
     44   USE ldfslp         ! lateral physics : slopes 
     45   USE bdy_oce  , ONLY: ln_bdy 
     46   USE bdytra         ! open boundary condition (bdy_tra routine) 
    4247   ! 
    43    USE in_out_manager  ! I/O manager 
    44    USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    45    USE prtctl          ! Print control 
    46    USE timing          ! Timing 
     48   USE in_out_manager ! I/O manager 
     49   USE lbclnk         ! ocean lateral boundary conditions (or mpp link) 
     50   USE prtctl         ! Print control 
     51   USE timing         ! Timing 
    4752 
    4853   IMPLICIT NONE 
     
    365370      ! 
    366371   END SUBROUTINE tra_atf_qco_lf 
    367  
     372#endif 
     373    
    368374   !!====================================================================== 
    369375END MODULE traatf_qco 
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/TRA/traqsr.F90

    r14215 r14618  
    7777CONTAINS 
    7878 
    79    SUBROUTINE tra_qsr( kt, Kmm, pts, Krhs ) 
     79   SUBROUTINE tra_qsr( kt, Kmm, pts, Krhs, kstg ) 
    8080      !!---------------------------------------------------------------------- 
    8181      !!                  ***  ROUTINE tra_qsr  *** 
     
    103103      !!              Morel, A. et Berthon, JF, 1989, Limnol Oceanogr 34(8), 1545-1562 
    104104      !!---------------------------------------------------------------------- 
    105       INTEGER,                                   INTENT(in   ) :: kt            ! ocean time-step 
    106       INTEGER,                                   INTENT(in   ) :: Kmm, Krhs     ! time level indices 
    107       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts           ! active tracers and RHS of tracer equation 
     105      INTEGER,                                   INTENT(in   ) ::   kt, Kmm, Krhs   ! ocean time-step and time level indices 
     106      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) ::   pts             ! active tracers and RHS of tracer equation 
     107      INTEGER , OPTIONAL                       , INTENT(in   ) ::   kstg            ! RK3 stage index 
    108108      ! 
    109109      INTEGER  ::   ji, jj, jk               ! dummy loop indices 
    110110      INTEGER  ::   irgb, isi, iei, isj, iej ! local integers 
     111      INTEGER  ::   istg_1, istg_3           !   -       - 
    111112      REAL(wp) ::   zchl, zcoef, z1_2        ! local scalars 
    112       REAL(wp) ::   zc0 , zc1 , zc2 , zc3    !    -         - 
    113       REAL(wp) ::   zzc0, zzc1, zzc2, zzc3   !    -         - 
    114       REAL(wp) ::   zz0 , zz1 , ze3t, zlui   !    -         - 
     113      REAL(wp) ::   zc0 , zc1 , zc2 , zc3    !   -      - 
     114      REAL(wp) ::   zzc0, zzc1, zzc2, zzc3   !   -      - 
     115      REAL(wp) ::   zz0 , zz1 , ze3t, zlui   !   -      - 
    115116      REAL(wp) ::   zCb, zCmax, zpsi, zpsimax, zrdpsi, zCze 
    116117      REAL(wp) ::   zlogc, zlogze, zlogCtot, zlogCze 
    117       REAL(wp), ALLOCATABLE, DIMENSION(:,:)   :: ze0, ze1, ze2, ze3 
    118       REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: ztrdt, zetot, ztmp3d 
     118      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   ze0, ze1, ze2, ze3 
     119      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ztrdt, zetot, ztmp3d 
    119120      !!---------------------------------------------------------------------- 
    120121      ! 
    121122      IF( ln_timing )   CALL timing_start('tra_qsr') 
     123      ! 
     124      IF( PRESENT( kstg ) ) THEN      ! RK3 : a few things have to be done at only a specific stage 
     125         istg_1 = kstg   ;   istg_3 = kstg 
     126      ELSE                            ! MLF : only one call by time step 
     127         istg_1 =   1    ;   istg_3 =   3 
     128      ENDIF 
    122129      ! 
    123130      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     
    142149      IF( ntej == Nje0 ) THEN ; iej = nn_hls ; ELSE ; iej = 0 ; ENDIF 
    143150 
    144       IF( kt == nit000 ) THEN          !==  1st time step  ==! 
     151      IF( kt == nit000 .AND. istg_1 == 1 ) THEN   !==  1st time step  ==! (RK3: only at stage 1) 
    145152         IF( ln_rstart .AND. .NOT.l_1st_euler ) THEN    ! read in restart 
    146153            z1_2 = 0.5_wp 
     
    155162            END_3D 
    156163         ENDIF 
    157       ELSE                             !==  Swap of qsr heat content  ==! 
     164      ELSEIF( istg_3 == 3 ) THEN                  !==  Swap of qsr heat content  ==! 
    158165         z1_2 = 0.5_wp 
    159166         DO_3D( isi, iei, isj, iej, 1, jpk ) 
     
    207214               ze3(ji,jj) = EXP( - zc3 )                                          ! ze3 = 1/zze 
    208215            END_2D 
    209  
    210 ! 
     216            ! 
    211217            DO_3D( isi, iei, isj, iej, 1, nksr + 1 ) 
    212218               ! zchl    = ALOG( ze0(ji,jj) ) 
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/TRA/trasbc.F90

    r14215 r14618  
    5151CONTAINS 
    5252 
    53    SUBROUTINE tra_sbc ( kt, Kmm, pts, Krhs ) 
     53   SUBROUTINE tra_sbc ( kt, Kmm, pts, Krhs, kstg ) 
    5454      !!---------------------------------------------------------------------- 
    5555      !!                  ***  ROUTINE tra_sbc  *** 
     
    7272      !!              - send trends to trdtra module for further diagnostics(l_trdtra=T) 
    7373      !!---------------------------------------------------------------------- 
    74       INTEGER,                                   INTENT(in   ) ::   kt         ! ocean time-step index 
    75       INTEGER,                                   INTENT(in   ) ::   Kmm, Krhs  ! time level indices 
    76       REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) ::   pts        ! active tracers and RHS of tracer Eq. 
     74      INTEGER,                                   INTENT(in   ) ::   kt, Kmm, Krhs   ! ocean time-step and time-level indices 
     75      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) ::   pts             ! active tracers and RHS of tracer Eq. 
     76      INTEGER , OPTIONAL                       , INTENT(in   ) ::   kstg            ! RK3 stage index 
    7777      ! 
    7878      INTEGER  ::   ji, jj, jk, jn               ! dummy loop indices 
    79       INTEGER  ::   ikt, ikb, isi, iei, isj, iej ! local integers 
     79      INTEGER  ::   istg_1, istg_3               ! local integers 
     80      INTEGER  ::   ikt, ikb, isi, iei, isj, iej !   -       - 
    8081      REAL(wp) ::   zfact, z1_e3t, zdep, ztim    ! local scalar 
    8182      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::  ztrdt, ztrds 
     
    8384      ! 
    8485      IF( ln_timing )   CALL timing_start('tra_sbc') 
     86      ! 
     87      IF( PRESENT( kstg ) ) THEN      ! RK3 : a few things have to be done at only a specific stage 
     88         istg_1 = kstg   ;   istg_3 = kstg 
     89      ELSE                            ! MLF : only one call by time step 
     90         istg_1 =   1    ;   istg_3 =   3 
     91      ENDIF 
    8592      ! 
    8693      IF( ntile == 0 .OR. ntile == 1 )  THEN                       ! Do only on the first tile 
     
    104111 
    105112!!gm  This should be moved into sbcmod.F90 module ? (especially now that ln_traqsr is read in namsbc namelist) 
    106       IF( .NOT.ln_traqsr ) THEN     ! no solar radiation penetration 
     113      IF( .NOT.ln_traqsr .AND. istg_1 == 1 ) THEN     ! no solar radiation penetration (RK3: only at stage 1) 
    107114         DO_2D( isi, iei, isj, iej ) 
    108             qns(ji,jj) = qns(ji,jj) + qsr(ji,jj)      ! total heat flux in qns 
    109             qsr(ji,jj) = 0._wp                        ! qsr set to zero 
     115            qns(ji,jj) = qns(ji,jj) + qsr(ji,jj)         ! total heat flux in qns 
     116            qsr(ji,jj) = 0._wp                           ! qsr set to zero 
    110117         END_2D 
    111118      ENDIF 
     
    115122      !---------------------------------------- 
    116123      !                             !==  Set before sbc tracer content fields  ==! 
    117       IF( kt == nit000 ) THEN             !* 1st time-step 
     124      IF( kt == nit000 .AND. istg_1 == 1 ) THEN             !* 1st time-step 
    118125         IF( ln_rstart .AND. .NOT.l_1st_euler ) THEN      ! Restart: read in restart file 
    119126            zfact = 0.5_wp 
     
    131138            END_2D 
    132139         ENDIF 
    133       ELSE                                !* other time-steps: swap of forcing fields 
     140      ELSEIF( istg_3 == 3 ) THEN          !* other time-steps: swap of forcing fields (RK3: only at stage 3) 
    134141         zfact = 0.5_wp 
    135142         DO_2D( isi, iei, isj, iej ) 
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/stpmlf.F90

    r14547 r14618  
    3535   !!            4.x  !  2020-08  (S. Techene, G. Madec)  quasi eulerian coordinate time stepping 
    3636   !!---------------------------------------------------------------------- 
    37 #if defined key_qco   ||   defined key_linssh 
     37#if ! defined key_RK3 
     38# if defined key_qco   ||   defined key_linssh 
    3839   !!---------------------------------------------------------------------- 
    3940   !!   'key_qco'                        Quasi-Eulerian vertical coordinate 
     
    527528   END SUBROUTINE finalize_lbc 
    528529 
    529 #else 
     530# else 
    530531   !!---------------------------------------------------------------------- 
    531532   !!   default option             EMPTY MODULE           qco not activated 
    532533   !!---------------------------------------------------------------------- 
     534# endif 
    533535#endif 
    534536    
  • NEMO/branches/2021/dev_r14318_RK3_stage1/src/OCE/stprk3.F90

    r14548 r14618  
    210210       
    211211!==>>>  at Nbb  no more Nnn  
    212        
     212      
    213213      IF( ln_diacfl  )   CALL dia_cfl   ( kstp,      Nbb )      ! Courant number diagnostics 
    214214                         CALL dia_hth   ( kstp,      Nbb )      ! Thermocline depth (20 degres isotherm depth) 
    215215      IF( ln_diadct  )   CALL dia_dct   ( kstp,      Nbb )      ! Transports 
    216                          CALL dia_ar5   ( kstp,      Nbb )      ! ar5 diag 
     216!!st                         CALL dia_ar5   ( kstp,      Nbb )      ! ar5 diag 
    217217                         CALL dia_ptr   ( kstp,      Nbb )      ! Poleward adv/ldf TRansports diagnostics 
    218218                         CALL dia_wri   ( kstp,      Nbb )      ! ocean model: outputs 
     
    245245      ! Control 
    246246      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    247                          CALL stp_ctl      ( kstp, Nnn ) 
    248  
     247                         CALL stp_ctl      ( kstp, Nbb ) 
    249248#if defined key_agrif 
    250249      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    282281      ENDIF 
    283282#endif 
    284       ! 
    285       IF( l_1st_euler ) THEN         ! recover Leap-frog timestep 
    286          rDt   = 2._wp * rn_Dt 
    287          r1_Dt = 1._wp / rDt 
    288          l_1st_euler = .FALSE. 
    289       ENDIF 
    290283      ! 
    291284      IF( ln_timing )   CALL timing_stop('stp_RK3') 
  • 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.