source: trunk/SOURCES/Draggings_modules/dragging_prescr_beta_mod.f90

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

Addition of a non-linear Coulomb friction law + iterations on velocity computation for convergence / for m=1 gives results very similar to the standard linear law in dragging_param_beta

File size: 18.7 KB
Line 
1!> \file dragging_prescr_beta_mod.f90
2!! Module pour calculer le beta a partir d'une carte sur les noeuds centres.
3!<
4
5!> \namespace dragging_prescr_beta
6!! Calcule le beta a partir de vitesses de bilan
7!! @note Il faut partir d'un cptr
8!! \author ...
9!! \date ...
10!! @note Used module
11!! @note   - use module3D_phy
12!<
13
14module dragging_prescr_beta
15
16
17  use module3d_phy
18  use interface_input
19  use fake_beta_iter_vitbil_mod
20 
21  implicit none
22  real, dimension(nx,ny) :: Vcol_x           !< uniquement pour compatibilite
23  real, dimension(nx,ny) :: Vcol_y           !<
24  real, dimension(nx,ny) :: Vsl_x            !<
25  real, dimension(nx,ny) :: Vsl_y            !<
26
27
28  real, dimension(nx,ny) :: drag_centre      !< beta on major node (average)
29  real, dimension(nx,ny) :: betamoy_x        !< pour faire les moyennes glissantes
30  real, dimension(nx,ny) :: betamoy_y        !<
31  integer, dimension(nx,ny) :: nbvoisx       !<
32  integer, dimension(nx,ny) :: nbvoisy       !<
33
34
35  real :: beta_limgz                         !< when beta gt beta_limgz -> not gzmx
36  real :: beta_min                           !< minimum value of beta
37  real :: beta_mult                          !< coefficient of beta field
38  integer  :: ill,jll,nmoy
39
40  real :: maxi                               !< calcul correction dS a la grounding line
41  real :: mini
42
43  logical :: corr_def = .False.               !< for deformation correction (compatibility)
44
45contains
46
47  !-----------------------------------------------------------------------------------------
48  !> SUBROUTINE: init_dragging
49  !! Cette routine fait l'initialisation du dragging.
50  !<
51  subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
52                                ! nouvelle version : lit les fichiers nc
53    implicit none
54    character(len=100) :: beta_c_file               !  beta on centered grid
55    character(len=100) :: betamx_file               !  beta on staggered grid mx
56    character(len=100) :: betamy_file               !  beta on staggered grid mx
57
58    namelist/beta_prescr/beta_c_file,beta_limgz,beta_min,beta_mult
59
60    if (itracebug.eq.1)  call tracebug(' Prescr_beta       subroutine init_dragging')
61
62    inv_beta=0
63
64    ! lecture des parametres du run                   
65    !--------------------------------------------------------------------
66
67    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
68428 format(A)
69    read(num_param,beta_prescr)
70
71    write(num_rep_42,428) '!___________________________________________________________' 
72    write(num_rep_42,428) '!  read beta on centered grid'
73    write(num_rep_42,beta_prescr)                   
74    write(num_rep_42,428) '!  beta_file : nom des fichier qui contiennent les beta_c'
75    write(num_rep_42,428) '!  above beta_limgz, gzmx is false'
76    write(num_rep_42,428) '!  if grounded, beta_min minimum value '
77    write(num_rep_42,428) '!  beta_mult : coefficient just after reading' 
78    write(num_rep_42,*)
79    write(num_rep_42,428) '!___________________________________________________________' 
80   
81
82    beta_c_file = trim(dirnameinp)//trim(beta_c_file)
83    betamx_file = trim(dirnameinp)//trim(betamx_file)
84    betamy_file = trim(dirnameinp)//trim(betamy_file)
85    betamax = beta_limgz
86    betamax_2d(:,:) = betamax
87
88!  read the beta value on centered and staggered grids
89!    call lect_input(1,'beta_c',1,drag_centre,beta_c_file,trim(dirnameinp)//trim(runname)//'.nc')
90!    call lect_input(1,'betamx',1,drag_mx,betamx_file,trim(dirnameinp)//trim(runname)//'.nc')
91!    call lect_input(1,'betamy',1,drag_my,betamy_file,trim(dirnameinp)//trim(runname)//'.nc')
92
93
94    call lect_input(1,'beta_c',1,drag_centre,beta_c_file,'')
95
96! multiplication par le coefficient beta_mult
97
98    drag_centre(:,:) = drag_centre(:,:)*beta_mult
99
100    where (mk_init(:,:).eq.-2)                                             ! iles a probleme
101       drag_centre(:,:) = 2.* abs(beta_limgz)
102    end where
103
104
105    where ((H(:,:).lt.0.1).and.(Bsoc(:,:).gt.0))                            ! zones actuellement non couvertes
106       drag_centre(:,:) = 1000.                                             ! valeur relativement faible qui evite
107    end where                                                               ! les freinages
108
109
110
111!    call lect_input(1,'betamx',1,drag_mx,betamx_file,'')
112!    call lect_input(1,'betamy',1,drag_my,betamy_file,'')
113
114
115     call  beta_c2stag(nx,ny,drag_mx,drag_my,drag_centre)         ! redistribute on staggered grid
116
117     if  (beta_limgz.gt.0.) then
118        beta_centre(:,:) =  drag_centre(:,:)
119        betamx(:,:)      =  drag_mx(:,:)
120        betamy(:,:)      =  drag_my(:,:)
121
122     else if (beta_limgz.lt.0.) then                              ! attention sans doute pour plastic
123
124        beta_centre(:,:) = - beta_limgz
125        betamx(:,:)      = - beta_limgz
126        betamy(:,:)      = - beta_limgz
127        drag_centre(:,:) = - beta_limgz 
128        drag_mx(:,:)     = - beta_limgz
129        drag_my(:,:)     = - beta_limgz 
130        beta_limgz = 1.
131     end if
132
133
134! divers traitements qui peuvent etre appliques a beta : mis sous forme de routine
135
136    ! call beta_plateau_slope   !  give value beta_limgz to places where slope is >1.e-2
137    ! call beta_stag2c          !  average the staggered values on central ones
138    ! call beta_plateau_S       !  give value beta_limgz to places where S > 2500 m
139
140   call beta_min_value(beta_min)      !  valeur mini que peut avoir beta (en Pa an m-1)
141                                ! on peut envisager une valeur ~ 10
142                                ! rappel : beta doit être positif.
143
144    ! call  beta_c2stag(nx,ny,drag_mx,drag_my,drag_centre)         ! redistribute on staggered grid
145
146
147    !  mstream attribution
148    !________________________
149    ! mstream -> regions where streaming is allowed
150    ! drag_mx -> value of beta if an ice stream is allowed
151    ! betamx  -> value of dragging (may depend on other conditions)
152
153    ! to ensure having some "Dirichlet" points, there must be at least
154    ! a mstream=0 in every island
155
156
157! un tour de lissage
158!    call beta_stag2c(nx,ny,drag_mx,drag_my,drag_centre)                  ! average  staggered
159!                                                                         ! on central
160
161!    call beta_c2stag(nx,ny,drag_mx,drag_my,drag_centre)                  ! centered distributed
162!                                                                         ! on staggered
163!
164
165
166    beta_centre(:,:) =  drag_centre(:,:)                                 ! surtout pour sorties
167
168
169
170    where (drag_mx(:,:).ge.beta_limgz)
171       mstream_mx(:,:) = 0
172       betamx(:,:)     = beta_limgz
173       drag_mx(:,:)    = beta_limgz
174    elsewhere
175       mstream_mx(:,:) = 1
176       betamx(:,:)     = drag_mx(:,:)
177    endwhere
178
179    where (drag_my(:,:).ge.beta_limgz)
180       mstream_my(:,:) = 0
181       betamy(:,:)     = beta_limgz
182       drag_my(:,:)    = beta_limgz
183    elsewhere
184       mstream_my(:,:) = 1
185       betamy(:,:)     = drag_my(:,:)
186    endwhere
187
188    if (itracebug.eq.1)  call tracebug(' Prescr_beta  apres lecture')
189    mstream_mx(:,:) = 0
190    mstream_my(:,:) = 0
191
192
193    call gzm_beta_prescr
194!    call write_datfile2(nx,ny,betamx,betamy,'beta-LBq15-08.dat')
195
196!    call make_gaussienne(n_g,sigma_beta,M_gauss)
197
198    niter_nolin = 1
199
200    return
201
202  end subroutine init_dragging
203  !________________________________________________________________________________
204
205  !> SUBROUTINE: dragging
206  !! Defini la localisation des streams et le frottement basal
207  !<
208  !-------------------------------------------------------------------------
209  subroutine dragging   ! defini la localisation des streams et le frottement basal
210
211    implicit none
212
213
214
215    if (itracebug.eq.1)  call tracebug(' beta_prescr         subroutine dragging')
216
217
218!   tentative d'iteration
219
220!    call beta_stag2c(nx,ny,drag_mx,drag_my,drag_centre)                  ! average  staggered
221                                                                          ! on central
222!   call iter_gauss(n_g,nx,ny,M_gauss,-1.,drag_centre,corrstream_centre,hdot)  ! correct central
223!   call beta_c2stag(nx,ny,drag_mx,drag_my,corrstream_centre)             ! centered distributed
224                                                                          ! on staggered
225!    drag_centre(:,:) = corrstream_centre(:,:)
226!    beta_centre(:,:) =  drag_centre(:,:)                                 ! surtout pour sorties
227
228    ! pour l'instant pas d'autre condition que mstream
229    where (mstream_mx(:,:).eq.1)
230       betamx(:,:) = drag_mx(:,:)
231    elsewhere 
232       betamx(:,:)= beta_limgz
233    end where
234
235    where (mstream_my(:,:).eq.1)
236       betamy(:,:) = drag_my(:,:)
237    elsewhere 
238       betamy(:,:)= beta_limgz
239    end where
240
241    betamx(:,:)=drag_mx(:,:)
242    betamy(:,:)=drag_my(:,:)
243
244
245    return
246  end subroutine dragging
247
248  !----------------------------------------------------------------------------------
249  !>SUBROUTINE: gzm_beta_prescr
250  !! Calcul de gzmx
251  !<
252
253  subroutine mstream_dragging
254
255    implicit none
256
257    call  gzm_beta_prescr        ! determine gz et flgz et met a 0 le beta des points flottants
258   
259    return
260  end subroutine mstream_dragging
261
262  subroutine gzm_beta_prescr     ! attribue flgzmx et gzmx (et my)
263   
264    logical ::  test_Tx          ! test if there is a basal melting point in the surrounding
265    logical ::  test_Ty          ! test if there is a basal melting point in the surrounding
266    logical ::  test_beta_x      ! test if there is a low drag point in the surrounding
267    logical ::  test_beta_y      ! test if there is a low drag point in the surrounding
268
269    real    ::  beta_forc_fleuve ! below this value -> ice stream
270
271!    beta_forc_fleuve = 5.e3
272    beta_forc_fleuve = beta_limgz
273
274    ! calcul de gzmx
275
276
277    gzmx(:,:)=.true.
278    gzmy(:,:)=.true.
279    flgzmx(:,:)=.false.
280    flgzmy(:,:)=.false.
281
282
283    ! points flottants : flgzmx mais pas gzmx
284    do j=2,ny
285       do i=2,nx
286          if (flot(i,j).and.(flot(i-1,j))) then
287             flotmx(i,j)=.true.
288             flgzmx(i,j)=.true.
289             gzmx(i,j)=.false.
290             betamx(i,j)=0.
291          end if
292          if (flot(i,j).and.(flot(i,j-1))) then
293             flotmy(i,j)=.true.
294             flgzmy(i,j)=.true.
295             gzmy(i,j)=.false.
296             betamy(i,j)=0.
297          end if
298       end do
299    end do
300
301    if (itracebug.eq.1)  call tracebug(' apres criteres flot')
302
303    if (itracebug.eq.1) write(num_tracebug,*) 'gzmx', gzmx(127,127) 
304
305    !--------- autres criteres
306
307    ! points poses
308    gzmx(:,:)=gzmx(:,:).and.(betamx(:,:).lt.beta_limgz)   !  Pas de calcul pour les points
309    gzmy(:,:)=gzmy(:,:).and.(betamy(:,:).lt.beta_limgz)   !  au fort beta
310
311
312
313    ! critere (lache) sur la temperature
314    do j = 2, ny-1
315       do i =2, nx-1 
316
317!   test s'il y a un point tempere dans les environs
318          test_Tx = any( ibase( i-1:i   , j-1:j+1 ).gt.1)
319          test_Ty = any( ibase( i-1:i+1 , j-1:j   ).gt.1)
320
321!   test s'il y a un point low drag dans les environs
322          test_beta_x = any( drag_centre( i-1:i   , j-1:j+1 ) .lt. beta_forc_fleuve )
323          test_beta_y = any( drag_centre( i-1:i+1 , j-1:j   ) .lt. beta_forc_fleuve )
324
325!   critere pour gzmx
326
327! Attention : les deux lignes suivants sont en commentaire pour voir l'effet
328
329!          gzmx(i,j) =  gzmx(i,j) .and. (test_Tx.or.test_beta_x)
330!          gzmy(i,j) =  gzmy(i,j) .and. (test_Ty.or.test_beta_y)
331
332       end do
333    end do
334
335
336    if (itracebug.eq.1)  call tracebug(' apres autres criteres ')
337    if (itracebug.eq.1) write(num_tracebug,*) 'gzmx', gzmx(127,127) 
338
339    flgzmx(:,:) = flgzmx(:,:) .or. gzmx(:,:)
340    flgzmy(:,:) = flgzmy(:,:) .or. gzmy(:,:)
341
342    where ((.not.flotmx(:,:)).and.(.not.gzmx(:,:)))
343       betamx(:,:) = beta_limgz
344    end where
345
346    where ((.not.flotmy(:,:)).and.(.not.gzmy(:,:)))
347       betamy(:,:) = beta_limgz
348    end where
349
350        fleuvemx(:,:)=gzmx(:,:)
351        fleuvemy(:,:)=gzmy(:,:)
352   
353  end subroutine gzm_beta_prescr
354
355
356  !----------------------------------------------------------------------------------
357  ! Routines qui permettent des manip simples sur beta.
358
359
360!______________________________________________________________________________________
361!>SUBROUTINE:  beta_plateau_slope
362!! give value beta_limgz to places where slope is >1.e-2
363!<
364
365  subroutine beta_plateau_slope      ! give value beta_limgz to places where slope is >1.e-2
366
367    call slope_surf
368
369    where (abs(sdx(:,:)).gt.1.e-2)
370       drag_mx(:,:) = beta_limgz
371    end where
372
373    where (abs(sdy(:,:)).gt.1.e-2)
374       drag_my(:,:) = beta_limgz
375    end where
376  end subroutine beta_plateau_slope
377!______________________________________________________________________________________
378!> SUBROUTINE beta_plateau_S
379!! give value beta_limgz to places where S > 2500 m
380!<
381  subroutine beta_plateau_S     !  give value beta_limgz to places where S > 2500 m
382
383    where (S(:,:).gt.2500.)     ! moyen simple de nettoyer l'interieur. WAIS pas touche
384       drag_centre(:,:)=beta_limgz
385       drag_mx(:,:)=beta_limgz
386       drag_my(:,:)=beta_limgz
387    end where
388  end subroutine beta_plateau_S
389!______________________________________________________________________________________
390!>SUBROUTINE:   beta_min_value
391!! valeur mini que peut avoir beta (en Pa an m-1)
392!<
393
394  subroutine beta_min_value(val_betamin)
395
396  real :: val_betamin     ! valeur mini que peut avoir beta (en Pa an m-1)
397
398    drag_centre(:,:) = max(drag_centre(:,:),val_betamin)
399    drag_mx(:,:) = max(drag_mx(:,:),val_betamin)
400    drag_my(:,:) = max(drag_my(:,:),val_betamin)
401    where(flot(:,:))
402       drag_centre(:,:) = 0.
403    end where
404
405  end subroutine beta_min_value
406!______________________________________________________________________________________
407!>SUBROUTINE:   beta_stag2c
408!! average the staggered values on central ones
409!<
410
411  subroutine beta_stag2c(nxx,nyy,R_mx,R_my,R_centre)      ! average the staggered values on central ones
412
413    implicit none
414    integer                                ::  nxx,nyy                ! dimension des tableaux
415    real, dimension(nxx,nyy), intent(in)   ::  R_mx                   ! valeur du beta sur maille mx
416    real, dimension(nxx,nyy), intent(in)   ::  R_my                   ! valeur du beta sur maille my
417    real, dimension(nxx,nyy), intent(out)  ::  R_centre               ! valeur du beta sur maille centree   
418   
419    R_centre(:,:)=0.
420    do j=2,ny-1
421       do i=2,nx-1
422          R_centre(i,j)= ((R_mx(i,j)+R_mx(i+1,j)) &
423               + (R_my(i,j)+R_my(i,j+1)))*0.25
424       end do
425    end do
426  end subroutine beta_stag2c
427!______________________________________________________________________________________
428!>SUBROUTINE:   beta_c2stag_min 
429!! redistribute central values on staggered grid
430!<
431
432  subroutine beta_c2stag(nxx,nyy,R_mx,R_my,R_centre)              ! redistribute central values on staggered grid
433
434    implicit none
435    integer                                ::  nxx,nyy                ! dimension des tableaux
436    real, dimension(nxx,nyy), intent(out)  ::  R_mx                   ! valeur du beta sur maille mx
437    real, dimension(nxx,nyy), intent(out)  ::  R_my                   ! valeur du beta sur maille my
438    real, dimension(nxx,nyy), intent(in)   ::  R_centre               ! valeur du beta sur maille centree
439
440    do j=2,nyy 
441       do i=2,nxx
442          R_mx(i,j)= (R_centre(i-1,j)+R_centre(i,j))*0.5
443          R_my(i,j)= (R_centre(i,j-1)+R_centre(i,j))*0.5
444       end do
445    end do
446
447  end subroutine beta_c2stag
448!______________________________________________________________________________________
449!>SUBROUTINE:   beta_c2stag_min 
450!! redistribute central values on staggered grid
451!<
452
453  subroutine beta_c2stag_min(nxx,nyy,R_mx,R_my,R_centre)              ! redistribute central values on staggered grid
454
455    implicit none
456    integer                                ::  nxx,nyy                ! dimension des tableaux
457    real, dimension(nxx,nyy), intent(inout)::  R_mx                   ! valeur du beta sur maille mx
458    real, dimension(nxx,nyy), intent(inout)::  R_my                   ! valeur du beta sur maille my
459    real, dimension(nxx,nyy), intent(in)   ::  R_centre               ! valeur du beta sur maille centree
460
461    do j=2,nyy 
462       do i=2,nxx
463          R_mx(i,j)= min(R_mx(i,j),(R_centre(i-1,j)+R_centre(i,j))*0.5)
464          R_my(i,j)= min(R_my(i,j),(R_centre(i,j-1)+R_centre(i,j))*0.5)
465       end do
466    end do
467
468  end subroutine beta_c2stag_min
469!______________________________________________________________________________________
470!>SUBROUTINE:   make_gaussienne
471!! cree une matrice gaussienne normalisee a 1
472!<
473
474  subroutine  make_gaussienne(n,sigma,M)
475
476    implicit none
477    integer,intent(in)                       :: n          ! demi taille du tableau gaussienne
478    real, intent(in)                         :: sigma      ! sigma gaussienne exprime en dx
479    real,dimension(-n:n,-n:n),intent(inout)  :: M          ! tableau gaussienne
480
481
482!   variables locales
483    real                                     :: dist2      ! distance au centre (au carre)
484    real                                     :: som        ! pour calibration
485    real                                     :: sigma2     ! pour gaussienne
486
487    som = 0.
488    sigma2 = sigma * sigma * 2
489    do j = -n, n
490       do i = -n, n
491          dist2 = (i*i+j*j)
492          M(i,j) = exp(-dist2/sigma2)
493          som=som+M(i,j)
494       end do
495    end do
496
497    M(:,:) = M(:,:) / som
498
499!    do j = -n, n
500!       write(6,'(5(f0.4,1x))') (M(i,j),i=-n,n)
501!    end do
502
503  end subroutine make_gaussienne
504
505!______________________________________________________________________________________
506!>SUBROUTINE:  iter_gauss
507!! utilise hdot pour faire les iterations sur beta_centre
508!<
509  subroutine iter_gauss(n,nxx,nyy,M,coef,beta_in,beta_out,Ecart)
510
511    integer,                  intent(in)    :: n          ! taille matrice gaussienne
512    integer,                  intent(in)    :: nxx,nyy    ! taille du domaine
513    real,dimension(-n:n,-n:n),intent(in)    :: M          ! matrice gaussienne
514    real,                     intent(in)    :: coef       ! coefficient multiplicateur
515    real,dimension(nx,ny),    intent(in)    :: beta_in    ! la valeur initiale
516    real,dimension(nx,ny),    intent(out)   :: beta_out   ! la valeur apres iteration
517    real,dimension(nx,ny),    intent(in)    :: Ecart      ! l'ecart a l'obs
518
519! variables locales
520    integer                                 :: ii,jj
521
522    do j= 1 + n, nyy - n
523       do i = 1 + n, nxx - n
524          beta_out(i,j) = beta_in(i,j)
525
526          do jj = -n, n
527             do ii = -n, n 
528                beta_out(i,j) = beta_out(i,j) + coef * M(ii,jj) * ecart(i+ii,j+jj)
529             end do
530          end do
531
532       end do
533    end do
534
535    beta_out(:,:) = max (beta_out(:,:),0.)
536  end subroutine iter_gauss
537
538end module dragging_prescr_beta
539
540!!$! add the nodes of topographical instability (on major nodes)
541! a garder eventuellement pour runs instabilite
542!!$beta_file = trim(dirnameinp)//'sensib_H_Gael-15km.dat'
543!!$call lect_datfile(nx,ny,drag_instab_topo,1,beta_file)   
544!!$
545!!$where (drag_instab_topo(:,:).gt.0.)
546!!$   drag_centre(:,:)=min(drag_centre(:,:),400.)
547!!$endwhere
Note: See TracBrowser for help on using the repository browser.