Ignore:
Timestamp:
06/24/16 09:19:59 (8 years ago)
Author:
dumas
Message:

OpenMP parallelization in remplimat-shelves-tabTu, litho and taubed

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/SOURCES/litho-0.4.f90

    r65 r74  
    3737 
    3838  
    39       subroutine litho 
     39subroutine litho 
     40  !$ USE OMP_LIB 
     41  USE module3D_phy 
     42  USE param_phy_mod 
     43  USE ISO_DECLAR ! module de declaration des variables specifiques a l'isostasie 
    4044 
    41       USE module3D_phy 
    42       USE param_phy_mod 
    43       USE ISO_DECLAR ! module de declaration des variables specifiques a l'isostasie 
     45  implicit none 
    4446 
    45       implicit none 
    46  
    47       INTEGER :: II,SOM1,SOM2 
    48       REAL, dimension(:,:), allocatable :: WLOC 
     47  INTEGER :: II,SOM1,SOM2 
     48  REAL, dimension(:,:), allocatable :: WLOC 
    4949 
    5050 
    51 !----- allocation de WLOC  et de croix ----------- 
     51  !----- allocation de WLOC  et de croix ----------- 
    5252 
    53       if (.not.allocated(WLOC)) THEN 
    54          allocate(WLOC(-LBLOC:LBLOC,-LBLOC:LBLOC),stat=err) 
    55          if (err/=0) then 
    56             print *,"Erreur a l'allocation du tableau WLOC",err 
    57             stop 4 
    58          end if 
    59       end if 
     53  if (.not.allocated(WLOC)) THEN 
     54     allocate(WLOC(-LBLOC:LBLOC,-LBLOC:LBLOC),stat=err) 
     55     if (err/=0) then 
     56        print *,"Erreur a l'allocation du tableau WLOC",err 
     57        stop 4 
     58     end if 
     59  end if 
    6060 
    61 !----- fin de l'allocation -------------- 
     61  !----- fin de l'allocation -------------- 
    6262 
    63 !     calcul de la deflexion par sommation des voisins appartenant 
    64 !     au bloc de taille 2*LBLOC 
    65        som1=0. 
    66        som2=0. 
     63  !     calcul de la deflexion par sommation des voisins appartenant 
     64  !     au bloc de taille 2*LBLOC 
     65  som1=0. 
     66  som2=0. 
    6767 
    68 !      On somme aussi les contributions des points exterieurs au domaine 
    69 !      lorsque la charge est due a l'ocean. On suppose alors  que 
    70 !      ces points ont la meme charge que les limites 
    71  
    72 boucleij : do J=1,NY 
    73            do I=1,NX 
     68  !      On somme aussi les contributions des points exterieurs au domaine 
     69  !      lorsque la charge est due a l'ocean. On suppose alors  que 
     70  !      ces points ont la meme charge que les limites 
     71  !$OMP PARALLEL 
     72  !$OMP DO PRIVATE(Wloc) 
     73  boucleij : do J=1,NY 
     74     do I=1,NX 
    7475 
    7576 
    76               W1(I,J)=0.  
    77               ii=0 
     77        W1(I,J)=0.  
     78 !   ii=0 
    7879 
    79 ! Wloc : enfoncement centre sur ij du a chaque charges autour 
    80       Wloc(:,:) = We(:,:) * charge(i-lbloc:i+lbloc,j-lbloc:j+lbloc) 
     80 ! Wloc : enfoncement centre sur ij du a chaque charges autour 
     81        Wloc(:,:) = We(:,:) * charge(i-lbloc:i+lbloc,j-lbloc:j+lbloc) 
    8182 
    8283 
    83       W1(i,j) = sum (Wloc(:,:))   ! effet en ij de toutes les charges autour (distance < lbloc) 
     84        W1(i,j) = sum (Wloc(:,:))   ! effet en ij de toutes les charges autour (distance < lbloc) 
    8485 
    85 !    ------------------------ FIN DE L'INTEGRATION SUR LE PAVE LBLOC 
    86               som1=som1+W1(I,J) 
    87               som2=som2-charge(I,J)/ROMG 
     86 !    ------------------------ FIN DE L'INTEGRATION SUR LE PAVE LBLOC 
     87 !   som1=som1+W1(I,J) 
     88 !   som2=som2-charge(I,J)/ROMG 
    8889 
    89       end do 
    90       end do  boucleij 
    91         
    92 !     write(6,*)'enfoncement total avec rigidite de la lithosphere:', 
    93 !    &  som1 
    94 !     write(6,*)'enfoncement total avec isostasie locale :', 
    95 !    &  som2 
     90     end do 
     91  end do  boucleij 
     92  !$OMP END DO 
     93  !$OMP END PARALLEL 
     94 
     95  !     write(6,*)'enfoncement total avec rigidite de la lithosphere:', 
     96  !    &  som1 
     97  !     write(6,*)'enfoncement total avec isostasie locale :', 
     98  !    &  som2 
    9699 
    97100 
    98     end subroutine litho 
    99        
     101end subroutine litho 
     102 
Note: See TracChangeset for help on using the changeset viewer.