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 6748 for branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF – NEMO

Ignore:
Timestamp:
2016-06-28T11:53:56+02:00 (8 years ago)
Author:
mocavero
Message:

GYRE hybrid parallelization

Location:
branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfbfr.F90

    r6140 r6748  
    112112         IF ( ln_loglayer.AND. .NOT.ln_linssh ) THEN ! "log layer" bottom friction coefficient 
    113113 
     114!$OMP PARALLEL DO private(jj,ji,ikbt,ztmp) 
    114115            DO jj = 1, jpj 
    115116               DO ji = 1, jpi 
     
    123124! (ISF) 
    124125            IF ( ln_isfcav ) THEN 
     126!$OMP PARALLEL DO private(jj,ji,ikbt,ztmp) 
    125127               DO jj = 1, jpj 
    126128                  DO ji = 1, jpi 
     
    135137            !    
    136138         ELSE 
     139!$OMP PARALLEL WORKSHARE 
    137140            zbfrt(:,:) = bfrcoef2d(:,:) 
    138141            ztfrt(:,:) = tfrcoef2d(:,:) 
    139          ENDIF 
    140  
     142!$OMP END PARALLEL WORKSHARE 
     143         ENDIF 
     144 
     145!$OMP PARALLEL DO private(jj,ji,ikbu,ikbv,zvu,zuv,zecu,zecv) 
    141146         DO jj = 2, jpjm1 
    142147            DO ji = 2, jpim1 
     
    173178 
    174179         IF( ln_isfcav ) THEN 
     180!$OMP PARALLEL DO private(jj,ji,ikbu,ikbv,zvu,zuv,zecu,zecv) 
    175181            DO jj = 2, jpjm1 
    176182               DO ji = 2, jpim1 
     
    266272      CASE( 0 ) 
    267273         IF(lwp) WRITE(numout,*) '      free-slip ' 
     274!$OMP PARALLEL WORKSHARE 
    268275         bfrua(:,:) = 0.e0 
    269276         bfrva(:,:) = 0.e0 
    270277         tfrua(:,:) = 0.e0 
    271278         tfrva(:,:) = 0.e0 
     279!$OMP END PARALLEL WORKSHARE 
    272280         ! 
    273281      CASE( 1 ) 
     
    296304         ENDIF 
    297305         ! 
     306!$OMP PARALLEL WORKSHARE 
    298307         bfrua(:,:) = - bfrcoef2d(:,:) 
    299308         bfrva(:,:) = - bfrcoef2d(:,:) 
     309!$OMP END PARALLEL WORKSHARE 
    300310         ! 
    301311         IF ( ln_isfcav ) THEN 
     
    310320            ENDIF 
    311321            ! 
     322!$OMP PARALLEL WORKSHARE 
    312323            tfrua(:,:) = - tfrcoef2d(:,:) 
    313324            tfrva(:,:) = - tfrcoef2d(:,:) 
     325!$OMP END PARALLEL WORKSHARE 
    314326         END IF 
    315327         ! 
     
    371383         ! 
    372384         IF( ln_loglayer.AND. ln_linssh ) THEN ! set "log layer" bottom friction once for all 
     385!$OMP PARALLEL DO private(jj,ji,ikbt,ztmp) 
    373386            DO jj = 1, jpj 
    374387               DO ji = 1, jpi 
     
    380393            END DO 
    381394            IF ( ln_isfcav ) THEN 
     395!$OMP PARALLEL DO private(jj,ji,ikbt,ztmp) 
    382396               DO jj = 1, jpj 
    383397                  DO ji = 1, jpi 
     
    419433      zmaxtfr = -1.e10_wp    ! initialise tracker for maximum of bottom friction coefficient 
    420434      ! 
     435!$OMP PARALLEL DO private(jj,ji,ikbu,ikbv,zfru,zfrv,ictu,ictv,zminbfr,zmaxbfr,zmintfr,zmaxtfr) 
    421436      DO jj = 2, jpjm1 
    422437         DO ji = 2, jpim1 
  • branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfevd.F90

    r6140 r6748  
    7676         zavm_evd(:,:,:) = avm(:,:,:)           ! set avm prior to evd application 
    7777         ! 
     78!$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 
    7879         DO jk = 1, jpkm1  
    7980            DO jj = 2, jpj             ! no vector opt. 
     
    9798         ! 
    9899      CASE DEFAULT         ! enhance vertical eddy diffusivity only (if rn2<-1.e-12)  
     100!$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 
    99101         DO jk = 1, jpkm1 
    100102!!!         WHERE( rn2(:,:,jk) <= -1.e-12 ) avt(:,:,jk) = tmask(:,:,jk) * avevd   ! agissant sur T SEUL!  
  • branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfmxl.F90

    r6352 r6748  
    9696 
    9797      ! w-level of the mixing and mixed layers 
     98!$OMP PARALLEL WORKSHARE 
    9899      nmln(:,:)  = nlb10               ! Initialization to the number of w ocean point 
    99100      hmlp(:,:)  = 0._wp               ! here hmlp used as a dummy variable, integrating vertically N^2 
     101!$OMP END PARALLEL WORKSHARE 
    100102      zN2_c = grav * rho_c * r1_rau0   ! convert density criteria into N^2 criteria 
    101103      DO jk = nlb10, jpkm1 
     
    110112      ! 
    111113      ! w-level of the turbocline and mixing layer (iom_use) 
     114!$OMP PARALLEL WORKSHARE 
    112115      imld(:,:) = mbkt(:,:) + 1        ! Initialization to the number of w ocean point 
     116!$OMP END PARALLEL WORKSHARE 
     117 
    113118      DO jk = jpkm1, nlb10, -1         ! from the bottom to nlb10  
    114119         DO jj = 1, jpj 
     
    119124      END DO 
    120125      ! depth of the mixing and mixed layers 
     126!$OMP PARALLEL DO schedule(static) private(jj, ji, iiki, iikn) 
    121127      DO jj = 1, jpj 
    122128         DO ji = 1, jpi 
  • branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/ZDF/zdftke.F90

    r6497 r6748  
    179179      ! 
    180180      IF( kt /= nit000 ) THEN   ! restore before value to compute tke 
     181!$OMP PARALLEL WORKSHARE 
    181182         avt (:,:,:) = avt_k (:,:,:)  
    182183         avm (:,:,:) = avm_k (:,:,:)  
    183184         avmu(:,:,:) = avmu_k(:,:,:)  
    184185         avmv(:,:,:) = avmv_k(:,:,:)  
     186!$OMP END PARALLEL WORKSHARE 
    185187      ENDIF  
    186188      ! 
     
    189191      CALL tke_avn      ! now avt, avm, avmu, avmv 
    190192      ! 
     193!$OMP PARALLEL WORKSHARE 
    191194      avt_k (:,:,:) = avt (:,:,:)  
    192195      avm_k (:,:,:) = avm (:,:,:)  
    193196      avmu_k(:,:,:) = avmu(:,:,:)  
    194197      avmv_k(:,:,:) = avmv(:,:,:)  
     198!$OMP END PARALLEL WORKSHARE 
    195199      ! 
    196200#if defined key_agrif 
     
    253257      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    254258      IF ( ln_isfcav ) THEN 
     259!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    255260         DO jj = 2, jpjm1            ! en(mikt(ji,jj))   = rn_emin 
    256261            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    259264         END DO 
    260265      END IF 
     266!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    261267      DO jj = 2, jpjm1            ! en(1)   = rn_ebb taum / rau0  (min value rn_emin0) 
    262268         DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    295301         zpelc(:,:,1) =  MAX( rn2b(:,:,1), 0._wp ) * gdepw_n(:,:,1) * e3w_n(:,:,1) 
    296302         DO jk = 2, jpk 
    297             zpelc(:,:,jk)  = zpelc(:,:,jk-1) + MAX( rn2b(:,:,jk), 0._wp ) * gdepw_n(:,:,jk) * e3w_n(:,:,jk) 
     303!$OMP PARALLEL DO schedule(static) private(jj, ji) 
     304            DO jj =1, jpj 
     305               DO ji=1, jpi 
     306                  zpelc(ji,jj,jk)  = zpelc(ji,jj,jk-1) + MAX( rn2b(ji,jj,jk), 0._wp ) * gdepw_n(ji,jj,jk) * e3w_n(ji,jj,jk) 
     307               END DO 
     308            END DO 
    298309         END DO 
    299310         !                        !* finite Langmuir Circulation depth 
     
    309320         END DO 
    310321         !                               ! finite LC depth 
     322!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    311323         DO jj = 1, jpj  
    312324            DO ji = 1, jpi 
     
    315327         END DO 
    316328         zcof = 0.016 / SQRT( zrhoa * zcdrag ) 
     329!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zus, zind, zwlc) 
    317330         DO jk = 2, jpkm1         !* TKE Langmuir circulation source term added to en 
    318331            DO jj = 2, jpjm1 
     
    338351      !                     ! zdiag : diagonal zd_up : upper diagonal zd_lw : lower diagonal 
    339352      ! 
     353!$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 
    340354      DO jk = 2, jpkm1           !* Shear production at uw- and vw-points (energy conserving form) 
    341355         DO jj = 1, jpjm1 
     
    356370         ! Note that zesh2 is also computed in the next loop. 
    357371         ! We decided to compute it twice to keep code readability and avoid an IF case in the DO loops 
     372!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zesh2, zri) 
    358373         DO jk = 2, jpkm1 
    359374            DO jj = 2, jpjm1 
     
    372387      ENDIF 
    373388      !          
     389!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zcof, zzd_up, zzd_lw, zesh2) 
    374390      DO jk = 2, jpkm1           !* Matrix and right hand side in en 
    375391         DO jj = 2, jpjm1 
     
    405421      !                          !* Matrix inversion from level 2 (tke prescribed at level 1) 
    406422      DO jk = 3, jpkm1                             ! First recurrence : Dk = Dk - Lk * Uk-1 / Dk-1 
     423!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    407424         DO jj = 2, jpjm1 
    408425            DO ji = fs_2, fs_jpim1    ! vector opt. 
     
    411428         END DO 
    412429      END DO 
     430!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    413431      DO jj = 2, jpjm1                             ! Second recurrence : Lk = RHSk - Lk / Dk-1 * Lk-1 
    414432         DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    417435      END DO 
    418436      DO jk = 3, jpkm1 
     437!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    419438         DO jj = 2, jpjm1 
    420439            DO ji = fs_2, fs_jpim1    ! vector opt. 
     
    423442         END DO 
    424443      END DO 
     444!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    425445      DO jj = 2, jpjm1                             ! thrid recurrence : Ek = ( Lk - Uk * Ek+1 ) / Dk 
    426446         DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    429449      END DO 
    430450      DO jk = jpk-2, 2, -1 
     451!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    431452         DO jj = 2, jpjm1 
    432453            DO ji = fs_2, fs_jpim1    ! vector opt. 
     
    435456         END DO 
    436457      END DO 
     458!$OMP PARALLEL DO schedule(static) private(jk,jj, ji) 
    437459      DO jk = 2, jpkm1                             ! set the minimum value of tke 
    438460         DO jj = 2, jpjm1 
     
    450472       
    451473      IF( nn_etau == 1 ) THEN           !* penetration below the mixed layer (rn_efr fraction) 
     474!$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 
    452475         DO jk = 2, jpkm1 
    453476            DO jj = 2, jpjm1 
     
    459482         END DO 
    460483      ELSEIF( nn_etau == 2 ) THEN       !* act only at the base of the mixed layer (jk=nmln)  (rn_efr fraction) 
     484!$OMP PARALLEL DO schedule(static) private(jj, ji, jk) 
    461485         DO jj = 2, jpjm1 
    462486            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    467491         END DO 
    468492      ELSEIF( nn_etau == 3 ) THEN       !* penetration belox the mixed layer (HF variability) 
     493!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, ztx2, zty2, ztau, zdif) 
    469494         DO jk = 2, jpkm1 
    470495            DO jj = 2, jpjm1 
     
    545570      ! 
    546571      ! initialisation of interior minimum value (avoid a 2d loop with mikt) 
     572!$OMP PARALLEL WORKSHARE 
    547573      zmxlm(:,:,:)  = rmxl_min     
    548574      zmxld(:,:,:)  = rmxl_min 
     575!$OMP END PARALLEL WORKSHARE 
    549576      ! 
    550577      IF( ln_mxl0 ) THEN            ! surface mixing length = F(stress) : l=vkarmn*2.e5*taum/(rau0*g) 
     578!$OMP PARALLEL DO schedule(static) private(jj, ji, zraug) 
    551579         DO jj = 2, jpjm1 
    552580            DO ji = fs_2, fs_jpim1 
     
    556584         END DO 
    557585      ELSE  
     586!$OMP PARALLEL WORKSHARE 
    558587         zmxlm(:,:,1) = rn_mxl0 
     588!$OMP END PARALLEL WORKSHARE 
    559589      ENDIF 
    560590      ! 
     591!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zrn2) 
    561592      DO jk = 2, jpkm1              ! interior value : l=sqrt(2*e/n^2) 
    562593         DO jj = 2, jpjm1 
     
    570601      !                     !* Physical limits for the mixing length 
    571602      ! 
     603!$OMP PARALLEL WORKSHARE 
    572604      zmxld(:,:, 1 ) = zmxlm(:,:,1)   ! surface set to the minimum value  
    573605      zmxld(:,:,jpk) = rmxl_min       ! last level  set to the minimum value 
     606!$OMP END PARALLEL WORKSHARE 
    574607      ! 
    575608      SELECT CASE ( nn_mxl ) 
     
    578611      ! where wmask = 0 set zmxlm == e3w_n 
    579612      CASE ( 0 )           ! bounded by the distance to surface and bottom 
     613!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zemxl) 
    580614         DO jk = 2, jpkm1 
    581615            DO jj = 2, jpjm1 
     
    591625         ! 
    592626      CASE ( 1 )           ! bounded by the vertical scale factor 
     627!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zemxl) 
    593628         DO jk = 2, jpkm1 
    594629            DO jj = 2, jpjm1 
     
    603638      CASE ( 2 )           ! |dk[xml]| bounded by e3t : 
    604639         DO jk = 2, jpkm1         ! from the surface to the bottom : 
     640!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    605641            DO jj = 2, jpjm1 
    606642               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    610646         END DO 
    611647         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : 
     648!$OMP PARALLEL DO schedule(static) private(jj, ji, zemxl) 
    612649            DO jj = 2, jpjm1 
    613650               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    621658      CASE ( 3 )           ! lup and ldown, |dk[xml]| bounded by e3t : 
    622659         DO jk = 2, jpkm1         ! from the surface to the bottom : lup 
     660!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    623661            DO jj = 2, jpjm1 
    624662               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    628666         END DO 
    629667         DO jk = jpkm1, 2, -1     ! from the bottom to the surface : ldown 
     668!$OMP PARALLEL DO schedule(static) private(jj, ji) 
    630669            DO jj = 2, jpjm1 
    631670               DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    634673            END DO 
    635674         END DO 
     675!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zemlm, zemlp) 
    636676         DO jk = 2, jpkm1 
    637677            DO jj = 2, jpjm1 
     
    648688      ! 
    649689# if defined key_c1d 
     690!$OMP PARALLEL WORKSHARE 
    650691      e_dis(:,:,:) = zmxld(:,:,:)      ! c1d configuration : save mixing and dissipation turbulent length scales 
    651692      e_mix(:,:,:) = zmxlm(:,:,:) 
     693!$OMP END PARALLEL WORKSHARE 
    652694# endif 
    653695 
     
    655697      !                     !  Vertical eddy viscosity and diffusivity  (avmu, avmv, avt) 
    656698      !                     !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     699!$OMP PARALLEL DO schedule(static) private(jk, jj, ji, zsqen, zav) 
    657700      DO jk = 1, jpkm1            !* vertical eddy viscosity & diffivity at w-points 
    658701         DO jj = 2, jpjm1 
     
    668711      CALL lbc_lnk( avm, 'W', 1. )      ! Lateral boundary conditions (sign unchanged) 
    669712      ! 
     713!$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 
    670714      DO jk = 2, jpkm1            !* vertical eddy viscosity at wu- and wv-points 
    671715         DO jj = 2, jpjm1 
     
    679723      ! 
    680724      IF( nn_pdl == 1 ) THEN      !* Prandtl number case: update avt 
     725!$OMP PARALLEL DO schedule(static) private(jk, jj, ji) 
    681726         DO jk = 2, jpkm1 
    682727            DO jj = 2, jpjm1 
     
    804849      ENDIF 
    805850      !                               !* set vertical eddy coef. to the background value 
     851!$OMP PARALLEL DO schedule(static) private(jk) 
    806852      DO jk = 1, jpk 
    807853         avt (:,:,jk) = avtb(jk) * wmask (:,:,jk) 
     
    857903           ELSE                                     ! No TKE array found: initialisation 
    858904              IF(lwp) WRITE(numout,*) ' ===>>>> : previous run without tke scheme, en computed by iterative loop' 
     905!$OMP PARALLEL WORKSHARE 
    859906              en (:,:,:) = rn_emin * tmask(:,:,:) 
     907!$OMP END PARALLEL WORKSHARE 
    860908              CALL tke_avn                               ! recompute avt, avm, avmu, avmv and dissl (approximation) 
    861909              ! 
     910!$OMP PARALLEL WORKSHARE 
    862911              avt_k (:,:,:) = avt (:,:,:) 
    863912              avm_k (:,:,:) = avm (:,:,:) 
    864913              avmu_k(:,:,:) = avmu(:,:,:) 
    865914              avmv_k(:,:,:) = avmv(:,:,:) 
     915!$OMP END PARALLEL WORKSHARE 
    866916              ! 
    867917              DO jit = nit000 + 1, nit000 + 10   ;   CALL zdf_tke( jit )   ;   END DO 
    868918           ENDIF 
    869919        ELSE                                   !* Start from rest 
     920!$OMP PARALLEL WORKSHARE 
    870921           en(:,:,:) = rn_emin * tmask(:,:,:) 
     922!$OMP END PARALLEL WORKSHARE 
     923!$OMP PARALLEL DO schedule(static) private(jk) 
    871924           DO jk = 1, jpk                           ! set the Kz to the background value 
    872925              avt (:,:,jk) = avtb(jk) * wmask (:,:,jk) 
Note: See TracChangeset for help on using the changeset viewer.