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 4427 for branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

Ignore:
Timestamp:
2014-02-04T13:14:00+01:00 (10 years ago)
Author:
trackstand2
Message:

First files changed on last FINISS work package. Stephen's work although
commited by Andy P.

Location:
branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/dynkeg.F90

    r3432 r4427  
    8787      DO jj = 2, jpj            ! Horizontal kinetic energy at T-point 
    8888         DO ji = 2, jpi 
    89             DO jk = 1, jpkm1 
     89            DO jk = 1, mbkmax(ji,jj)-1 
    9090               zhke(ji,jj,jk) = 0.25 * (   un(ji-1,jj  ,jk) * un(ji-1,jj  ,jk)   & 
    9191                  &                      + un(ji  ,jj  ,jk) * un(ji  ,jj  ,jk)   & 
    92                                          + vn(ji  ,jj-1,jk) * vn(ji  ,jj-1,jk)   & 
     92                  &                      + vn(ji  ,jj-1,jk) * vn(ji  ,jj-1,jk)   & 
    9393                  &                      + vn(ji  ,jj  ,jk) * vn(ji  ,jj  ,jk)   ) 
    9494            END DO   
     
    9797      DO jj = 2, jpjm1          ! add the gradient of kinetic energy to the general momentum trends 
    9898         DO ji = 2, jpim1 
    99             DO jk = 1, jpkm1 
     99            DO jk = 1, mbkmax(ji,jj)-1 
    100100               ua(ji,jj,jk) = ua(ji,jj,jk) - ( zhke(ji+1,jj  ,jk) - zhke(ji,jj,jk) ) / e1u(ji,jj) 
    101101               va(ji,jj,jk) = va(ji,jj,jk) - ( zhke(ji  ,jj+1,jk) - zhke(ji,jj,jk) ) / e2v(ji,jj) 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/dynldf_lap.F90

    r3211 r4427  
    8383      DO jj = 2, jpjm1 
    8484         DO ji = 2, jpim1 
    85             DO jk = 1, jpkm1 
     85            DO jk = 1, mbkmax(ji,jj)-1 
    8686#else 
    8787      !                                                ! =============== 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/dynzad.F90

    r3837 r4427  
    9999            zwuw(ji,jj, 1 ) = 0.e0     ! Surface values set to zero 
    100100            zwvw(ji,jj, 1 ) = 0.e0 
    101             DO jk = 2, jpkm1 
     101            DO jk = 2, mbkmax(ji,jj)-1 
    102102               zwuw(ji,jj,jk) =   ( zww(ji+1,jj  )*wn(ji+1,jj  ,jk) + zww(ji,jj)*wn(ji,jj,jk) )   & 
    103103                  &             * ( un(ji,jj,jk-1)-un(ji,jj,jk) )  
     
    105105                  &             * ( vn(ji,jj,jk-1)-vn(ji,jj,jk) ) 
    106106            END DO   
    107             zwuw(ji,jj,jpk) = 0.e0     ! Bottom values set to zero 
    108             zwvw(ji,jj,jpk) = 0.e0 
     107            zwuw(ji,jj,mbkmax(ji,jj)) = 0.e0     ! Bottom values set to zero 
     108            zwvw(ji,jj,mbkmax(ji,jj)) = 0.e0 
    109109         END DO    
    110110      END DO 
     
    136136      DO jj = 2, jpjm1              ! Vertical momentum advection at u- and v-points 
    137137         DO ji = 2, jpim1 
    138             DO jk = 1, jpkm1 
     138            DO jk = 1, mbkmax(ji,jj)-1 
    139139#else 
    140140      DO jk = 1, jpkm1              ! Vertical momentum advection at u- and v-points 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90

    r4415 r4427  
    9494      ! 
    9595      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
     96#if defined key_z_first 
     97      INTEGER  ::   klim         ! upper bound on k loop 
     98#endif 
    9699      REAL(wp) ::   zcoefu, zcoefv, zcoeff, z2dt, z1_2dt, z1_rau0   ! local scalars 
    97100      !!---------------------------------------------------------------------- 
     
    107110         IF(lwp) WRITE(numout,*) '~~~~~~~ ' 
    108111         ! 
     112#if defined key_z_first 
     113         DO jj=1,jpj 
     114            DO ji=1,jpi 
     115               DO jk=mbkmax(ji,jj), jpk 
     116                  wn(ji,jj,jk) = 0._wp 
     117               END DO 
     118            END DO 
     119         END DO 
     120#else 
    109121         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all) 
     122#endif 
    110123         ! 
    111124         IF( lk_vvl ) THEN                    ! before and now Sea SSH at u-, v-, f-points (vvl case only) 
     
    161174            ! gdept_1(1:jpkm1,:,:) = (gdept(1:jpkm1,:,:)*(1.+sshn(:,:)*mut(1:jpkm1,:,:)))  
    162175            ! which contains non-conforming array expressions. 
    163          DO jj=1,jpj,1 
    164             DO ji=1,jpi,1 
    165                DO jk=1,jpk,1 
    166                   fsdept(ji,jj,jk) = fsdept_n(ji,jj,jk)   ! now local depths stored in fsdep. arrays 
    167                END DO 
    168             END DO 
    169          END DO 
    170          DO jj=1,jpj,1 
    171             DO ji=1,jpi,1 
    172                DO jk=1,jpk,1 
    173                   fsdepw(ji,jj,jk) = fsdepw_n(ji,jj,jk) 
    174                END DO 
    175             END DO 
    176          END DO 
    177          DO jj=1,jpj,1 
    178             DO ji=1,jpi,1 
    179                DO jk=1,jpk,1 
    180                   fsde3w(ji,jj,jk) = fsde3w_n(ji,jj,jk) 
    181                END DO 
    182             END DO 
    183          END DO 
    184             ! 
    185          DO jj=1,jpj,1 
    186             DO ji=1,jpi,1 
    187                DO jk=1,jpk,1 
    188                   fse3t (ji,jj,jk) = fse3t_n (ji,jj,jk)   ! vertical scale factors stored in fse3. arrays 
    189                END DO 
    190             END DO 
    191          END DO 
    192          DO jj=1,jpj,1 
    193             DO ji=1,jpi,1 
    194                DO jk=1,jpk,1 
    195                   fse3u (ji,jj,jk) = fse3u_n (ji,jj,jk) 
    196                END DO 
    197             END DO 
    198          END DO 
    199          DO jj=1,jpj,1 
    200             DO ji=1,jpi,1 
    201                DO jk=1,jpk,1 
    202                   fse3v (ji,jj,jk) = fse3v_n (ji,jj,jk) 
    203                END DO 
    204             END DO 
    205          END DO 
    206          DO jj=1,jpj,1 
    207             DO ji=1,jpi,1 
    208                DO jk=1,jpk,1 
    209                   fse3f (ji,jj,jk) = fse3f_n (ji,jj,jk) 
    210                END DO 
    211             END DO 
    212          END DO 
    213          DO jj=1,jpj,1 
    214             DO ji=1,jpi,1 
    215                DO jk=1,jpk,1 
    216                   fse3w (ji,jj,jk) = fse3w_n (ji,jj,jk) 
    217                END DO 
    218             END DO 
    219          END DO 
    220  
    221  
    222          DO jj=1,jpj,1 
    223             DO ji=1,jpi,1 
    224                DO jk=1,jpk,1 
    225                   fse3uw(ji,jj,jk) = fse3uw_n(ji,jj,jk) 
    226                END DO 
    227             END DO 
    228          END DO 
    229  
    230          DO jj=1,jpj,1 
    231             DO ji=1,jpi,1 
    232                DO jk=1,jpk,1 
    233                   fse3vw(ji,jj,jk) = fse3vw_n(ji,jj,jk) 
    234                END DO 
     176         DO jj=1,jpj 
     177            DO ji=1,jpi 
     178               klim=mbkmax(ji,jj) 
     179               ! now local depths stored in fsdep. arrays 
     180               fsdept(ji,jj,1:klim) = fsdept_n(ji,jj,1:klim) 
     181               fsdepw(ji,jj,1:klim) = fsdepw_n(ji,jj,1:klim) 
     182               fsde3w(ji,jj,1:klim) = fsde3w_n(ji,jj,1:klim) 
     183               ! vertical scale factors stored in fse3. arrays 
     184               fse3t (ji,jj,1:klim) = fse3t_n (ji,jj,1:klim) 
     185               fse3u (ji,jj,1:klim) = fse3u_n (ji,jj,1:klim) 
     186               fse3v (ji,jj,1:klim) = fse3v_n (ji,jj,1:klim) 
     187               fse3f (ji,jj,1:klim) = fse3f_n (ji,jj,1:klim) 
     188               fse3w (ji,jj,1:klim) = fse3w_n (ji,jj,1:klim) 
     189               fse3uw(ji,jj,1:klim) = fse3uw_n(ji,jj,1:klim) 
     190               fse3vw(ji,jj,1:klim) = fse3vw_n(ji,jj,1:klim) 
    235191            END DO 
    236192         END DO 
     
    279235      DO jj = 1, jpj 
    280236         DO ji = 1, jpi 
    281             DO jk = 1, jpkm1                           ! Horizontal divergence of barotropic transports 
     237            DO jk = 1, mbkmax(ji,jj)-1                 ! Horizontal divergence of barotropic transports 
    282238               zhdiv(ji,jj) = zhdiv(ji,jj) + fse3t(ji,jj,jk) * hdivn(ji,jj,jk) 
    283239            END DO 
     
    355311      DO jj = 1, jpj 
    356312         DO ji = 1, jpi 
    357             DO jk = jpkm1, 1, -1                      ! integrate from the bottom the hor. divergence 
     313            DO jk = mbkmax(ji,jj)-1, 1, -1                      ! integrate from the bottom the hor. divergence 
    358314                wn(ji,jj,jk) = wn(ji,jj,jk+1)                               & 
    359315                   &         -   fse3t_n(ji,jj,jk) * hdivn(ji,jj,jk)        & 
     
    390346         DO jj = 1, jpj 
    391347            DO ji = 1, jpi 
    392                DO jk = 1, jpk 
     348               DO jk = 1, mbkmax(ji,jj) 
    393349                  z3d(ji,jj,jk) = wn(ji,jj,jk) * z2d(ji,jj) 
    394350               END DO 
Note: See TracChangeset for help on using the changeset viewer.