source: trunk/SOURCES/dragging_param_beta_nolin_mod.f90 @ 193

Last change on this file since 193 was 193, checked in by aquiquet, 6 years ago

Correcting a missing variable initialisation: toblim for parameterised beta

File size: 7.2 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_nolin
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,dimension(nx,ny) :: uslmag ! vitesse basale - afq
47
48real :: coef_ile      ! coef frottement zones iles   (0.1)
49
50real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
51real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
52real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
53real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
54logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
55
56contains
57!-------------------------------------------------------------------------------
58
59!> SUBROUTINE: init_dragging
60!! Cette routine fait l'initialisation du dragging.
61!>
62subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
63
64implicit none
65
66namelist/drag_param_beta/beta_slope,beta_expo,betamax,betamin,coef_ile
67
68if (itracebug.eq.1)  call tracebug(' dragging avec neff et slope')
69
70
71! formats pour les ecritures dans 42
72428 format(A)
73
74! lecture des parametres du run                      block drag neff
75!--------------------------------------------------------------------
76
77rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
78read(num_param,drag_param_beta)
79
80write(num_rep_42,428)'!___________________________________________________________' 
81write(num_rep_42,428) '&drag_param_beta_nolin    ! nom du bloc dragging param beta'
82write(num_rep_42,*)
83write(num_rep_42,*) 'beta_slope            = ', beta_slope
84write(num_rep_42,*) 'beta_expo             = ', beta_expo
85write(num_rep_42,*) 'betamax               = ', betamax
86write(num_rep_42,*) 'betamin               = ', betamin
87write(num_rep_42,*)'/'                           
88write(num_rep_42,428) '! Slope & expo of beta = (- slope x Neff ** expo ) / Uslmag'
89
90inv_beta=0
91!-------------------------------------------------------------------
92! masque stream
93mstream_mx(:,:)=1
94mstream_my(:,:)=1
95
96
97! coefficient permettant de modifier le basal drag.
98drag_mx(:,:)=1.
99drag_my(:,:)=1.
100     
101betamax_2d(:,:) = betamax
102
103toblim = 0.7e5 ! afq -- pour les iles, mais est-ce vraiment utile?
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 la localisation des streams et 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!$OMP DO
121do j=2,ny
122   do i=2,nx
123      ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
124      ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
125   end do
126end do
127!$OMP END DO
128
129
130! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
131
132
133!$OMP WORKSHARE
134fleuvemx(:,:)=.false.
135fleuvemy(:,:)=.false.
136cote(:,:)=.false.
137gzmx(:,:)=.true.
138gzmy(:,:)=.true.
139flgzmx(:,:)=.false.
140flgzmy(:,:)=.false.
141!$OMP END WORKSHARE
142
143! detection des cotes
144!$OMP DO
145do  j=2,ny-1 
146   do i=2,nx-1 
147      if ((.not.flot(i,j)).and. & 
148      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
149         cote(i,j)=.true.
150      endif
151   end do
152end do
153!$OMP END DO
154
155! calcul de neff (pression effective sur noeuds majeurs)
156if (sum(neffmx(:,:)).le.0.) neffmx(:,:) =1.e8
157if (sum(neffmy(:,:)).le.0.) neffmy(:,:) =1.e8
158
159!$OMP DO
160do j=1,ny-1
161   do i=1,nx-1
162        neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4
163   enddo
164enddo
165!$OMP END DO
166!aurel, for the borders:
167!$OMP WORKSHARE
168neff(:,ny)=1.e5
169neff(nx,:)=1.e5
170! calcul de hwat (staggered)
171hwatmx(:,:)=0.
172hwatmy(:,:)=0.
173!$OMP END WORKSHARE
174!$OMP DO
175do i=2,nx
176   hwatmx(i,:)=(hwater(i-1,:)+hwater(i,:))/2.
177enddo
178!$OMP END DO
179!$OMP DO
180do j=2,ny
181   hwatmy(:,j)=(hwater(:,j-1)+hwater(:,j))/2.
182enddo
183!$OMP END DO
184
185!$OMP WORKSHARE
186
187! new parametrisation of beta on Neff:
188betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo)
189betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo)
190
191where (ilemx(:,:)) betamx(:,:) = betamx(:,:) * coef_ile
192where (ilemy(:,:)) betamy(:,:) = betamy(:,:) * coef_ile
193
194!uslmag(:,:) = (ux(:,:,nz)**2+uy(:,:,nz)**2)**(0.5)
195!where (uslmag(:,:).gt.0)
196!   betamx(:,:) = betamx(:,:) / uslmag(:,:)
197!   betamy(:,:) = betamy(:,:) / uslmag(:,:)
198!endwhere
199
200where (abs(ux(:,:,nz)).gt.0)
201   betamx(:,:) = betamx(:,:) / min(abs(ux(:,:,nz)),100.)
202endwhere
203where (abs(uy(:,:,nz)).gt.0)
204   betamy(:,:) = betamy(:,:) / min(abs(uy(:,:,nz)),100.)
205endwhere
206
207betamx(:,:)=max(betamx(:,:),betamin)
208betamy(:,:)=max(betamy(:,:),betamin)
209betamx(:,:)=min(betamx(:,:),betamax)
210betamy(:,:)=min(betamy(:,:),betamax)
211
212where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
213where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
214
215!$OMP END WORKSHARE
216
217
218! calcul de gzmx
219
220! points flottants : flgzmx mais pas gzmx
221!$OMP DO
222do j=2,ny
223   do i=2,nx
224      if (flot(i,j).and.(flot(i-1,j))) then
225         flotmx(i,j)=.true.
226         flgzmx(i,j)=.true.
227         gzmx(i,j)=.false.
228         betamx(i,j)=0.
229      end if
230      if (flot(i,j).and.(flot(i,j-1))) then
231         flotmy(i,j)=.true.
232         flgzmy(i,j)=.true.
233         gzmy(i,j)=.false.
234         betamy(i,j)=0.
235      end if
236   end do
237end do
238!$OMP END DO
239
240!--------- autres criteres
241
242! points poses
243!$OMP WORKSHARE
244gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
245gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
246
247flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
248flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
249
250fleuvemx(:,:)=gzmx(:,:)
251fleuvemy(:,:)=gzmy(:,:)
252!$OMP END WORKSHARE
253
254!$OMP DO
255do j=2,ny-1
256   do i=2,nx-1
257      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
258          + (betamy(i,j)+betamy(i,j+1)))*0.25
259   end do
260end do
261!$OMP END DO
262!$OMP END PARALLEL
263
264return
265end subroutine dragging
266
267end module dragging_param_beta_nolin
268
Note: See TracBrowser for help on using the repository browser.