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

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

Cleaning branch: use only statement in module3D

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