source: trunk/SOURCES/dragging_param_beta_sedim_mod.f90 @ 332

Last change on this file since 332 was 332, checked in by aquiquet, 3 years ago

Dragging: distinction between drag computation and mask update / flottab only updates masks

File size: 9.1 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_sedim
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) :: h_sedimmx   ! sediment thickness (or rebounded bedrock)
49real,dimension(nx,ny) :: h_sedimmy   ! sediment thickness (or rebounded bedrock)
50real :: seuil_sedim   ! drag reduced beyond this threshold
51real :: coef_sedim    ! drag reduction factor
52real :: neffmin = 1.e5 ! as in neffect, but this should be the exact same variable
53
54real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite avec spinup cat
55real, dimension(nx,ny) :: Vcol_y           !< uniquement pour compatibilite avec spinup cat
56real, dimension(nx,ny) :: Vsl_x            !< uniquement pour compatibilite avec spinup cat
57real, dimension(nx,ny) :: Vsl_y            !< uniquement pour compatibilite avec spinup cat
58logical :: corr_def = .false.              !< for deformation correction, pour compatibilite beta
59
60contains
61!-------------------------------------------------------------------------------
62
63!> SUBROUTINE: init_dragging
64!! Cette routine fait l'initialisation du dragging.
65!>
66subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
67
68implicit none
69
70real,dimension(nx,ny) :: h_sedim   ! sediment thickness (or rebounded bedrock)
71character(len=150) :: file_sedim   ! file for sediment thickness for HN or rebounded bsoc for Antar
72character(len=150) :: file_ncdf    ! fichier netcdf issue des fichiers .dat
73
74namelist/drag_param_beta_sedim/beta_slope,beta_expo,betamax,betamin,coef_ile,file_sedim,seuil_sedim,coef_sedim
75
76if (itracebug.eq.1)  call tracebug(' dragging param avec sedim')
77
78
79! lecture des parametres du run                      block drag neff
80!--------------------------------------------------------------------
81
82rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
83read(num_param,drag_param_beta_sedim)
84
85write(num_rep_42,'(A)')'!___________________________________________________________' 
86write(num_rep_42,'(A)') '&drag_param_beta_sedim          ! nom du bloc dragging param beta'
87write(num_rep_42,*)
88write(num_rep_42,*) 'beta_slope            = ', beta_slope
89write(num_rep_42,*) 'beta_expo             = ', beta_expo
90write(num_rep_42,*) 'betamax               = ', betamax
91write(num_rep_42,*) 'betamin               = ', betamin
92write(num_rep_42,'(A,A,A)') 'file_sedim     = "',trim(file_sedim),'"'
93write(num_rep_42,*) 'seuil_sedim    = ', seuil_sedim
94write(num_rep_42,*) 'coef_sedim     = ', coef_sedim
95write(num_rep_42,*)'/'                           
96write(num_rep_42,'(A)') '! Slope & expo of beta = - slope x Neff ** expo'
97write(num_rep_42,'(A)') '! Where h_sedim > seuil_sedim, beta*coef_sedim'
98
99!-------------------------------------------------------------------
100! masque stream
101      mstream_mx(:,:)=1
102      mstream_my(:,:)=1
103
104
105! coefficient permettant de modifier le basal drag.
106      drag_mx(:,:)=1.
107      drag_my(:,:)=1.
108     
109      betamax_2d(:,:) = betamax
110
111! modification of basal drag depends on where we have sediments
112! for Antarctica, we assume that what is present-day below sea level has sediment
113      file_sedim=trim(dirnameinp)//trim(file_sedim)
114      call lect_input(1,'h_sedim',1,h_sedim(:,:),file_sedim,file_ncdf)
115 
116      h_sedimmx(1,:)=h_sedim(1,:)
117      h_sedimmy(:,1)=h_sedim(:,1)
118      do i=2,nx
119         h_sedimmx(i,:)=(h_sedim(i-1,:)+h_sedim(i,:))/2.
120      enddo
121      do j=2,ny
122         h_sedimmy(:,j)=(h_sedim(:,j-1)+h_sedim(:,j))/2.
123      enddo
124 
125      !afq -- toblim = 0. !0.7e5 ! afq -- pour les iles, mais est-ce vraiment utile?
126
127      return
128      end subroutine init_dragging
129!________________________________________________________________________________
130
131!> SUBROUTINE: dragging
132!! Defini la localisation des streams et le frottement basal
133!>
134!-------------------------------------------------------------------------
135subroutine dragging   ! defini le frottement basal
136!$ USE OMP_LIB
137
138!         les iles n'ont pas de condition neff mais ont la condition toblim
139!         (ce bloc etait avant dans flottab)
140
141!$OMP PARALLEL
142
143!    -- !$OMP DO
144!    -- do j=2,ny
145!    --    do i=2,nx
146!    --       ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
147!    --       ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
148!    --    end do
149!    -- end do
150!    -- !$OMP END DO
151
152
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)=neffmin
168neff(nx,:)=neffmin
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
186where (h_sedimmx(:,:).le.seuil_sedim)
187  betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo)
188elsewhere
189  betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo) * coef_sedim
190endwhere
191
192where (h_sedimmy(:,:).le.seuil_sedim)
193  betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo)
194elsewhere
195  betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo) * coef_sedim
196endwhere
197!where (h_sedimmx(:,:).gt.seuil_sedim) betamx(:,:) = beta_slope*(neffmin**beta_expo)/neffmin
198!where (h_sedimmy(:,:).gt.seuil_sedim) betamy(:,:) = beta_slope*(neffmin**beta_expo)/neffmin
199
200where (ilemx(:,:)) betamx(:,:) = betamx(:,:) * coef_ile
201where (ilemy(:,:)) betamy(:,:) = betamy(:,:) * coef_ile
202
203betamx(:,:)=max(betamx(:,:),betamin)
204betamy(:,:)=max(betamy(:,:),betamin)
205betamx(:,:)=min(betamx(:,:),betamax)
206betamy(:,:)=min(betamy(:,:),betamax)
207
208where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
209where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
210
211!$OMP END WORKSHARE
212
213
214! calcul de gzmx
215
216! points flottants
217!$OMP DO
218do j=2,ny
219   do i=2,nx
220      if (flot(i,j).and.(flot(i-1,j))) then
221         betamx(i,j)=0.
222      end if
223      if (flot(i,j).and.(flot(i,j-1))) then
224         betamy(i,j)=0.
225      end if
226   end do
227end do
228!$OMP END DO
229
230!$OMP DO
231do j=2,ny-1
232   do i=2,nx-1
233      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
234          + (betamy(i,j)+betamy(i,j+1)))*0.25
235   end do
236end do
237!$OMP END DO
238!$OMP END PARALLEL
239
240return
241end subroutine dragging
242!________________________________________________________________________________
243
244!> SUBROUTINE: mstream_dragging
245!! Defini la localisation des streams
246!>
247!-------------------------------------------------------------------------
248subroutine mstream_dragging   ! defini la localisation des streams
249
250!$OMP PARALLEL
251
252!$OMP WORKSHARE
253fleuvemx(:,:)=.false.
254fleuvemy(:,:)=.false.
255cote(:,:)=.false.
256gzmx(:,:)=.true.
257gzmy(:,:)=.true.
258flgzmx(:,:)=.false.
259flgzmy(:,:)=.false.
260!$OMP END WORKSHARE
261
262! detection des cotes
263!$OMP DO
264do  j=2,ny-1
265   do i=2,nx-1
266      if ((.not.flot(i,j)).and. &
267      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
268         cote(i,j)=.true.
269      endif
270   end do
271end do
272!$OMP END DO
273
274
275! calcul de gzmx
276
277! points flottants : flgzmx mais pas gzmx
278!$OMP DO
279do j=2,ny
280   do i=2,nx
281      if (flot(i,j).and.(flot(i-1,j))) then
282         flotmx(i,j)=.true.
283         flgzmx(i,j)=.true.
284         gzmx(i,j)=.false.
285         betamx(i,j)=0.
286      end if
287      if (flot(i,j).and.(flot(i,j-1))) then
288         flotmy(i,j)=.true.
289         flgzmy(i,j)=.true.
290         gzmy(i,j)=.false.
291         betamy(i,j)=0.
292      end if
293   end do
294end do
295!$OMP END DO
296
297!--------- autres criteres
298
299! points poses
300!$OMP WORKSHARE
301gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
302gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
303
304flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
305flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
306
307fleuvemx(:,:)=gzmx(:,:)
308fleuvemy(:,:)=gzmy(:,:)
309!$OMP END WORKSHARE
310
311!$OMP END PARALLEL
312
313return
314end subroutine mstream_dragging
315
316end module dragging_param_beta_sedim
317
Note: See TracBrowser for help on using the repository browser.