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

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

Grisli-iloveclim branch merged to trunk at revision 196

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