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 3211 for branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/traadv_qck.F90 – NEMO

Ignore:
Timestamp:
2011-12-11T16:00:26+01:00 (12 years ago)
Author:
spickles2
Message:

Stephen Pickles, 11 Dec 2011

Commit to bring the rest of the DCSE NEMO development branch
in line with the latest development version. This includes
array index re-ordering of all OPA_SRC/.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/traadv_qck.F90

    r2715 r3211  
    3535   REAL(wp) :: r1_6 = 1./ 6.   ! 1/6 ratio 
    3636 
     37   !! * Control permutation of array indices 
     38#  include "oce_ftrans.h90" 
     39#  include "dom_oce_ftrans.h90" 
     40#  include "trc_oce_ftrans.h90" 
     41 
    3742   !! * Substitutions 
    3843#  include "domzgr_substitute.h90" 
     
    8590      INTEGER                              , INTENT(in   ) ::   kjpt            ! number of tracers 
    8691      REAL(wp), DIMENSION(        jpk     ), INTENT(in   ) ::   p2dt            ! vertical profile of tracer time-step 
    87       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun, pvn, pwn   ! 3 ocean velocity components 
    88       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn        ! before and now tracer fields 
    89       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta             ! tracer trend  
     92 
     93      !! DCSE_NEMO: This style defeats ftrans 
     94!     REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun, pvn, pwn   ! 3 ocean velocity components 
     95!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn        ! before and now tracer fields 
     96!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta             ! tracer trend  
     97 
     98!FTRANS pun pvn pwn :I :I :z 
     99!FTRANS ptb ptn :I :I :z : 
     100!FTRANS pta :I :I :z : 
     101      REAL(wp), INTENT(in   ) ::   pun(jpi,jpj,jpk)        ! ocean velocity component (u) 
     102      REAL(wp), INTENT(in   ) ::   pvn(jpi,jpj,jpk)        ! ocean velocity component (v) 
     103      REAL(wp), INTENT(in   ) ::   pwn(jpi,jpj,jpk)        ! ocean velocity component (w) 
     104      REAL(wp), INTENT(in   ) ::   ptb(jpi,jpj,jpk,kjpt)   ! tracer fields (before) 
     105      REAL(wp), INTENT(in   ) ::   ptn(jpi,jpj,jpk,kjpt)   ! tracer fields (now) 
     106      REAL(wp), INTENT(inout) ::   pta(jpi,jpj,jpk,kjpt)   ! tracer trend  
     107 
    90108      !!---------------------------------------------------------------------- 
    91109 
     
    107125      CALL tra_adv_cen2_k( kt, cdtype, pwn,         ptn, pta, kjpt ) 
    108126      ! 
     127 
     128      !! * Reset control of array index permutation 
     129!FTRANS CLEAR 
     130#  include "oce_ftrans.h90" 
     131#  include "dom_oce_ftrans.h90" 
     132#  include "trc_oce_ftrans.h90" 
     133 
    109134   END SUBROUTINE tra_adv_qck 
    110135 
     
    118143      USE oce     , ONLY:   zwx => ua       ! ua used as workspace 
    119144      USE wrk_nemo, ONLY:   zfu => wrk_3d_1 , zfc => wrk_3d_2, zfd => wrk_3d_3   ! 3D workspace 
     145 
     146      !! DCSE_NEMO: need additional directives for renamed module variables 
     147!FTRANS zwx zfu zfc zfd :I :I :z 
     148 
    120149      ! 
    121150      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index 
     
    123152      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers 
    124153      REAL(wp), DIMENSION(        jpk     ), INTENT(in   ) ::   p2dt       ! vertical profile of tracer time-step 
    125       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun        ! i-velocity components 
    126       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn   ! before and now tracer fields 
    127       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend  
     154 
     155      !! DCSE_NEMO: This style defeats ftrans 
     156!     REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun        ! i-velocity components 
     157!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn   ! before and now tracer fields 
     158!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend  
     159 
     160!FTRANS pun :I :I :z 
     161!FTRANS ptb ptn pta :I :I :z : 
     162      REAL(wp), INTENT(in   ) ::   pun(jpi,jpj,jpk)         ! i-velocity component 
     163      REAL(wp), INTENT(in   ) ::   ptb(jpi,jpj,jpk,kjpt)    ! tracer field (before) 
     164      REAL(wp), INTENT(in   ) ::   ptn(jpi,jpj,jpk,kjpt)    ! tracer field (now) 
     165      REAL(wp), INTENT(inout) ::   pta(jpi,jpj,jpk,kjpt)    ! tracer trend  
     166 
    128167      !! 
    129168      INTEGER  :: ji, jj, jk, jn   ! dummy loop indices 
     
    140179         zfd(:,:,:) = 0.0     ;   zwx(:,:,:) = 0.0      
    141180         !                                                   
     181#if defined key_z_first 
     182         !--- Computation of the upstream and downstream value of the tracer and the mask 
     183         DO jj = 2, jpjm1 
     184            DO ji = 2, jpim1 
     185               DO jk = 1, jpkm1                                 
     186#else 
    142187         DO jk = 1, jpkm1                                 
    143188            !                                              
    144             !--- Computation of the ustream and downstream value of the tracer and the mask 
     189            !--- Computation of the upstream and downstream value of the tracer and the mask 
    145190            DO jj = 2, jpjm1 
    146191               DO ji = fs_2, fs_jpim1   ! vector opt. 
     192#endif 
    147193                  ! Upstream in the x-direction for the tracer 
    148194                  zfc(ji,jj,jk) = ptb(ji-1,jj,jk,jn) 
     
    158204         ! --------------------------- 
    159205         ! 
     206#if defined key_z_first 
     207         DO jj = 2, jpjm1 
     208            DO ji = 2, jpim1 
     209               DO jk = 1, jpkm1                              
     210#else 
    160211         DO jk = 1, jpkm1                              
    161212            DO jj = 2, jpjm1 
    162213               DO ji = fs_2, fs_jpim1   ! vector opt.          
     214#endif 
    163215                  zdir = 0.5 + SIGN( 0.5, pun(ji,jj,jk) )   ! if pun > 0 : zdir = 1 otherwise zdir = 0  
    164216                  zfu(ji,jj,jk) = zdir * zfc(ji,jj,jk ) + ( 1. - zdir ) * zfd(ji+1,jj,jk)  ! FU in the x-direction for T  
     
    167219         END DO 
    168220         ! 
     221#if defined key_z_first 
     222         DO jj = 2, jpjm1 
     223            DO ji = 2, jpim1 
     224               DO jk = 1, jpkm1   
     225                  zdt =  p2dt(jk) 
     226#else 
    169227         DO jk = 1, jpkm1   
    170228            zdt =  p2dt(jk) 
    171229            DO jj = 2, jpjm1 
    172230               DO ji = fs_2, fs_jpim1   ! vector opt.    
     231#endif 
    173232                  zdir = 0.5 + SIGN( 0.5, pun(ji,jj,jk) )   ! if pun > 0 : zdir = 1 otherwise zdir = 0  
    174233                  zdx  = ( zdir * e1t(ji,jj) + ( 1. - zdir ) * e1t(ji+1,jj) ) * e2u(ji,jj) * fse3u(ji,jj,jk) 
     
    187246         ! 
    188247         ! Mask at the T-points in the x-direction (mask=0 or mask=1) 
     248#if defined key_z_first 
     249         DO jj = 2, jpjm1 
     250            DO ji = 2, jpim1 
     251               DO jk = 1, jpkm1   
     252#else 
    189253         DO jk = 1, jpkm1   
    190254            DO jj = 2, jpjm1 
    191255               DO ji = fs_2, fs_jpim1   ! vector opt.                
     256#endif 
    192257                  zfu(ji,jj,jk) = tmask(ji-1,jj,jk) + tmask(ji,jj,jk) + tmask(ji+1,jj,jk) - 2. 
    193258               END DO 
     
    198263         ! 
    199264         ! Tracer flux on the x-direction 
     265#if defined key_z_first 
     266         DO jj = 2, jpjm1 
     267            DO ji = 2, jpim1 
     268               DO jk = 1, jpkm1   
     269#else 
    200270         DO jk = 1, jpkm1   
    201             ! 
    202271            DO jj = 2, jpjm1 
    203272               DO ji = fs_2, fs_jpim1   ! vector opt.                
     273#endif 
    204274                  zdir = 0.5 + SIGN( 0.5, pun(ji,jj,jk) )   ! if pun > 0 : zdir = 1 otherwise zdir = 0  
    205275                  !--- If the second ustream point is a land point 
     
    210280               END DO 
    211281            END DO 
     282#if defined key_z_first 
     283         END DO 
     284         ! Computation of the trend 
     285         DO jj = 2, jpjm1 
     286            DO ji = 2, jpim1 
     287               DO jk = 1, jpkm1 
     288#else 
    212289            ! 
    213290            ! Computation of the trend 
    214291            DO jj = 2, jpjm1 
    215292               DO ji = fs_2, fs_jpim1   ! vector opt.   
     293#endif 
    216294                  zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
    217295                  ! horizontal advective trends 
     
    230308      IF( wrk_not_released(3, 1,2,3) )   CALL ctl_stop('tra_adv_qck_i: failed to release workspace arrays') 
    231309      ! 
     310 
     311      !! * Reset control of array index permutation 
     312!FTRANS CLEAR 
     313#  include "oce_ftrans.h90" 
     314#  include "dom_oce_ftrans.h90" 
     315#  include "trc_oce_ftrans.h90" 
     316 
    232317   END SUBROUTINE tra_adv_qck_i 
    233318 
     
    241326      USE oce     , ONLY:   zwy => ua       ! ua used as workspace 
    242327      USE wrk_nemo, ONLY:   zfu => wrk_3d_1 , zfc => wrk_3d_2, zfd => wrk_3d_3   ! 3D workspace 
     328 
     329      !! DCSE_NEMO: need additional directives for renamed module variables 
     330!FTRANS zwy zfu zfc zfd :I :I :z 
     331 
    243332      ! 
    244333      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index 
     
    246335      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers 
    247336      REAL(wp), DIMENSION(        jpk     ), INTENT(in   ) ::   p2dt       ! vertical profile of tracer time-step 
    248       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pvn        ! j-velocity components 
    249       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn   ! before and now tracer fields 
    250       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend  
     337 
     338      !! DCSE_NEMO: This style defeats ftrans 
     339!     REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pvn        ! j-velocity components 
     340!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn   ! before and now tracer fields 
     341!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend  
     342 
     343!FTRANS pvn :I :I :z 
     344!FTRANS ptb ptn pta :I :I :z : 
     345      REAL(wp), INTENT(in   ) ::   pvn(jpi,jpj,jpk)         ! j-velocity component 
     346      REAL(wp), INTENT(in   ) ::   ptb(jpi,jpj,jpk,kjpt)    ! tracer field (before) 
     347      REAL(wp), INTENT(in   ) ::   ptn(jpi,jpj,jpk,kjpt)    ! tracer field (now) 
     348      REAL(wp), INTENT(inout) ::   pta(jpi,jpj,jpk,kjpt)    ! tracer trend  
     349 
    251350      !! 
    252351      INTEGER  :: ji, jj, jk, jn   ! dummy loop indices 
     
    264363         zfd(:,:,:) = 0.0     ;   zwy(:,:,:) = 0.0      
    265364         !                                                   
     365#if defined key_z_first 
     366         !--- Computation of the ustream and downstream value of the tracer and the mask 
     367         DO jj = 2, jpjm1 
     368            DO ji = 2, jpim1 
     369               DO jk = 1, jpkm1                                 
     370#else 
    266371         DO jk = 1, jpkm1                                 
    267372            !                                              
     
    269374            DO jj = 2, jpjm1 
    270375               DO ji = fs_2, fs_jpim1   ! vector opt. 
     376#endif 
    271377                  ! Upstream in the x-direction for the tracer 
    272378                  zfc(ji,jj,jk) = ptb(ji,jj-1,jk,jn) 
     
    283389         ! --------------------------- 
    284390         ! 
     391#if defined key_z_first 
     392         DO jj = 2, jpjm1 
     393            DO ji = 2, jpim1 
     394               DO jk = 1, jpkm1                              
     395#else 
    285396         DO jk = 1, jpkm1                              
    286397            DO jj = 2, jpjm1 
    287398               DO ji = fs_2, fs_jpim1   ! vector opt.          
     399#endif 
    288400                  zdir = 0.5 + SIGN( 0.5, pvn(ji,jj,jk) )   ! if pun > 0 : zdir = 1 otherwise zdir = 0  
    289401                  zfu(ji,jj,jk) = zdir * zfc(ji,jj,jk ) + ( 1. - zdir ) * zfd(ji,jj+1,jk)  ! FU in the x-direction for T  
     
    292404         END DO 
    293405         ! 
     406#if defined key_z_first 
     407         DO jj = 2, jpjm1 
     408            DO ji = 2, jpim1 
     409               DO jk = 1, jpkm1   
     410                  zdt =  p2dt(jk) 
     411#else 
    294412         DO jk = 1, jpkm1   
    295413            zdt =  p2dt(jk) 
    296414            DO jj = 2, jpjm1 
    297415               DO ji = fs_2, fs_jpim1   ! vector opt.    
     416#endif 
    298417                  zdir = 0.5 + SIGN( 0.5, pvn(ji,jj,jk) )   ! if pun > 0 : zdir = 1 otherwise zdir = 0  
    299418                  zdx  = ( zdir * e2t(ji,jj) + ( 1. - zdir ) * e2t(ji,jj+1) ) * e1v(ji,jj) * fse3v(ji,jj,jk) 
     
    313432         ! 
    314433         ! Mask at the T-points in the x-direction (mask=0 or mask=1) 
     434#if defined key_z_first 
     435         DO jj = 2, jpjm1 
     436            DO ji = 2, jpim1 
     437               DO jk = 1, jpkm1   
     438#else 
    315439         DO jk = 1, jpkm1   
    316440            DO jj = 2, jpjm1 
    317441               DO ji = fs_2, fs_jpim1   ! vector opt.                
     442#endif 
    318443                  zfu(ji,jj,jk) = tmask(ji,jj-1,jk) + tmask(ji,jj,jk) + tmask(ji,jj+1,jk) - 2. 
    319444               END DO 
     
    324449         ! 
    325450         ! Tracer flux on the x-direction 
     451#if defined key_z_first 
     452         DO jj = 2, jpjm1 
     453            DO ji = 2, jpim1 
     454               DO jk = 1, jpkm1   
     455#else 
    326456         DO jk = 1, jpkm1   
    327457            ! 
    328458            DO jj = 2, jpjm1 
    329459               DO ji = fs_2, fs_jpim1   ! vector opt.                
     460#endif 
    330461                  zdir = 0.5 + SIGN( 0.5, pvn(ji,jj,jk) )   ! if pun > 0 : zdir = 1 otherwise zdir = 0  
    331462                  !--- If the second ustream point is a land point 
     
    336467               END DO 
    337468            END DO 
     469#if defined key_z_first 
     470         END DO 
     471         ! Computation of the trend 
     472         DO jj = 2, jpjm1 
     473            DO ji = 2, jpim1 
     474               DO jk = 1, jpkm1 
     475#else 
    338476            ! 
    339477            ! Computation of the trend 
    340478            DO jj = 2, jpjm1 
    341479               DO ji = fs_2, fs_jpim1   ! vector opt.   
     480#endif 
    342481                  zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
    343482                  ! horizontal advective trends 
     
    361500      IF( wrk_not_released(3, 1,2,3) )   CALL ctl_stop('tra_adv_qck_j: failed to release workspace arrays') 
    362501      ! 
     502 
     503      !! * Reset control of array index permutation 
     504!FTRANS CLEAR 
     505#  include "oce_ftrans.h90" 
     506#  include "dom_oce_ftrans.h90" 
     507#  include "trc_oce_ftrans.h90" 
     508 
    363509   END SUBROUTINE tra_adv_qck_j 
    364510 
     
    370516      !!---------------------------------------------------------------------- 
    371517      USE oce, ONLY:   zwz => ua   ! ua used as workspace 
     518 
     519      !! DCSE_NEMO: need additional directives for renamed module variables 
     520!FTRANS zwz :I :I :z 
     521 
    372522      ! 
    373523      INTEGER                              , INTENT(in   ) ::   kt       ! ocean time-step index 
    374524      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype   ! =TRA or TRC (tracer indicator) 
    375525      INTEGER                              , INTENT(in   ) ::   kjpt     ! number of tracers 
    376       REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pwn      ! vertical velocity  
    377       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptn      ! before and now tracer fields 
    378       REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta      ! tracer trend  
     526 
     527      !! DCSE_NEMO: This style defeats ftrans 
     528!     REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pwn      ! vertical velocity  
     529!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptn      ! tracer fields (now) 
     530!     REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta      ! tracer trend  
     531 
     532!FTRANS pwn :I :I :z 
     533!FTRANS ptn pta :I :I :z : 
     534      REAL(wp), INTENT(in   ) ::   pwn(jpi,jpj,jpk)         ! vertical velocity 
     535      REAL(wp), INTENT(in   ) ::   ptn(jpi,jpj,jpk,kjpt)    ! tracer fields (now) 
     536      REAL(wp), INTENT(inout) ::   pta(jpi,jpj,jpk,kjpt)    ! tracer trend  
     537 
    379538      ! 
    380539      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices 
     
    393552         ENDIF 
    394553         ! 
     554#if defined key_z_first 
     555         DO jj = 2, jpjm1 
     556            DO ji = 2, jpim1 
     557               DO jk = 2, jpkm1            ! Interior point: second order centered tracer flux at w-point 
     558#else 
    395559         DO jk = 2, jpkm1                  ! Interior point: second order centered tracer flux at w-point 
    396560            DO jj = 2, jpjm1 
    397561               DO ji = fs_2, fs_jpim1   ! vector opt. 
     562#endif 
    398563                  zwz(ji,jj,jk) = 0.5 * pwn(ji,jj,jk) * ( ptn(ji,jj,jk-1,jn) + ptn(ji,jj,jk,jn) ) 
    399564               END DO 
     
    401566         END DO 
    402567         ! 
     568#if defined key_z_first 
     569         DO jj = 2, jpjm1 
     570            DO ji = 2, jpim1 
     571               DO jk = 1, jpkm1    !==  Tracer flux divergence added to the general trend  ==! 
     572#else 
    403573         DO jk = 1, jpkm1          !==  Tracer flux divergence added to the general trend  ==! 
    404574            DO jj = 2, jpjm1 
    405575               DO ji = fs_2, fs_jpim1   ! vector opt. 
     576#endif 
    406577                  zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 
    407578                  ! k- vertical advective trends  
     
    417588      END DO 
    418589      ! 
     590 
     591      !! * Reset control of array index permutation 
     592!FTRANS CLEAR 
     593#  include "oce_ftrans.h90" 
     594#  include "dom_oce_ftrans.h90" 
     595#  include "trc_oce_ftrans.h90" 
     596 
    419597   END SUBROUTINE tra_adv_cen2_k 
    420598 
     
    427605      !! ** Method :    
    428606      !!---------------------------------------------------------------------- 
    429       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfu   ! second upwind point 
    430       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfd   ! first douwning point 
    431       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfc   ! the central point (or the first upwind point) 
    432       REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   puc   ! input as Courant number ; output as flux 
     607 
     608      !! DCSE_NEMO: This style defeats ftrans 
     609 
     610!     REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfu   ! second upwind point 
     611!     REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfd   ! first douwning point 
     612!     REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pfc   ! the central point (or the first upwind point) 
     613!     REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::   puc   ! input as Courant number ; output as flux 
     614 
     615!FTRANS pfu pfd pfc puc :I :I :z 
     616      REAL(wp), INTENT(in   ) ::   pfu(jpi,jpj,jpk)   ! second upwind point 
     617      REAL(wp), INTENT(in   ) ::   pfd(jpi,jpj,jpk)   ! first douwning point 
     618      REAL(wp), INTENT(in   ) ::   pfc(jpi,jpj,jpk)   ! the central point (or the first upwind point) 
     619      REAL(wp), INTENT(inout) ::   puc(jpi,jpj,jpk)   ! input as Courant number ; output as flux 
     620 
    433621      !! 
    434622      INTEGER  ::  ji, jj, jk               ! dummy loop indices  
     
    437625      !---------------------------------------------------------------------- 
    438626 
     627#if defined key_z_first 
     628      DO jj = 1, jpj 
     629         DO ji = 1, jpi 
     630            DO jk = 1, jpkm1 
     631#else 
    439632      DO jk = 1, jpkm1 
    440633         DO jj = 1, jpj 
    441634            DO ji = 1, jpi 
     635#endif 
    442636               zc     = puc(ji,jj,jk)                         ! Courant number 
    443637               zcurv  = pfd(ji,jj,jk) + pfu(ji,jj,jk) - 2. * pfc(ji,jj,jk) 
Note: See TracChangeset for help on using the changeset viewer.