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