source: branches/GRISLIv3/SOURCES/isostasie_mod-0.3.f90 @ 376

Last change on this file since 376 was 376, checked in by aquiquet, 16 months ago

Cleaning branch: use only statements for isostasy

File size: 4.2 KB
Line 
1!> \file isostasie_mod-0.3.f90
2!!Package isostasie
3!<
4
5!> \namespace isostasie_mod
6!! Module pour le calcul de l'isostasie
7!! \author ...
8!! \date ...
9!! @note Used modules
10!! @note    - use iso_declar
11!! @note    - use module3D_phy
12!! @note    - use param_phy_mod
13!<
14
15
16
17module isostasie_mod
18 
19   use iso_declar,only: nbed 
20
21contains
22
23!> SUBROUTINE: init_iso
24!! Routine qui permet d'initialiser les variables
25!>
26  subroutine init_iso ! routine qui permet d'initialiser les variables
27
28    use iso_declar,only: nlith,dt_iso,tausoc,dl,rl,lbloc,we,charge
29    use module3D_phy, only: geoplace,icouple,marine,err,dx,nx,ny,h0,bsoc0,sealevel_2d,w0,w1,i,j
30    use param_phy_mod, only: ro,row,rog,rowg,romg
31
32    implicit none
33
34    !     nbed=0 pas d'isostasie
35    !     nbed=1 temps de reaction
36
37
38    if (GEOPLACE.eq.'eismint') then
39       NBED=0
40       NLITH=0
41    else if (GEOPLACE(1:6).eq.'marine') then
42       NBED=1
43       NLITH=1
44    else if ((geoplace.eq.'anteis1').or.(geoplace.eq.'ant20km')) then
45
46       if (icouple.eq.2) then
47          NBED=0
48          NLITH=0
49          if (marine) then
50             NBED=1
51             NLITH=1
52          endif
53       else
54          NBED=1
55          NLITH=1
56       endif
57    else
58       NBED=1
59       !       switch sur la lithosphere
60       NLITH=1
61    endif
62
63    !      temps de reaction en annees
64    if (NBED.eq.1) tausoc=3000
65
66    if (NLITH.eq.1) then
67       !        DL lithosphere flexural rigidity (N.m)
68       DL=9.87E24
69       !        radius of relative stiffness  (metre)
70       RL=131910.
71       !        LBLOC, 400 km de part et d'autre
72       LBLOC=int((400000.-0.1)/DX)+1
73       !         LBLOC=int((480000.-0.1)/DX)+1
74       !         LBLOC=int(400000./DX)+1
75
76
77       !******** allocation des tableaux en fonction de la valeur de LBLOC *****
78
79       if (.not.allocated(WE)) then
80          allocate(WE(-LBLOC:LBLOC,-LBLOC:LBLOC),stat=err)
81          if (err/=0) then
82             print *,"Erreur a l'allocation du tableau WE",err
83             stop 4
84          end if
85       end if
86       !    PRINT*,'shape(we)',SHAPE(we),'lbloc',lbloc
87
88       if (.not.allocated(CHARGE)) then
89          allocate(CHARGE(1-LBLOC:NX+LBLOC,1-LBLOC:NY+LBLOC),stat=err)
90          if (err/=0) then
91             print *,"Erreur a l'allocation du tableau CHARGE",err
92             stop 4
93          end if
94       end if
95       !********* FIN DE L'ALLOCATION DES TABLEAUX ***********
96       call tab_litho
97    end if
98
99    !******** initialisation de CHARGE ***********
100
101    !       pour calcul du socle initial on calcule la charge
102    !       avec l'etat initial suppose en equilibre S0, H0, Bsoc0
103
104    do i=1,nx
105       do j=1,ny
106          if (ro*h0(i,j).ge.-row*(BSOC0(i,j)-sealevel_2d(i,j))) then
107
108             !        glace ou terre
109             charge(i,j)=rog*h0(i,j)
110
111          else
112             !        ocean
113             charge(i,j)=-rowg*(Bsoc0(i,j)-sealevel_2d(i,j))
114          endif
115       end do
116    end do
117
118    do j=1,ny  ! parties de charge a l'exterieure de la grille
119       charge(1-lbloc:0,j)=charge(1,j)
120       charge(nx+1:nx+lbloc,j)=charge(nx,j)
121    end do
122    do i=1,nx
123       charge(i,1-lbloc:0)=charge(i,1)
124       charge(i,ny+1:ny+lbloc)=charge(i,ny)
125    end do
126
127    do j=1,ny  ! parties de charge a l'exterieure de la grille
128       charge(1-lbloc:0,j)=charge(1,j)
129       charge(nx+1:nx+lbloc,j)=charge(nx,j)
130    end do
131    do i=1,nx
132       charge(i,1-lbloc:0)=charge(i,1)
133       charge(i,ny+1:ny+lbloc)=charge(i,ny)
134    end do
135    charge(1-lbloc:0,1-lbloc:0)=charge(1,1) !valeurs aux quatres coins
136    charge(1-lbloc:0,ny+1:ny+lbloc)=charge(1,ny) ! exterieurs au domaine
137    charge(nx+1:nx+lbloc,1-lbloc:0)=charge(nx,1)
138    charge(nx+1:nx+lbloc,ny+1:ny+lbloc)=charge(nx,ny)
139
140    !     deflection de la lithosphere
141    if (nlith.eq.1) then
142
143       call litho
144
145       do j=1,ny
146          do i=1,nx
147             w0(i,j)=w1(i,j)
148          end do
149       end do
150    else
151       do j=1,ny
152          do i=1,nx
153             w0(i,j)=charge(i,j)/romg
154          end do
155       end do
156    end if
157
158    dt_iso = 50.
159
160  end subroutine init_iso
161
162
163
164!> SUBROUTINE: bedrock
165!! Routine qui calcule l'isostasie
166!>
167
168subroutine bedrock
169  call taubed  ! routine qui calcule l'isostasie
170end subroutine bedrock
171
172end module isostasie_mod
Note: See TracBrowser for help on using the repository browser.