source: trunk/SOURCES/diffusiv-polyn-0.6.f90 @ 68

Last change on this file since 68 was 68, checked in by dumas, 8 years ago

No sliding in module_choix | Antarctic 15km configuration usable | Use depth instade of dist_talus in bmelt-ant-regions_mod.f90

File size: 6.1 KB
Line 
1!> \file diffusiv-polyn-0.5.f90
2!! Contient la subroutine de calcul diffusivites
3!<
4
5!> SUBROUTINE: diffusiv
6!! Calcule de  diffusivites
7!! \author Catherine
8!! \date   18 mars 2006
9!! @note Used modules:
10!! @note    - use module3D_phy
11!! @note    - use module_choix
12!!
13!<
14      subroutine diffusiv()
15
16!     ===============================================================
17!     **     DIFFUSIVITIES                    18 mars 2006   ***
18!     glissement avec reserve d'eau 29 Avril 97
19!     choix sur le type de discretisation (locale iloc=1 )   
20!
21!     Modification Vince --> 2 fev 95
22!     pour couplage avec Shelf
23!
24!     Modification C. Dumas Dec 99
25!     DDX, SA,S2A,... ont une dimension en plus (loi de deformation)
26!
27!     Modifications  de Cat du 11 juin 2000
28!     Modif Christophe 24 janv 2002 : passage au fortran 90
29!     pour optimisation de la routine
30!     Modif Vincent dec 2003 : on utilise plus frontmx/y > commentaires
31!
32!     Modif Cat mars 2003
33!                  modif mars 2007 : moyennes de S2a avec S2a_mx et S2a_my
34
35!     ===============================================================
36
37USE module3D_phy
38
39use module_choix
40
41implicit none
42
43REAL :: ulgliss,ultot,uldef,glenexp
44REAL :: INV_4DX ! inverse de dx pour eviter les divisions =1/(4*dx)
45REAL :: INV_4DY ! inverse de dy pour eviter les divisions =1/(4*dy)
46
47if (itracebug.eq.1)  call tracebug(' Entree dans routine diffusiv')
48
49!     la valeur de ff est donnee dans initial
50!     limitation de la vitesse de glissement de déformation et totale
51ulgliss = V_limit
52ultot   = V_limit
53uldef   = 2000.
54
55
56INV_4DX=1/(4*dx)
57INV_4DY=1/(4*dy)
58
59! initialisation de la diffusion
60diffmx(:,:)=0.
61diffmy(:,:)=0.
62
63
64! calcul de SDX, SDY et de la pente au carre en mx et en my :
65
66call slope_surf
67
68call sliding     ! au sens vitesse de glissement
69
70!  le glissement est maintenant dans un module a part choisi dans le module choix
71!  pour l'instant seules les lois Heino (loigliss=4) et Bindshadler(loigliss=2)
72!  sont programmees.
73
74!      ddbx et ddby termes dus au glissement
75!      relation avec la vitesse de glissement ubx et uby
76!      ubx=-ddbx*sdx et uby=-ddby*sdy
77!      -------------------------------------------------
78
79! Calcul de la vitesse de glissement, qu'elle soit calculee par diagno ou par sliding
80!------------------------------------------------------------------------------------
81!
82where (flgzmx(:,:))
83    ubx(:,:) = uxflgz(:,:)
84elsewhere
85    ubx(:,:) = ddbx(:,:)* (-sdx(:,:))     ! on pourrait mettre une limitation
86endwhere
87
88where (flgzmy(:,:))
89    uby(:,:) = uyflgz(:,:)
90elsewhere
91    uby(:,:) = ddby(:,:)* (-sdy(:,:))
92endwhere
93
94
95if (itracebug.eq.1) call tracebug(' diffusiv  apres calcul glissement')
96
97
98
99! Deformation SIA : Calcul de ddy  et ddx pour chaque valeur de iglen
100!------------------------------------------------------------------
101
102do iglen=n1poly,n2poly
103   glenexp=max(0.,(glen(iglen)-1.)/2.)
104
105
106   do j=1,ny
107      do i=1,nx
108
109         if (.not.flotmy(i,j)) then !  On calcule quand la deformation même pour les ice streams
110            ddy(i,j,iglen)=((slope2my(i,j)**  &  ! slope2my calcule en debut de routine
111                 glenexp)*(rog)**glen(iglen)) &
112                 *hmy(i,j)**(glen(iglen)+1)
113
114
115         endif
116      end do
117   end do
118end do
119
120do iglen=n1poly,n2poly
121   glenexp=max(0.,(glen(iglen)-1.)/2.)
122
123
124   do j=1,ny 
125      do i=1,nx
126
127         if (.not.flotmx(i,j)) then      ! bug y->x corrige le 5 avril 2008 (
128                                   
129            ddx(i,j,iglen)=((slope2mx(i,j)**  & ! slope2mx calcule en debut de routine
130                 glenexp)*(rog)**glen(iglen)) &
131                 *hmx(i,j)**(glen(iglen)+1)
132         endif
133      end do
134   end do
135end do
136
137! vitesse due a la déformation-----------------------------------------------
138
139
140ydef: do j=2,ny
141       do i=1,nx 
142
143        if (.not.flotmy(i,j)) then !  On calcule la deformation même pour les ice streams
144
145         uydef(i,j)=0.
146         diffmy(i,j)=0.
147
148         do iglen=n1poly,n2poly
149           diffmy(i,j)=diffmy(i,j)+ddy(i,j,iglen)*(s2a_my(i,j,1,iglen))
150         end do
151
152         uydef(i,j)=diffmy(i,j)*(-sdy(i,j))    ! partie deformation
153         uydef(i,j)=min(uydef(i,j),uldef)
154         uydef(i,j)=max(uydef(i,j),-uldef)
155
156         diffmy(i,j)=diffmy(i,j)+ddby(i,j)     ! pour utilisation comme diffusion, a multiplier par H
157
158         uybar(i,j)=uby(i,j)+uydef(i,j)
159
160
161         else    ! points flottants
162            uydef(i,j)=0.   ! normalement uxbed est attribue         
163         endif
164    end do
165   end do ydef
166
167xdef: do j=1,ny
168       do i=2,nx 
169
170        if (.not.flotmx(i,j)) then !  On calcule la deformation même pour les ice streams
171
172         uxdef(i,j)=0.
173         diffmx(i,j)=0.
174
175         do iglen=n1poly,n2poly
176           diffmx(i,j)=diffmx(i,j)+ddx(i,j,iglen)*(s2a_mx(i,j,1,iglen))
177         end do
178
179         uxdef(i,j)=diffmx(i,j)*(-sdx(i,j))    ! partie deformation
180         uxdef(i,j)=min(uxdef(i,j),uldef)
181         uxdef(i,j)=max(uxdef(i,j),-uldef)
182
183
184         diffmx(i,j)=diffmx(i,j)+ddbx(i,j)     ! pour utilisation comme diffusion,
185                                               ! a multiplier par H
186         uxbar(i,j)=ubx(i,j)+uxdef(i,j)
187
188
189         else    ! points flottants
190            uxdef(i,j)=0.   ! normalement uxbed est attribue         
191         endif
192    end do
193   end do xdef
194
195if (itracebug.eq.1) call tracebug(' diffusiv  avant limit')
196
197
198! limitation par ultot----------------------------------------------------------
199
200! la limitation selon x
201where(.not.flgzmx(:,:))
202   uxbar(:,:)=min(uxbar(:,:),ultot)
203   uxbar(:,:)=max(uxbar(:,:),-ultot)
204end where
205
206! la limitation selon y
207where(.not.flgzmy(:,:))
208   uybar(:,:)=min(uybar(:,:),ultot)
209   uybar(:,:)=max(uybar(:,:),-ultot)
210end where
211
212!    cas particulier des sommets ou il ne reste plus de neige.
213do j=2,ny-1
214   do i=2,nx-1
215      if ((H(i,j).le.1.).and.(.not.flot(i,j))) then
216         uxbar(i+1,j)=min(uxbar(i+1,j),ultot)  ! n'agit que si ux ->
217         uxbar(i,j)=max(uxbar(i,j),-ultot)     ! n'agit que si ux <-
218         uybar(i,j+1)=min(uybar(i,j+1),ultot)
219         uybar(i,j)=max(uybar(i,j),-ultot)
220      endif
221   end do
222end do
223
224if (itracebug.eq.1)  call tracebug(' fin diffusiv ')
225
226return
227end subroutine diffusiv
228
Note: See TracBrowser for help on using the repository browser.