source: trunk/SOURCES/Antarctique_general_files/dragging_prescr_beta_mod_Ant.f90 @ 4

Last change on this file since 4 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

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