source: trunk/SOURCES/Antarctique_general_files/dragging-vit_bil_CISM_gen_mod.f90 @ 23

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

initial import GRISLI trunk

File size: 9.2 KB
Line 
1!> \file dragging-vit_bil_CISM_gen_mod.f90
2!! Module qui defini les zones de stream d'apres les vitesses de bilan
3!<
4
5!> \namespace dragging_vit_bil_CISM_gen
6!! Module qui defini les zones de stream
7!! \author ...
8!! \date ...
9!! @note Defini les zones de stream avec :
10!! @note   - un critere sur la hauteur d'eau
11!! @note   - un critere de fleuves actuels
12!! @note Attention il faut fonctionner en additionnant SIA et Stream pour éviter les pb
13!! de fleuves en pointille
14!! @note Used module
15!! @note   - use module3D_phy
16!! @note   - use param_phy_mod
17!<
18
19
20module dragging_vit_bil_CISM_gen
21
22  ! Defini les zones de stream avec :
23  ! * un critere sur la hauteur d'eau
24  ! * un critere de fleuves actuels
25  ! attention il faut fonctionner en additionnant SIA et Stream pour éviter les pb
26  ! de fleuves en pointille
27
28
29  use module3d_phy
30  use param_phy_mod
31  use interface_input
32  implicit none
33  logical,dimension(nx,ny) :: fleuvemx         !< fleuves sont les tableaux courants (dep. time)   
34  logical,dimension(nx,ny) :: fleuvemy
35  logical,dimension(nx,ny) :: fleuve     
36  logical,dimension(nx,ny) :: cote   
37
38  character(len=100) :: balance_vel_file        !<  vitesses de bilan
39  real, dimension(nx,ny) :: vitbil_max          !<  vitesse de bilan max sur la maille
40  real ::  seuil_vel                            !<  seuil sur les vitesses pour definir le masque stream
41
42  real :: valmax
43  integer :: imax,jmax
44  integer :: i_moins1,i_plus1,j_moins1,j_plus1
45  integer :: lmax=20
46  integer :: idep,jdep,iloc,jloc
47
48
49  real :: tostick      ! pour la glace posee
50  real :: tob_ile      ! pour les iles
51
52contains
53  !-------------------------------------------------------------------------------
54
55
56  !> SUBROUTINE: init_dragging
57  !! Cette routine fait l'initialisation du dragging.
58  !<   
59
60  subroutine init_dragging      ! Cette routine fait l'initialisation du dragging.
61
62    implicit none
63
64    namelist/drag_vit_bil_CISM_gen/hwatstream,cf,betamax,toblim,seuil_vel,balance_vel_file   
65
66    if (itracebug.eq.1)  call tracebug(' Init_dragging avec dragging_vit_bil_CISM_gen')
67
68
69    ! formats pour les ecritures dans 42
70428 format(A)
71
72    ! lecture des parametres du run                   
73    !--------------------------------------------------------------------
74
75    rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
76    read(num_param,drag_vit_bil_CISM_gen)
77
78    write(num_rep_42,428) '!___________________________________________________________' 
79    write(num_rep_42,428) '! module dragging_vit_bil_CISM_gen '
80    write(num_rep_42,drag_vit_bil_CISM_gen)                   
81    write(num_rep_42,428) '! hwatstream (m) :  critere de passage en stream vitesses de bilan'
82    write(num_rep_42,428) '!  si hwater > hwatstream '
83    write(num_rep_42,428) '! cf coefficient de la loi de frottement fonction Neff'
84    write(num_rep_42,428) '! seulement pour les points cotiers'
85    write(num_rep_42,428) '! betamax : (Pa) frottement maxi sous les streams '
86    write(num_rep_42,428) '! toblim : (Pa)  pour les iles '
87    write(num_rep_42,428) '! seuil_vel (m/an) : seuil sur les vitesses pour definir le masque stream '
88    write(num_rep_42,428) '! balance_vel_file : nom du fichier qui contient les vitesse de bilan'
89    write(num_rep_42,*)
90
91    tostick=1.e5   ! valeurs pour les points non flgzmx
92    tob_ile=betamax/2.
93    moteurmax=toblim
94    !moteurmax=0.
95
96    !-------------------------------------------------------------------
97    ! masque stream : lecture des vitesses de bilan
98    ! attention on prend la colonne 3 qui est la vitesse max dans la maille
99
100    balance_vel_file  = trim(dirnameinp)//trim(balance_vel_file)
101
102    call lect_input(3,'vitbil',1,vitbil_max,balance_vel_file,trim(dirnameinp)//trim(runname)//'.nc')
103    !call lect_datfile(nx,ny,vitbil_max,1,balance_vel_file)
104    !call lect_datfile(nx,ny,vitbil_max,3,balance_vel_file)    colonne 3 (valeur max) donne des fleuves trop larges) 
105
106    mstream_mx(:,:)=0
107    mstream_my(:,:)=0
108
109    do j=1,ny
110       do i=1,nx
111          if (vitbil_max(i,j).gt.seuil_vel) then
112             mstream_mx(i:i+1,j) = 1
113             mstream_my(i,j:j+1) = 1
114             mstream(i,j) = 1                             ! allowed ice streams
115          end if
116       end do
117    end do
118
119    debug_3D(:,:,25) = mstream(:,:)
120
121
122    ! coefficient permettant de modifier le basal drag.
123    drag_mx(:,:)=1.
124    drag_my(:,:)=1.
125
126
127    return
128  end subroutine init_dragging
129  !________________________________________________________________________________
130
131  !> SUBROUTINE: dragging
132  !! defini la localisation des streams et le frottement basal
133  !<
134
135  !-------------------------------------------------------------------------
136  subroutine dragging   ! defini la localisation des streams et le frottement basal
137
138
139    !         les iles n'ont pas de condition neff mais ont la condition toblim
140    !         (ce bloc etait avant dans flottab)
141
142
143    do j=2,ny
144       do i=2,nx
145          ilemx(i,j)=ilemx(i,j).and.(abs(rog*Hmx(i,j)*sdx(i,j)).lt.toblim)
146          ilemy(i,j)=ilemy(i,j).and.(abs(rog*Hmy(i,j)*sdy(i,j)).lt.toblim)
147       end do
148    end do
149
150    fleuvemx(:,:)=.false.
151    fleuvemy(:,:)=.false.
152    fleuve(:,:)=.false.
153    cote(:,:)=.false.
154    cotemx(:,:)=.false.
155    cotemy(:,:)=.false.
156
157    ! detection des cotes
158    do  j=2,ny-1 
159       do i=2,nx-1 
160          if ((.not.flot(i,j)).and. & 
161               ((flot(i+1,j)).or.(flot(i,j+1)).or.(flot(i-1,j)).or.(flot(i,j-1)))) then
162             cote(i,j)=.true.
163             cotemx(i:i+1,j) = .true.
164             cotemy(i,j:j+1) = .true.
165          endif
166       end do
167    end do
168
169    ! le calcul des fleuves se fait sur les mailles majeures
170    ! normalement, les points cotiers sont deja tagges gzmx dans la routine flottab
171
172    do j=2,ny-1
173       do i=2,nx-1
174
175          if ((hwater(i,j).gt.hwatstream).and.(mstream(i,j).eq.1).and.(.not.flot(i,j))) then
176             fleuve(i,j)=.true.
177             fleuvemx(i:i+1,j)=.true.
178             fleuvemy(i,j:j+1)=.true.
179             gzmx(i:i+1,j)=.true.
180             gzmy(i,j:j+1)=.true.
181          endif
182       end do
183    end do
184
185
186    ! calcul du frottement basal (ce bloc etait avant dans neffect)
187
188    drag_mx : do j=2,ny-1
189       do i=2,nx-1
190
191          if ((.not.ilemx(i,j)).and.(fleuvemx(i,j))) gzmx(i,j)=.true. 
192
193          if (cotemx(i,j)) then                        ! point cotier peut etre diminue si pression d'eau
194             betamx(i,j)=cf*neffmx(i,j) 
195             betamx(i,j)=min(betamx(i,j),betamax)
196          else if ((gzmx(i,j)).and.(.not.cotemx(i,j))) then  ! stream -> betamax
197             betamx(i,j)=betamax
198             betamx(i,j)=min(betamx(i,j),betamax)
199             betamx(i,j)=max(betamx(i,j),10.)
200
201             !         if (cf*neffmx(i,j).gt.1500.) then
202             !            gzmx(i,j)=.false.
203             !            betamx(i,j)=tostick
204             !         endif
205
206          else if (ilemx(i,j)) then                     ! ile : tob_ile=betamax/2 + peut etre diminue si pression d'eau
207             betamx(i,j)=cf*neffmx(i,j) 
208             betamx(i,j)=min(betamx(i,j),tob_ile)       
209          else if (flotmx(i,j)) then                    ! flottant ou ile
210             betamx(i,j)=0.
211          else                                         ! grounded, SIA
212             betamx(i,j)=tostick                       ! frottement glace posee (1 bar)
213          endif
214
215       end do
216    end do drag_mx
217
218    drag_my : do j=2,ny-1
219       do i=2,nx-1
220
221          if ((.not.ilemy(i,j)).and.(fleuvemy(i,j))) gzmy(i,j)=.true. 
222
223          if (cotemy(i,j)) then                        ! point cotier peut etre diminue si pression d'eau
224             betamy(i,j)=cf*neffmy(i,j) 
225             betamy(i,j)=min(betamy(i,j),betamax)
226          else if ((gzmy(i,j)).and.(.not.cotemy(i,j))) then  ! stream -> betamax
227             betamy(i,j)=betamax
228             betamy(i,j)=min(betamy(i,j),betamax)
229             betamy(i,j)=max(betamy(i,j),10.)
230
231             !         if (cf*neffmy(i,j).gt.1500.) then
232             !            gzmy(i,j)=.false.
233             !            betamy(i,j)=tostick
234             !         endif
235
236          else if (ilemy(i,j)) then                     ! ile : tob_ile=betamax/2 + peut etre diminue si pression d'eau
237             betamy(i,j)=cf*neffmy(i,j) 
238             betamy(i,j)=min(betamy(i,j),tob_ile)       
239          else if (flotmy(i,j)) then                    ! flottant ou ile
240             betamy(i,j)=0.
241          else                                         ! grounded, SIA
242             betamy(i,j)=tostick                       ! frottement glace posee (1 bar)
243          endif
244
245       end do
246    end do drag_my
247
248
249    where (fleuve(:,:))
250       debug_3D(:,:,1)=1
251    elsewhere (flot(:,:))
252       debug_3D(:,:,1)=2
253    elsewhere
254       debug_3D(:,:,1)=0
255    endwhere
256
257    where (cote(:,:))
258       debug_3D(:,:,2)=1
259    elsewhere
260       debug_3D(:,:,2)=0
261    endwhere
262
263    where (fleuvemx(:,:))
264       debug_3D(:,:,3)=1
265    elsewhere
266       debug_3D(:,:,3)=0
267    endwhere
268
269    where (flotmx(:,:))
270       debug_3D(:,:,3)=-1
271    endwhere
272
273    !_____________________
274
275
276    where (fleuvemy(:,:))
277       debug_3D(:,:,4)=1
278    elsewhere
279       debug_3D(:,:,4)=0
280    endwhere
281
282    where (flotmy(:,:))
283       debug_3D(:,:,4)=-1
284    end where
285
286    debug_3D(:,:,23)= abs(RO*G*HMX(:,:)*sdx(:,:)*1.e-5) 
287    debug_3D(:,:,24)= abs(RO*G*HMY(:,:)*sdy(:,:)*1.e-5) 
288
289    return
290  end subroutine dragging
291
292end module dragging_vit_bil_CISM_gen
293
Note: See TracBrowser for help on using the repository browser.