source: trunk/SOURCES/dragging_param_beta_nolin_mod.f90 @ 161

Last change on this file since 161 was 161, checked in by dumas, 7 years ago

Attempt to debug the classical crash in diagnoshelf : determin_tache now called by steps_time_loop | isostasy after diagnoshelf | flgzmx & flgzmy .false. if there is no ice

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