source: branches/iLoveclim/SOURCES/dragging_param_beta_nolin_mod.f90 @ 187

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

Grisli-iloveclim branch merged to trunk at revision 185

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