source: branches/iLoveclim/SOURCES/dragging_param_beta_mod.f90 @ 146

Last change on this file since 146 was 146, checked in by aquiquet, 7 years ago

Grisli-iLoveclim branch: merged to trunk at revision 145

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