source: trunk/SOURCES/dragging_param_beta_mod.f90 @ 334

Last change on this file since 334 was 333, checked in by aquiquet, 3 years ago

Addition of a non-linear Coulomb friction law + iterations on velocity computation for convergence / for m=1 gives results very similar to the standard linear law in dragging_param_beta

File size: 7.3 KB
Line 
1!> \file dragging_neff_contmaj_mod.f90
2!! Defini les zones de stream avec 3 criteres
3!<
4
5!> \namespace  dragging_neff_contmaj
6!! Defini les zones de stream avec 3 criteres
7!! \author ...
8!! \date ...
9!! @note Les trois criteres sont:
10!! @note   * un critere sur la pression effective
11!! @note   * un critere de continuite depuis la cote
12!! @note   * un critere sur la courbure du socle (si negatif, vallees)
13!! @note Used module
14!! @note   - use module3D_phy
15!! @note   - use param_phy_mod
16!<
17
18module dragging_param_beta
19
20! Defini les zones de stream avec :
21! * un critere sur la pression effective
22! * un critere de continuite depuis la cote
23! * un critere sur la courbure du socle (si negatif, vallees)
24
25use module3d_phy
26use param_phy_mod
27use interface_input
28use io_netcdf_grisli
29use fake_beta_iter_vitbil_mod
30
31implicit none
32
33logical,dimension(nx,ny) :: cote
34
35real,dimension(nx,ny) :: beta_param ! local variable
36
37real :: betamin   ! betamin : (Pa) frottement mini sous les streams
38
39real :: beta_slope        ! = A in:     B = A x Neff ** K
40real :: beta_expo         ! = K in:     B = A x Neff ** K
41
42real,dimension(nx,ny) :: neff   ! pression effective noeuds majeurs
43real,dimension(nx,ny) :: hwatmx ! hauteur d'eau staggered grille - afq
44real,dimension(nx,ny) :: hwatmy ! hauteur d'eau staggered grille - afq
45
46real :: coef_ile      ! coef frottement zones iles   (0.1)
47
48real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
49real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
50real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
51real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
52logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
53
54contains
55!-------------------------------------------------------------------------------
56
57!> SUBROUTINE: init_dragging
58!! Cette routine fait l'initialisation du dragging.
59!>
60subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
61
62implicit none
63
64namelist/drag_param_beta/beta_slope,beta_expo,betamax,betamin,coef_ile
65
66if (itracebug.eq.1)  call tracebug(' dragging avec neff et slope')
67
68
69! formats pour les ecritures dans 42
70428 format(A)
71
72! lecture des parametres du run                      block drag neff
73!--------------------------------------------------------------------
74
75rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
76read(num_param,drag_param_beta)
77
78write(num_rep_42,428)'!___________________________________________________________' 
79write(num_rep_42,428) '&drag_param_beta          ! nom du bloc dragging param beta'
80write(num_rep_42,*)
81write(num_rep_42,*) 'beta_slope            = ', beta_slope
82write(num_rep_42,*) 'beta_expo             = ', beta_expo
83write(num_rep_42,*) 'betamax               = ', betamax
84write(num_rep_42,*) 'betamin               = ', betamin
85write(num_rep_42,*)'/'                           
86write(num_rep_42,428) '! Slope & expo of beta = - slope x Neff ** expo'
87
88inv_beta=0
89!-------------------------------------------------------------------
90! masque stream
91mstream_mx(:,:)=1
92mstream_my(:,:)=1
93
94
95! coefficient permettant de modifier le basal drag.
96drag_mx(:,:)=1.
97drag_my(:,:)=1.
98     
99betamax_2d(:,:) = betamax
100
101!afq -- toblim = 0. !0.7e5 ! afq -- pour les iles, mais est-ce vraiment utile?
102
103niter_nolin = 1
104
105return
106end subroutine init_dragging
107!________________________________________________________________________________
108
109!> SUBROUTINE: dragging
110!! Defini la localisation des streams et le frottement basal
111!>
112!-------------------------------------------------------------------------
113subroutine dragging   ! defini le frottement basal
114!$ USE OMP_LIB
115
116!         les iles n'ont pas de condition neff mais ont la condition toblim
117!         (ce bloc etait avant dans flottab)
118
119!$OMP PARALLEL
120
121!    -- !$OMP DO
122!    -- do j=2,ny
123!    --    do i=2,nx
124!    --       ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
125!    --       ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
126!    --    end do
127!    -- end do
128!    -- !$OMP END DO
129
130
131
132! calcul de neff (pression effective sur noeuds majeurs)
133if (sum(neffmx(:,:)).le.0.) neffmx(:,:) =1.e8
134if (sum(neffmy(:,:)).le.0.) neffmy(:,:) =1.e8
135
136!$OMP DO
137do j=1,ny-1
138   do i=1,nx-1
139        neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4
140   enddo
141enddo
142!$OMP END DO
143!aurel, for the borders:
144!$OMP WORKSHARE
145neff(:,ny)=1.e5
146neff(nx,:)=1.e5
147! calcul de hwat (staggered)
148hwatmx(:,:)=0.
149hwatmy(:,:)=0.
150!$OMP END WORKSHARE
151!$OMP DO
152do i=2,nx
153   hwatmx(i,:)=(hwater(i-1,:)+hwater(i,:))/2.
154enddo
155!$OMP END DO
156!$OMP DO
157do j=2,ny
158   hwatmy(:,j)=(hwater(:,j-1)+hwater(:,j))/2.
159enddo
160!$OMP END DO
161
162!$OMP WORKSHARE
163
164! new parametrisation of beta on Neff:
165betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo)
166betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo)
167
168where (ilemx(:,:)) betamx(:,:) = betamx(:,:) * coef_ile
169where (ilemy(:,:)) betamy(:,:) = betamy(:,:) * coef_ile
170
171
172betamx(:,:)=max(betamx(:,:),betamin)
173betamy(:,:)=max(betamy(:,:),betamin)
174betamx(:,:)=min(betamx(:,:),betamax)
175betamy(:,:)=min(betamy(:,:),betamax)
176
177where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
178where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
179
180!$OMP END WORKSHARE
181
182! points flottants
183!$OMP DO
184do j=2,ny
185   do i=2,nx
186      if (flot(i,j).and.(flot(i-1,j))) then
187         betamx(i,j)=0.
188      end if
189      if (flot(i,j).and.(flot(i,j-1))) then
190         betamy(i,j)=0.
191      end if
192   end do
193end do
194!$OMP END DO
195
196
197!$OMP DO
198do j=2,ny-1
199   do i=2,nx-1
200      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
201          + (betamy(i,j)+betamy(i,j+1)))*0.25
202   end do
203end do
204!$OMP END DO
205!$OMP END PARALLEL
206
207return
208end subroutine dragging
209!________________________________________________________________________________
210
211!> SUBROUTINE: mstream_dragging
212!! Defini la localisation des streams
213!>
214!-------------------------------------------------------------------------
215subroutine mstream_dragging   ! defini la localisation des streams
216
217!$OMP PARALLEL
218
219!$OMP WORKSHARE
220fleuvemx(:,:)=.false.
221fleuvemy(:,:)=.false.
222cote(:,:)=.false.
223gzmx(:,:)=.true.
224gzmy(:,:)=.true.
225flgzmx(:,:)=.false.
226flgzmy(:,:)=.false.
227!$OMP END WORKSHARE
228
229! detection des cotes
230!$OMP DO
231do  j=2,ny-1 
232   do i=2,nx-1 
233      if ((.not.flot(i,j)).and. & 
234      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
235         cote(i,j)=.true.
236      endif
237   end do
238end do
239!$OMP END DO
240
241
242! calcul de gzmx
243
244! points flottants : flgzmx mais pas gzmx
245!$OMP DO
246do j=2,ny
247   do i=2,nx
248      if (flot(i,j).and.(flot(i-1,j))) then
249         flotmx(i,j)=.true.
250         flgzmx(i,j)=.true.
251         gzmx(i,j)=.false.
252         betamx(i,j)=0.
253      end if
254      if (flot(i,j).and.(flot(i,j-1))) then
255         flotmy(i,j)=.true.
256         flgzmy(i,j)=.true.
257         gzmy(i,j)=.false.
258         betamy(i,j)=0.
259      end if
260   end do
261end do
262!$OMP END DO
263
264!--------- autres criteres
265
266! points poses
267!$OMP WORKSHARE
268gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
269gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
270
271flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
272flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
273
274fleuvemx(:,:)=gzmx(:,:)
275fleuvemy(:,:)=gzmy(:,:)
276!$OMP END WORKSHARE
277
278!$OMP END PARALLEL
279
280return
281end subroutine mstream_dragging
282
283end module dragging_param_beta
284
Note: See TracBrowser for help on using the repository browser.