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/LDF/ldfslp.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/LDF/ldfslp.F90

    r2715 r3211  
    4646   !                                                                !! Griffies operator 
    4747   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   wslp2                !: wslp**2 from Griffies quarter cells 
    48    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi_g, triadj_g   !: skew flux  slopes relative to geopotentials  
    49    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi  , triadj     !: isoneutral slopes relative to model-coordinate 
     48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   & 
     49             &   triadi_g, triadj_g   !: skew flux  slopes relative to geopotentials  
     50   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   & 
     51             &   triadi  , triadj     !: isoneutral slopes relative to model-coordinate 
    5052 
    5153   !                                                              !! Madec operator 
     
    6264   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   zdzrho , zdyrho, zdxrho     ! Horizontal and vertical density gradients 
    6365   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   zti_mlb, ztj_mlb            ! for Griffies operator only 
     66 
     67   !! * Control permutation of array indices 
     68#  include "ldfslp_ftrans.h90" 
     69!FTRANS zdxrho :I :I :z : 
     70!FTRANS zdyrho :I :I :z : 
     71!FTRANS zdzrho :I :I :z : 
     72#  include "oce_ftrans.h90" 
     73#  include "dom_oce_ftrans.h90" 
     74#  include "ldftra_oce_ftrans.h90" 
     75#  include "ldfdyn_oce_ftrans.h90" 
    6476 
    6577   !! * Substitutions 
     
    119131      USE oce     , ONLY:   zgrv => ta       , zwz => sa   ! (ta,sa) used as workspace 
    120132      USE wrk_nemo, ONLY:   zdzr => wrk_3d_1               ! 3D workspace 
     133      !! DCSE_NEMO: need additional directives for renamed module variables 
     134!FTRANS zgru :I :I :z 
     135!FTRANS zww  :I :I :z 
     136!FTRANS zgrv :I :I :z 
     137!FTRANS zwz  :I :I :z 
     138!FTRANS zdzr :I :I :z 
     139 
    121140      !! 
    122141      INTEGER , INTENT(in)                   ::   kt    ! ocean time-step index 
     142!FTRANS prd :I :I :z 
     143!FTRANS pn2 :I :I :z 
    123144      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   prd   ! in situ density 
    124145      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   pn2   ! Brunt-Vaisala frequency (locally ref.) 
     
    145166      zwz(:,:,:) = 0._wp 
    146167      ! 
     168#if defined key_z_first 
     169      DO jj = 1, jpjm1           !==   i- & j-gradient of density   ==! 
     170         DO ji = 1, jpim1 
     171            DO jk = 1, jpk 
     172#else 
    147173      DO jk = 1, jpk             !==   i- & j-gradient of density   ==! 
    148174         DO jj = 1, jpjm1 
    149175            DO ji = 1, fs_jpim1   ! vector opt. 
     176#endif 
    150177               zgru(ji,jj,jk) = umask(ji,jj,jk) * ( prd(ji+1,jj  ,jk) - prd(ji,jj,jk) )  
    151178               zgrv(ji,jj,jk) = vmask(ji,jj,jk) * ( prd(ji  ,jj+1,jk) - prd(ji,jj,jk) )  
     
    154181      END DO 
    155182      IF( ln_zps ) THEN                           ! partial steps correction at the bottom ocean level 
    156 # if defined key_vectopt_loop   
     183!! DCSE_NEMO: Attention! key_vectopt_loop will break key_z_first 
     184# if ( defined key_vectopt_loop ) && ! ( defined key_z_first )  
    157185         DO jj = 1, 1 
    158186            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling) 
     
    167195      ENDIF 
    168196      ! 
     197#if defined key__first 
     198      DO jj = 1, jpj 
     199         DO ji = 1, jpi 
     200            zdzr(ji,jj,1) = 0._wp   !==   Local vertical density gradient at T-point   == !   (evaluated from N^2) 
     201            DO jk = 2, jpkm1 
     202               zdzr(ji,jj,jk) = zm1_g * ( prd(ji,jj,jk) + 1._wp )              & 
     203                  &                   * ( pn2(ji,jj,jk) + pn2(ji,jj,jk+1) ) * ( 1._wp - 0.5_wp * tmask(ji,jj,jk+1) ) 
     204            END DO 
     205         END DO 
     206      END DO 
     207#else 
    169208      zdzr(:,:,1) = 0._wp        !==   Local vertical density gradient at T-point   == !   (evaluated from N^2) 
    170209      DO jk = 2, jpkm1 
     
    177216            &                 * ( pn2(:,:,jk) + pn2(:,:,jk+1) ) * ( 1._wp - 0.5_wp * tmask(:,:,jk+1) ) 
    178217      END DO 
     218#endif 
    179219      ! 
    180220      !                          !==   Slopes just below the mixed layer   ==! 
     
    185225      ! ===========================      | vslp = d/dj( prd ) / d/dz( prd ) 
    186226      !                
     227#if defined key_z_first 
     228      DO jj = 2, jpjm1               !* Slopes at u and v points 
     229         DO ji = 2, jpim1 
     230            DO jk = 2, jpkm1 
     231#else 
    187232      DO jk = 2, jpkm1                            !* Slopes at u and v points 
    188233         DO jj = 2, jpjm1 
    189234            DO ji = fs_2, fs_jpim1   ! vector opt. 
     235#endif 
    190236               !                                      ! horizontal and vertical density gradient at u- and v-points 
    191237               zau = zgru(ji,jj,jk) / e1u(ji,jj) 
     
    223269      ! 
    224270      !                                            !* horizontal Shapiro filter 
     271#if defined key_z_first 
     272      DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only 
     273         DO ji = 2, jpim1   
     274            DO jk = 2, jpkm1 
     275#else 
    225276      DO jk = 2, jpkm1 
    226277         DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only 
    227278            DO ji = 2, jpim1   
     279#endif 
    228280               uslp(ji,jj,jk) = z1_16 * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      & 
    229281                  &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      & 
     
    238290            END DO 
    239291         END DO 
     292#if defined key_z_first 
     293      END DO 
     294      DO jj = 3, jpj-2                               ! other rows 
     295         DO ji = 2, jpim1 
     296            DO jk = 2, jpkm1 
     297#else 
    240298         DO jj = 3, jpj-2                               ! other rows 
    241299            DO ji = fs_2, fs_jpim1   ! vector opt. 
     300#endif 
    242301               uslp(ji,jj,jk) = z1_16 * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      & 
    243302                  &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      & 
     
    252311            END DO 
    253312         END DO 
     313#if defined key_z_first 
     314      END DO 
     315      !                                           !* decrease along coastal boundaries 
     316      DO jj = 2, jpjm1 
     317         DO ji = 2, jpim1 
     318            DO jk = 2, jpkm1 
     319#else 
    254320         !                                        !* decrease along coastal boundaries 
    255321         DO jj = 2, jpjm1 
    256322            DO ji = fs_2, fs_jpim1   ! vector opt. 
     323#endif 
    257324               uslp(ji,jj,jk) = uslp(ji,jj,jk) * ( umask(ji,jj+1,jk) + umask(ji,jj-1,jk  ) ) * 0.5_wp   & 
    258325                  &                            * ( umask(ji,jj  ,jk) + umask(ji,jj  ,jk+1) ) * 0.5_wp 
     
    267334      ! ===========================      | wslpj = mij( d/dj( prd ) / d/dz( prd ) 
    268335      !                
     336#if defined key_z_first 
     337      DO jj = 2, jpjm1 
     338         DO ji = 2, jpim1 
     339            DO jk = 2, jpkm1 
     340#else 
    269341      DO jk = 2, jpkm1 
    270342         DO jj = 2, jpjm1 
    271343            DO ji = fs_2, fs_jpim1   ! vector opt. 
     344#endif 
    272345               !                                  !* Local vertical density gradient evaluated from N^2 
    273346               zbw = zm1_2g * pn2 (ji,jj,jk) * ( prd (ji,jj,jk) + prd (ji,jj,jk-1) + 2. ) 
     
    305378      ! 
    306379      !                                           !* horizontal Shapiro filter 
     380#if defined key_z_first 
     381      DO jj = 2, jpjm1, MAX(1, jpj-3)                           ! rows jj=2 and =jpjm1 only 
     382         DO ji = 2, jpim1 
     383            DO jk = 2, jpkm1 
     384#else 
    307385      DO jk = 2, jpkm1 
    308386         DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only 
    309387            DO ji = 2, jpim1 
     388#endif 
    310389               wslpi(ji,jj,jk) = (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     & 
    311390                  &                +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     & 
     
    321400            END DO 
    322401         END DO   
     402#if defined key_z_first 
     403      END DO 
     404      DO jj = 3, jpj-2                                  ! other rows 
     405         DO ji = 2, jpim1 
     406            DO jk = 2, jpkm1 
     407#else 
    323408         DO jj = 3, jpj-2                               ! other rows 
    324409            DO ji = fs_2, fs_jpim1   ! vector opt. 
     410#endif 
    325411               wslpi(ji,jj,jk) = (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     & 
    326412                  &                +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     & 
     
    336422            END DO 
    337423         END DO 
     424#if defined key_z_first 
     425      END DO 
     426      !                                           !* decrease along coastal boundaries 
     427      DO jj = 2, jpjm1 
     428         DO ji = 2, jpim1 
     429            DO jk = 2, jpkm1 
     430#else 
    338431         !                                        !* decrease along coastal boundaries 
    339432         DO jj = 2, jpjm1 
    340433            DO ji = fs_2, fs_jpim1   ! vector opt. 
     434#endif 
    341435               zck =   ( umask(ji,jj,jk) + umask(ji-1,jj,jk) )   & 
    342436                  &  * ( vmask(ji,jj,jk) + vmask(ji,jj-1,jk) ) * 0.25 
     
    383477   END SUBROUTINE ldf_slp 
    384478    
     479   !! * Reset control of array index permutation 
     480!FTRANS CLEAR 
     481#  include "ldfslp_ftrans.h90" 
     482!FTRANS zdxrho :I :I :z : 
     483!FTRANS zdyrho :I :I :z : 
     484!FTRANS zdzrho :I :I :z : 
     485#  include "oce_ftrans.h90" 
     486#  include "dom_oce_ftrans.h90" 
     487#  include "ldftra_oce_ftrans.h90" 
     488#  include "ldfdyn_oce_ftrans.h90" 
    385489 
    386490   SUBROUTINE ldf_slp_grif ( kt ) 
     
    404508      USE wrk_nemo, ONLY:   zalpha  => wrk_3d_4 , zbeta => wrk_3d_5    ! alpha, beta at T points, at depth fsgdept 
    405509      USE wrk_nemo, ONLY:   z1_mlbw => wrk_2d_1 
     510      !! DCSE_NEMO: need additional directives for renamed module variables 
     511!FTRANS zdit   :I :I :z 
     512!FTRANS zdis   :I :I :z 
     513!FTRANS zdjt   :I :I :z 
     514!FTRANS zdjs   :I :I :z 
     515!FTRANS zdkt   :I :I :z 
     516!FTRANS zdks   :I :I :z 
     517!FTRANS zalpha :I :I :z 
     518!FTRANS zbeta  :I :I :z 
    406519      ! 
    407520      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
     
    426539      CALL eos_alpbet( tsb, zalpha, zbeta )     !==  before thermal and haline expension coeff. at T-points  ==! 
    427540      ! 
     541#if defined key_z_first 
     542      DO jj = 1, jpjm1 
     543         DO ji = 1, jpim1 
     544            DO jk = 1, jpkm1                    !==  before lateral T & S gradients at T-level jk  ==! 
     545#else 
    428546      DO jk = 1, jpkm1                          !==  before lateral T & S gradients at T-level jk  ==! 
    429547         DO jj = 1, jpjm1 
    430548            DO ji = 1, fs_jpim1   ! vector opt. 
     549#endif 
    431550               zdit(ji,jj,jk) = ( tb(ji+1,jj,jk) - tb(ji,jj,jk) ) * umask(ji,jj,jk)   ! i-gradient of T and S at jj 
    432551               zdis(ji,jj,jk) = ( sb(ji+1,jj,jk) - sb(ji,jj,jk) ) * umask(ji,jj,jk) 
     
    437556      END DO 
    438557      IF( ln_zps ) THEN                               ! partial steps: correction at the last level 
    439 # if defined key_vectopt_loop 
     558# if ( defined key_vectopt_loop ) && ! ( defined key_z_first ) 
    440559         DO jj = 1, 1 
    441560            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling) 
     
    452571      ENDIF 
    453572      ! 
     573#if defined key_z_first 
     574      DO jj = 1, jpj 
     575         DO ji = 1, jpi 
     576            zdkt(ji,jj,1) = 0._wp            !==  before vertical T & S gradient at w-level  ==! 
     577            zdks(ji,jj,1) = 0._wp 
     578            DO jk = 2, jpk 
     579               zdkt(ji,jj,jk) = ( tb(ji,jj,jk-1) - tb(ji,jj,jk) ) * tmask(ji,jj,jk) 
     580               zdks(ji,jj,jk) = ( sb(ji,jj,jk-1) - sb(ji,jj,jk) ) * tmask(ji,jj,jk) 
     581            END DO 
     582         END DO 
     583      END DO 
     584#else 
    454585      zdkt(:,:,1) = 0._wp                    !==  before vertical T & S gradient at w-level  ==! 
    455586      zdks(:,:,1) = 0._wp 
     
    458589         zdks(:,:,jk) = ( sb(:,:,jk-1) - sb(:,:,jk) ) * tmask(:,:,jk) 
    459590      END DO 
     591#endif 
    460592      ! 
    461593      ! 
    462594      DO jl = 0, 1                           !==  density i-, j-, and k-gradients  ==! 
    463595         ip = jl   ;   jp = jl         ! guaranteed nonzero gradients ( absolute value larger than repsln) 
     596#if defined key_z_first 
     597         DO jj = 1, jpjm1                          ! NB: not masked due to the minimum value set 
     598            DO ji = 1, jpim1 
     599               DO jk = 1, jpkm1                    ! done each pair of triad 
     600#else 
    464601         DO jk = 1, jpkm1                          ! done each pair of triad 
    465602            DO jj = 1, jpjm1                       ! NB: not masked due to the minimum value set 
    466603               DO ji = 1, fs_jpim1   ! vector opt.  
     604#endif 
    467605                  zdxrho_raw = ( zalpha(ji+ip,jj   ,jk) * zdit(ji,jj,jk) + zbeta(ji+ip,jj   ,jk) * zdis(ji,jj,jk) ) / e1u(ji,jj) 
    468606                  zdyrho_raw = ( zalpha(ji   ,jj+jp,jk) * zdjt(ji,jj,jk) + zbeta(ji   ,jj+jp,jk) * zdjs(ji,jj,jk) ) / e2v(ji,jj) 
     
    474612      END DO 
    475613     DO kp = 0, 1                           !==  density i-, j-, and k-gradients  ==! 
     614#if defined key_z_first 
     615         DO jj = 1, jpj                         ! NB: not masked due to the minimum value set 
     616            DO ji = 1, jpi 
     617               DO jk = 1, jpkm1                    ! done each pair of triad 
     618#else 
    476619         DO jk = 1, jpkm1                          ! done each pair of triad 
    477620            DO jj = 1, jpj                       ! NB: not masked due to the minimum value set 
    478621               DO ji = 1, jpi   ! vector opt.  
     622#endif 
    479623                  zdzrho_raw = ( zalpha(ji,jj,jk) * zdkt(ji,jj,jk+kp) + zbeta(ji,jj,jk) * zdks(ji,jj,jk+kp) )   & 
    480624                     &       / fse3w(ji,jj,jk+kp) 
     
    530674         DO jl = 0, 1 
    531675            ip = jl   ;   jp = jl         ! i- and j-indices of triads (i-k and j-k planes) 
     676#if defined key_z_first 
     677            DO jj = 1, jpjm1 
     678               DO ji = 1, jpim1 
     679                  DO jk = 1, jpkm1 
     680#else 
    532681            DO jk = 1, jpkm1 
    533682               DO jj = 1, jpjm1 
    534683                  DO ji = 1, fs_jpim1   ! vector opt. 
     684#endif 
    535685                     ! 
    536686                     ! Calculate slope relative to geopotentials used for GM skew fluxes 
     
    605755   END SUBROUTINE ldf_slp_grif 
    606756 
     757   !! * Reset control of array index permutation 
     758!FTRANS CLEAR 
     759#  include "ldfslp_ftrans.h90" 
     760!FTRANS zdxrho :I :I :z : 
     761!FTRANS zdyrho :I :I :z : 
     762!FTRANS zdzrho :I :I :z : 
     763#  include "oce_ftrans.h90" 
     764#  include "dom_oce_ftrans.h90" 
     765#  include "ldftra_oce_ftrans.h90" 
     766#  include "ldfdyn_oce_ftrans.h90" 
    607767 
    608768   SUBROUTINE ldf_slp_mxl( prd, pn2, p_gru, p_grv, p_dzr ) 
     
    622782      !!                omlmask         :  mixed layer mask 
    623783      !!---------------------------------------------------------------------- 
     784!FTRANS prd   :I :I :z 
     785!FTRANS pn2   :I :I :z 
     786!FTRANS p_gru :I :I :z 
     787!FTRANS p_grv :I :I :z 
     788!FTRANS p_dzr :I :I :z 
    624789      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   prd            ! in situ density 
    625790      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   pn2            ! Brunt-Vaisala frequency (locally ref.) 
     
    646811      !                          !==   surface mixed layer mask   ! 
    647812      DO jk = 1, jpk                      ! =1 inside the mixed layer, =0 otherwise 
    648 # if defined key_vectopt_loop 
     813# if ( defined key_vectopt_loop ) && ! ( defined key_z_first ) 
    649814         DO jj = 1, 1 
    650815            DO ji = 1, jpij   ! vector opt. (forced unrolling) 
     
    672837      !----------------------------------------------------------------------- 
    673838      ! 
    674 # if defined key_vectopt_loop 
     839# if ( defined key_vectopt_loop ) && ! ( defined key_z_first ) 
    675840      DO jj = 1, 1 
    676841         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling) 
     
    727892   END SUBROUTINE ldf_slp_mxl 
    728893 
     894   !! * Reset control of array index permutation 
     895!FTRANS CLEAR 
     896#  include "ldfslp_ftrans.h90" 
     897!FTRANS zdxrho :I :I :z : 
     898!FTRANS zdyrho :I :I :z : 
     899!FTRANS zdzrho :I :I :z : 
     900#  include "oce_ftrans.h90" 
     901#  include "dom_oce_ftrans.h90" 
     902#  include "ldftra_oce_ftrans.h90" 
     903#  include "ldfdyn_oce_ftrans.h90" 
    729904 
    730905   SUBROUTINE ldf_slp_init 
     
    780955            ! set the slope of diffusion to the slope of s-surfaces 
    781956            !      ( c a u t i o n : minus sign as fsdep has positive value ) 
     957#if defined key_z_first 
     958            DO jj = 2, jpjm1 
     959               DO ji = 2, jpim1 
     960                  DO jk = 1, jpk 
     961#else 
    782962            DO jk = 1, jpk 
    783963               DO jj = 2, jpjm1 
    784964                  DO ji = fs_2, fs_jpim1   ! vector opt. 
     965#endif 
    785966                     uslp (ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept(ji+1,jj,jk) - fsdept(ji ,jj ,jk) ) * umask(ji,jj,jk) 
    786967                     vslp (ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept(ji,jj+1,jk) - fsdept(ji ,jj ,jk) ) * vmask(ji,jj,jk) 
Note: See TracChangeset for help on using the changeset viewer.