source: branches/iLoveclim/SOURCES/branche-Cat-spinup-dec2011/dragging_prescr_beta_mod.f90 @ 254

Last change on this file since 254 was 77, checked in by dumas, 8 years ago

Merge branche iLOVECLIM sur rev 76

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