source: branches/GRISLIv3/SOURCES/litho-0.4.f90 @ 483

Last change on this file since 483 was 467, checked in by aquiquet, 5 months ago

Cleaning branch: continuing module3D cleaning

File size: 3.1 KB
Line 
1!> \file litho-0.4.f90     
2!! Donne la repartition des enfoncements
3!! en fonction de la rigidite de la lithosphere.
4!<
5
6!> SUBROUTINE: litho()
7!!Petit routine qui donne la repartition des enfoncements
8!!en fonction de la rigidite de la lithosphere.
9!! \author ....
10!! \date  10 septembre 2009     
11!! @note Cette version utilise la syntaxe f90 pour calculer W1
12!!     - avantage : plus facile a lire
13!!         pas de pb d'execution avac ifor 10
14!!     - inconvenient : pourrait perdre la symetrie.
15!! @note  - En entree
16!! @note      WE(-LBLOC:LBLOC,-LBLOC:LBLOC)  : deflection due a une charge unitaire
17!! defini dans tab-litho
18!! @note      LBLOC : relie a la distance : distance en noeud autour de laquelle
19!! la flexure de la lithosphere est calculee
20!! @note      CHARGE(1-LBLOC:NX+LBLOC,1-LBLOC:NY+LBLOC) : poids par unite de surface
21!! (unite ?) au temps time, calcule avant  'appel a litho
22!! dans taubed ou initial2
23!! @note  - En sortie
24!! @note      W1(nx,ny)   : enfoncement courant
25!! @note  - tableaux de travail.
26!! @note      WLOC(-LBLOC:LBLOC,-LBLOC:LBLOC) : enfoncement au point i,j du a la charge
27!! au point IP,JP
28!! @note      LPX et LPY sont les coordonnees entre (-LBLOC:LBLOC,-LBLOC:LBLOC) autour
29!! du point I,J pour calculer tous les WLOC
30!! @note      IP et JP sont les coordonnees de CHARGE dans un repere :
31!! (1-LBLOC:NX+LBLOC,1-LBLOC:NY+LBLOC)
32!! @note Used modules:
33!! @note    - use module3D_phy
34!! @note    - use ISO_DECLAR
35!! @note    - use param_phy_mod
36!<
37
38 
39subroutine litho
40  !$ USE OMP_LIB
41  USE module3D_phy, only: err
42  use geography, only: nx,ny
43  use isostasie_mod, only: w1
44  USE iso_declar, only: lbloc,we,charge ! module de declaration des variables specifiques a l'isostasie
45
46  implicit none
47
48  INTEGER :: SOM1,SOM2
49  REAL, dimension(:,:), allocatable :: WLOC
50
51  integer :: i,j
52
53  !----- allocation de WLOC  et de croix -----------
54
55  if (.not.allocated(WLOC)) THEN
56     allocate(WLOC(-LBLOC:LBLOC,-LBLOC:LBLOC),stat=err)
57     if (err/=0) then
58        print *,"Erreur a l'allocation du tableau WLOC",err
59        stop 4
60     end if
61  end if
62
63  !----- fin de l'allocation --------------
64
65  !     calcul de la deflexion par sommation des voisins appartenant
66  !     au bloc de taille 2*LBLOC
67  som1=0.
68  som2=0.
69
70  !      On somme aussi les contributions des points exterieurs au domaine
71  !      lorsque la charge est due a l'ocean. On suppose alors  que
72  !      ces points ont la meme charge que les limites
73  !$OMP PARALLEL
74  !$OMP DO PRIVATE(Wloc)
75  boucleij : do J=1,NY
76     do I=1,NX
77
78
79        W1(I,J)=0. 
80 !   ii=0
81
82 ! Wloc : enfoncement centre sur ij du a chaque charges autour
83        Wloc(:,:) = We(:,:) * charge(i-lbloc:i+lbloc,j-lbloc:j+lbloc)
84
85
86        W1(i,j) = sum (Wloc(:,:))   ! effet en ij de toutes les charges autour (distance < lbloc)
87
88 !    ------------------------ FIN DE L'INTEGRATION SUR LE PAVE LBLOC
89 !   som1=som1+W1(I,J)
90 !   som2=som2-charge(I,J)/ROMG
91
92     end do
93  end do  boucleij
94  !$OMP END DO
95  !$OMP END PARALLEL
96
97  !     write(6,*)'enfoncement total avec rigidite de la lithosphere:',
98  !    &  som1
99  !     write(6,*)'enfoncement total avec isostasie locale :',
100  !    &  som2
101
102
103end subroutine litho
104
Note: See TracBrowser for help on using the repository browser.