Ignore:
Timestamp:
06/15/16 17:13:33 (8 years ago)
Author:
dumas
Message:

First parallelization instructions with OpenMP tested in debug : exactly the same results | isostasy called only every 50 years

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/SOURCES/moy_mxmy.f90

    r4 r71  
    1818subroutine moy_mxmy(n1,n2,X2D,X_mx,X_my) 
    1919! fait la moyenne ponderee d'un tableau X2D sur les mailles staggered 
    20  
     20 !$ USE OMP_LIB 
    2121implicit none 
    2222integer, intent(in) :: n1,n2   !< dimension des tableaux 
     
    2929 
    3030 
     31  !$OMP PARALLEL PRIVATE(i_moins1,j_moins1,i_plus1,j_plus1) 
     32  !$OMP DO 
     33        do j=1,n2 
     34                do i=1,n1 
     35                        i_moins1=max(1,i-1) 
     36                        j_moins1=max(1,j-1) 
     37                        i_plus1=min(n1,i+1) 
     38                        j_plus1=min(n2,j+1) 
    3139 
    32 do j=1,n2 
    33    do i=1,n1 
    34       i_moins1=max(1,i-1) 
    35       j_moins1=max(1,j-1) 
    36       i_plus1=min(n1,i+1) 
    37       j_plus1=min(n2,j+1) 
    38  
    39       X_mx(i,j)=0.25*(X2D(i,j)+X2D(i_moins1,j)) & 
     40                        X_mx(i,j)=0.25*(X2D(i,j)+X2D(i_moins1,j)) & 
    4041           + 0.125*((X2D(i_moins1,j_plus1)+X2D(i,j_plus1))   & 
    4142           +       (X2D(i_moins1,j_moins1)+X2D(i,j_moins1))) 
    4243 
    43       X_my(i,j)=0.25*(X2D(i,j)+X2D(i,j_moins1)) & 
     44                        X_my(i,j)=0.25*(X2D(i,j)+X2D(i,j_moins1)) & 
    4445           + 0.125*((X2D(i_plus1,j_moins1)+X2D(i_plus1,j))   & 
    4546           +       (X2D(i_moins1,j_moins1)+X2D(i_moins1,j))) 
    46  
    47  
    48    end do 
    49 end do 
     47                end do 
     48        end do 
     49  !$OMP END DO 
     50  !$OMP END PARALLEL 
    5051 
    5152return 
Note: See TracChangeset for help on using the changeset viewer.