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

Changeset 10743


Ignore:
Timestamp:
2019-03-12T12:06:29+01:00 (6 years ago)
Author:
davestorkey
Message:

First block of changes:

  1. New state variables uu, vv, e3X and corresponding temporary pointers.
  2. Updated code for time-filtering and time-level-swapping uu, vv, e3X.
  3. DYN routines not updated to use new variables at this stage.
Location:
NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DOM/dom_oce.F90

    r10068 r10743  
    6262   INTEGER , PUBLIC ::   neuler         !: restart euler forward option (0=Euler) 
    6363   REAL(wp), PUBLIC ::   r2dt           !: = 2*rdt except at nit000 (=rdt) if neuler=0 
     64 
     65   !!---------------------------------------------------------------------- 
     66   !! time level indices 
     67   !!---------------------------------------------------------------------- 
     68   INTEGER, PUBLIC :: Nm1, Nnn, Np1, Nrhs, Nm1_2lev, Nnn_2lev 
    6469 
    6570   !!---------------------------------------------------------------------- 
     
    129134   LOGICAL, PUBLIC ::   ln_sco       !: s-coordinate or hybrid z-s coordinate 
    130135   LOGICAL, PUBLIC ::   ln_isfcav    !: presence of ISF  
    131    !                                                        !  ref.   ! before  !   now   ! after  ! 
    132    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3t_0 ,   e3t_b ,   e3t_n ,  e3t_a   !: t- vert. scale factor [m] 
    133    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3u_0 ,   e3u_b ,   e3u_n ,  e3u_a   !: u- vert. scale factor [m] 
    134    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3v_0 ,   e3v_b ,   e3v_n ,  e3v_a   !: v- vert. scale factor [m] 
    135    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3f_0           ,   e3f_n            !: f- vert. scale factor [m] 
    136    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3w_0 ,   e3w_b ,   e3w_n            !: w- vert. scale factor [m] 
    137    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::    e3uw_0 ,  e3uw_b ,  e3uw_n            !: uw-vert. scale factor [m] 
    138    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::    e3vw_0 ,  e3vw_b ,  e3vw_n            !: vw-vert. scale factor [m] 
     136   !                                                        !  ref. scale factors 
     137   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3t_0  !: t- vert. scale factor [m] 
     138   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3u_0  !: u- vert. scale factor [m] 
     139   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3v_0  !: v- vert. scale factor [m] 
     140   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3f_0  !: f- vert. scale factor [m] 
     141   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::     e3w_0  !: w- vert. scale factor [m] 
     142   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::    e3uw_0  !: uw-vert. scale factor [m] 
     143   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::    e3vw_0  !: vw-vert. scale factor [m] 
     144   !                                                        !  time-dependent scale factors 
     145   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:,:,:) ::   e3t, e3u, e3v, e3w, e3uw, e3vw  !: vert. scale factor [m] 
     146   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, TARGET, DIMENSION(:,:,:) ::     e3f                             !: F-point vert. scale factor [m] 
    139147 
    140148   !                                                        !  ref.   ! before  !   now   ! 
     
    149157   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::           r1_hu_b , r1_hu_n , r1_hu_a   !: inverse of u-depth [1/m] 
    150158   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::           r1_hv_b , r1_hv_n , r1_hv_a   !: inverse of v-depth [1/m] 
     159 
     160   !! TEMPORARY POINTERS FOR DEVELOPMENT ONLY 
     161   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::     e3t_b ,   e3t_n ,  e3t_a   !: t- vert. scale factor [m] 
     162   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::     e3u_b ,   e3u_n ,  e3u_a   !: u- vert. scale factor [m] 
     163   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::     e3v_b ,   e3v_n ,  e3v_a   !: v- vert. scale factor [m] 
     164   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::               e3f_n            !: f- vert. scale factor [m] 
     165   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::     e3w_b ,   e3w_n            !: w- vert. scale factor [m] 
     166   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::    e3uw_b ,  e3uw_n            !: uw-vert. scale factor [m] 
     167   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:) ::    e3vw_b ,  e3vw_n            !: vw-vert. scale factor [m] 
     168   !! TEMPORARY POINTERS FOR DEVELOPMENT ONLY 
    151169 
    152170   INTEGER, PUBLIC ::   nla10              !: deepest    W level Above  ~10m (nlb10 - 1) 
     
    260278         &      gdept_n(jpi,jpj,jpk) , gdepw_n(jpi,jpj,jpk) , gde3w_n(jpi,jpj,jpk) , STAT=ierr(4) ) 
    261279         ! 
    262       ALLOCATE( e3t_0(jpi,jpj,jpk) , e3u_0(jpi,jpj,jpk) , e3v_0(jpi,jpj,jpk) , e3f_0(jpi,jpj,jpk) , e3w_0(jpi,jpj,jpk) ,   & 
    263          &      e3t_b(jpi,jpj,jpk) , e3u_b(jpi,jpj,jpk) , e3v_b(jpi,jpj,jpk) ,                      e3w_b(jpi,jpj,jpk) ,   &  
    264          &      e3t_n(jpi,jpj,jpk) , e3u_n(jpi,jpj,jpk) , e3v_n(jpi,jpj,jpk) , e3f_n(jpi,jpj,jpk) , e3w_n(jpi,jpj,jpk) ,   &  
    265          &      e3t_a(jpi,jpj,jpk) , e3u_a(jpi,jpj,jpk) , e3v_a(jpi,jpj,jpk) ,                                             & 
    266          !                                                          ! 
    267          &      e3uw_0(jpi,jpj,jpk) , e3vw_0(jpi,jpj,jpk) ,         & 
    268          &      e3uw_b(jpi,jpj,jpk) , e3vw_b(jpi,jpj,jpk) ,         &                
    269          &      e3uw_n(jpi,jpj,jpk) , e3vw_n(jpi,jpj,jpk) ,     STAT=ierr(5) )                        
     280      ALLOCATE( e3t_0(jpi,jpj,jpk)       , e3u_0(jpi,jpj,jpk)       , e3v_0(jpi,jpj,jpk)       , e3f_0(jpi,jpj,jpk) ,    & 
     281         &      e3w_0(jpi,jpj,jpk)       , e3uw_0(jpi,jpj,jpk)      , e3vw_0(jpi,jpj,jpk)      ,                         & 
     282         &      e3t  (jpi,jpj,jpk,jpt)   , e3u  (jpi,jpj,jpk,jpt)   , e3v  (jpi,jpj,jpk,jpt)   ,                         &  
     283         &      e3w  (jpi,jpj,jpk,jpt-1) , e3uw (jpi,jpj,jpk,jpt-1) , e3vw (jpi,jpj,jpk,jpt-1) ,                         &  
     284         &      e3f  (jpi,jpj,jpk),     STAT=ierr(5) )                        
    270285         ! 
    271286      ALLOCATE( ht_0(jpi,jpj) , hu_0(jpi,jpj) , hv_0(jpi,jpj) ,                                           & 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DOM/domvvl.F90

    r10425 r10743  
    563563 
    564564 
    565    SUBROUTINE dom_vvl_sf_swp( kt ) 
     565   SUBROUTINE dom_vvl_sf_swp( kt, kNm1, kNnn, kNm1_2lev, kNnn_2lev ) 
    566566      !!---------------------------------------------------------------------- 
    567567      !!                ***  ROUTINE dom_vvl_sf_swp  *** 
     
    588588      !!---------------------------------------------------------------------- 
    589589      INTEGER, INTENT( in ) ::   kt   ! time step 
     590      INTEGER, INTENT( in ) ::   kNm1, kNnn, kNm1_2lev, kNnn_2lev ! time level indices 
    590591      ! 
    591592      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
     
    605606      ! Time filter and swap of scale factors 
    606607      ! ===================================== 
    607       ! - ML - e3(t/u/v)_b are allready computed in dynnxt. 
     608      ! - ML - e3(t/u/v)_b are already computed in dynnxt. 
    608609      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
    609610         IF( neuler == 0 .AND. kt == nit000 ) THEN 
     
    618619      gdepw_b(:,:,:) = gdepw_n(:,:,:) 
    619620 
    620       e3t_n(:,:,:) = e3t_a(:,:,:) 
    621       e3u_n(:,:,:) = e3u_a(:,:,:) 
    622       e3v_n(:,:,:) = e3v_a(:,:,:) 
     621!! TO BE DELETED 
     622!!$      e3t_n(:,:,:) = e3t_a(:,:,:) 
     623!!$      e3u_n(:,:,:) = e3u_a(:,:,:) 
     624!!$      e3v_n(:,:,:) = e3v_a(:,:,:) 
     625!! TO BE DELETED 
    623626 
    624627      ! Compute all missing vertical scale factor and depths 
     
    629632      ! - JC - hu_b, hv_b, hur_b, hvr_b also 
    630633       
    631       CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F'  ) 
     634      CALL dom_vvl_interpol( e3u(:,:,:,kNnn), e3f(:,:,:), 'F'  ) 
    632635       
    633636      ! Vertical scale factor interpolations 
    634       CALL dom_vvl_interpol( e3t_n(:,:,:),  e3w_n(:,:,:), 'W'  ) 
    635       CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' ) 
    636       CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' ) 
    637       CALL dom_vvl_interpol( e3t_b(:,:,:),  e3w_b(:,:,:), 'W'  ) 
    638       CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' ) 
    639       CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' ) 
     637      CALL dom_vvl_interpol( e3t(:,:,:,kNnn),  e3w(:,:,:,kNnn_2lev), 'W'  ) 
     638      CALL dom_vvl_interpol( e3u(:,:,:,kNnn), e3uw(:,:,:,kNnn_2lev), 'UW' ) 
     639      CALL dom_vvl_interpol( e3v(:,:,:,kNnn), e3vw(:,:,:,kNnn_2lev), 'VW' ) 
     640      CALL dom_vvl_interpol( e3t(:,:,:,kNm1),  e3w(:,:,:,kNm1_2lev), 'W'  ) 
     641      CALL dom_vvl_interpol( e3u(:,:,:,kNm1), e3uw(:,:,:,kNm1_2lev), 'UW' ) 
     642      CALL dom_vvl_interpol( e3v(:,:,:,kNm1), e3vw(:,:,:,kNm1_2lev), 'VW' ) 
    640643 
    641644      ! t- and w- points depth (set the isf depth as it is in the initial step) 
    642       gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1) 
     645      gdept_n(:,:,1) = 0.5_wp * e3w(:,:,1,kNnn_2lev) 
    643646      gdepw_n(:,:,1) = 0.0_wp 
    644647      gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:) 
     
    649652                                                                 ! 1 for jk = mikt 
    650653               zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk)) 
    651                gdepw_n(ji,jj,jk) = gdepw_n(ji,jj,jk-1) + e3t_n(ji,jj,jk-1) 
    652                gdept_n(ji,jj,jk) =    zcoef  * ( gdepw_n(ji,jj,jk  ) + 0.5 * e3w_n(ji,jj,jk) )  & 
    653                    &             + (1-zcoef) * ( gdept_n(ji,jj,jk-1) +       e3w_n(ji,jj,jk) )  
     654               gdepw_n(ji,jj,jk) = gdepw_n(ji,jj,jk-1) + e3t(ji,jj,jk-1,kNnn) 
     655               gdept_n(ji,jj,jk) =    zcoef  * ( gdepw_n(ji,jj,jk  ) + 0.5 * e3w(ji,jj,jk,kNnn_2lev) )  & 
     656                   &             + (1-zcoef) * ( gdept_n(ji,jj,jk-1) +       e3w(ji,jj,jk,kNnn_2lev) )  
    654657               gde3w_n(ji,jj,jk) = gdept_n(ji,jj,jk) - sshn(ji,jj) 
    655658            END DO 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DYN/dynnxt.F90

    r10425 r10743  
    6464CONTAINS 
    6565 
    66    SUBROUTINE dyn_nxt ( kt ) 
     66   SUBROUTINE dyn_nxt ( kt, kNm1, kNp1, puu, pvv, pe3t, pe3u, pe3v ) 
    6767      !!---------------------------------------------------------------------- 
    6868      !!                  ***  ROUTINE dyn_nxt  *** 
     
    9292      !!               un,vn   now horizontal velocity of next time-step 
    9393      !!---------------------------------------------------------------------- 
    94       INTEGER, INTENT( in ) ::   kt      ! ocean time-step index 
     94      INTEGER, INTENT( in ) ::   kt         ! ocean time-step index 
     95      INTEGER, INTENT( in ) ::   kNm1, kNp1 ! before and after time level indices 
     96      REAL(wp), INTENT( inout), DIMENSION(jpi,jpj,jpk) :: puu, pvv ! now velocities to be time filtered 
     97      REAL(wp), INTENT( inout), DIMENSION(jpi,jpj,jpk) :: pe3t, pe3u, pe3v ! now scale factors to be time filtered 
    9598      ! 
    9699      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
    97100      INTEGER  ::   ikt          ! local integers 
    98       REAL(wp) ::   zue3a, zue3n, zue3b, zuf, zcoef    ! local scalars 
    99       REAL(wp) ::   zve3a, zve3n, zve3b, zvf, z1_2dt   !   -      - 
     101      REAL(wp) ::   zue3a, zue3n, zue3b, zcoef    ! local scalars 
     102      REAL(wp) ::   zve3a, zve3n, zve3b, z1_2dt   !   -      - 
    100103      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   zue, zve 
    101       REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ze3u_f, ze3v_f, zua, zva  
     104      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ze3t_f, ze3u_f, ze3v_f, zua, zva  
    102105      !!---------------------------------------------------------------------- 
    103106      ! 
     
    115118         ! Ensure below that barotropic velocities match time splitting estimate 
    116119         ! Compute actual transport and replace it with ts estimate at "after" time step 
    117          zue(:,:) = e3u_a(:,:,1) * ua(:,:,1) * umask(:,:,1) 
    118          zve(:,:) = e3v_a(:,:,1) * va(:,:,1) * vmask(:,:,1) 
     120         zue(:,:) = e3u(:,:,1,kNp1) * uu(:,:,1,kNp1) * umask(:,:,1) 
     121         zve(:,:) = e3v(:,:,1,kNp1) * vv(:,:,1,kNp1) * vmask(:,:,1) 
    119122         DO jk = 2, jpkm1 
    120             zue(:,:) = zue(:,:) + e3u_a(:,:,jk) * ua(:,:,jk) * umask(:,:,jk) 
    121             zve(:,:) = zve(:,:) + e3v_a(:,:,jk) * va(:,:,jk) * vmask(:,:,jk) 
     123            zue(:,:) = zue(:,:) + e3u(:,:,jk,kNp1) * uu(:,:,jk,kNp1) * umask(:,:,jk) 
     124            zve(:,:) = zve(:,:) + e3v(:,:,jk,kNp1) * vv(:,:,jk,kNp1) * vmask(:,:,jk) 
    122125         END DO 
    123126         DO jk = 1, jpkm1 
    124             ua(:,:,jk) = ( ua(:,:,jk) - zue(:,:) * r1_hu_a(:,:) + ua_b(:,:) ) * umask(:,:,jk) 
    125             va(:,:,jk) = ( va(:,:,jk) - zve(:,:) * r1_hv_a(:,:) + va_b(:,:) ) * vmask(:,:,jk) 
     127            uu(:,:,jk,kNp1) = ( uu(:,:,jk,kNp1) - zue(:,:) * r1_hu_a(:,:) + ua_b(:,:) ) * umask(:,:,jk) 
     128            vv(:,:,jk,kNp1) = ( vv(:,:,jk,kNp1) - zve(:,:) * r1_hv_a(:,:) + va_b(:,:) ) * vmask(:,:,jk) 
    126129         END DO 
    127130         ! 
     
    132135            ! so that asselin contribution is removed at the same time  
    133136            DO jk = 1, jpkm1 
    134                un(:,:,jk) = ( un(:,:,jk) - un_adv(:,:)*r1_hu_n(:,:) + un_b(:,:) )*umask(:,:,jk) 
    135                vn(:,:,jk) = ( vn(:,:,jk) - vn_adv(:,:)*r1_hv_n(:,:) + vn_b(:,:) )*vmask(:,:,jk) 
     137               puu(:,:,jk) = ( puu(:,:,jk) - un_adv(:,:)*r1_hu_n(:,:) + un_b(:,:) )*umask(:,:,jk) 
     138               pvv(:,:,jk) = ( pvv(:,:,jk) - vn_adv(:,:)*r1_hv_n(:,:) + vn_b(:,:) )*vmask(:,:,jk) 
    136139            END DO   
    137140         ENDIF 
     
    144147# endif 
    145148      ! 
    146       CALL lbc_lnk_multi( 'dynnxt', ua, 'U', -1., va, 'V', -1. )     !* local domain boundaries 
     149      CALL lbc_lnk_multi( 'dynnxt', uu(:,:,:,kNp1), 'U', -1., vv(:,:,:,kNp1), 'V', -1. )     !* local domain boundaries 
    147150      ! 
    148151      !                                !* BDY open boundaries 
     
    157160         ! 
    158161         !                                  ! Kinetic energy and Conversion 
    159          IF( ln_KE_trd  )   CALL trd_dyn( ua, va, jpdyn_ken, kt ) 
     162         IF( ln_KE_trd  )   CALL trd_dyn( uu(:,:,:,kNp1), vv(:,:,:,kNp1), jpdyn_ken, kt ) 
    160163         ! 
    161164         IF( ln_dyn_trd ) THEN              ! 3D output: total momentum trends 
    162             zua(:,:,:) = ( ua(:,:,:) - ub(:,:,:) ) * z1_2dt 
    163             zva(:,:,:) = ( va(:,:,:) - vb(:,:,:) ) * z1_2dt 
     165            zua(:,:,:) = ( uu(:,:,:,kNp1) - uu(:,:,:,kNm1) ) * z1_2dt 
     166            zva(:,:,:) = ( vv(:,:,:,kNp1) - vv(:,:,:,kNm1) ) * z1_2dt 
    164167            CALL iom_put( "utrd_tot", zua )        ! total momentum trends, except the asselin time filter 
    165168            CALL iom_put( "vtrd_tot", zva ) 
    166169         ENDIF 
    167170         ! 
    168          zua(:,:,:) = un(:,:,:)             ! save the now velocity before the asselin filter 
    169          zva(:,:,:) = vn(:,:,:)             ! (caution: there will be a shift by 1 timestep in the 
     171         zua(:,:,:) = puu(:,:,:)             ! save the now velocity before the asselin filter 
     172         zva(:,:,:) = pvv(:,:,:)             ! (caution: there will be a shift by 1 timestep in the 
    170173         !                                  !  computation of the asselin filter trends) 
    171174      ENDIF 
     
    173176      ! Time filter and swap of dynamics arrays 
    174177      ! ------------------------------------------ 
    175       IF( neuler == 0 .AND. kt == nit000 ) THEN        !* Euler at first time-step: only swap 
    176          DO jk = 1, jpkm1 
    177             un(:,:,jk) = ua(:,:,jk)                         ! un <-- ua 
    178             vn(:,:,jk) = va(:,:,jk) 
    179          END DO 
    180          IF( .NOT.ln_linssh ) THEN                          ! e3._b <-- e3._n 
     178!! TO BE DELETED 
     179!!$      IF( neuler == 0 .AND. kt == nit000 ) THEN        !* Euler at first time-step: only swap 
     180!!$         DO jk = 1, jpkm1 
     181!!$            un(:,:,jk) = ua(:,:,jk)                         ! un <-- ua 
     182!!$            vn(:,:,jk) = va(:,:,jk) 
     183!!$         END DO  
     184!!$         IF( .NOT.ln_linssh ) THEN                          ! e3._b <-- e3._n 
    181185!!gm BUG ????    I don't understand why it is not : e3._n <-- e3._a   
    182             DO jk = 1, jpkm1 
     186!!$            DO jk = 1, jpkm1 
    183187!               e3t_b(:,:,jk) = e3t_n(:,:,jk) 
    184188!               e3u_b(:,:,jk) = e3u_n(:,:,jk) 
    185189!               e3v_b(:,:,jk) = e3v_n(:,:,jk) 
    186190               ! 
    187                e3t_n(:,:,jk) = e3t_a(:,:,jk) 
    188                e3u_n(:,:,jk) = e3u_a(:,:,jk) 
    189                e3v_n(:,:,jk) = e3v_a(:,:,jk) 
    190             END DO 
     191!!$               e3t_n(:,:,jk) = e3t_a(:,:,jk) 
     192!!$               e3u_n(:,:,jk) = e3u_a(:,:,jk) 
     193!!$               e3v_n(:,:,jk) = e3v_a(:,:,jk) 
     194!!$            END DO 
    191195!!gm BUG end 
    192          ENDIF 
     196!!$         ENDIF 
     197!! TO BE DELETED 
    193198                                                            !  
    194199          
    195       ELSE                                             !* Leap-Frog : Asselin filter and swap 
     200      IF( .NOT.( neuler == 0 .AND. kt == nit000 ) ) THEN    !* Leap-Frog : Asselin filter and swap 
    196201         !                                ! =============! 
    197202         IF( ln_linssh ) THEN             ! Fixed volume ! 
     
    200205               DO jj = 1, jpj 
    201206                  DO ji = 1, jpi     
    202                      zuf = un(ji,jj,jk) + atfp * ( ub(ji,jj,jk) - 2._wp * un(ji,jj,jk) + ua(ji,jj,jk) ) 
    203                      zvf = vn(ji,jj,jk) + atfp * ( vb(ji,jj,jk) - 2._wp * vn(ji,jj,jk) + va(ji,jj,jk) ) 
     207                     puu(ji,jj,jk) = puu(ji,jj,jk) + atfp * ( uu(ji,jj,jk,kNm1) - 2._wp * puu(ji,jj,jk) + uu(ji,jj,jk,kNp1) ) 
     208                     pvv(ji,jj,jk) = pvv(ji,jj,jk) + atfp * ( vv(ji,jj,jk,kNm1) - 2._wp * pvv(ji,jj,jk) + vv(ji,jj,jk,kNp1) ) 
    204209                     ! 
    205                      ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity 
    206                      vb(ji,jj,jk) = zvf 
    207                      un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua 
    208                      vn(ji,jj,jk) = va(ji,jj,jk) 
     210!! TO BE DELETED 
     211!!$                     ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity 
     212!!$                     vb(ji,jj,jk) = zvf 
     213!!$                     un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua 
     214!!$                     vn(ji,jj,jk) = va(ji,jj,jk) 
     215!! TO BE DELETED 
    209216                  END DO 
    210217               END DO 
     
    213220         ELSE                             ! Variable volume ! 
    214221            !                             ! ================! 
    215             ! Before scale factor at t-points 
    216             ! (used as a now filtered scale factor until the swap) 
     222            ! Time-filtered scale factor at t-points 
    217223            ! ---------------------------------------------------- 
     224            ALLOCATE( ze3t_f(jpi,jpj,jpk) ) 
    218225            DO jk = 1, jpkm1 
    219                e3t_b(:,:,jk) = e3t_n(:,:,jk) + atfp * ( e3t_b(:,:,jk) - 2._wp * e3t_n(:,:,jk) + e3t_a(:,:,jk) ) 
     226               ze3t_f(:,:,jk) = pe3t(:,:,jk) + atfp * ( e3t(:,:,jk,kNm1) - 2._wp * pe3t(:,:,jk) + e3t(:,:,jk,kNp1) ) 
    220227            END DO 
    221228            ! Add volume filter correction: compatibility with tracer advection scheme 
     
    223230            zcoef = atfp * rdt * r1_rau0 
    224231 
    225             e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) ) * tmask(:,:,1) 
     232            ze3t_f(:,:,1) = ze3t_f(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) ) * tmask(:,:,1) 
    226233 
    227234            IF ( ln_rnf ) THEN 
     
    231238                        DO ji = 1, jpi 
    232239                           IF( jk <=  nk_rnf(ji,jj)  ) THEN 
    233                                e3t_b(ji,jj,jk) =   e3t_b(ji,jj,jk) - zcoef *  ( - rnf_b(ji,jj) + rnf(ji,jj) ) & 
    234                                       &          * ( e3t_n(ji,jj,jk) / h_rnf(ji,jj) ) * tmask(ji,jj,jk) 
     240                               ze3t_f(ji,jj,jk) =   ze3t_f(ji,jj,jk) - zcoef *  ( - rnf_b(ji,jj) + rnf(ji,jj) ) & 
     241                                      &          * ( pe3t(ji,jj,jk) / h_rnf(ji,jj) ) * tmask(ji,jj,jk) 
    235242                           ENDIF 
    236243                        ENDDO 
     
    238245                  ENDDO 
    239246               ELSE 
    240                   e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef *  ( -rnf_b(:,:) + rnf(:,:))*tmask(:,:,1) 
     247                  ze3t_f(:,:,1) = ze3t_f(:,:,1) - zcoef *  ( -rnf_b(:,:) + rnf(:,:))*tmask(:,:,1) 
    241248               ENDIF 
    242249            END IF 
     
    247254                     DO ji = 1, jpi 
    248255                        IF( misfkt(ji,jj) <=jk .and. jk < misfkb(ji,jj)  ) THEN 
    249                            e3t_b(ji,jj,jk) = e3t_b(ji,jj,jk) - zcoef * ( fwfisf_b(ji,jj) - fwfisf(ji,jj) ) & 
    250                                 &          * ( e3t_n(ji,jj,jk) * r1_hisf_tbl(ji,jj) ) * tmask(ji,jj,jk) 
     256                           ze3t_f(ji,jj,jk) = ze3t_f(ji,jj,jk) - zcoef * ( fwfisf_b(ji,jj) - fwfisf(ji,jj) ) & 
     257                                &          * ( pe3t(ji,jj,jk) * r1_hisf_tbl(ji,jj) ) * tmask(ji,jj,jk) 
    251258                        ELSEIF ( jk==misfkb(ji,jj) ) THEN 
    252                            e3t_b(ji,jj,jk) = e3t_b(ji,jj,jk) - zcoef * ( fwfisf_b(ji,jj) - fwfisf(ji,jj) ) & 
    253                                 &          * ( e3t_n(ji,jj,jk) * r1_hisf_tbl(ji,jj) ) * ralpha(ji,jj) * tmask(ji,jj,jk) 
     259                           ze3t_f(ji,jj,jk) = ze3t_f(ji,jj,jk) - zcoef * ( fwfisf_b(ji,jj) - fwfisf(ji,jj) ) & 
     260                                &          * ( pe3t(ji,jj,jk) * r1_hisf_tbl(ji,jj) ) * ralpha(ji,jj) * tmask(ji,jj,jk) 
    254261                        ENDIF 
    255262                     END DO 
     
    258265            END IF 
    259266            ! 
     267            pe3t(:,:,1:jpkm1) = ze3t_f(:,:,1:jpkm1)        ! filtered scale factor at T-points 
     268            ! 
    260269            IF( ln_dynadv_vec ) THEN      ! Asselin filter applied on velocity 
    261270               ! Before filtered scale factor at (u/v)-points 
    262                CALL dom_vvl_interpol( e3t_b(:,:,:), e3u_b(:,:,:), 'U' ) 
    263                CALL dom_vvl_interpol( e3t_b(:,:,:), e3v_b(:,:,:), 'V' ) 
     271               CALL dom_vvl_interpol( pe3t(:,:,:), pe3u(:,:,:), 'U' ) 
     272               CALL dom_vvl_interpol( pe3t(:,:,:), pe3v(:,:,:), 'V' ) 
    264273               DO jk = 1, jpkm1 
    265274                  DO jj = 1, jpj 
    266275                     DO ji = 1, jpi 
    267                         zuf = un(ji,jj,jk) + atfp * ( ub(ji,jj,jk) - 2._wp * un(ji,jj,jk) + ua(ji,jj,jk) ) 
    268                         zvf = vn(ji,jj,jk) + atfp * ( vb(ji,jj,jk) - 2._wp * vn(ji,jj,jk) + va(ji,jj,jk) ) 
     276                        puu(ji,jj,jk) = puu(ji,jj,jk) + atfp * ( uu(ji,jj,jk,kNm1) - 2._wp * puu(ji,jj,jk) + uu(ji,jj,jk,kNp1) ) 
     277                        pvv(ji,jj,jk) = pvv(ji,jj,jk) + atfp * ( vv(ji,jj,jk,kNm1) - 2._wp * pvv(ji,jj,jk) + vv(ji,jj,jk,kNp1) ) 
    269278                        ! 
    270                         ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity 
    271                         vb(ji,jj,jk) = zvf 
    272                         un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua 
    273                         vn(ji,jj,jk) = va(ji,jj,jk) 
     279!! TO BE DELETED 
     280!!$                        ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity 
     281!!$                        vb(ji,jj,jk) = zvf 
     282!!$                        un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua 
     283!!$                        vn(ji,jj,jk) = va(ji,jj,jk) 
     284!! TO BE DELETED 
    274285                     END DO 
    275286                  END DO 
     
    279290               ! 
    280291               ALLOCATE( ze3u_f(jpi,jpj,jpk) , ze3v_f(jpi,jpj,jpk) ) 
    281                ! Before filtered scale factor at (u/v)-points stored in ze3u_f, ze3v_f 
    282                CALL dom_vvl_interpol( e3t_b(:,:,:), ze3u_f, 'U' ) 
    283                CALL dom_vvl_interpol( e3t_b(:,:,:), ze3v_f, 'V' ) 
     292               ! Now filtered scale factor at (u/v)-points stored in ze3u_f, ze3v_f 
     293               CALL dom_vvl_interpol( pe3t(:,:,:), ze3u_f, 'U' ) 
     294               CALL dom_vvl_interpol( pe3t(:,:,:), ze3v_f, 'V' ) 
    284295               DO jk = 1, jpkm1 
    285296                  DO jj = 1, jpj 
    286297                     DO ji = 1, jpi                   
    287                         zue3a = e3u_a(ji,jj,jk) * ua(ji,jj,jk) 
    288                         zve3a = e3v_a(ji,jj,jk) * va(ji,jj,jk) 
    289                         zue3n = e3u_n(ji,jj,jk) * un(ji,jj,jk) 
    290                         zve3n = e3v_n(ji,jj,jk) * vn(ji,jj,jk) 
    291                         zue3b = e3u_b(ji,jj,jk) * ub(ji,jj,jk) 
    292                         zve3b = e3v_b(ji,jj,jk) * vb(ji,jj,jk) 
     298                        zue3a = e3u(ji,jj,jk,kNp1) * uu(ji,jj,jk,kNp1) 
     299                        zve3a = e3v(ji,jj,jk,kNp1) * vv(ji,jj,jk,kNp1) 
     300                        zue3n = pe3u(ji,jj,jk)     * puu(ji,jj,jk) 
     301                        zve3n = pe3v(ji,jj,jk)     * pvv(ji,jj,jk) 
     302                        zue3b = e3u(ji,jj,jk,kNm1) * uu(ji,jj,jk,kNm1) 
     303                        zve3b = e3v(ji,jj,jk,kNm1) * vv(ji,jj,jk,kNm1) 
    293304                        ! 
    294                         zuf = ( zue3n + atfp * ( zue3b - 2._wp * zue3n  + zue3a ) ) / ze3u_f(ji,jj,jk) 
    295                         zvf = ( zve3n + atfp * ( zve3b - 2._wp * zve3n  + zve3a ) ) / ze3v_f(ji,jj,jk) 
     305                        puu(ji,jj,jk) = ( zue3n + atfp * ( zue3b - 2._wp * zue3n  + zue3a ) ) / ze3u_f(ji,jj,jk) 
     306                        pvv(ji,jj,jk) = ( zve3n + atfp * ( zve3b - 2._wp * zve3n  + zve3a ) ) / ze3v_f(ji,jj,jk) 
    296307                        ! 
    297                         ub(ji,jj,jk) = zuf                     ! ub <-- filtered velocity 
    298                         vb(ji,jj,jk) = zvf 
    299                         un(ji,jj,jk) = ua(ji,jj,jk)            ! un <-- ua 
    300                         vn(ji,jj,jk) = va(ji,jj,jk) 
     308!! TO BE DELETED 
     309!!$                        ub(ji,jj,jk) = zuf                     ! ub <-- filtered velocity 
     310!!$                        vb(ji,jj,jk) = zvf 
     311!!$                        un(ji,jj,jk) = ua(ji,jj,jk)            ! un <-- ua 
     312!!$                        vn(ji,jj,jk) = va(ji,jj,jk) 
     313!! TO BE DELETED 
    301314                     END DO 
    302315                  END DO 
    303316               END DO 
    304                e3u_b(:,:,1:jpkm1) = ze3u_f(:,:,1:jpkm1)        ! e3u_b <-- filtered scale factor 
    305                e3v_b(:,:,1:jpkm1) = ze3v_f(:,:,1:jpkm1) 
     317               pe3u(:,:,1:jpkm1) = ze3u_f(:,:,1:jpkm1)   
     318               pe3v(:,:,1:jpkm1) = ze3v_f(:,:,1:jpkm1) 
    306319               ! 
    307320               DEALLOCATE( ze3u_f , ze3v_f ) 
    308321            ENDIF 
    309322            ! 
     323            DEALLOCATE( ze3t_f ) 
    310324         ENDIF 
    311325         ! 
    312326         IF( ln_dynspg_ts .AND. ln_bt_fw ) THEN 
    313             ! Revert "before" velocities to time split estimate 
     327            ! Revert filtered "now" velocities to time split estimate 
    314328            ! Doing it here also means that asselin filter contribution is removed   
    315             zue(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1) 
    316             zve(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1)     
     329            zue(:,:) = pe3u(:,:,1) * puu(:,:,1) * umask(:,:,1) 
     330            zve(:,:) = pe3v(:,:,1) * pvv(:,:,1) * vmask(:,:,1)     
    317331            DO jk = 2, jpkm1 
    318                zue(:,:) = zue(:,:) + e3u_b(:,:,jk) * ub(:,:,jk) * umask(:,:,jk) 
    319                zve(:,:) = zve(:,:) + e3v_b(:,:,jk) * vb(:,:,jk) * vmask(:,:,jk)     
     332               zue(:,:) = zue(:,:) + pe3u(:,:,jk) * puu(:,:,jk) * umask(:,:,jk) 
     333               zve(:,:) = zve(:,:) + pe3v(:,:,jk) * pvv(:,:,jk) * vmask(:,:,jk)     
    320334            END DO 
    321335            DO jk = 1, jpkm1 
    322                ub(:,:,jk) = ub(:,:,jk) - (zue(:,:) * r1_hu_n(:,:) - un_b(:,:)) * umask(:,:,jk) 
    323                vb(:,:,jk) = vb(:,:,jk) - (zve(:,:) * r1_hv_n(:,:) - vn_b(:,:)) * vmask(:,:,jk) 
     336               puu(:,:,jk) = puu(:,:,jk) - (zue(:,:) * r1_hu_n(:,:) - un_b(:,:)) * umask(:,:,jk) 
     337               pvv(:,:,jk) = pvv(:,:,jk) - (zve(:,:) * r1_hv_n(:,:) - vn_b(:,:)) * vmask(:,:,jk) 
    324338            END DO 
    325339         ENDIF 
     
    331345      ! integration 
    332346      ! 
     347      ! DS IMMERSE :: This is very confusing as it stands. But should  
     348      ! hopefully be simpler when we do the time-level swapping for the 
     349      ! external mode solution.  
    333350      ! 
    334351      IF(.NOT.ln_linssh ) THEN 
    335          hu_b(:,:) = e3u_b(:,:,1) * umask(:,:,1) 
    336          hv_b(:,:) = e3v_b(:,:,1) * vmask(:,:,1) 
     352         hu_b(:,:) = pe3u(:,:,1) * umask(:,:,1) 
     353         hv_b(:,:) = pe3v(:,:,1) * vmask(:,:,1) 
    337354         DO jk = 2, jpkm1 
    338             hu_b(:,:) = hu_b(:,:) + e3u_b(:,:,jk) * umask(:,:,jk) 
    339             hv_b(:,:) = hv_b(:,:) + e3v_b(:,:,jk) * vmask(:,:,jk) 
     355            hu_b(:,:) = hu_b(:,:) + pe3u(:,:,jk) * umask(:,:,jk) 
     356            hv_b(:,:) = hv_b(:,:) + pe3v(:,:,jk) * vmask(:,:,jk) 
    340357         END DO 
    341358         r1_hu_b(:,:) = ssumask(:,:) / ( hu_b(:,:) + 1._wp - ssumask(:,:) ) 
     
    343360      ENDIF 
    344361      ! 
    345       un_b(:,:) = e3u_a(:,:,1) * un(:,:,1) * umask(:,:,1) 
    346       ub_b(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1) 
    347       vn_b(:,:) = e3v_a(:,:,1) * vn(:,:,1) * vmask(:,:,1) 
    348       vb_b(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1) 
     362      un_b(:,:) = e3u(:,:,1,kNp1) * uu(:,:,1,kNp1) * umask(:,:,1) 
     363      ub_b(:,:) = pe3u(:,:,1)     * puu(:,:,1) * umask(:,:,1) 
     364      vn_b(:,:) = e3v(:,:,1,kNp1) * vv(:,:,1,kNp1) * vmask(:,:,1) 
     365      vb_b(:,:) = pe3v(:,:,1)     * pvv(:,:,1) * vmask(:,:,1) 
    349366      DO jk = 2, jpkm1 
    350          un_b(:,:) = un_b(:,:) + e3u_a(:,:,jk) * un(:,:,jk) * umask(:,:,jk) 
    351          ub_b(:,:) = ub_b(:,:) + e3u_b(:,:,jk) * ub(:,:,jk) * umask(:,:,jk) 
    352          vn_b(:,:) = vn_b(:,:) + e3v_a(:,:,jk) * vn(:,:,jk) * vmask(:,:,jk) 
    353          vb_b(:,:) = vb_b(:,:) + e3v_b(:,:,jk) * vb(:,:,jk) * vmask(:,:,jk) 
     367         un_b(:,:) = un_b(:,:) + e3u(:,:,jk,kNp1) * uu(:,:,jk,kNp1) * umask(:,:,jk) 
     368         ub_b(:,:) = ub_b(:,:) + pe3u(:,:,jk)     * puu(:,:,jk) * umask(:,:,jk) 
     369         vn_b(:,:) = vn_b(:,:) + e3v(:,:,jk,kNp1) * vv(:,:,jk,kNp1) * vmask(:,:,jk) 
     370         vb_b(:,:) = vb_b(:,:) + pe3v(:,:,jk)     * pvv(:,:,jk) * vmask(:,:,jk) 
    354371      END DO 
    355372      un_b(:,:) = un_b(:,:) * r1_hu_a(:,:) 
     
    363380      ENDIF 
    364381      IF( l_trddyn ) THEN                ! 3D output: asselin filter trends on momentum 
    365          zua(:,:,:) = ( ub(:,:,:) - zua(:,:,:) ) * z1_2dt 
    366          zva(:,:,:) = ( vb(:,:,:) - zva(:,:,:) ) * z1_2dt 
     382         zua(:,:,:) = ( puu(:,:,:) - zua(:,:,:) ) * z1_2dt 
     383         zva(:,:,:) = ( pvv(:,:,:) - zva(:,:,:) ) * z1_2dt 
    367384         CALL trd_dyn( zua, zva, jpdyn_atf, kt ) 
    368385      ENDIF 
    369386      ! 
    370       IF(ln_ctl)   CALL prt_ctl( tab3d_1=un, clinfo1=' nxt  - Un: ', mask1=umask,   & 
    371          &                       tab3d_2=vn, clinfo2=' Vn: '       , mask2=vmask ) 
     387      IF(ln_ctl)   CALL prt_ctl( tab3d_1=uu(:,:,:,kNp1), clinfo1=' nxt  - uu(:,:,:,kNp1): ', mask1=umask,   & 
     388         &                       tab3d_2=vv(:,:,:,kNp1), clinfo2=' vv(:,:,:,kNp1): '       , mask2=vmask ) 
    372389      !  
    373390      IF( ln_dynspg_ts )   DEALLOCATE( zue, zve ) 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/nemogcm.F90

    r10588 r10743  
    401401      CALL nemo_alloc() 
    402402 
     403      ! Initialise time level indices 
     404      Nm1 = 1; Nnn = 2; Np1 = 3; Nrhs = Np1 
     405      Nm1_2lev = 1; Nnn_2lev = 2 
     406 
     407      ! Initialisation of temporary pointers (to be deleted after development finished) 
     408      CALL update_pointers() 
     409 
    403410      !                             !-------------------------------! 
    404411      !                             !  NEMO general initialization  ! 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/oce.F90

    r10425 r10743  
    1717   PUBLIC oce_alloc ! routine called by nemo_init in nemogcm.F90 
    1818 
    19    !! dynamics and tracer fields                            ! before ! now    ! after  ! the after trends becomes the fields 
    20    !! --------------------------                            ! fields ! fields ! trends ! only after tra_zdf and dyn_spg 
    21    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   ub   ,  un    , ua     !: i-horizontal velocity        [m/s] 
    22    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::   vb   ,  vn    , va     !: j-horizontal velocity        [m/s] 
     19   !! dynamics and tracer fields  NOTE THAT "TARGET" ATTRIBUTE CAN BE REMOVED AFTER IMMERSE DEVELOPMENT FINISHED                             
     20   !! --------------------------                             
     21   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:), TARGET ::   uu   ,  vv     !: horizontal velocities        [m/s] 
    2322   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::           wn             !: vertical velocity            [m/s] 
    2423   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)   ::           wi             !: vertical vel. (adaptive-implicit) [m/s] 
     
    6968   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   fraqsr_1lev        !: fraction of solar net radiation absorbed in the first ocean level [-] 
    7069 
     70   !! TEMPORARY POINTERS FOR DEVELOPMENT ONLY 
     71   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:)   ::   ub   ,  un    , ua     !: i-horizontal velocity        [m/s] 
     72   REAL(wp), PUBLIC, POINTER, SAVE, DIMENSION(:,:,:)   ::   vb   ,  vn    , va     !: j-horizontal velocity        [m/s] 
     73 
    7174   !!---------------------------------------------------------------------- 
    7275   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    8487      ! 
    8588      ierr(:) = 0  
    86       ALLOCATE( ub   (jpi,jpj,jpk)      , un   (jpi,jpj,jpk)      , ua(jpi,jpj,jpk)       ,     & 
    87          &      vb   (jpi,jpj,jpk)      , vn   (jpi,jpj,jpk)      , va(jpi,jpj,jpk)       ,     &           
     89      ALLOCATE( uu   (jpi,jpj,jpk, jpt) , vv   (jpi,jpj,jpk,jpt)  ,                             & 
    8890         &      wn   (jpi,jpj,jpk)      , hdivn(jpi,jpj,jpk)      ,                             & 
    8991         &      tsb  (jpi,jpj,jpk,jpts) , tsn  (jpi,jpj,jpk,jpts) , tsa(jpi,jpj,jpk,jpts) ,     & 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/par_oce.F90

    r10068 r10743  
    5151   INTEGER, PUBLIC ::   jpj   !                                                    !: second dimension 
    5252   INTEGER, PUBLIC ::   jpk   ! = jpkglo                                           !: third  dimension 
     53   INTEGER, PUBLIC, PARAMETER :: jpt = 3                                           !: fourth dimension (time level) 
    5354   INTEGER, PUBLIC ::   jpim1 ! = jpi-1                                            !: inner domain indices 
    5455   INTEGER, PUBLIC ::   jpjm1 ! = jpj-1                                            !:   -     -      - 
  • NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/step.F90

    r10364 r10743  
    4444 
    4545   PUBLIC   stp   ! called by nemogcm.F90 
     46   PUBLIC   update_pointers ! called by nemo_init 
    4647 
    4748   !!---------------------------------------------------------------------- 
     
    276277!!jc2: dynnxt must be the latest call. e3t_b are indeed updated in that routine 
    277278                         CALL tra_nxt       ( kstp )  ! finalize (bcs) tracer fields at next time step and swap 
    278                          CALL dyn_nxt       ( kstp )  ! finalize (bcs) velocities at next time step and swap (always called after tra_nxt) 
     279                         CALL dyn_nxt       ( kstp, Nm1, Np1, uu(:,:,:,Nnn), vv(:,:,:,Nnn), e3t(:,:,:,Nnn), e3u(:,:,:,Nnn), e3v(:,:,:,Nnn) )  !  
    279280                         CALL ssh_swp       ( kstp )  ! swap of sea surface height 
    280       IF(.NOT.ln_linssh) CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors 
     281      ! 
     282      ! Swap time levels 
     283      Nrhs = Nm1 
     284      Nm1 = Nnn 
     285      Nnn = Np1 
     286      Np1 = Nrhs 
     287      ! 
     288      ! Update temporary pointers 
     289      CALL update_pointers() 
     290      ! 
     291      ! Note that 2-time-level indices don't need to be swapped because both "before" and "now" fields are derived in dom_vvl_sf_swp 
     292      IF(.NOT.ln_linssh) CALL dom_vvl_sf_swp( kstp, Nm1, Nnn, Nm1_2lev, Nnn_2lev )  ! interpolate vertical scale factors for Nnn time level 
    281293      ! 
    282294      IF( ln_diahsb  )   CALL dia_hsb       ( kstp )  ! - ML - global conservation diagnostics 
     
    329341   END SUBROUTINE stp 
    330342    
     343   SUBROUTINE update_pointers 
     344      !!---------------------------------------------------------------------- 
     345      !!                     ***  ROUTINE update_pointers  *** 
     346      !! 
     347      !! ** Purpose :   Associate temporary pointer arrays. 
     348      !!                For IMMERSE development phase only - to be deleted 
     349      !! 
     350      !! ** Method  : 
     351      !!---------------------------------------------------------------------- 
     352 
     353      ub => uu(:,:,:,Nm1); un => uu(:,:,:,Nnn); ua => uu(:,:,:,Np1) 
     354      vb => vv(:,:,:,Nm1); vn => vv(:,:,:,Nnn); va => vv(:,:,:,Np1) 
     355 
     356      e3t_b => e3t(:,:,:,Nm1); e3t_n => e3t(:,:,:,Nnn); e3t_a => e3t(:,:,:,Np1) 
     357      e3u_b => e3u(:,:,:,Nm1); e3u_n => e3u(:,:,:,Nnn); e3u_a => e3u(:,:,:,Np1) 
     358      e3v_b => e3v(:,:,:,Nm1); e3v_n => e3v(:,:,:,Nnn); e3v_a => e3v(:,:,:,Np1) 
     359 
     360      e3f_n => e3f(:,:,:) 
     361 
     362      e3w_b  => e3w (:,:,:,Nm1_2lev); e3w_n  => e3w (:,:,:,Nnn_2lev) 
     363      e3uw_b => e3uw(:,:,:,Nm1_2lev); e3uw_n => e3uw(:,:,:,Nnn_2lev) 
     364      e3vw_b => e3vw(:,:,:,Nm1_2lev); e3vw_n => e3vw(:,:,:,Nnn_2lev) 
     365 
     366   END SUBROUTINE update_pointers 
     367 
    331368   !!====================================================================== 
    332369END MODULE step 
Note: See TracChangeset for help on using the changeset viewer.