source: trunk/SOURCES/Hudson_files/sliding-hudson_mod.f90 @ 334

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

initial import GRISLI trunk

File size: 10.2 KB
Line 
1module sliding_hudson
2
3! Défini la loi de glissement dans les expériences Heino
4! modification par rapport a la version de Catherine : ice-shelves + streams autorises
5
6  USE module3D_phy
7  USE sedim_declar
8
9  implicit none
10
11  real :: Cr                            ! coefficient glissement pour hard rock
12  real :: Cs                            ! coefficient glissement pour sediment
13  real :: coefmax
14  real :: coefslid
15  real :: coefdrag
16!logical,dimension(nx,ny) :: gzmx_heino     ! pour transporter la valeur dans flottab
17!logical,dimension(nx,ny) :: gzmy_heino
18  real :: longslope
19  real :: seuil_gliss                   ! seuil sur hw_mx pour avoir du glissement
20  real, dimension(nx,ny) :: hw_mx
21  real, dimension(nx,ny) :: hw_my
22
23!character(len=80) :: filin
24
25
26contains
27!-------------------------------------------------------------------------------
28  subroutine init_sliding
29
30! formats pour les ecritures dans 42
31428 format(A)
32
33    namelist/sliding_Heino/seuil_gliss,coefmax,Cr,Cs   
34! lecture des parametres du run                      block sliding_heino
35!--------------------------------------------------------------------
36
37    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
38    read(num_param,sliding_Heino)
39
40    write(num_rep_42,428)'!___________________________________________________________' 
41    write(num_rep_42,428) '&sliding_Heino                  ! nom du bloc sliding Heino'
42    write(num_rep_42,*)
43    write(num_rep_42,*) 'seuil_gliss    = ', seuil_gliss 
44    write(num_rep_42,*) 'coefmax        = ', coefmax
45    write(num_rep_42,*) 'Cr             = ', Cr
46    write(num_rep_42,*) 'Cs             = ', Cs
47    write(num_rep_42,*)'/'                           
48    write(num_rep_42,428) '! seuil_gliss (m) : glissement si eau basale des 6 voisins > seuil_gliss'
49    write(num_rep_42,428) '! coefmax : hauteur d eau max dans le sediment (=hwatermax)'
50    write(num_rep_42,428) '! Cr : coeficient glissemnt sur hard rock (loi en puissance 3)'
51    write(num_rep_42,428) '! Cs : coefficient de la loi de frottement sediment (loi lineaire)'
52    write(num_rep_42,*)
53
54!cdc if (itracebug.eq.1)  call tracebug(' Heino: entree dans routine init_sliding')
55
56! Cette routine fait l'initialisation du glissement et est appelée juste apres
57! la lecture du masque
58
59! Cette routine fait l'initialisation du glissement
60
61!loigliss=4       ! 4 -> que glissement Heino (pas de mode dragging)
62!if (loigliss.eq.4) then
63!   flgzmx(:,:)=.false.
64!   flgzmy(:,:)=.false.
65!endif
66!    coefmax=hwatermax
67!    seuil_gliss=0.0001
68                      ! 5 -> glissement Heino + dragging
69
70!    ecriture dans le fichier parametres
71    write(42,*)
72    write(42,*)'glissement Heino -  initialisation masque'
73!    write(42,*)'loi glissement=',loigliss
74!     write(num_rep_42,*)'coef=1, gliss si les 2 voisins temperes'
75    write(42,*)'glissement si eau basale (6 voisins) >=seuil_gliss=',seuil_gliss
76    write(42,*)'coef dependant de la hauteur d eau    coefmax=',coefmax
77    write(42,*)' moyenne de la hauteur d eau a partir de 6 voisins'
78    write(42,*)'----------------------------------------------------------'
79
80
81! fichier mask fourni par l'intercomparaison Heino
82! on refait la lecture
83
84!      filin='mask-shelf-40.dat'
85
86!      filin = TRIM(DIRNAMEINP)//TRIM(filin)
87
88!       write(42,*) 'type de socle '
89!       write(42,*) 'mksedim donne par ',filin
90
91!open(11,file=filin)
92!do i=1,6
93!      read(11,*)  ! 6 lignes de commentaires
94!end do
95
96
97! lecture proprement dite de la carte sediment
98!do j=ny, 1, -1
99!    read(11,'(87i1)') (mksedim(i,j), i=1,nx)
100!end do
101
102
103!close(11)
104
105
106! test pour une calotte axi-symetrique avec glissement Cr partout
107!write(num_rep_42,*) 'masque sediment mis a 1 '
108!mksedim(:,:)=min(mksedim(:,:),1)
109
110! a la frontiere le demi noeud est sediment.
111! calcul du masque sur les demi mailles : frontiere -> sedim
112!mkxsedim(:,:)=max(mksedim(:,:),eoshift(mksedim(:,:),shift=-1,boundary=0,dim=1))
113!mkysedim(:,:)=max(mksedim(:,:),eoshift(mksedim(:,:),shift=-1,boundary=0,dim=2))
114! write(42,*) 'sur les bords zone sediment, masque exterieur '
115
116
117!mkxsedim(:,:)=min(mksedim(:,:),eoshift(mksedim(:,:),shift=-1,boundary=0,dim=1))
118!mkysedim(:,:)=min(mksedim(:,:),eoshift(mksedim(:,:),shift=-1,boundary=0,dim=2))
119
120!write(num_rep_42,*) 'sur les bords zone sediment, demi-maille sediment si les 2 mailles le sont'
121
122    mslid_mx(:,:)=mkxsedim(:,:)
123    mslid_my(:,:)=mkysedim(:,:)
124
125!    Cr=1.e5    ! en a-1  coefficient hard rock (loi en puissance 3)
126!    Cs=500.    ! en a-1  coefficient sediment (lineaire)
127
128
129    write(42,*)'coefficients:   Cr=',Cr,'    Cs=',Cs,' (a-1)'
130!  write(num_rep_42,*)'glissement seulement si les 2 voisins sont au point de fusion'
131
132
133
134    return
135  end subroutine init_sliding
136
137
138  subroutine sliding
139
140
141! cette routine calcule le terme ddbx et ddby dans le cas Heino
142!      -------------------------------------------------
143!      GLISSEMENT
144!      DDBX et DDBY termes dus au glissement
145!      relation avec la vitesse de glissement UXB et UYB
146!      UXB=-DDBX*SDX et UYB=-DDBY*SDY
147!      -------------------------------------------------
148
149    implicit none
150
151    integer :: i_moins1,j_moins1,i_plus1,j_plus1
152
153!!$real, dimension(nx,ny) :: hw_mx
154!!$real, dimension(nx,ny) :: hw_my
155
156!cdc if (itracebug.eq.1)  call tracebug(' Heino: entree dans routine sliding')
157
158
159    call moy_mxmy(nx,ny,Hwater,hw_mx,hw_my)
160
161
162
163    slidy: do j=1,NY
164       do I=1,NX
165
166          pose_y:  if (.not.flgzmy(i,j)) then
167
168             if (HMY(I,J).lt.1.) then        ! glace tres fine
169                ddby(i,j)=0.
170                mslid_my(i,j)=0
171
172!   else if ((T(i,j,nz).lt.Tpmp(i,j,nz)).and.(T(i,j-1,nz).lt.Tpmp(i,j-1,nz))) then
173!      ddby(i,j)=0.                 ! pas de glissement seulement
174!                                   ! si les 2 voisins sont froids
175
176!   else if ((T(i,j,nz).lt.Tpmp(i,j,nz)).or.(T(i,j-1,nz).lt.Tpmp(i,j-1,nz))) then
177
178!      else if ((ibase(i,j).eq.1).or.(ibase(i,j-1).eq.1)) then
179!      else if ((ibase(i,j).eq.1).and.(ibase(i,j-1).eq.1)) then
180
181             else if (hw_my(i,j).lt.seuil_gliss) then
182                ddby(i,j)=0.                 ! pas de glissement si tous les voisins
183                                   ! sont a sec
184                mslid_my(i,j)=0
185!      flgzmy(i,j)=.false. ! on ne bloque pas les ice-shelves
186
187             else                            ! glissement
188!      coefslid=(hwater(i,j)+hwater(i,j-1))*0.5
189                coefslid=hw_my(i,j)
190                coefslid=max(coefslid,0.)
191                coefslid=min(coefslid/coefmax,1.)
192!      coefslid=1.
193
194
195!jalv le mksedim maintenant dans hudson: mksedim=0 ==> hard rock (=1 ==> ocean)
196!                if (mkysedim(i,j).eq.1) then  ! loi hard rock
197!                if (mkysedim(i,j).eq.0) then  ! loi hard rock
198
199!jalv: pas besoin de if mksedim... si l'on enleve le sliding de type Calov
200! (y a pas besoin de specifier dans quel type de sediments on est)
201!                if ((mkysedim(i,j).eq.0).or.(mkysedim(i,j).eq.2)) then  ! loi hard rock
202                ddby(i,j)=Cr*slope2my(i,j)*Hmy(i,j)*coefslid          ! loi hard rock
203
204
205                mslid_my(i,j)=1
206
207 
208!jalv:si on met en comentaire les lignes suivantes:
209!     test pas de sliding de type Calov dans la zone sediment
210!     (sliding normal avec la dependence en Cr et non en Cs)
211
212!               else if ((mkysedim(i,j).eq.2)) then  ! loi sediment
213!                   ddby(i,j)=Cs*Hmy(i,j)*coefslid
214
215!                   mslid_my(i,j)=2
216
217!jalv: test longslope
218
219!         if (abs(sdy(i,j)).gt.1.e-8) then
220!            j_1=max(j-1,1)
221!            j1=min(j+1,ny)
222!            longslope=(sdy(i,j_1)+sdy(i,j)+sdy(i,j1))/3.
223!            ddby(i,j)=ddby(i,j)/sdy(i,j)*longslope
224!         endif
225
226!jalv: pas besoin maintenant si c'est bien fait
227!                else
228!                   ddby(i,j)=0.
229!                   mslid_my(i,j)=0
230
231!                endif
232             endif
233
234          endif pose_y
235       end do
236    end do slidy
237
238
239    slidx: do j=1,NY
240       do i=1,NY
241          pose_x:  if (.not.flgzmx(i,j)) then
242
243             if (HMX(I,J).lt.1.) then        ! glace tres fine
244                ddbx(i,j)=0.
245                mslid_mx(i,j)=0
246
247!   else if ((T(i,j,nz).lt.Tpmp(i,j,nz)).and.(T(i-1,j,nz).lt.Tpmp(i-1,j,nz))) then
248!      ddbx(i,j)=0.                 ! pas de glissement seulement
249!                                   ! si les 2 voisins sont froids
250
251!   else if ((T(i,j,nz).lt.Tpmp(i,j,nz)).or.(T(i-1,j,nz).lt.Tpmp(i-1,j,nz))) then
252!     else if ((ibase(i,j).eq.1).or.(ibase(i-1,j).eq.1)) then
253
254!      else if ((ibase(i,j).eq.1).and.(ibase(i-1,j).eq.1)) then
255             else if (hw_mx(i,j).lt.seuil_gliss) then
256               
257                ddbx(i,j)=0.                 ! pas de glissement si au moins un seul
258                                   ! des 2 voisins est froid
259                mslid_mx(i,j)=0
260!      flgzmx(i,j)=.false. ! on ne bloque pas les ice-shelves
261
262
263             else                            ! glissement
264
265!      coefslid=(hwater(i,j)+hwater(i-1,j))*0.5
266                coefslid=hw_mx(i,j)
267                coefslid=max(coefslid,0.)
268                coefslid=min(coefslid/coefmax,1.)
269!      coefslid=1.
270
271!jalv le mksedim maintenant dans hudson: mksedim=0 ==> hard rock (=1 ==> ocean)
272!                if (mkxsedim(i,j).eq.0) then  ! loi hard rock
273
274!jalv: pas besoin de if mksedim... si l'on enleve le sliding de type Calov
275! (y a pas besoin de specifier dans quel type de sediments on est)
276!                if ((mkxsedim(i,j).eq.0).or.(mkysedim(i,j).eq.2)) then  ! loi hard rock
277                ddbx(i,j)=Cr*slope2mx(i,j)*Hmx(i,j)*coefslid
278                mslid_mx(i,j)=1
279
280
281!jalv: si on met en comentaire les lignes suivantes:
282!      test pas de sliding de type Calov dans la zone sediment
283!      (sliding normal avec la dependence en Cr et non en Cs)
284!                else if ((mkxsedim(i,j).eq.2)) then  ! loi sediment
285!                   ddbx(i,j)=Cs*Hmx(i,j)*coefslid
286!                   mslid_mx(i,j)=2
287
288!jalv: test longslope
289
290!         if (abs(sdx(i,j)).gt.1.e-8) then
291!            i_1=max(i-1,1)
292!            i1=min(i+1,nx)
293!            longslope=(sdx(i_1,j)+sdx(i,j)+sdx(i1,j))/3.
294!            ddbx(i,j)=ddbx(i,j)/sdx(i,j)*longslope
295!         endif
296
297!jalv: pas besoin maintenant si c'est bien fait
298!                else
299!                   ddbx(i,j)=0.
300!                   mslid_mx(i,j)=0
301
302!                endif
303             endif
304
305          endif pose_x
306       end do
307    end do slidx
308
309    return
310  end subroutine  sliding
311
312
313END MODULE sliding_hudson
314
315
Note: See TracBrowser for help on using the repository browser.