Ignore:
Timestamp:
06/22/16 15:43:40 (8 years ago)
Author:
dumas
Message:

OpenMP parallelization in Temperature-routines, deformation_mod_2lois and resol_adv_diff_2D-sept2009

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/SOURCES/resol_adv_diff_2D-sept2009.f90

    r65 r73  
    101101subroutine resol_adv_diff_2D_vect(Dfx,Dfy,advx,advy,H_presc,i_Hpresc,bil,vieuxH,newH) 
    102102 
     103!$ USE OMP_LIB 
    103104implicit none 
    104105 
     
    148149if (itracebug.eq.1)  call tracebug(' Entree dans routine resolution_diffusion') 
    149150 
    150  
     151!$OMP PARALLEL 
     152!$OMP WORKSHARE 
    151153! calcul de bdx et hdx 
    152154hdx(:,:)=dx1*(vieuxH(:,:)-eoshift(vieuxH(:,:),shift=-1,boundary=0.,dim=1)) 
     
    164166Frelax(:,:)=0. 
    165167DeltaH(:,:)=0. 
    166  
     168!$OMP END WORKSHARE 
    167169 
    168170! schema spatial 
    169171 
    170172if (upwind.eq.1.) then                 !schema amont 
    171  
     173!$OMP WORKSHARE 
    172174   where (Advx(:,:).ge.0.) 
    173175      c_west(:,:)=1. 
     
    185187      c_north(:,:)=1. 
    186188   end where 
    187  
     189!$OMP END WORKSHARE 
    188190else if (upwind.lt.1.) then             ! schema centre 
     191!$OMP WORKSHARE 
    189192      c_west(:,:)=0.5 
    190193      c_east(:,:)=0.5 
    191194      c_south(:,:)=0.5 
    192195      c_north(:,:)=0.5 
     196!$OMP END WORKSHARE       
    193197end if 
    194198 
     199 
    195200! attribution des elements des diagonales 
     201!$OMP DO PRIVATE(frdx,frdy,fraxw,fraxe,frays,frayn) 
    196202do j=2,ny-1 
    197203  do i=2,nx-1 
     
    254260     end do 
    255261  end do 
    256  
    257  
     262!$OMP END DO 
     263 
     264!$OMP WORKSHARE 
    258265where (i_hpresc(:,:) .eq.1)          ! thickness prescribed 
    259266   frelax(:,:) = H_presc(:,:) 
     
    264271   erelax(:,:) = 0. 
    265272end where 
     273!$OMP END WORKSHARE 
     274!$OMP END PARALLEL 
    266275 
    267276stopp = .false. 
    268277ntour=0 
    269278 
    270  
    271  
    272279relax_loop: do while(.not.stopp) 
    273280ntour=ntour+1 
    274  
    275  
     281!$OMP PARALLEL 
     282!$OMP DO PRIVATE(reste) 
    276283do j=2,ny-1 
    277284   do i=2,nx-1 
     
    281288           + crelax(i,j)*newH(i,j))- frelax(i,j) 
    282289 
    283 if (ntour.eq.1)  debug_3D(i,j,49)=(((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & 
    284            + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & 
    285            + crelax(i,j)*newH(i,j)) 
     290!if (ntour.eq.1)  debug_3D(i,j,49)=(((arelax(i,j)*newH(i-1,j) +drelax(i,j)*newH(i,j-1)) & 
     291!           + (brelax(i,j)*newH(i+1,j) + erelax(i,j)*newH(i,j+1))) & 
     292!           + crelax(i,j)*newH(i,j)) 
    286293 
    287294 
     
    290297   end do 
    291298end do 
    292  
    293 debug_3D(:,:,50)=arelax(:,:) 
    294 debug_3D(:,:,51)=brelax(:,:) 
    295 debug_3D(:,:,52)=crelax(:,:) 
    296 debug_3D(:,:,53)=drelax(:,:) 
    297 debug_3D(:,:,54)=erelax(:,:) 
    298 debug_3D(:,:,55)=frelax(:,:) 
    299  
    300  
    301  
     299!$OMP END DO 
     300 
     301!debug_3D(:,:,50)=arelax(:,:) 
     302!debug_3D(:,:,51)=brelax(:,:) 
     303!debug_3D(:,:,52)=crelax(:,:) 
     304!debug_3D(:,:,53)=drelax(:,:) 
     305!debug_3D(:,:,54)=erelax(:,:) 
     306!debug_3D(:,:,55)=frelax(:,:) 
     307 
     308 
     309!$OMP WORKSHARE 
    302310newH(:,:)=newH(:,:)-deltaH(:,:) 
    303  
    304  
     311!$OMP END WORKSHARE 
     312!$OMP END PARALLEL 
    305313 
    306314 
     
    310318delh=0 
    311319 
    312  
     320!$OMP PARALLEL 
     321!$OMP DO REDUCTION(+:delh) 
    313322do j=2,ny-1 
    314323   do i=2,nx-1 
     
    316325   end do 
    317326end do 
     327!$OMP END DO 
     328!$OMP END PARALLEL 
    318329 
    319330if (delh.gt.0.) then 
     
    330341 
    331342! thickness at the upwind node 
    332 do j = 1, ny 
    333    do i = 2, nx 
    334       debug_3D(i,j,92) = c_west(i,j) * newH(i-1,j) + c_east(i,j) * newH(i,j) 
    335    end do 
    336 end do 
    337 do j = 2, ny 
    338    do i = 1, nx 
    339       debug_3D(i,j,93) = c_south(i,j) * newH(i,j-1) + c_north(i,j) * newH(i,j) 
    340    end do 
    341 end do 
    342  
    343  if (itracebug.eq.1)  call tracebug(' fin routine resolution_diffusion') 
     343!do j = 1, ny 
     344!   do i = 2, nx 
     345!      debug_3D(i,j,92) = c_west(i,j) * newH(i-1,j) + c_east(i,j) * newH(i,j) 
     346!   end do 
     347!end do 
     348!do j = 2, ny 
     349!   do i = 1, nx 
     350!      debug_3D(i,j,93) = c_south(i,j) * newH(i,j-1) + c_north(i,j) * newH(i,j) 
     351!   end do 
     352!end do 
     353 
     354if (itracebug.eq.1)  call tracebug(' fin routine resolution_diffusion') 
     355 
     356  
    344357return 
    345358 
Note: See TracChangeset for help on using the changeset viewer.