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 4479 for branches/2011/DEV_r2739_STFC_dCSE – NEMO

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

Remove jpkf as un-needed now we just reset jpk instead

Location:
branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC
Files:
9 edited

Legend:

Unmodified
Added
Removed
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DOM/domzgr.F90

    r4475 r4479  
    772772      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released 
    773773      USE wrk_nemo, ONLY:   zmbk => wrk_2d_1 
    774       USE par_oce,  ONLY:   jpkf, jpkfm1 
    775774      USE arpdebugging, ONLY: dump_array 
    776775      !! 
     
    824823               mbkmax(ji,jj) = MIN(jpk, MAX(mbkt(ji,jj)+1, mbku(ji,jj), mbkv(ji,jj))) 
    825824               ! write(*,*) narea, ': SMPDBG: ', ji, jj, mbkt(ji,jj), mbku(ji,jj), mbkv(ji,jj), mbkmax(ji,jj) 
    826 !!$               IF(mbkmax(ji,jj) /= jpkf)THEN 
    827 !!$                  WRITE (*,*) narea,': ARPDBG: mbkmax at ',ji,jj,' in {mpp_init3,zgr_bot_level} = ',mbkmax(ji,jj),jpkf 
    828 !!$               END IF 
    829825            END DO 
    830826         END DO 
     
    834830         ! mbkmax computed in mpp_init3 because needed before ANY halo  
    835831         ! swaps can be performed 
    836          jpkf = MAXVAL( mbkmax(:,:) ) 
    837          WRITE(*,*) narea,': ARPDBG: shallowest pt and jpkf = ', & 
    838                     MINVAL(mbkmax(:,:)), jpkf 
    839832         ! Play rather fast and loose and just change the value of jpk 
    840833         ! here... 
    841          jpk = jpkf 
     834         jpk = MAXVAL( mbkmax(:,:) ) 
     835         WRITE(*,*) narea,': ARPDBG: shallowest pt and new, sub-domain jpk = ', & 
     836                    MINVAL(mbkmax(:,:)), jpk 
    842837         jpkm1 = jpk - 1 
    843838      ELSE 
    844839         WRITE(*,*) narea,': ARPDBG: NOT trimming domain in z' 
    845840         mbkmax(:,:) = jpk 
    846          jpkf = jpk 
    847841      END IF 
    848       jpkfm1 = jpkf - 1 
    849842 
    850843      ! 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/dynzdf_imp.F90

    r4436 r4479  
    105105      DO jj = 2, jpjm1  
    106106         DO ji = 2, jpim1 
    107             DO jk = 1, jpkfm1 
     107            DO jk = 1, jpkm1 
    108108               zcoef = - p2dt / fse3u(ji,jj,jk) 
    109109               zzwi          = zcoef * avmu (ji,jj,jk  ) / fse3uw(ji,jj,jk  ) 
     
    119119      END DO 
    120120#else 
    121       DO jk = 1, jpkfm1        ! Matrix 
     121      DO jk = 1, jpkm1        ! Matrix 
    122122         DO jj = 2, jpjm1  
    123123            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    160160            ua(ji,jj,1) = ub(ji,jj,1) + p2dt * (  ua(ji,jj,1) + 0.5_wp * ( utau_b(ji,jj) + utau(ji,jj) )   & 
    161161               &                                                       / ( fse3u(ji,jj,1) * rau0       )  ) 
    162             DO jk = 2, jpkfm1 
     162            DO jk = 2, jpkm1 
    163163               zzwibd = zwi(ji,jj,jk) / zwd(ji,jj,jk-1) 
    164164               !==  First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1   (increasing k)  == 
     
    169169            END DO 
    170170            !==  third recurrence : SOLk = ( Lk - Uk * Ek+1 ) / Dk  == 
    171             ua(ji,jj,jpkfm1) = ua(ji,jj,jpkfm1) / zwd(ji,jj,jpkfm1) 
    172             DO jk = jpkf-2, 1, -1 
     171            ua(ji,jj,jpkm1) = ua(ji,jj,jpkm1) / zwd(ji,jj,jpkm1) 
     172            DO jk = jpk-2, 1, -1 
    173173               ua(ji,jj,jk) = ( ua(ji,jj,jk) - zws(ji,jj,jk) * ua(ji,jj,jk+1) ) / zwd(ji,jj,jk) 
    174174            END DO 
    175175            ! Normalization to obtain the general momentum trend ua 
    176             DO jk = 1, jpkfm1 
     176            DO jk = 1, jpkm1 
    177177               ua(ji,jj,jk) = ( ua(ji,jj,jk) - ub(ji,jj,jk) ) * z1_p2dt 
    178178            END DO 
     
    180180      END DO 
    181181#else 
    182       DO jk = 2, jpkfm1        !==  First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1   (increasing k)  == 
     182      DO jk = 2, jpkm1        !==  First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1   (increasing k)  == 
    183183         DO jj = 2, jpjm1    
    184184            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    194194         END DO 
    195195      END DO 
    196       DO jk = 2, jpkfm1 
     196      DO jk = 2, jpkm1 
    197197         DO jj = 2, jpjm1    
    198198            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    205205      DO jj = 2, jpjm1        !==  third recurrence : SOLk = ( Lk - Uk * Ek+1 ) / Dk  == 
    206206         DO ji = fs_2, fs_jpim1   ! vector opt. 
    207             ua(ji,jj,jpkfm1) = ua(ji,jj,jpkfm1) / zwd(ji,jj,jpkfm1) 
    208          END DO 
    209       END DO 
    210       DO jk = jpkf-2, 1, -1 
     207            ua(ji,jj,jpkm1) = ua(ji,jj,jpkm1) / zwd(ji,jj,jpkm1) 
     208         END DO 
     209      END DO 
     210      DO jk = jpk-2, 1, -1 
    211211         DO jj = 2, jpjm1    
    212212            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    216216      END DO 
    217217      ! Normalization to obtain the general momentum trend ua 
    218       DO jk = 1, jpkfm1 
     218      DO jk = 1, jpkm1 
    219219         DO jj = 2, jpjm1    
    220220            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    237237      DO jj = 2, jpjm1    
    238238         DO ji = 2, jpim1 
    239             DO jk = 1, jpkfm1        ! Matrix 
     239            DO jk = 1, jpkm1        ! Matrix 
    240240               zcoef = -p2dt / fse3v(ji,jj,jk) 
    241241               zzwi          = zcoef * avmv (ji,jj,jk  ) / fse3vw(ji,jj,jk  ) 
     
    251251      END DO 
    252252#else 
    253       DO jk = 1, jpkfm1        ! Matrix 
     253      DO jk = 1, jpkm1        ! Matrix 
    254254         DO jj = 2, jpjm1    
    255255            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    292292            va(ji,jj,1) = vb(ji,jj,1) + p2dt * (  va(ji,jj,1) + 0.5_wp * ( vtau_b(ji,jj) + vtau(ji,jj) )   & 
    293293               &                                                       / ( fse3v(ji,jj,1) * rau0       )  ) 
    294             DO jk = 2, jpkfm1 
     294            DO jk = 2, jpkm1 
    295295               zzwibd = zwi(ji,jj,jk) / zwd(ji,jj,jk-1) 
    296296               !==  First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1   (increasing k)  == 
     
    301301            END DO 
    302302            !==  third recurrence : SOLk = ( Lk - Uk * SOLk+1 ) / Dk  == 
    303             va(ji,jj,jpkfm1) = va(ji,jj,jpkfm1) / zwd(ji,jj,jpkfm1) 
    304             DO jk = jpkf-2, 1, -1 
     303            va(ji,jj,jpkm1) = va(ji,jj,jpkm1) / zwd(ji,jj,jpkm1) 
     304            DO jk = jpk-2, 1, -1 
    305305               va(ji,jj,jk) = ( va(ji,jj,jk) - zws(ji,jj,jk) * va(ji,jj,jk+1) ) / zwd(ji,jj,jk) 
    306306            END DO 
    307307            ! Normalization to obtain the general momentum trend va 
    308             DO jk = 1, jpkfm1 
     308            DO jk = 1, jpkm1 
    309309               va(ji,jj,jk) = ( va(ji,jj,jk) - vb(ji,jj,jk) ) * z1_p2dt 
    310310            END DO 
     
    312312      END DO 
    313313#else 
    314       DO jk = 2, jpkfm1        !==  First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1   (increasing k)  == 
     314      DO jk = 2, jpkm1        !==  First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1   (increasing k)  == 
    315315         DO jj = 2, jpjm1    
    316316            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    326326         END DO 
    327327      END DO 
    328       DO jk = 2, jpkfm1 
     328      DO jk = 2, jpkm1 
    329329         DO jj = 2, jpjm1 
    330330            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    337337      DO jj = 2, jpjm1        !==  third recurrence : SOLk = ( Lk - Uk * SOLk+1 ) / Dk  == 
    338338         DO ji = fs_2, fs_jpim1   ! vector opt. 
    339             va(ji,jj,jpkfm1) = va(ji,jj,jpkfm1) / zwd(ji,jj,jpkfm1) 
    340          END DO 
    341       END DO 
    342       DO jk = jpkf-2, 1, -1 
     339            va(ji,jj,jpkm1) = va(ji,jj,jpkm1) / zwd(ji,jj,jpkm1) 
     340         END DO 
     341      END DO 
     342      DO jk = jpk-2, 1, -1 
    343343         DO jj = 2, jpjm1    
    344344            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    349349 
    350350      ! Normalization to obtain the general momentum trend va 
    351       DO jk = 1, jpkfm1 
     351      DO jk = 1, jpkm1 
    352352         DO jj = 2, jpjm1    
    353353            DO ji = fs_2, fs_jpim1   ! vector opt. 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/LBC/exchmod.F90

    r4476 r4479  
    300300    USE par_oce, ONLY: wp, jpreci, jprecj, jpim1 
    301301    USE dom_oce, ONLY: nlci, nlcj, nldi, nlei, nldj, nlej, & 
    302                        nperio, nbondi, npolj, narea, jpkf 
     302                       nperio, nbondi, npolj, narea 
    303303    USE mapcomm_mod, ONLY: Iminus, Iplus, NONE, ilbext, iubext, cyclic_bc 
    304304    USE mapcomm_mod, ONLY: trimmed, eidx, widx 
     
    359359       ! we can limit the length of our z loops to the 
    360360       ! no. of levels above the ocean floor. 
    361        IF(kdim1 == jpkorig)kdim1 = jpkf 
     361       IF(kdim1 == jpkorig)kdim1 = jpk 
    362362    ELSEIF ( PRESENT(ib3) ) THEN 
    363363#if defined key_z_first 
     
    369369       ! we can limit the length of our z loops to the 
    370370       ! no. of levels above the ocean floor. 
    371        IF(kdim1 == jpk)kdim1 = jpkf 
     371       IF(kdim1 == jpk)kdim1 = jpk 
    372372    ELSEIF ( PRESENT(b2) ) THEN 
    373373       kdim1 = SIZE(b2,dim=2) 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90

    r4458 r4479  
    172172      DO jj = 1, jpjm1           !==   i- & j-gradient of density   ==! 
    173173         DO ji = 1, jpim1 
    174             DO jk = 1, mbkmax(ji,jj) ! jpkf 
    175 #else 
    176       DO jk = 1, jpkf            !==   i- & j-gradient of density   ==! 
     174            DO jk = 1, mbkmax(ji,jj) ! jpk 
     175#else 
     176      DO jk = 1, jpk            !==   i- & j-gradient of density   ==! 
    177177         DO jj = 1, jpjm1 
    178178            DO ji = 1, fs_jpim1   ! vector opt. 
     
    202202         DO ji = 1, jpi 
    203203            zdzr(ji,jj,1) = 0._wp   !==   Local vertical density gradient at T-point   == !   (evaluated from N^2) 
    204             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
     204            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
    205205               zdzr(ji,jj,jk) = zm1_g * ( prd(ji,jj,jk) + 1._wp )              & 
    206206                  &                   * ( pn2(ji,jj,jk) + pn2(ji,jj,jk+1) ) * ( 1._wp - 0.5_wp * tmask(ji,jj,jk+1) ) 
     
    210210#else 
    211211      zdzr(:,:,1) = 0._wp        !==   Local vertical density gradient at T-point   == !   (evaluated from N^2) 
    212       DO jk = 2, jpkfm1 
     212      DO jk = 2, jpkm1 
    213213         !                                ! zdzr = d/dz(prd)= - ( prd ) / grav * mk(pn2) -- at t point 
    214214         !                                !   trick: tmask(ik  )  = 0   =>   all pn2   = 0   =>   zdzr = 0 
     
    231231      DO jj = 2, jpjm1               !* Slopes at u and v points 
    232232         DO ji = 2, jpim1 
    233             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
    234 #else 
    235       DO jk = 2, jpkfm1                            !* Slopes at u and v points 
     233            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
     234#else 
     235      DO jk = 2, jpkm1                            !* Slopes at u and v points 
    236236         DO jj = 2, jpjm1 
    237237            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    275275      DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only 
    276276         DO ji = 2, jpim1   
    277             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
    278 #else 
    279       DO jk = 2, jpkfm1 
     277            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
     278#else 
     279      DO jk = 2, jpkm1 
    280280         DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only 
    281281            DO ji = 2, jpim1   
     
    297297      DO jj = 3, jpj-2                               ! other rows 
    298298         DO ji = 2, jpim1 
    299             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
     299            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
    300300#else 
    301301         DO jj = 3, jpj-2                               ! other rows 
     
    319319      DO jj = 2, jpjm1 
    320320         DO ji = 2, jpim1 
    321             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
     321            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
    322322#else 
    323323         !                                        !* decrease along coastal boundaries 
     
    340340      DO jj = 2, jpjm1 
    341341         DO ji = 2, jpim1 
    342             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
    343 #else 
    344       DO jk = 2, jpkfm1 
     342            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
     343#else 
     344      DO jk = 2, jpkm1 
    345345         DO jj = 2, jpjm1 
    346346            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    384384      DO jj = 2, jpjm1, MAX(1, jpj-3)                           ! rows jj=2 and =jpjm1 only 
    385385         DO ji = 2, jpim1 
    386             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
    387 #else 
    388       DO jk = 2, jpkfm1 
     386            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
     387#else 
     388      DO jk = 2, jpkm1 
    389389         DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only 
    390390            DO ji = 2, jpim1 
     
    407407      DO jj = 3, jpj-2                                  ! other rows 
    408408         DO ji = 2, jpim1 
    409             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
     409            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
    410410#else 
    411411         DO jj = 3, jpj-2                               ! other rows 
     
    430430      DO jj = 2, jpjm1 
    431431         DO ji = 2, jpim1 
    432             DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
     432            DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
    433433#else 
    434434         !                                        !* decrease along coastal boundaries 
     
    450450         !                                                    ! Gibraltar Strait 
    451451         ij0 =  50   ;   ij1 =  53 
    452          ii0 =  69   ;   ii1 =  71   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
     452         ii0 =  69   ;   ii1 =  71   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
    453453         ij0 =  51   ;   ij1 =  53 
    454          ii0 =  68   ;   ii1 =  71   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
    455          ii0 =  69   ;   ii1 =  71   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
    456          ii0 =  69   ;   ii1 =  71   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
     454         ii0 =  68   ;   ii1 =  71   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
     455         ii0 =  69   ;   ii1 =  71   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
     456         ii0 =  69   ;   ii1 =  71   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
    457457         ! 
    458458         !                                                    ! Mediterranean Sea 
    459459         ij0 =  49   ;   ij1 =  56 
    460          ii0 =  71   ;   ii1 =  90   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
     460         ii0 =  71   ;   ii1 =  90   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
    461461         ij0 =  50   ;   ij1 =  56 
    462          ii0 =  70   ;   ii1 =  90   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
    463          ii0 =  71   ;   ii1 =  90   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
    464          ii0 =  71   ;   ii1 =  90   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpkf ) = 0._wp 
     462         ii0 =  70   ;   ii1 =  90   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
     463         ii0 =  71   ;   ii1 =  90   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
     464         ii0 =  71   ;   ii1 =  90   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , :jpk ) = 0._wp 
    465465      ENDIF 
    466466 
     
    547547      DO jj = 1, jpjm1 
    548548         DO ji = 1, jpim1 
    549             DO jk = 1, jpkfm1                    !==  before lateral T & S gradients at T-level jk  ==! 
    550 #else 
    551       DO jk = 1, jpkfm1                          !==  before lateral T & S gradients at T-level jk  ==! 
     549            DO jk = 1, jpkm1                    !==  before lateral T & S gradients at T-level jk  ==! 
     550#else 
     551      DO jk = 1, jpkm1                          !==  before lateral T & S gradients at T-level jk  ==! 
    552552         DO jj = 1, jpjm1 
    553553            DO ji = 1, fs_jpim1   ! vector opt. 
     
    581581            zdkt(ji,jj,1) = 0._wp            !==  before vertical T & S gradient at w-level  ==! 
    582582            zdks(ji,jj,1) = 0._wp 
    583             DO jk = 2, jpkf 
     583            DO jk = 2, jpk 
    584584               zdkt(ji,jj,jk) = ( tb(ji,jj,jk-1) - tb(ji,jj,jk) ) * tmask(ji,jj,jk) 
    585585               zdks(ji,jj,jk) = ( sb(ji,jj,jk-1) - sb(ji,jj,jk) ) * tmask(ji,jj,jk) 
     
    590590      zdkt(:,:,1) = 0._wp                    !==  before vertical T & S gradient at w-level  ==! 
    591591      zdks(:,:,1) = 0._wp 
    592       DO jk = 2, jpkf 
     592      DO jk = 2, jpk 
    593593         zdkt(:,:,jk) = ( tb(:,:,jk-1) - tb(:,:,jk) ) * tmask(:,:,jk) 
    594594         zdks(:,:,jk) = ( sb(:,:,jk-1) - sb(:,:,jk) ) * tmask(:,:,jk) 
     
    602602         DO jj = 1, jpjm1                          ! NB: not masked due to the minimum value set 
    603603            DO ji = 1, jpim1 
    604                DO jk = 1, jpkfm1                    ! done each pair of triad 
    605 #else 
    606          DO jk = 1, jpkfm1                          ! done each pair of triad 
     604               DO jk = 1, jpkm1                    ! done each pair of triad 
     605#else 
     606         DO jk = 1, jpkm1                          ! done each pair of triad 
    607607            DO jj = 1, jpjm1                       ! NB: not masked due to the minimum value set 
    608608               DO ji = 1, fs_jpim1   ! vector opt.  
     
    620620         DO jj = 1, jpj                         ! NB: not masked due to the minimum value set 
    621621            DO ji = 1, jpi 
    622                DO jk = 1, jpkfm1                    ! done each pair of triad 
    623 #else 
    624          DO jk = 1, jpkfm1                          ! done each pair of triad 
     622               DO jk = 1, jpkm1                    ! done each pair of triad 
     623#else 
     624         DO jk = 1, jpkm1                          ! done each pair of triad 
    625625            DO jj = 1, jpj                       ! NB: not masked due to the minimum value set 
    626626               DO ji = 1, jpi   ! vector opt.  
     
    682682            DO jj = 1, jpjm1 
    683683               DO ji = 1, jpim1 
    684                   DO jk = 1, jpkfm1 
    685 #else 
    686             DO jk = 1, jpkfm1 
     684                  DO jk = 1, jpkm1 
     685#else 
     686            DO jk = 1, jpkm1 
    687687               DO jj = 1, jpjm1 
    688688                  DO ji = 1, fs_jpim1   ! vector opt. 
     
    818818      DO jj = 1, jpj 
    819819         DO ji = 1, jpi 
    820             DO jk = 1, mbkmax(ji,jj) ! jpkf ! =1 inside the mixed layer, =0 otherwise 
    821 #else 
    822       DO jk = 1, jpkf                ! =1 inside the mixed layer, =0 otherwise 
     820            DO jk = 1, mbkmax(ji,jj) ! jpk ! =1 inside the mixed layer, =0 otherwise 
     821#else 
     822      DO jk = 1, jpk                ! =1 inside the mixed layer, =0 otherwise 
    823823# if defined key_vectopt_loop  
    824824         DO jj = 1, 1 
     
    987987                  DO jk = 1, mbkmax(ji,jj) 
    988988#else 
    989             DO jk = 1, jpkf 
     989            DO jk = 1, jpk 
    990990               DO jj = 2, jpjm1 
    991991                  DO ji = fs_2, fs_jpim1   ! vector opt. 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/traadv_tvd.F90

    r4452 r4479  
    157157         END DO 
    158158#else 
    159          zwx(:,:,jpkf) = 0.e0    ;    zwz(:,:,jpkf) = 0.e0 
    160          zwy(:,:,jpkf) = 0.e0    ;    zwi(:,:,jpkf) = 0.e0 
     159         zwx(:,:,jpk) = 0.e0    ;    zwz(:,:,jpk) = 0.e0 
     160         zwy(:,:,jpk) = 0.e0    ;    zwi(:,:,jpk) = 0.e0 
    161161#endif 
    162162         ! 2. upstream advection with initial mass fluxes & intermediate update 
     
    167167         DO jj = 1, jpjm1 
    168168            DO ji = 1, jpim1 
    169                DO jk = 1, mbkmax(ji,jj)-1 ! jpkfm1 
    170 #else 
    171          DO jk = 1, jpkfm1 
     169               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
     170#else 
     171         DO jk = 1, jpkm1 
    172172            DO jj = 1, jpjm1 
    173173               DO ji = 1, fs_jpim1   ! vector opt. 
     
    195195         DO jj = 1, jpj 
    196196            DO ji = 1, jpi 
    197                DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
    198 #else 
    199          DO jk = 2, jpkfm1 
     197               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
     198#else 
     199         DO jk = 2, jpkm1 
    200200            DO jj = 1, jpj 
    201201               DO ji = 1, jpi 
     
    214214         DO jj = 2, jpjm1 
    215215            DO ji = 2, jpim1 
    216                DO jk = 1, mbkmax(ji,jj)-1 ! jpkfm1 
     216               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
    217217                  z2dtt = p2dt(jk) 
    218218#else 
    219          DO jk = 1, jpkfm1 
     219         DO jk = 1, jpkm1 
    220220            z2dtt = p2dt(jk) 
    221221            DO jj = 2, jpjm1 
     
    243243         IF( l_trd )  THEN  
    244244            ! store intermediate advective trends 
    245             ztrdx(:,:,1:jpkf) = zwx(:,:,1:jpkf)   ;    ztrdy(:,:,1:jpkf) = zwy(:,:,1:jpkf)  ;   ztrdz(:,:,1:jpkf) = zwz(:,:,1:jpkf) 
     245            ztrdx(:,:,1:jpk) = zwx(:,:,1:jpk)   ;    ztrdy(:,:,1:jpk) = zwy(:,:,1:jpk)  ;   ztrdz(:,:,1:jpk) = zwz(:,:,1:jpk) 
    246246         END IF 
    247247         !                                 ! "Poleward" heat and salt transports (contribution of upstream fluxes) 
     
    258258         DO jj = 1, jpjm1 
    259259            DO ji = 1, jpim1 
    260                DO jk = 1, mbkmax(ji,jj)-1 ! jpkfm1 
    261 #else 
    262          DO jk = 1, jpkfm1 
     260               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
     261#else 
     262         DO jk = 1, jpkm1 
    263263            DO jj = 1, jpjm1 
    264264               DO ji = 1, fs_jpim1   ! vector opt. 
     
    277277            DO ji = 1, jpi 
    278278               zwz(ji,jj,1) = 0.e0   ! Surface value 
    279                DO jk = 2, mbkmax(ji,jj)-1 ! jpkfm1 
     279               DO jk = 2, mbkmax(ji,jj)-1 ! jpkm1 
    280280#else 
    281281         zwz(:,:,1) = 0.e0           ! Surface value 
    282282         ! 
    283          DO jk = 2, jpkfm1            ! Interior value 
     283         DO jk = 2, jpkm1            ! Interior value 
    284284            DO jj = 1, jpj 
    285285               DO ji = 1, jpi 
     
    309309         DO jj = 2, jpjm1 
    310310            DO ji = 2, jpim1 
    311                DO jk = 1, mbkmax(ji,jj)-1 ! jpkfm1 
    312 #else 
    313          DO jk = 1, jpkfm1 
     311               DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
     312#else 
     313         DO jk = 1, jpkm1 
    314314            DO jj = 2, jpjm1 
    315315               DO ji = fs_2, fs_jpim1   ! vector opt.   
     
    419419      END DO 
    420420#else 
    421       zbetup(:,:,jpkf) = 0._wp   ;   zbetdo(:,:,jpkf) = 0._wp 
     421      zbetup(:,:,jpk) = 0._wp   ;   zbetdo(:,:,jpk) = 0._wp 
    422422#endif 
    423423 
     
    433433      DO jj = 2, jpjm1 
    434434         DO ji = 2, jpim1 
    435             DO jk = 1, mbkmax(ji,jj)-1 ! jpkfm1 
     435            DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
    436436               ikm1 = MAX(jk-1,1) 
    437437               z2dtt = p2dt(jk) 
    438438#else 
    439       DO jk = 1, jpkfm1 
     439      DO jk = 1, jpkm1 
    440440         ikm1 = MAX(jk-1,1) 
    441441         z2dtt = p2dt(jk) 
     
    484484      DO jj = 2, jpjm1 
    485485         DO ji = 2, jpim1 
    486             DO jk = 1, mbkmax(ji,jj)-1 ! jpkfm1 
    487 #else 
    488       DO jk = 1, jpkfm1 
     486            DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
     487#else 
     488      DO jk = 1, jpkm1 
    489489         DO jj = 2, jpjm1 
    490490            DO ji = fs_2, fs_jpim1   ! vector opt. 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso.F90

    r4427 r4479  
    189189#else 
    190190         !!bug ajout.... why?   ( 1,jpj,:) and (jpi,1,:) should be sufficient.... 
    191          zdit (1,:,1:jpkf) = 0.e0     ;     zdit (jpi,:,1:jpkf) = 0.e0 
    192          zdjt (1,:,1:jpkf) = 0.e0     ;     zdjt (jpi,:,1:jpkf) = 0.e0 
     191         zdit (1,:,1:jpk) = 0.e0     ;     zdit (jpi,:,1:jpk) = 0.e0 
     192         zdjt (1,:,1:jpk) = 0.e0     ;     zdjt (jpi,:,1:jpk) = 0.e0 
    193193         !!end 
    194          DO jk = 1, jpkfm1 ! jpkm1 
     194         DO jk = 1, jpkm1 ! jpkm1 
    195195            DO jj = 1, jpjm1 
    196196               DO ji = 1, fs_jpim1   ! vector opt. 
     
    326326!CDIR PARALLEL DO PRIVATE( zdk1t )  
    327327         !                                                ! =============== 
    328          DO jk = 1, jpkfm1 ! jpkm1                        ! Horizontal slab 
     328         DO jk = 1, jpkm1 ! jpkm1                        ! Horizontal slab 
    329329            !                                             ! =============== 
    330330            ! 1. Vertical tracer gradient at level jk and jk+1 
     
    395395                  DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
    396396#else 
    397             DO jk = 1, jpkfm1 ! jpkm1 
     397            DO jk = 1, jpkm1 ! jpkm1 
    398398               DO jj = 2, jpjm1 
    399399                  DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    412412                  DO jk = 1, mbkmax(ji,jj)-1 ! jpkm1 
    413413#else 
    414             DO jk = 1, jpkfm1 ! jpkm1 
     414            DO jk = 1, jpkm1 ! jpkm1 
    415415               DO jj = 2, jpjm1 
    416416                  DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    441441         ztfw(:,:,:) = 0.0_wp 
    442442#else 
    443          ztfw(1,:,1:jpkf) = 0.e0     ;     ztfw(jpi,:,1:jpkf) = 0.e0 
    444          ztfw(:,:, 1 ) = 0.e0      ;      ztfw(:,:,jpkf) = 0.e0 
     443         ztfw(1,:,1:jpk) = 0.e0     ;     ztfw(jpi,:,1:jpk) = 0.e0 
     444         ztfw(:,:, 1 ) = 0.e0      ;      ztfw(:,:,jpk) = 0.e0 
    445445#endif 
    446446 
     
    451451               DO jk = 2, mbkmax(ji,jj)-1 
    452452#else 
    453          DO jk = 2, jpkfm1 
     453         DO jk = 2, jpkm1 
    454454            DO jj = 2, jpjm1 
    455455               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    481481               DO jk = 1, mbkmax(ji,jj)-1 
    482482#else 
    483          DO jk = 1, jpkfm1 
     483         DO jk = 1, jpkm1 
    484484            DO jj = 2, jpjm1 
    485485               DO ji = fs_2, fs_jpim1   ! vector opt. 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90

    r4425 r4479  
    287287            DO ji = 1, jpi 
    288288               zpelc(ji,jj,1) =  MAX( rn2b(ji,jj,1), 0._wp ) * fsdepw(ji,jj,1) * fse3w(ji,jj,1) 
    289                DO jk = 2, jpkf 
     289               DO jk = 2, jpk 
    290290                  zpelc(ji,jj,jk)  = zpelc(ji,jj,jk-1) + MAX( rn2b(ji,jj,jk), 0._wp ) * fsdepw(ji,jj,jk) * fse3w(ji,jj,jk) 
    291291               END DO 
     
    294294#else 
    295295         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * fsdepw(:,:,1) * fse3w(:,:,1) 
    296          DO jk = 2, jpkf 
     296         DO jk = 2, jpk 
    297297            zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * fsdepw(:,:,jk) * fse3w(:,:,jk) 
    298298         END DO 
     
    305305            DO ji = 1, jpi               !      with us=0.016*wind(starting from jpk-1) 
    306306               zus  = zcof * taum(ji,jj) 
    307                DO jk = jpkfm1, 2, -1 
    308 #else 
    309          DO jk = jpkfm1, 2, -1 
     307               DO jk = jpkm1, 2, -1 
     308#else 
     309         DO jk = jpkm1, 2, -1 
    310310            DO jj = 1, jpj               ! Last w-level at which zpelc>=0.5*us*us  
    311311               DO ji = 1, jpi            !      with us=0.016*wind(starting from jpk-1) 
     
    332332            DO ji = 2, jpim1 
    333333               zus  = zcof * SQRT( taum(ji,jj) )           ! Stokes drift 
    334                DO jk = 2, jpkfm1 
    335 #else 
    336 !CDIR NOVERRCHK 
    337          DO jk = 2, jpkfm1         !* TKE Langmuir circulation source term added to en 
     334               DO jk = 2, jpkm1 
     335#else 
     336!CDIR NOVERRCHK 
     337         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en 
    338338!CDIR NOVERRCHK 
    339339            DO jj = 2, jpjm1 
     
    365365      DO jj = 1, jpj 
    366366         DO ji = 1, jpi 
    367             DO jk = 2, jpkfm1 
    368 #else 
    369       DO jk = 2, jpkfm1           !* Shear production at uw- and vw-points (energy conserving form) 
     367            DO jk = 2, jpkm1 
     368#else 
     369      DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form) 
    370370         DO jj = 1, jpj                 ! here avmu, avmv used as workspace 
    371371            DO ji = 1, jpi 
     
    387387      DO jj = 2, jpjm1 
    388388         DO ji = 2, jpim1 
    389             DO jk = 2, jpkfm1     !* Matrix and right hand side in en 
    390 #else 
    391       DO jk = 2, jpkfm1           !* Matrix and right hand side in en 
     389            DO jk = 2, jpkm1     !* Matrix and right hand side in en 
     390#else 
     391      DO jk = 2, jpkm1           !* Matrix and right hand side in en 
    392392         DO jj = 2, jpjm1 
    393393            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    417417         DO ji = 2, jpim1 
    418418            ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 
    419             DO jk = 3, jpkfm1 
     419            DO jk = 3, jpkm1 
    420420               zdiag(ji,jj,jk) = zdiag(ji,jj,jk) - zd_lw(ji,jj,jk) * zd_up(ji,jj,jk-1) / zdiag(ji,jj,jk-1) 
    421421            END DO 
    422422            ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 
    423423            zd_lw(ji,jj,2) = en(ji,jj,2) - zd_lw(ji,jj,2) * en(ji,jj,1)    ! Surface boudary conditions on tke 
    424             DO jk = 3, jpkfm1 
     424            DO jk = 3, jpkm1 
    425425               zd_lw(ji,jj,jk) = en(ji,jj,jk) - zd_lw(ji,jj,jk) / zdiag(ji,jj,jk-1) *zd_lw(ji,jj,jk-1) 
    426426            END DO 
    427427            ! Third recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk 
    428             en(ji,jj,jpkfm1) = zd_lw(ji,jj,jpkfm1) / zdiag(ji,jj,jpkfm1) 
    429             DO jk = jpkf-2, 2, -1 
     428            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1) 
     429            DO jk = jpk-2, 2, -1 
    430430               en(ji,jj,jk) = ( zd_lw(ji,jj,jk) - zd_up(ji,jj,jk) * en(ji,jj,jk+1) ) / zdiag(ji,jj,jk) 
    431431            END DO 
    432             DO jk = 2, jpkfm1                       ! set the minimum value of tke 
     432            DO jk = 2, jpkm1                       ! set the minimum value of tke 
    433433               en(ji,jj,jk) = MAX( en(ji,jj,jk), rn_emin ) * tmask(ji,jj,jk) 
    434434            END DO 
     
    436436      END DO 
    437437#else 
    438       DO jk = 3, jpkfm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 
     438      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 
    439439         DO jj = 2, jpjm1 
    440440            DO ji = fs_2, fs_jpim1    ! vector opt. 
     
    448448         END DO 
    449449      END DO 
    450       DO jk = 3, jpkfm1 
     450      DO jk = 3, jpkm1 
    451451         DO jj = 2, jpjm1 
    452452            DO ji = fs_2, fs_jpim1    ! vector opt. 
     
    457457      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk 
    458458         DO ji = fs_2, fs_jpim1    ! vector opt. 
    459             en(ji,jj,jpkfm1) = zd_lw(ji,jj,jpkfm1) / zdiag(ji,jj,jpkfm1) 
    460          END DO 
    461       END DO 
    462       DO jk = jpkf-2, 2, -1 
     459            en(ji,jj,jpkm1) = zd_lw(ji,jj,jpkm1) / zdiag(ji,jj,jpkm1) 
     460         END DO 
     461      END DO 
     462      DO jk = jpk-2, 2, -1 
    463463         DO jj = 2, jpjm1 
    464464            DO ji = fs_2, fs_jpim1    ! vector opt. 
     
    467467         END DO 
    468468      END DO 
    469       DO jk = 2, jpkfm1                             ! set the minimum value of tke 
     469      DO jk = 2, jpkm1                             ! set the minimum value of tke 
    470470         DO jj = 2, jpjm1 
    471471            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    483483         DO jj = 2, jpjm1 
    484484            DO ji = 2, jpim1 
    485                DO jk = 2, jpkfm1 
    486 #else 
    487          DO jk = 2, jpkfm1 
     485               DO jk = 2, jpkm1 
     486#else 
     487         DO jk = 2, jpkm1 
    488488            DO jj = 2, jpjm1 
    489489               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    507507         !!            unless we also make zdif a 2-d (jpi,jpj) array 
    508508!CDIR NOVERRCHK 
    509          DO jk = 2, jpkfm1 
     509         DO jk = 2, jpkm1 
    510510!CDIR NOVERRCHK 
    511511            DO jj = 2, jpjm1 
     
    601601      DO jj = 2, jpjm1 
    602602         DO ji = 2, jpim1 
    603             zmxlm(ji,jj,jpkf) = rmxl_min     ! last level set to the interior minium value 
    604             DO jk = 2, jpkfm1        ! interior value : l=sqrt(2*e/n^2) 
    605 #else 
    606       zmxlm(:,:,jpkf)  = rmxl_min    ! last level set to the interior minium value 
    607       ! 
    608 !CDIR NOVERRCHK 
    609       DO jk = 2, jpkfm1              ! interior value : l=sqrt(2*e/n^2) 
     603            zmxlm(ji,jj,jpk) = rmxl_min     ! last level set to the interior minium value 
     604            DO jk = 2, jpkm1        ! interior value : l=sqrt(2*e/n^2) 
     605#else 
     606      zmxlm(:,:,jpk)  = rmxl_min    ! last level set to the interior minium value 
     607      ! 
     608!CDIR NOVERRCHK 
     609      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2) 
    610610!CDIR NOVERRCHK 
    611611         DO jj = 2, jpjm1 
     
    622622      ! 
    623623      zmxld(:,:, 1  ) = zmxlm(:,:,1)   ! surface set to the zmxlm   value 
    624       zmxld(:,:,jpkf) = rmxl_min       ! last level  set to the minimum value 
     624      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value 
    625625      ! 
    626626      SELECT CASE ( nn_mxl ) 
     
    630630         DO jj = 2, jpjm1 
    631631            DO ji = 2, jpim1 
    632                DO jk = 2, jpkfm1 
    633 #else 
    634          DO jk = 2, jpkfm1 
     632               DO jk = 2, jpkm1 
     633#else 
     634         DO jk = 2, jpkm1 
    635635            DO jj = 2, jpjm1 
    636636               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    648648         DO jj = 2, jpjm1 
    649649            DO ji = 2, jpim1 
    650                DO jk = 2, jpkfm1 
    651 #else 
    652          DO jk = 2, jpkfm1 
     650               DO jk = 2, jpkm1 
     651#else 
     652         DO jk = 2, jpkm1 
    653653            DO jj = 2, jpjm1 
    654654               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    665665         DO jj = 2, jpjm1 
    666666            DO ji = 2, jpim1 
    667                DO jk = 2, jpkfm1   ! from the surface to the bottom : 
     667               DO jk = 2, jpkm1   ! from the surface to the bottom : 
    668668                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) ) 
    669669               END DO 
    670                DO jk = jpkfm1, 2, -1     ! from the bottom to the surface : 
     670               DO jk = jpkm1, 2, -1     ! from the bottom to the surface : 
    671671                  zemxl = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) ) 
    672672                  zmxlm(ji,jj,jk) = zemxl 
     
    676676         END DO 
    677677#else 
    678          DO jk = 2, jpkfm1         ! from the surface to the bottom : 
     678         DO jk = 2, jpkm1         ! from the surface to the bottom : 
    679679            DO jj = 2, jpjm1 
    680680               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    683683            END DO 
    684684         END DO 
    685          DO jk = jpkfm1, 2, -1     ! from the bottom to the surface : 
     685         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : 
    686686            DO jj = 2, jpjm1 
    687687               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    698698         DO jj = 2, jpjm1 
    699699            DO ji = 2, jpim1 
    700                DO jk = 2, jpkfm1         ! from the surface to the bottom : lup 
     700               DO jk = 2, jpkm1         ! from the surface to the bottom : lup 
    701701                  zmxld(ji,jj,jk) = MIN( zmxld(ji,jj,jk-1) + fse3t(ji,jj,jk-1), zmxlm(ji,jj,jk) ) 
    702702               END DO 
    703                DO jk = jpkfm1, 2, -1     ! from the bottom to the surface : ldown 
     703               DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown 
    704704                  zmxlm(ji,jj,jk) = MIN( zmxlm(ji,jj,jk+1) + fse3t(ji,jj,jk+1), zmxlm(ji,jj,jk) ) 
    705705               END DO 
    706                DO jk = 2, jpkfm1 
     706               DO jk = 2, jpkm1 
    707707                  zemlm = MIN ( zmxld(ji,jj,jk),  zmxlm(ji,jj,jk) ) 
    708708                  zemlp = SQRT( zmxld(ji,jj,jk) * zmxlm(ji,jj,jk) ) 
     
    713713         END DO 
    714714#else 
    715          DO jk = 2, jpkfm1         ! from the surface to the bottom : lup 
     715         DO jk = 2, jpkm1         ! from the surface to the bottom : lup 
    716716            DO jj = 2, jpjm1 
    717717               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    720720            END DO 
    721721         END DO 
    722          DO jk = jpkfm1, 2, -1     ! from the bottom to the surface : ldown 
     722         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown 
    723723            DO jj = 2, jpjm1 
    724724               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    728728         END DO 
    729729!CDIR NOVERRCHK 
    730          DO jk = 2, jpkfm1 
     730         DO jk = 2, jpkm1 
    731731!CDIR NOVERRCHK 
    732732            DO jj = 2, jpjm1 
     
    755755      DO jj = 2, jpjm1 
    756756         DO ji = 2, jpim1 
    757             DO jk = 1, jpkfm1      !* vertical eddy viscosity & diffivity at w-points 
    758 #else 
    759 !CDIR NOVERRCHK 
    760       DO jk = 1, jpkfm1            !* vertical eddy viscosity & diffivity at w-points 
     757            DO jk = 1, jpkm1      !* vertical eddy viscosity & diffivity at w-points 
     758#else 
     759!CDIR NOVERRCHK 
     760      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points 
    761761!CDIR NOVERRCHK 
    762762         DO jj = 2, jpjm1 
     
    777777      DO jj = 2, jpjm1 
    778778         DO ji = 2, jpim1 
    779             DO jk = 2, jpkfm1      !* vertical eddy viscosity at u- and v-points 
    780 #else 
    781       DO jk = 2, jpkfm1            !* vertical eddy viscosity at u- and v-points 
     779            DO jk = 2, jpkm1      !* vertical eddy viscosity at u- and v-points 
     780#else 
     781      DO jk = 2, jpkm1            !* vertical eddy viscosity at u- and v-points 
    782782         DO jj = 2, jpjm1 
    783783            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    795795         DO jj = 2, jpjm1 
    796796            DO ji = 2, jpim1 
    797                DO jk = 2, jpkfm1 
    798 #else 
    799          DO jk = 2, jpkfm1 
     797               DO jk = 2, jpkm1 
     798#else 
     799         DO jk = 2, jpkm1 
    800800            DO jj = 2, jpjm1 
    801801               DO ji = fs_2, fs_jpim1   ! vector opt. 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/nemogcm.F90

    r4423 r4479  
    242242         jpkm1 = jpk-1                                            ! inner domain indices 
    243243         jpkorig = jpk                        ! Copy of jpk that is NOT modified 
    244          jpkf    = jpk                        ! Max depth of this sub-domain. Initially set to jpk here 
    245                                               ! but altered later in domzgr 
    246244      ENDIF 
    247245 
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/par_oce.F90

    r4421 r4479  
    177177   INTEGER, PUBLIC  ::   jpj   ! = ( jpjglo-2*jprecj + (jpnj-1) ) / jpnj + 2*jprecj   !: second dimension 
    178178   INTEGER, PUBLIC  ::   jpk   ! = jpkdta                                             !: third dimension 
    179    INTEGER, PUBLIC  ::   jpkf  ! <= jpk                                             !: Max wet level of MPP subdomain 
    180    INTEGER, PUBLIC  ::   jpkorig ! = jpk before it is reset to jpkf 
     179   INTEGER, PUBLIC  ::   jpkorig ! = jpk before it is reduced to a sub-domain-local value 
    181180   INTEGER, PUBLIC  ::   jpim1 ! = jpi-1                                            !: inner domain indices 
    182181   INTEGER, PUBLIC  ::   jpjm1 ! = jpj-1                                            !:   -     -      - 
    183182   INTEGER, PUBLIC  ::   jpkm1 ! = jpk-1                                            !:   -     -      - 
    184    INTEGER, PUBLIC  ::   jpkfm1 ! = jpkf-1                                          !:   -     -      - 
    185183   INTEGER, PUBLIC  ::   jpij  ! = jpi*jpj                                          !:  jpi x jpj 
    186184 
Note: See TracChangeset for help on using the changeset viewer.