source: trunk/SOURCES/dragging_param_beta_mod.f90 @ 193

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

Correcting a missing variable initialisation: toblim for parameterised beta

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
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
101toblim = 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!$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
193betamx(:,:)=max(betamx(:,:),betamin)
194betamy(:,:)=max(betamy(:,:),betamin)
195betamx(:,:)=min(betamx(:,:),betamax)
196betamy(:,:)=min(betamy(:,:),betamax)
197
198where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
199where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
200
201!$OMP END WORKSHARE
202
203
204! calcul de gzmx
205
206! points flottants : flgzmx mais pas gzmx
207!$OMP DO
208do j=2,ny
209   do i=2,nx
210      if (flot(i,j).and.(flot(i-1,j))) then
211         flotmx(i,j)=.true.
212         flgzmx(i,j)=.true.
213         gzmx(i,j)=.false.
214         betamx(i,j)=0.
215      end if
216      if (flot(i,j).and.(flot(i,j-1))) then
217         flotmy(i,j)=.true.
218         flgzmy(i,j)=.true.
219         gzmy(i,j)=.false.
220         betamy(i,j)=0.
221      end if
222   end do
223end do
224!$OMP END DO
225
226!--------- autres criteres
227
228! points poses
229!$OMP WORKSHARE
230gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
231gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
232
233flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
234flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
235
236fleuvemx(:,:)=gzmx(:,:)
237fleuvemy(:,:)=gzmy(:,:)
238!$OMP END WORKSHARE
239
240!$OMP DO
241do j=2,ny-1
242   do i=2,nx-1
243      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
244          + (betamy(i,j)+betamy(i,j+1)))*0.25
245   end do
246end do
247!$OMP END DO
248!$OMP END PARALLEL
249
250return
251end subroutine dragging
252
253end module dragging_param_beta
254
Note: See TracBrowser for help on using the repository browser.