source: trunk/SOURCES/dragging_param_beta_mod.f90 @ 332

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

Dragging: distinction between drag computation and mask update / flottab only updates masks

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 
103return
104end subroutine init_dragging
105!________________________________________________________________________________
106
107!> SUBROUTINE: dragging
108!! Defini la localisation des streams et le frottement basal
109!>
110!-------------------------------------------------------------------------
111subroutine dragging   ! defini le frottement basal
112!$ USE OMP_LIB
113
114!         les iles n'ont pas de condition neff mais ont la condition toblim
115!         (ce bloc etait avant dans flottab)
116
117!$OMP PARALLEL
118
119!    -- !$OMP DO
120!    -- do j=2,ny
121!    --    do i=2,nx
122!    --       ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
123!    --       ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
124!    --    end do
125!    -- end do
126!    -- !$OMP END DO
127
128
129
130! calcul de neff (pression effective sur noeuds majeurs)
131if (sum(neffmx(:,:)).le.0.) neffmx(:,:) =1.e8
132if (sum(neffmy(:,:)).le.0.) neffmy(:,:) =1.e8
133
134!$OMP DO
135do j=1,ny-1
136   do i=1,nx-1
137        neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4
138   enddo
139enddo
140!$OMP END DO
141!aurel, for the borders:
142!$OMP WORKSHARE
143neff(:,ny)=1.e5
144neff(nx,:)=1.e5
145! calcul de hwat (staggered)
146hwatmx(:,:)=0.
147hwatmy(:,:)=0.
148!$OMP END WORKSHARE
149!$OMP DO
150do i=2,nx
151   hwatmx(i,:)=(hwater(i-1,:)+hwater(i,:))/2.
152enddo
153!$OMP END DO
154!$OMP DO
155do j=2,ny
156   hwatmy(:,j)=(hwater(:,j-1)+hwater(:,j))/2.
157enddo
158!$OMP END DO
159
160!$OMP WORKSHARE
161
162! new parametrisation of beta on Neff:
163betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo)
164betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo)
165
166where (ilemx(:,:)) betamx(:,:) = betamx(:,:) * coef_ile
167where (ilemy(:,:)) betamy(:,:) = betamy(:,:) * coef_ile
168
169
170betamx(:,:)=max(betamx(:,:),betamin)
171betamy(:,:)=max(betamy(:,:),betamin)
172betamx(:,:)=min(betamx(:,:),betamax)
173betamy(:,:)=min(betamy(:,:),betamax)
174
175where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
176where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
177
178!$OMP END WORKSHARE
179
180! points flottants
181!$OMP DO
182do j=2,ny
183   do i=2,nx
184      if (flot(i,j).and.(flot(i-1,j))) then
185         betamx(i,j)=0.
186      end if
187      if (flot(i,j).and.(flot(i,j-1))) then
188         betamy(i,j)=0.
189      end if
190   end do
191end do
192!$OMP END DO
193
194
195!$OMP DO
196do j=2,ny-1
197   do i=2,nx-1
198      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
199          + (betamy(i,j)+betamy(i,j+1)))*0.25
200   end do
201end do
202!$OMP END DO
203!$OMP END PARALLEL
204
205return
206end subroutine dragging
207!________________________________________________________________________________
208
209!> SUBROUTINE: mstream_dragging
210!! Defini la localisation des streams
211!>
212!-------------------------------------------------------------------------
213subroutine mstream_dragging   ! defini la localisation des streams
214
215!$OMP PARALLEL
216
217!$OMP WORKSHARE
218fleuvemx(:,:)=.false.
219fleuvemy(:,:)=.false.
220cote(:,:)=.false.
221gzmx(:,:)=.true.
222gzmy(:,:)=.true.
223flgzmx(:,:)=.false.
224flgzmy(:,:)=.false.
225!$OMP END WORKSHARE
226
227! detection des cotes
228!$OMP DO
229do  j=2,ny-1 
230   do i=2,nx-1 
231      if ((.not.flot(i,j)).and. & 
232      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
233         cote(i,j)=.true.
234      endif
235   end do
236end do
237!$OMP END DO
238
239
240! calcul de gzmx
241
242! points flottants : flgzmx mais pas gzmx
243!$OMP DO
244do j=2,ny
245   do i=2,nx
246      if (flot(i,j).and.(flot(i-1,j))) then
247         flotmx(i,j)=.true.
248         flgzmx(i,j)=.true.
249         gzmx(i,j)=.false.
250         betamx(i,j)=0.
251      end if
252      if (flot(i,j).and.(flot(i,j-1))) then
253         flotmy(i,j)=.true.
254         flgzmy(i,j)=.true.
255         gzmy(i,j)=.false.
256         betamy(i,j)=0.
257      end if
258   end do
259end do
260!$OMP END DO
261
262!--------- autres criteres
263
264! points poses
265!$OMP WORKSHARE
266gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
267gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
268
269flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
270flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
271
272fleuvemx(:,:)=gzmx(:,:)
273fleuvemy(:,:)=gzmy(:,:)
274!$OMP END WORKSHARE
275
276!$OMP END PARALLEL
277
278return
279end subroutine mstream_dragging
280
281end module dragging_param_beta
282
Note: See TracBrowser for help on using the repository browser.