source: trunk/SOURCES/dragging_param_beta_mod.f90 @ 223

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

We apply schoof flux even on island tagged points. Dragging param is corrected so that we do not alter islands as computed in flottab.

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