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 12424 for NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src – NEMO

Ignore:
Timestamp:
2020-02-20T16:06:41+01:00 (4 years ago)
Author:
davestorkey
Message:
  1. Rename r2dt -> rDt
  2. Rename r1_2dt -> r1_Dt
  3. Reorganise management of initial Euler timestep for leapfrogging.

This version passes all SETTE tests and bit-compares with the trunk @ 12377

Location:
NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE
Files:
16 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DIA/diacfl.F90

    r12397 r12424  
    6060      ! 
    6161      DO_3D_11_11( 1, jpk ) 
    62          zCu_cfl(ji,jj,jk) = ABS( uu(ji,jj,jk,Kmm) ) * r2dt / e1u  (ji,jj)      ! for i-direction 
    63          zCv_cfl(ji,jj,jk) = ABS( vv(ji,jj,jk,Kmm) ) * r2dt / e2v  (ji,jj)      ! for j-direction 
    64          zCw_cfl(ji,jj,jk) = ABS( ww(ji,jj,jk) ) * r2dt / e3w(ji,jj,jk,Kmm)   ! for k-direction 
     62         zCu_cfl(ji,jj,jk) = ABS( uu(ji,jj,jk,Kmm) ) * rDt / e1u  (ji,jj)      ! for i-direction 
     63         zCv_cfl(ji,jj,jk) = ABS( vv(ji,jj,jk,Kmm) ) * rDt / e2v  (ji,jj)      ! for j-direction 
     64         zCw_cfl(ji,jj,jk) = ABS( ww(ji,jj,jk) ) * rDt / e3w(ji,jj,jk,Kmm)   ! for k-direction 
    6565      END_3D 
    6666      ! 
     
    112112         WRITE(numcfl,*) '******************************************' 
    113113         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cu', rCu_max, nCu_loc(1), nCu_loc(2), nCu_loc(3) 
    114          WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', r2dt/rCu_max 
     114         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', rDt/rCu_max 
    115115         WRITE(numcfl,*) '******************************************' 
    116116         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cv', rCv_max, nCv_loc(1), nCv_loc(2), nCv_loc(3) 
    117          WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', r2dt/rCv_max 
     117         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', rDt/rCv_max 
    118118         WRITE(numcfl,*) '******************************************' 
    119119         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cw', rCw_max, nCw_loc(1), nCw_loc(2), nCw_loc(3) 
    120          WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', r2dt/rCw_max 
     120         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', rDt/rCw_max 
    121121         CLOSE( numcfl )  
    122122         ! 
     
    125125         WRITE(numout,*) 'dia_cfl : Maximum Courant number information for the run ' 
    126126         WRITE(numout,*) '~~~~~~~' 
    127          WRITE(numout,*) '   Max Cu = ', rCu_max, ' at (i,j,k) = (',nCu_loc(1),nCu_loc(2),nCu_loc(3),') => dt/C = ', r2dt/rCu_max 
    128          WRITE(numout,*) '   Max Cv = ', rCv_max, ' at (i,j,k) = (',nCv_loc(1),nCv_loc(2),nCv_loc(3),') => dt/C = ', r2dt/rCv_max 
    129          WRITE(numout,*) '   Max Cw = ', rCw_max, ' at (i,j,k) = (',nCw_loc(1),nCw_loc(2),nCw_loc(3),') => dt/C = ', r2dt/rCw_max 
     127         WRITE(numout,*) '   Max Cu = ', rCu_max, ' at (i,j,k) = (',nCu_loc(1),nCu_loc(2),nCu_loc(3),') => dt/C = ', rDt/rCu_max 
     128         WRITE(numout,*) '   Max Cv = ', rCv_max, ' at (i,j,k) = (',nCv_loc(1),nCv_loc(2),nCv_loc(3),') => dt/C = ', rDt/rCv_max 
     129         WRITE(numout,*) '   Max Cw = ', rCw_max, ' at (i,j,k) = (',nCw_loc(1),nCw_loc(2),nCw_loc(3),') => dt/C = ', rDt/rCw_max 
    130130      ENDIF 
    131131      ! 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DOM/dom_oce.F90

    r12406 r12424  
    5959   !                                   !!! associated variables 
    6060   LOGICAL , PUBLIC ::   l_1st_euler    !: Euler 1st time-step flag (=T if ln_restart=F or ln_1st_euler=T) 
    61    REAL(wp), PUBLIC ::   r2dt, r1_2dt   !: = 2*rdt and 1/(2*rdt) except at nit000 (=rdt and 1/rdt) if l_1st_euler=.true. 
     61   REAL(wp), PUBLIC ::   rDt, r1_Dt     !: Current model timestep and reciprocal 
     62                                        !:  = 2 * rn_Dt if leapfrog and l_1st_euler = F 
     63                                        !:  =     rn_Dt if leapfrog and l_1st_euler = T 
     64                                        !:  =     rn_Dt if RK3 
    6265 
    6366   !!---------------------------------------------------------------------- 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DOM/domain.F90

    r12406 r12424  
    415415      !          ! conversion DOCTOR names into model names (this should disappear soon) 
    416416      atfp = rn_atfp 
    417 !!      Initialise current model timestep rDt = 2*rn_Dt if MLF or rDt = rn_Dt if RK3 
    418 !!      rDt  = rn_Dt 
     417      !! Initialise current model timestep rDt = 2*rn_Dt if MLF or rDt = rn_Dt if RK3 
     418      rDt  = 2._wp * rn_Dt 
     419      r1_Dt = 1._wp / rDt 
    419420 
    420421      IF( TRIM(Agrif_CFixed()) == '0' ) THEN 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DOM/domvvl.F90

    r12406 r12424  
    431431         ! --------------------------------------------- 
    432432         CALL lbc_lnk( 'domvvl', tilde_e3t_a(:,:,:), 'T', 1._wp ) 
    433          tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + r2dt * tmask(:,:,:) * tilde_e3t_a(:,:,:) 
     433         tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + rDt * tmask(:,:,:) * tilde_e3t_a(:,:,:) 
    434434 
    435435         ! Maximum deformation control 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DYN/dynatf.F90

    r12406 r12424  
    162162         ! 
    163163         IF( ln_dyn_trd ) THEN              ! 3D output: total momentum trends 
    164             zua(:,:,:) = ( puu(:,:,:,Kaa) - puu(:,:,:,Kbb) ) * r1_2dt 
    165             zva(:,:,:) = ( pvv(:,:,:,Kaa) - pvv(:,:,:,Kbb) ) * r1_2dt 
     164            zua(:,:,:) = ( puu(:,:,:,Kaa) - puu(:,:,:,Kbb) ) * r1_Dt 
     165            zva(:,:,:) = ( pvv(:,:,:,Kaa) - pvv(:,:,:,Kbb) ) * r1_Dt 
    166166            CALL iom_put( "utrd_tot", zua )        ! total momentum trends, except the asselin time filter 
    167167            CALL iom_put( "vtrd_tot", zva ) 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DYN/dynspg_ts.F90

    r12406 r12424  
    150150      LOGICAL  ::   ll_init               ! =T : special startup of 2d equations 
    151151      INTEGER  ::   noffset               ! local integers  : time offset for bdy update 
    152       REAL(wp) ::   r1_2dt_b, z1_hu, z1_hv          ! local scalars 
     152      REAL(wp) ::   r1_Dt_b, z1_hu, z1_hv          ! local scalars 
    153153      REAL(wp) ::   za0, za1, za2, za3              !   -      - 
    154154      REAL(wp) ::   zztmp, zldg               !   -      - 
     
    180180!     zwdramp = 1._wp / (rn_wdmin2 - rn_wdmin1) ! more general ramp 
    181181      !                                         ! inverse of baroclinic time step  
    182       r1_2dt_b = 1._wp / r2dt  
     182      r1_Dt_b = 1._wp / rDt  
    183183      ! 
    184184      ll_init     = ln_bt_av                    ! if no time averaging, then no specific restart  
     
    729729      IF( ln_dynadv_vec .OR. ln_linssh ) THEN 
    730730         DO jk=1,jpkm1 
    731             puu(:,:,jk,Krhs) = puu(:,:,jk,Krhs) + ( puu_b(:,:,Kaa) - puu_b(:,:,Kbb) ) * r1_2dt_b 
    732             pvv(:,:,jk,Krhs) = pvv(:,:,jk,Krhs) + ( pvv_b(:,:,Kaa) - pvv_b(:,:,Kbb) ) * r1_2dt_b 
     731            puu(:,:,jk,Krhs) = puu(:,:,jk,Krhs) + ( puu_b(:,:,Kaa) - puu_b(:,:,Kbb) ) * r1_Dt_b 
     732            pvv(:,:,jk,Krhs) = pvv(:,:,jk,Krhs) + ( pvv_b(:,:,Kaa) - pvv_b(:,:,Kbb) ) * r1_Dt_b 
    733733         END DO 
    734734      ELSE 
     
    745745         ! 
    746746         DO jk=1,jpkm1 
    747             puu(:,:,jk,Krhs) = puu(:,:,jk,Krhs) + r1_hu(:,:,Kmm) * ( puu_b(:,:,Kaa) - puu_b(:,:,Kbb) * hu(:,:,Kbb) ) * r1_2dt_b 
    748             pvv(:,:,jk,Krhs) = pvv(:,:,jk,Krhs) + r1_hv(:,:,Kmm) * ( pvv_b(:,:,Kaa) - pvv_b(:,:,Kbb) * hv(:,:,Kbb) ) * r1_2dt_b 
     747            puu(:,:,jk,Krhs) = puu(:,:,jk,Krhs) + r1_hu(:,:,Kmm) * ( puu_b(:,:,Kaa) - puu_b(:,:,Kbb) * hu(:,:,Kbb) ) * r1_Dt_b 
     748            pvv(:,:,jk,Krhs) = pvv(:,:,jk,Krhs) + r1_hv(:,:,Kmm) * ( pvv_b(:,:,Kaa) - pvv_b(:,:,Kbb) * hv(:,:,Kbb) ) * r1_Dt_b 
    749749         END DO 
    750750         ! Save barotropic velocities not transport: 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DYN/dynzdf.F90

    r12397 r12424  
    107107      IF( ln_dynadv_vec .OR. ln_linssh ) THEN   ! applied on velocity 
    108108         DO jk = 1, jpkm1 
    109             puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kbb) + r2dt * puu(:,:,jk,Krhs) ) * umask(:,:,jk) 
    110             pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kbb) + r2dt * pvv(:,:,jk,Krhs) ) * vmask(:,:,jk) 
     109            puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kbb) + rDt * puu(:,:,jk,Krhs) ) * umask(:,:,jk) 
     110            pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kbb) + rDt * pvv(:,:,jk,Krhs) ) * vmask(:,:,jk) 
    111111         END DO 
    112112      ELSE                                      ! applied on thickness weighted velocity 
    113113         DO jk = 1, jpkm1 
    114114            puu(:,:,jk,Kaa) = (         e3u(:,:,jk,Kbb) * puu(:,:,jk,Kbb)  & 
    115                &          + r2dt * e3u(:,:,jk,Kmm) * puu(:,:,jk,Krhs)  ) / e3u(:,:,jk,Kaa) * umask(:,:,jk) 
     115               &          + rDt * e3u(:,:,jk,Kmm) * puu(:,:,jk,Krhs)  ) / e3u(:,:,jk,Kaa) * umask(:,:,jk) 
    116116            pvv(:,:,jk,Kaa) = (         e3v(:,:,jk,Kbb) * pvv(:,:,jk,Kbb)  & 
    117                &          + r2dt * e3v(:,:,jk,Kmm) * pvv(:,:,jk,Krhs)  ) / e3v(:,:,jk,Kaa) * vmask(:,:,jk) 
     117               &          + rDt * e3v(:,:,jk,Kmm) * pvv(:,:,jk,Krhs)  ) / e3v(:,:,jk,Kaa) * vmask(:,:,jk) 
    118118         END DO 
    119119      ENDIF 
     
    133133            ze3ua =  ( 1._wp - r_vvl ) * e3u(ji,jj,iku,Kmm) + r_vvl * e3u(ji,jj,iku,Kaa) 
    134134            ze3va =  ( 1._wp - r_vvl ) * e3v(ji,jj,ikv,Kmm) + r_vvl * e3v(ji,jj,ikv,Kaa) 
    135             puu(ji,jj,iku,Kaa) = puu(ji,jj,iku,Kaa) + r2dt * 0.5*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) * uu_b(ji,jj,Kaa) / ze3ua 
    136             pvv(ji,jj,ikv,Kaa) = pvv(ji,jj,ikv,Kaa) + r2dt * 0.5*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) * vv_b(ji,jj,Kaa) / ze3va 
     135            puu(ji,jj,iku,Kaa) = puu(ji,jj,iku,Kaa) + rDt * 0.5*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) * uu_b(ji,jj,Kaa) / ze3ua 
     136            pvv(ji,jj,ikv,Kaa) = pvv(ji,jj,ikv,Kaa) + rDt * 0.5*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) * vv_b(ji,jj,Kaa) / ze3va 
    137137         END_2D 
    138138         IF( ln_isfcav ) THEN    ! Ocean cavities (ISF) 
     
    142142               ze3ua =  ( 1._wp - r_vvl ) * e3u(ji,jj,iku,Kmm) + r_vvl * e3u(ji,jj,iku,Kaa) 
    143143               ze3va =  ( 1._wp - r_vvl ) * e3v(ji,jj,ikv,Kmm) + r_vvl * e3v(ji,jj,ikv,Kaa) 
    144                puu(ji,jj,iku,Kaa) = puu(ji,jj,iku,Kaa) + r2dt * 0.5*( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) * uu_b(ji,jj,Kaa) / ze3ua 
    145                pvv(ji,jj,ikv,Kaa) = pvv(ji,jj,ikv,Kaa) + r2dt * 0.5*( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) * vv_b(ji,jj,Kaa) / ze3va 
     144               puu(ji,jj,iku,Kaa) = puu(ji,jj,iku,Kaa) + rDt * 0.5*( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) * uu_b(ji,jj,Kaa) / ze3ua 
     145               pvv(ji,jj,ikv,Kaa) = pvv(ji,jj,ikv,Kaa) + rDt * 0.5*( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) * vv_b(ji,jj,Kaa) / ze3va 
    146146            END_2D 
    147147         END IF 
     
    151151      ! 
    152152      !                    !* Matrix construction 
    153       zdt = r2dt * 0.5 
     153      zdt = rDt * 0.5 
    154154      IF( ln_zad_Aimp ) THEN   !! 
    155155         SELECT CASE( nldf_dyn ) 
     
    227227            iku = mbku(ji,jj)       ! ocean bottom level at u- and v-points 
    228228            ze3ua =  ( 1._wp - r_vvl ) * e3u(ji,jj,iku,Kmm) + r_vvl * e3u(ji,jj,iku,Kaa)   ! after scale factor at T-point 
    229             zwd(ji,jj,iku) = zwd(ji,jj,iku) - r2dt * 0.5*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) / ze3ua 
     229            zwd(ji,jj,iku) = zwd(ji,jj,iku) - rDt * 0.5*( rCdU_bot(ji+1,jj)+rCdU_bot(ji,jj) ) / ze3ua 
    230230         END_2D 
    231231         IF ( ln_isfcav ) THEN   ! top friction (always implicit) 
     
    234234               iku = miku(ji,jj)       ! ocean top level at u- and v-points  
    235235               ze3ua =  ( 1._wp - r_vvl ) * e3u(ji,jj,iku,Kmm) + r_vvl * e3u(ji,jj,iku,Kaa)   ! after scale factor at T-point 
    236                zwd(ji,jj,iku) = zwd(ji,jj,iku) - r2dt * 0.5*( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) / ze3ua 
     236               zwd(ji,jj,iku) = zwd(ji,jj,iku) - rDt * 0.5*( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) / ze3ua 
    237237            END_2D 
    238238         END IF 
     
    260260      DO_2D_00_00 
    261261         ze3ua =  ( 1._wp - r_vvl ) * e3u(ji,jj,1,Kmm) + r_vvl * e3u(ji,jj,1,Kaa)  
    262          puu(ji,jj,1,Kaa) = puu(ji,jj,1,Kaa) + r2dt * 0.5_wp * ( utau_b(ji,jj) + utau(ji,jj) )   & 
     262         puu(ji,jj,1,Kaa) = puu(ji,jj,1,Kaa) + rDt * 0.5_wp * ( utau_b(ji,jj) + utau(ji,jj) )   & 
    263263            &                                      / ( ze3ua * rau0 ) * umask(ji,jj,1)  
    264264      END_2D 
     
    277277      ! 
    278278      !                       !* Matrix construction 
    279       zdt = r2dt * 0.5 
     279      zdt = rDt * 0.5 
    280280      IF( ln_zad_Aimp ) THEN   !! 
    281281         SELECT CASE( nldf_dyn ) 
     
    352352            ikv = mbkv(ji,jj)       ! (deepest ocean u- and v-points) 
    353353            ze3va =  ( 1._wp - r_vvl ) * e3v(ji,jj,ikv,Kmm) + r_vvl * e3v(ji,jj,ikv,Kaa)   ! after scale factor at T-point 
    354             zwd(ji,jj,ikv) = zwd(ji,jj,ikv) - r2dt * 0.5*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) / ze3va            
     354            zwd(ji,jj,ikv) = zwd(ji,jj,ikv) - rDt * 0.5*( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj) ) / ze3va            
    355355         END_2D 
    356356         IF ( ln_isfcav ) THEN 
     
    358358               ikv = mikv(ji,jj)       ! (first wet ocean u- and v-points) 
    359359               ze3va =  ( 1._wp - r_vvl ) * e3v(ji,jj,ikv,Kmm) + r_vvl * e3v(ji,jj,ikv,Kaa)   ! after scale factor at T-point 
    360                zwd(ji,jj,ikv) = zwd(ji,jj,ikv) - r2dt * 0.5*( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) / ze3va 
     360               zwd(ji,jj,ikv) = zwd(ji,jj,ikv) - rDt * 0.5*( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) / ze3va 
    361361            END_2D 
    362362         ENDIF 
     
    384384      DO_2D_00_00 
    385385         ze3va =  ( 1._wp - r_vvl ) * e3v(ji,jj,1,Kmm) + r_vvl * e3v(ji,jj,1,Kaa)  
    386          pvv(ji,jj,1,Kaa) = pvv(ji,jj,1,Kaa) + r2dt * 0.5_wp * ( vtau_b(ji,jj) + vtau(ji,jj) )   & 
     386         pvv(ji,jj,1,Kaa) = pvv(ji,jj,1,Kaa) + rDt * 0.5_wp * ( vtau_b(ji,jj) + vtau(ji,jj) )   & 
    387387            &                                      / ( ze3va * rau0 ) * vmask(ji,jj,1)  
    388388      END_2D 
     
    399399      ! 
    400400      IF( l_trddyn )   THEN                      ! save the vertical diffusive trends for further diagnostics 
    401          ztrdu(:,:,:) = ( puu(:,:,:,Kaa) - puu(:,:,:,Kbb) ) / r2dt - ztrdu(:,:,:) 
    402          ztrdv(:,:,:) = ( pvv(:,:,:,Kaa) - pvv(:,:,:,Kbb) ) / r2dt - ztrdv(:,:,:) 
     401         ztrdu(:,:,:) = ( puu(:,:,:,Kaa) - puu(:,:,:,Kbb) ) / rDt - ztrdu(:,:,:) 
     402         ztrdv(:,:,:) = ( pvv(:,:,:,Kaa) - pvv(:,:,:,Kbb) ) / rDt - ztrdv(:,:,:) 
    403403         CALL trd_dyn( ztrdu, ztrdv, jpdyn_zdf, kt, Kmm ) 
    404404         DEALLOCATE( ztrdu, ztrdv )  
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/DYN/sshwzv.F90

    r12406 r12424  
    9494      !                                           !------------------------------! 
    9595      IF(ln_wd_il) THEN 
    96          CALL wad_lmt(pssh(:,:,Kbb), zcoef * (emp_b(:,:) + emp(:,:)), r2dt, Kmm, uu, vv ) 
     96         CALL wad_lmt(pssh(:,:,Kbb), zcoef * (emp_b(:,:) + emp(:,:)), rDt, Kmm, uu, vv ) 
    9797      ENDIF 
    9898 
     
    107107      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations. 
    108108      !  
    109       pssh(:,:,Kaa) = (  pssh(:,:,Kbb) - r2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:) 
     109      pssh(:,:,Kaa) = (  pssh(:,:,Kbb) - rDt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:) 
    110110      ! 
    111111#if defined key_agrif 
     
    182182            ! computation of w 
    183183            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk) + zhdiv(:,:,jk)    & 
    184                &                         + r1_2dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) )     ) * tmask(:,:,jk) 
     184               &                         + r1_Dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) )     ) * tmask(:,:,jk) 
    185185         END DO 
    186186         !          IF( ln_vvl_layer ) pww(:,:,:) = 0.e0 
     
    190190            ! computation of w 
    191191            pww(:,:,jk) = pww(:,:,jk+1) - (  e3t(:,:,jk,Kmm) * hdiv(:,:,jk)                 & 
    192                &                         + r1_2dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) )  ) * tmask(:,:,jk) 
     192               &                         + r1_Dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) )  ) * tmask(:,:,jk) 
    193193         END DO 
    194194      ENDIF 
     
    306306         DO_3D_00_00( 1, jpkm1 ) 
    307307            z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm) 
    308             ! 2*rn_Dt and not r2dt (for restartability) 
     308            ! 2*rn_Dt and not rDt (for restartability) 
    309309            Cu_adv(ji,jj,jk) = 2._wp * rn_Dt * ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) )                       &   
    310310               &                             + ( MAX( e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kmm)*uu(ji  ,jj,jk,Kmm) + un_td(ji  ,jj,jk), 0._wp ) -   & 
     
    319319         DO_3D_00_00( 1, jpkm1 ) 
    320320            z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm) 
    321             ! 2*rn_Dt and not r2dt (for restartability) 
     321            ! 2*rn_Dt and not rDt (for restartability) 
    322322            Cu_adv(ji,jj,jk) = 2._wp * rn_Dt * ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) )   &  
    323323               &                             + ( MAX( e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kmm)*uu(ji  ,jj,jk,Kmm), 0._wp ) -   & 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRA/traadv.F90

    r12397 r12424  
    144144         CALL tra_adv_cen    ( kt, nit000, 'TRA',         zuu, zvv, zww, Kmm, pts, jpts, Krhs, nn_cen_h, nn_cen_v ) 
    145145      CASE ( np_FCT )                                 ! FCT scheme      : 2nd / 4th order 
    146          CALL tra_adv_fct    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, nn_fct_h, nn_fct_v ) 
     146         CALL tra_adv_fct    ( kt, nit000, 'TRA', rDt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, nn_fct_h, nn_fct_v ) 
    147147      CASE ( np_MUS )                                 ! MUSCL 
    148          CALL tra_adv_mus    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, ln_mus_ups )  
     148         CALL tra_adv_mus    ( kt, nit000, 'TRA', rDt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, ln_mus_ups )  
    149149      CASE ( np_UBS )                                 ! UBS 
    150          CALL tra_adv_ubs    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, nn_ubs_v   ) 
     150         CALL tra_adv_ubs    ( kt, nit000, 'TRA', rDt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs, nn_ubs_v   ) 
    151151      CASE ( np_QCK )                                 ! QUICKEST 
    152          CALL tra_adv_qck    ( kt, nit000, 'TRA', r2dt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs ) 
     152         CALL tra_adv_qck    ( kt, nit000, 'TRA', rDt, zuu, zvv, zww, Kbb, Kmm, pts, jpts, Krhs ) 
    153153      ! 
    154154      END SELECT 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRA/traatf.F90

    r12406 r12424  
    162162      ! 
    163163      IF( l_trdtra .AND. ln_linssh ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt      
    164          zfact = 1._wp / r2dt              
     164         zfact = 1._wp / rDt              
    165165         DO jk = 1, jpkm1 
    166166            ztrdt(:,:,jk) = ( pts(:,:,jk,jp_tem,Kmm) - ztrdt(:,:,jk) ) * zfact 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRA/traldf_iso.F90

    r12397 r12424  
    173173               DO_3D_10_10( 2, jpkm1 ) 
    174174                  ze3w_2 = e3w(ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm) 
    175                   zcoef0 = r2dt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  ) 
    176                   akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * r1_2dt 
     175                  zcoef0 = rDt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  ) 
     176                  akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * r1_Dt 
    177177               END_3D 
    178178           ENDIF 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRA/traldf_triad.F90

    r12397 r12424  
    185185               DO_3D_10_10( 2, jpkm1 ) 
    186186                  ze3w_2 = e3w(ji,jj,jk,Kmm) * e3w(ji,jj,jk,Kmm) 
    187                   zcoef0 = r2dt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  ) 
    188                   akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * r1_2dt 
     187                  zcoef0 = rDt * (  akz(ji,jj,jk) + ah_wslp2(ji,jj,jk) / ze3w_2  ) 
     188                  akz(ji,jj,jk) = MAX( zcoef0 - 0.5_wp , 0._wp ) * ze3w_2 * r1_Dt 
    189189               END_3D 
    190190           ENDIF 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRA/tranpc.F90

    r12406 r12424  
    6767      LOGICAL  ::   l_bottom_reached, l_column_treated 
    6868      REAL(wp) ::   zta, zalfa, zsum_temp, zsum_alfa, zaw, zdz, zsum_z 
    69       REAL(wp) ::   zsa, zbeta, zsum_sali, zsum_beta, zbw, zrw, z1_r2dt 
     69      REAL(wp) ::   zsa, zbeta, zsum_sali, zsum_beta, zbw, zrw, z1_rDt 
    7070      REAL(wp), PARAMETER ::   zn2_zero = 1.e-14_wp             ! acceptance criteria for neutrality (N2==0) 
    7171      REAL(wp), DIMENSION(        jpk     )   ::   zvn2         ! vertical profile of N2 at 1 given point... 
     
    301301         ! 
    302302         IF( l_trdtra ) THEN         ! send the Non penetrative mixing trends for diagnostic 
    303             z1_r2dt = 1._wp / (2._wp * rn_Dt) 
    304             ztrdt(:,:,:) = ( pts(:,:,:,jp_tem,Kaa) - ztrdt(:,:,:) ) * z1_r2dt 
    305             ztrds(:,:,:) = ( pts(:,:,:,jp_sal,Kaa) - ztrds(:,:,:) ) * z1_r2dt 
     303            z1_rDt = 1._wp / (2._wp * rn_Dt) 
     304            ztrdt(:,:,:) = ( pts(:,:,:,jp_tem,Kaa) - ztrdt(:,:,:) ) * z1_rDt 
     305            ztrds(:,:,:) = ( pts(:,:,:,jp_sal,Kaa) - ztrds(:,:,:) ) * z1_rDt 
    306306            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_npc, ztrdt ) 
    307307            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_npc, ztrds ) 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRA/trazdf.F90

    r12397 r12424  
    7373      ! 
    7474      !                                      !* compute lateral mixing trend and add it to the general trend 
    75       CALL tra_zdf_imp( kt, nit000, 'TRA', r2dt, Kbb, Kmm, Krhs, pts, Kaa, jpts )  
     75      CALL tra_zdf_imp( kt, nit000, 'TRA', rDt, Kbb, Kmm, Krhs, pts, Kaa, jpts )  
    7676 
    7777!!gm WHY here !   and I don't like that ! 
     
    8585         DO jk = 1, jpkm1 
    8686            ztrdt(:,:,jk) = ( ( pts(:,:,jk,jp_tem,Kaa)*e3t(:,:,jk,Kaa) - pts(:,:,jk,jp_tem,Kbb)*e3t(:,:,jk,Kbb) ) & 
    87                &          / (e3t(:,:,jk,Kmm)*r2dt) ) - ztrdt(:,:,jk) 
     87               &          / (e3t(:,:,jk,Kmm)*rDt) ) - ztrdt(:,:,jk) 
    8888            ztrds(:,:,jk) = ( ( pts(:,:,jk,jp_sal,Kaa)*e3t(:,:,jk,Kaa) - pts(:,:,jk,jp_sal,Kbb)*e3t(:,:,jk,Kbb) ) & 
    89               &           / (e3t(:,:,jk,Kmm)*r2dt) ) - ztrds(:,:,jk) 
     89              &           / (e3t(:,:,jk,Kmm)*rDt) ) - ztrds(:,:,jk) 
    9090         END DO 
    9191!!gm this should be moved in trdtra.F90 and done on all trends 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/TRD/trdtra.F90

    r12397 r12424  
    244244 
    245245      !                   ! Potential ENergy trends 
    246       IF( ln_PE_trd  )   CALL trd_pen( ptrdx, ptrdy, ktrd, kt, r2dt, Kmm ) 
     246      IF( ln_PE_trd  )   CALL trd_pen( ptrdx, ptrdy, ktrd, kt, rDt, Kmm ) 
    247247 
    248248      !                   ! Mixed layer trends for active tracers 
     
    277277         CASE ( jptra_atf )        ;   CALL trd_mxl_zint( ptrdx, ptrdy, jpmxl_atf, '3D' )   ! asselin time filter (last trend) 
    278278                                   ! 
    279                                        CALL trd_mxl( kt, r2dt )                             ! trends: Mixed-layer (output) 
     279                                       CALL trd_mxl( kt, rDt )                             ! trends: Mixed-layer (output) 
    280280         END SELECT 
    281281         ! 
  • NEMO/branches/2020/KERNEL-03_Storkey_Coward_RK3_stage2/src/OCE/step.F90

    r12406 r12424  
    105105      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    106106      ! 
    107       IF( l_1st_euler ) THEN   ;   r2dt =         rn_Dt   ! start or restart with Euler 1st time-step 
    108       ELSE                     ;   r2dt = 2._wp * rn_Dt   ! restart with leapfrog 
    109       ENDIF 
    110       r1_2dt = 1._wp / r2dt 
     107      IF( l_1st_euler ) THEN   
     108         ! start or restart with Euler 1st time-step 
     109         rDt =  rn_Dt    
     110         r1_Dt = 1._wp / rDt 
     111      ENDIF 
    111112      ! 
    112113      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    344345#endif 
    345346      ! 
    346       IF( l_1st_euler ) l_1st_euler = .FALSE.      ! recover Leap-frog timestep 
     347      IF( l_1st_euler ) THEN         ! recover Leap-frog timestep 
     348         rDt = 2._wp * rn_Dt    
     349         r1_Dt = 1._wp / rDt 
     350         l_1st_euler = .FALSE.       
     351      ENDIF 
    347352      ! 
    348353      IF( ln_timing )   CALL timing_stop('stp') 
Note: See TracChangeset for help on using the changeset viewer.