source: trunk/SOURCES/dragging_param_beta_sedim_mod.f90 @ 194

Last change on this file since 194 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: 8.6 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! formats pour les ecritures dans 42
80428 format(A)
81
82! lecture des parametres du run                      block drag neff
83!--------------------------------------------------------------------
84
85rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
86read(num_param,drag_param_beta_sedim)
87
88write(num_rep_42,428)'!___________________________________________________________' 
89write(num_rep_42,428) '&drag_param_beta_sedim          ! nom du bloc dragging param beta'
90write(num_rep_42,*)
91write(num_rep_42,*) 'beta_slope            = ', beta_slope
92write(num_rep_42,*) 'beta_expo             = ', beta_expo
93write(num_rep_42,*) 'betamax               = ', betamax
94write(num_rep_42,*) 'betamin               = ', betamin
95write(num_rep_42,*) 'file_sedim     = ', file_sedim
96write(num_rep_42,*) 'seuil_sedim    = ', seuil_sedim
97write(num_rep_42,*) 'coef_sedim     = ', coef_sedim
98write(num_rep_42,*)'/'                           
99write(num_rep_42,428) '! Slope & expo of beta = - slope x Neff ** expo'
100write(num_rep_42,428) '! Where h_sedim > seuil_sedim, beta*coef_sedim'
101
102!-------------------------------------------------------------------
103! masque stream
104      mstream_mx(:,:)=1
105      mstream_my(:,:)=1
106
107
108! coefficient permettant de modifier le basal drag.
109      drag_mx(:,:)=1.
110      drag_my(:,:)=1.
111     
112      betamax_2d(:,:) = betamax
113
114! modification of basal drag depends on where we have sediments
115! for Antarctica, we assume that what is present-day below sea level has sediment
116      file_sedim=trim(dirnameinp)//trim(file_sedim)
117      call lect_input(1,'h_sedim',1,h_sedim(:,:),file_sedim,file_ncdf)
118 
119      h_sedimmx(1,:)=h_sedim(1,:)
120      h_sedimmy(:,1)=h_sedim(:,1)
121      do i=2,nx
122         h_sedimmx(i,:)=(h_sedim(i-1,:)+h_sedim(i,:))/2.
123      enddo
124      do j=2,ny
125         h_sedimmy(:,j)=(h_sedim(:,j-1)+h_sedim(:,j))/2.
126      enddo
127 
128      !afq -- toblim = 0. !0.7e5 ! afq -- pour les iles, mais est-ce vraiment utile?
129
130      return
131      end subroutine init_dragging
132!________________________________________________________________________________
133
134!> SUBROUTINE: dragging
135!! Defini la localisation des streams et le frottement basal
136!>
137!-------------------------------------------------------------------------
138subroutine dragging   ! defini la localisation des streams et le frottement basal
139!$ USE OMP_LIB
140
141!         les iles n'ont pas de condition neff mais ont la condition toblim
142!         (ce bloc etait avant dans flottab)
143
144!$OMP PARALLEL
145
146!    -- !$OMP DO
147!    -- do j=2,ny
148!    --    do i=2,nx
149!    --       ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
150!    --       ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
151!    --    end do
152!    -- end do
153!    -- !$OMP END DO
154
155
156! le gzmx initial a ete calculé dans flottab pour les points a cheval sur la grounding line
157
158
159!$OMP WORKSHARE
160fleuvemx(:,:)=.false.
161fleuvemy(:,:)=.false.
162cote(:,:)=.false.
163gzmx(:,:)=.true.
164gzmy(:,:)=.true.
165flgzmx(:,:)=.false.
166flgzmy(:,:)=.false.
167!$OMP END WORKSHARE
168
169! detection des cotes
170!$OMP DO
171do  j=2,ny-1 
172   do i=2,nx-1 
173      if ((.not.flot(i,j)).and. & 
174      ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
175         cote(i,j)=.true.
176      endif
177   end do
178end do
179!$OMP END DO
180
181! calcul de neff (pression effective sur noeuds majeurs)
182if (sum(neffmx(:,:)).le.0.) neffmx(:,:) =1.e8
183if (sum(neffmy(:,:)).le.0.) neffmy(:,:) =1.e8
184
185!$OMP DO
186do j=1,ny-1
187   do i=1,nx-1
188        neff(i,j)=(neffmx(i,j)+neffmx(i+1,j)+neffmy(i,j)+neffmy(i,j+1))/4
189   enddo
190enddo
191!$OMP END DO
192!aurel, for the borders:
193!$OMP WORKSHARE
194neff(:,ny)=neffmin
195neff(nx,:)=neffmin
196! calcul de hwat (staggered)
197hwatmx(:,:)=0.
198hwatmy(:,:)=0.
199!$OMP END WORKSHARE
200!$OMP DO
201do i=2,nx
202   hwatmx(i,:)=(hwater(i-1,:)+hwater(i,:))/2.
203enddo
204!$OMP END DO
205!$OMP DO
206do j=2,ny
207   hwatmy(:,j)=(hwater(:,j-1)+hwater(:,j))/2.
208enddo
209!$OMP END DO
210
211!$OMP WORKSHARE
212
213where (h_sedimmx(:,:).le.seuil_sedim)
214  betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo)
215elsewhere
216  betamx(:,:)= beta_slope*(neffmx(:,:)**beta_expo) * coef_sedim
217endwhere
218
219where (h_sedimmy(:,:).le.seuil_sedim)
220  betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo)
221elsewhere
222  betamy(:,:)= beta_slope*(neffmy(:,:)**beta_expo) * coef_sedim
223endwhere
224!where (h_sedimmx(:,:).gt.seuil_sedim) betamx(:,:) = beta_slope*(neffmin**beta_expo)/neffmin
225!where (h_sedimmy(:,:).gt.seuil_sedim) betamy(:,:) = beta_slope*(neffmin**beta_expo)/neffmin
226
227where (ilemx(:,:)) betamx(:,:) = betamx(:,:) * coef_ile
228where (ilemy(:,:)) betamy(:,:) = betamy(:,:) * coef_ile
229
230betamx(:,:)=max(betamx(:,:),betamin)
231betamy(:,:)=max(betamy(:,:),betamin)
232betamx(:,:)=min(betamx(:,:),betamax)
233betamy(:,:)=min(betamy(:,:),betamax)
234
235where ( hwatmx(:,:) .lt. 0.5 ) betamx(:,:) = betamax
236where ( hwatmy(:,:) .lt. 0.5 ) betamy(:,:) = betamax
237
238!$OMP END WORKSHARE
239
240
241! calcul de gzmx
242
243! points flottants : flgzmx mais pas gzmx
244!$OMP DO
245do j=2,ny
246   do i=2,nx
247      if (flot(i,j).and.(flot(i-1,j))) then
248         flotmx(i,j)=.true.
249         flgzmx(i,j)=.true.
250         gzmx(i,j)=.false.
251         betamx(i,j)=0.
252      end if
253      if (flot(i,j).and.(flot(i,j-1))) then
254         flotmy(i,j)=.true.
255         flgzmy(i,j)=.true.
256         gzmy(i,j)=.false.
257         betamy(i,j)=0.
258      end if
259   end do
260end do
261!$OMP END DO
262
263!--------- autres criteres
264
265! points poses
266!$OMP WORKSHARE
267gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.betamax)   !  Pas de calcul pour les points
268gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.betamax)   !  au fort beta
269
270flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
271flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
272
273fleuvemx(:,:)=gzmx(:,:)
274fleuvemy(:,:)=gzmy(:,:)
275!$OMP END WORKSHARE
276
277!$OMP DO
278do j=2,ny-1
279   do i=2,nx-1
280      beta_centre(i,j) = ((betamx(i,j)+betamx(i+1,j)) &
281          + (betamy(i,j)+betamy(i,j+1)))*0.25
282   end do
283end do
284!$OMP END DO
285!$OMP END PARALLEL
286
287return
288end subroutine dragging
289
290end module dragging_param_beta_sedim
291
Note: See TracBrowser for help on using the repository browser.