source: branches/GRISLIv3/SOURCES/taubed-0.3.f90 @ 446

Last change on this file since 446 was 446, checked in by aquiquet, 8 months ago

Cleaning branch: use only statement in module3D

File size: 4.1 KB
Line 
1!> \file taubed-0.3.f90
2!! Calcul la charge en chaque point de la grille puis appel la routine litho pour calculer la contribution
3!! de chaque point a la deflexion de la lithosphere   
4!<
5
6!> SUBROUTINE: taubed()
7!! \author ...
8!! \date  16 Novembre 1999
9!! @note  Routine qui calcul la charge en chaque point de la grille
10!! puis appel la routine litho pour calculer la contribution
11!! de chaque point a la deflexion de la lithosphere   
12!! @note En entree
13!!        - H,
14!!        - BSOC,
15!!        - SEALEVEL,
16!!        - RO,
17!!        - ROW
18!! @note En sortie
19!!        - CHARGE(1-LBLOC:NX+LBLOC,1-LBLOC:NY+LBLOC) : poids par unite de surface
20!!               (unite ?)   Elle est calculee initialement dans initial2
21!!               Poids de la colonne d'eau ou de la colonne de glace.
22!!               a l'exterieur du domaine : 1-LBLOC:1 et NX+1:NX+LBLOC
23!!               on donne les valeurs des bords de la grille
24!!               CHARGE est utilise par litho uniquement
25!!        - W1(NX,NY) est l'enfoncement courant, c'est le resultat
26!!               de la routine litho
27!!               W1 peut etre calcule de plusieurs facons selon le modele
28!!               d'isostasie utilise
29!! @note Used module:
30!! @note   - use module3D_phy
31!! @note   - use param_phy_mod
32!! @note   - use ISO_DECLAR
33!!
34!<
35
36!!     ****************************************************************
37!!     *           BEDROCK ADJUSTMENT avec temps de reaction          *
38!!     *  changement de nom B -> Bsoc                                 *
39!!     ****************************************************************
40
41
42
43subroutine taubed()
44
45  !$USE OMP_LIB
46  USE module3D_phy, only: h,bsoc0,bsoc,bdot,sealevel_2d,w0,w1
47  USE geography, only: nx,ny
48  USE param_phy_mod, only: ro,row,rog,rowg,rom
49  USE iso_declar,only: nlith,lbloc,charge,dt_iso,tausoc ! module de declaration des variables de l'isostasie
50
51  implicit none
52
53  integer :: i,j 
54
55  !    ********* calcul de W1 l'enfoncement d'equilibre au temps t
56  ! NLITH est defini dans isostasie et permet le choix du modele d'isostasie
57
58 
59  if (NLITH.eq.1) then
60     !       avec rigidite de la lithosphere
61     !$OMP PARALLEL
62     !$OMP DO
63     do J=1,NY 
64        do I=1,NX
65           if (RO*H(I,J).ge.-ROW*(BSOC(I,J)-SEALEVEL_2D(i,j))) then
66              !           glace ou terre
67              CHARGE(I,J)=ROG*H(I,J)
68           else
69              !           ocean
70              CHARGE(I,J)=-ROWG*(BSOC(I,J)-SEALEVEL_2D(i,j))
71           endif
72        end do
73     end do
74     !$OMP END DO
75
76     ! il faut remplir CHARGE dans les parties a l'exterieur de la grille :
77     ! a l'exterieur de la grille CHARGE est egale a la valeur sur le bord
78     !$OMP DO
79     do J=1,NY
80        CHARGE(1-LBLOC:0,J)=CHARGE(1,J)      ! bord W
81        CHARGE(NX+1:NX+LBLOC,J)=CHARGE(NX,J) ! bord E
82     end do
83     !$OMP END DO
84     !$OMP DO
85     do I=1,NX
86        CHARGE(I,1-LBLOC:0)=CHARGE(I,1)      ! bord S
87        CHARGE(I,NY+1:NY+LBLOC)=CHARGE(I,NY) ! bord N
88     end do
89     !$OMP END DO
90
91     ! valeur dans les quatres coins exterieurs au domaine       
92          !$OMP WORKSHARE
93     CHARGE(1-LBLOC:0,1-LBLOC:0)=CHARGE(1,1)           ! coin SW
94     CHARGE(1-LBLOC:0,NY+1:NY+LBLOC)=CHARGE(1,NY)      ! coin NW
95     CHARGE(NX+1:NX+LBLOC,1-LBLOC:0)=CHARGE(NX,1)      ! coin SE
96     CHARGE(NX+1:NX+LBLOC,NY+1:NY+LBLOC)=CHARGE(NX,NY) ! coin NE
97     !$OMP END WORKSHARE
98     !$OMP END PARALLEL
99     call litho
100
101  else
102     !     enfoncement local
103     !$OMP PARALLEL
104     !$OMP DO
105     do J=1,NY
106        do I=1,NX
107           if (RO*H(I,J).ge.-ROW*(BSOC(I,J)-SEALEVEL_2D(i,j))) then
108              !           glace ou terre
109              W1(I,J)=RO/ROM*H(I,J)
110           else
111              !           ocean
112              W1(I,J)=-ROW/ROM*(BSOC(I,J)-SEALEVEL_2D(i,j))
113           endif
114        end do
115     end do
116     !$OMP END DO
117     !$OMP END PARALLEL
118  endif
119
120  !     decroissance exponentielle de l'enfoncement
121  !$OMP PARALLEL
122  !$OMP DO
123  do J=1,NY
124     do I=1,NX
125        BDOT(I,J)=((Bsoc0(I,J)-BSOC(I,J))-(W1(I,J)-W0(I,J)))/TAUSOC
126        BSOC(I,J)=BSOC(I,J)+BDOT(I,J)*DT_ISO
127     end do
128  end do
129  !$OMP END DO
130  !$OMP END PARALLEL
131
132end subroutine taubed
Note: See TracBrowser for help on using the repository browser.