Ignore:
Timestamp:
06/29/16 16:21:13 (8 years ago)
Author:
dumas
Message:

OpenMP parallelization in conserv-mass-adv-diff_sept2009_mod.f90, diffusiv-polyn-0.6.f90, dragging_neff_slope_mod.f90, eaubasale-0.5_mod.f90 and relaxation_water_diffusion.f90

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/SOURCES/diffusiv-polyn-0.6.f90

    r68 r76  
    3434 
    3535!     =============================================================== 
    36  
     36!$ USE OMP_LIB 
    3737USE module3D_phy 
    3838 
     
    5858 
    5959! initialisation de la diffusion 
     60!$OMP PARALLEL 
     61!$OMP WORKSHARE 
    6062diffmx(:,:)=0. 
    6163diffmy(:,:)=0. 
    62  
     64!$OMP END WORKSHARE 
     65!$OMP END PARALLEL 
    6366 
    6467! calcul de SDX, SDY et de la pente au carre en mx et en my : 
     
    7982! Calcul de la vitesse de glissement, qu'elle soit calculee par diagno ou par sliding 
    8083!------------------------------------------------------------------------------------ 
    81 !  
     84!$OMP PARALLEL 
     85!$OMP WORKSHARE 
    8286where (flgzmx(:,:)) 
    8387    ubx(:,:) = uxflgz(:,:) 
     
    9195    uby(:,:) = ddby(:,:)* (-sdy(:,:)) 
    9296endwhere 
    93  
    94  
     97!$OMP END WORKSHARE 
     98!$OMP END PARALLEL 
    9599if (itracebug.eq.1) call tracebug(' diffusiv  apres calcul glissement') 
    96100 
     
    99103! Deformation SIA : Calcul de ddy  et ddx pour chaque valeur de iglen 
    100104!------------------------------------------------------------------ 
    101  
    102105do iglen=n1poly,n2poly 
    103106   glenexp=max(0.,(glen(iglen)-1.)/2.) 
    104  
    105  
     107  !$OMP PARALLEL  
     108  !$OMP DO 
    106109   do j=1,ny 
    107110      do i=1,nx 
    108  
    109111         if (.not.flotmy(i,j)) then !  On calcule quand la deformation même pour les ice streams 
    110112            ddy(i,j,iglen)=((slope2my(i,j)**  &  ! slope2my calcule en debut de routine 
    111113                 glenexp)*(rog)**glen(iglen)) & 
    112114                 *hmy(i,j)**(glen(iglen)+1) 
    113  
    114  
    115115         endif 
    116116      end do 
    117117   end do 
     118   !$OMP END DO 
     119   !$OMP END PARALLEL 
    118120end do 
     121 
     122 
    119123 
    120124do iglen=n1poly,n2poly 
    121125   glenexp=max(0.,(glen(iglen)-1.)/2.) 
    122  
    123  
     126   !$OMP PARALLEL 
     127   !$OMP DO 
    124128   do j=1,ny   
    125129      do i=1,nx 
    126  
    127130         if (.not.flotmx(i,j)) then      ! bug y->x corrige le 5 avril 2008 ( 
    128                                      
    129131            ddx(i,j,iglen)=((slope2mx(i,j)**  & ! slope2mx calcule en debut de routine 
    130132                 glenexp)*(rog)**glen(iglen)) & 
     
    133135      end do 
    134136   end do 
     137   !$OMP END DO 
     138   !$OMP END PARALLEL 
    135139end do 
    136140 
     141 
    137142! vitesse due a la déformation----------------------------------------------- 
    138  
    139  
     143!$OMP PARALLEL 
     144!$OMP DO 
    140145ydef: do j=2,ny 
    141146       do i=1,nx  
     
    164169    end do 
    165170   end do ydef 
    166  
     171!$OMP END DO 
     172 
     173!$OMP DO 
    167174xdef: do j=1,ny 
    168175       do i=2,nx  
     
    192199    end do 
    193200   end do xdef 
     201!$OMP END DO 
    194202 
    195203if (itracebug.eq.1) call tracebug(' diffusiv  avant limit') 
     
    199207 
    200208! la limitation selon x 
     209!$OMP WORKSHARE 
    201210where(.not.flgzmx(:,:)) 
    202211   uxbar(:,:)=min(uxbar(:,:),ultot) 
     
    209218   uybar(:,:)=max(uybar(:,:),-ultot) 
    210219end where 
     220!$OMP END WORKSHARE 
    211221 
    212222!    cas particulier des sommets ou il ne reste plus de neige. 
     223!$OMP DO 
    213224do j=2,ny-1 
    214225   do i=2,nx-1 
     
    221232   end do 
    222233end do 
    223  
     234!$OMP END DO 
     235!$OMP END PARALLEL 
    224236if (itracebug.eq.1)  call tracebug(' fin diffusiv ') 
    225237 
Note: See TracChangeset for help on using the changeset viewer.