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

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

OpenMP parallelization in conserv-mass-adv-diff_sept2009_mod.f90, diffusiv-polyn-0.6.f90, dragging_neff_slope_mod.f90, eaubasale-0.5_mod.f90 and relaxation_water_diffusion.f90

File size: 6.4 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!$ USE OMP_LIB
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
60!$OMP PARALLEL
61!$OMP WORKSHARE
62diffmx(:,:)=0.
63diffmy(:,:)=0.
64!$OMP END WORKSHARE
65!$OMP END PARALLEL
66
67! calcul de SDX, SDY et de la pente au carre en mx et en my :
68
69call slope_surf
70
71call sliding     ! au sens vitesse de glissement
72
73!  le glissement est maintenant dans un module a part choisi dans le module choix
74!  pour l'instant seules les lois Heino (loigliss=4) et Bindshadler(loigliss=2)
75!  sont programmees.
76
77!      ddbx et ddby termes dus au glissement
78!      relation avec la vitesse de glissement ubx et uby
79!      ubx=-ddbx*sdx et uby=-ddby*sdy
80!      -------------------------------------------------
81
82! Calcul de la vitesse de glissement, qu'elle soit calculee par diagno ou par sliding
83!------------------------------------------------------------------------------------
84!$OMP PARALLEL
85!$OMP WORKSHARE
86where (flgzmx(:,:))
87    ubx(:,:) = uxflgz(:,:)
88elsewhere
89    ubx(:,:) = ddbx(:,:)* (-sdx(:,:))     ! on pourrait mettre une limitation
90endwhere
91
92where (flgzmy(:,:))
93    uby(:,:) = uyflgz(:,:)
94elsewhere
95    uby(:,:) = ddby(:,:)* (-sdy(:,:))
96endwhere
97!$OMP END WORKSHARE
98!$OMP END PARALLEL
99if (itracebug.eq.1) call tracebug(' diffusiv  apres calcul glissement')
100
101
102
103! Deformation SIA : Calcul de ddy  et ddx pour chaque valeur de iglen
104!------------------------------------------------------------------
105do iglen=n1poly,n2poly
106   glenexp=max(0.,(glen(iglen)-1.)/2.)
107  !$OMP PARALLEL
108  !$OMP DO
109   do j=1,ny
110      do i=1,nx
111         if (.not.flotmy(i,j)) then !  On calcule quand la deformation même pour les ice streams
112            ddy(i,j,iglen)=((slope2my(i,j)**  &  ! slope2my calcule en debut de routine
113                 glenexp)*(rog)**glen(iglen)) &
114                 *hmy(i,j)**(glen(iglen)+1)
115         endif
116      end do
117   end do
118   !$OMP END DO
119   !$OMP END PARALLEL
120end do
121
122
123
124do iglen=n1poly,n2poly
125   glenexp=max(0.,(glen(iglen)-1.)/2.)
126   !$OMP PARALLEL
127   !$OMP DO
128   do j=1,ny 
129      do i=1,nx
130         if (.not.flotmx(i,j)) then      ! bug y->x corrige le 5 avril 2008 (
131            ddx(i,j,iglen)=((slope2mx(i,j)**  & ! slope2mx calcule en debut de routine
132                 glenexp)*(rog)**glen(iglen)) &
133                 *hmx(i,j)**(glen(iglen)+1)
134         endif
135      end do
136   end do
137   !$OMP END DO
138   !$OMP END PARALLEL
139end do
140
141
142! vitesse due a la déformation-----------------------------------------------
143!$OMP PARALLEL
144!$OMP DO
145ydef: do j=2,ny
146       do i=1,nx 
147
148        if (.not.flotmy(i,j)) then !  On calcule la deformation même pour les ice streams
149
150         uydef(i,j)=0.
151         diffmy(i,j)=0.
152
153         do iglen=n1poly,n2poly
154           diffmy(i,j)=diffmy(i,j)+ddy(i,j,iglen)*(s2a_my(i,j,1,iglen))
155         end do
156
157         uydef(i,j)=diffmy(i,j)*(-sdy(i,j))    ! partie deformation
158         uydef(i,j)=min(uydef(i,j),uldef)
159         uydef(i,j)=max(uydef(i,j),-uldef)
160
161         diffmy(i,j)=diffmy(i,j)+ddby(i,j)     ! pour utilisation comme diffusion, a multiplier par H
162
163         uybar(i,j)=uby(i,j)+uydef(i,j)
164
165
166         else    ! points flottants
167            uydef(i,j)=0.   ! normalement uxbed est attribue         
168         endif
169    end do
170   end do ydef
171!$OMP END DO
172
173!$OMP DO
174xdef: do j=1,ny
175       do i=2,nx 
176
177        if (.not.flotmx(i,j)) then !  On calcule la deformation même pour les ice streams
178
179         uxdef(i,j)=0.
180         diffmx(i,j)=0.
181
182         do iglen=n1poly,n2poly
183           diffmx(i,j)=diffmx(i,j)+ddx(i,j,iglen)*(s2a_mx(i,j,1,iglen))
184         end do
185
186         uxdef(i,j)=diffmx(i,j)*(-sdx(i,j))    ! partie deformation
187         uxdef(i,j)=min(uxdef(i,j),uldef)
188         uxdef(i,j)=max(uxdef(i,j),-uldef)
189
190
191         diffmx(i,j)=diffmx(i,j)+ddbx(i,j)     ! pour utilisation comme diffusion,
192                                               ! a multiplier par H
193         uxbar(i,j)=ubx(i,j)+uxdef(i,j)
194
195
196         else    ! points flottants
197            uxdef(i,j)=0.   ! normalement uxbed est attribue         
198         endif
199    end do
200   end do xdef
201!$OMP END DO
202
203if (itracebug.eq.1) call tracebug(' diffusiv  avant limit')
204
205
206! limitation par ultot----------------------------------------------------------
207
208! la limitation selon x
209!$OMP WORKSHARE
210where(.not.flgzmx(:,:))
211   uxbar(:,:)=min(uxbar(:,:),ultot)
212   uxbar(:,:)=max(uxbar(:,:),-ultot)
213end where
214
215! la limitation selon y
216where(.not.flgzmy(:,:))
217   uybar(:,:)=min(uybar(:,:),ultot)
218   uybar(:,:)=max(uybar(:,:),-ultot)
219end where
220!$OMP END WORKSHARE
221
222!    cas particulier des sommets ou il ne reste plus de neige.
223!$OMP DO
224do j=2,ny-1
225   do i=2,nx-1
226      if ((H(i,j).le.1.).and.(.not.flot(i,j))) then
227         uxbar(i+1,j)=min(uxbar(i+1,j),ultot)  ! n'agit que si ux ->
228         uxbar(i,j)=max(uxbar(i,j),-ultot)     ! n'agit que si ux <-
229         uybar(i,j+1)=min(uybar(i,j+1),ultot)
230         uybar(i,j)=max(uybar(i,j),-ultot)
231      endif
232   end do
233end do
234!$OMP END DO
235!$OMP END PARALLEL
236if (itracebug.eq.1)  call tracebug(' fin diffusiv ')
237
238return
239end subroutine diffusiv
240
Note: See TracBrowser for help on using the repository browser.