source: trunk/SOURCES/isostasie_mod-0.3.f90 @ 111

Last change on this file since 111 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

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