source: trunk/SOURCES/bilan_eau_mod.f90 @ 334

Last change on this file since 334 was 255, checked in by dumas, 5 years ago

Addition of output variables for ISMIP

File size: 5.7 KB
Line 
1!> \file bilan_eau.f90
2!! Clacul du bilan d'eau sur la calotte
3!<
4     
5
6!> MODULE: bilan_eau
7!! Calcule le bilan d'eau sur la grille GRISLI
8!! \author C. Dumas
9!! \date 20/02/2017
10!! @note Permet de vérifier que le bilan d'eau est ok et qu'il n'y a pas de fuite.
11!>
12
13! bilan d'eau :
14! ce qui tombe : Acc
15! ce qui sort : Bmelt, Abl, calving, ablbord
16! ce qui est stocké : masse de glace
17! le total de tous ces termes doit être constant (sur la grille, pas localement)
18
19module bilan_eau_mod
20
21
22use module3D_phy
23implicit none
24
25real :: sum_H
26
27
28real,dimension(nx,ny) :: tot_water    !< bilan d'eau
29real,dimension(nx,ny) :: H_beau_old, diff_H_2D, diff_H_water_bilan
30real,dimension(nx,ny) :: Bm_dtt       !< mass balance on ice points accumulated during dtt
31real,dimension(nx,ny) :: Bmelt_dtt    !< basal melting on ice points accumulated during dtt
32real,dimension(nx,ny) :: calv_dtt     !< calving sur dtt (pour calcul bilan d'eau)
33real,dimension(nx,ny) :: archimtab    ! point pose si > 0
34real,dimension(nx,ny) :: grline_dtt   ! grounding line flux during dtt
35
36real :: sum_H_old
37real :: diff_H
38real :: alpha_flot                      ! ro/row
39real :: dtt_flux   ! 1 an si dtt<1an sinon dtt_loc=dtt
40
41real,dimension(nx,ny) :: bm_dt,bmelt_dt
42
43
44contains
45
46subroutine init_bilan_eau
47 ! initialisation des variables
48  diff_H=0.
49  sum_H_old = sum(H(2:nx-1,2:ny-1),mask=ice(2:nx-1,2:ny-1)==1)
50  H_beau_old(:,:)=H(:,:)
51  tot_water(:,:)=0.
52  bm_dt(:,:)=0.
53  bmelt_dt(:,:)=0.
54  alpha_flot=ro/row
55  archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel_2d(:,:)
56  gr_line(:,:)=0
57  do j=1,ny
58    do i=1,nx
59      !afq    if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel_2d(i,j))) then ! grounded with ice
60      if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.)) then ! grounded with ice
61        if (archimtab(i-1,j).LT.0..and.Uxbar(i,j).LT.0..and..not.flot_marais(i-1,j)) gr_line(i,j)=1
62        if (archimtab(i+1,j).LT.0..and.Uxbar(i+1,j).GT.0..and..not.flot_marais(i+1,j)) gr_line(i,j)=1
63        if (archimtab(i,j-1).LT.0..and.Uybar(i,j).LT.0..and..not.flot_marais(i,j-1)) gr_line(i,j)=1
64        if (archimtab(i,j+1).LT.0..and.Uybar(i,j+1).GT.0..and..not.flot_marais(i,j+1)) gr_line(i,j)=1
65      endif
66    enddo
67  enddo
68  if (dtt.ge.1.) then
69    dtt_flux=dtt
70  else
71    dtt_flux=1.
72  endif
73
74end subroutine init_bilan_eau
75
76
77
78
79subroutine bilan_eau
80
81tot_water(:,:)=0.
82
83bm_dt(:,:)=0.
84bmelt_dt(:,:)=0.
85
86bm_dt(2:nx-1,2:ny-1)=bm(2:nx-1,2:ny-1)*dt
87bmelt_dt(2:nx-1,2:ny-1)=bmelt(2:nx-1,2:ny-1)*dt
88!~ where (Bm(:,:).lt.0.)
89!~      bm_dt(:,:)=max(Bm(:,:)*dt,-H(:,:))
90!~ elsewhere
91!~      bm_dt(:,:)=bm(:,:)*dt
92!~ endwhere
93
94!~ where (Bmelt(:,:).lt.0.)
95!~      bmelt_dt(:,:)=max(Bmelt(:,:)*dt,-H(:,:))
96!~ elsewhere
97!~      bmelt_dt(:,:)=bmelt(:,:)*dt
98!~ endwhere
99
100
101! on fait la somme sur la grille en excluant les bords :
102
103sum_H = sum(H(2:nx-1,2:ny-1),mask=ice(2:nx-1,2:ny-1)==1)
104
105diff_H = diff_H + (sum_H - sum_H_old) ! on calcul la variation de volume a chaque dt
106
107diff_H_2D(2:nx-1,2:ny-1)=H(2:nx-1,2:ny-1)-H_beau_old(2:nx-1,2:ny-1)
108
109where (ice(2:nx-1,2:ny-1).eq.1)
110        Bm_dtt(2:nx-1,2:ny-1) = Bm_dtt(2:nx-1,2:ny-1) + Bm_dt(2:nx-1,2:ny-1) !* dt ! somme Bm sur dt
111        bmelt_dtt(2:nx-1,2:ny-1) = bmelt_dtt(2:nx-1,2:ny-1) + bmelt_dt(2:nx-1,2:ny-1) ! * dt ! somme bmelt sur dt
112endwhere
113
114archimtab(:,:) = Bsoc(:,:)+H(:,:)*alpha_flot - sealevel_2d(:,:)
115gr_line(:,:)=0
116do j=1,ny
117  do i=1,nx
118     !afq    if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.).and.(Bsoc(i,j).LE.sealevel_2d(i,j))) then ! grounded with ice
119     if ((H(i,j).gt.0.).and.(archimtab(i,j).GE.0.)) then ! grounded with ice
120      if (archimtab(i-1,j).LT.0..and.Uxbar(i,j).LT.0..and..not.flot_marais(i-1,j)) gr_line(i,j)=1
121      if (archimtab(i+1,j).LT.0..and.Uxbar(i+1,j).GT.0..and..not.flot_marais(i+1,j)) gr_line(i,j)=1
122      if (archimtab(i,j-1).LT.0..and.Uybar(i,j).LT.0..and..not.flot_marais(i,j-1)) gr_line(i,j)=1
123      if (archimtab(i,j+1).LT.0..and.Uybar(i,j+1).GT.0..and..not.flot_marais(i,j+1)) gr_line(i,j)=1
124    endif
125  enddo
126enddo
127
128where (gr_line(:,:).eq.1)
129!~   grline_dtt(:,:)= (((uxbar(:,:)+eoshift(uxbar(:,:),shift=1,boundary=0.,dim=1))**2+ &
130!~                     (uybar(:,:)+eoshift(uybar(:,:),shift=1,boundary=0.,dim=2))**2)**0.5)*0.5 &
131!~                     *H(:,:) + grline_dtt(:,:)
132   grline_dtt(:,:)= - sqrt(                                                                        &
133                           ( (uxbar(:,:)+eoshift(uxbar(:,:),shift=1,boundary=0.,dim=1))*dy/2. )**2 +  &
134                           ( (uybar(:,:)+eoshift(uybar(:,:),shift=1,boundary=0.,dim=2))*dx/2. )**2 )  &
135                           * H(:,:) * dt + grline_dtt(:,:)
136endwhere
137
138if ((isynchro.eq.1.and.(mod(abs(TIME),1.).lt.dtmin)).or.(nt.eq.1)) then
139! on raisonne en bilan annuel pour simplifier :
140!~      where (ice(2:nx-1,2:ny-1).eq.1)
141!~              tot_water(2:nx-1,2:ny-1) = (Bm_dtt(2:nx-1,2:ny-1) - Bmelt_dtt(2:nx-1,2:ny-1) + calv_dtt(2:nx-1,2:ny-1) - ablbord_dtt(2:nx-1,2:ny-1)) / dtt
142!~              elsewhere
143!~              tot_water(2:nx-1,2:ny-1) = (calv_dtt(2:nx-1,2:ny-1) - ablbord_dtt(2:nx-1,2:ny-1)) / dtt
144!~      endwhere
145!cdc pas besoin du test sur ice ici, il a ete fait avant (et le masque ice varie a chaque dt)   
146  tot_water(2:nx-1,2:ny-1) = (Bm_dtt(2:nx-1,2:ny-1) - Bmelt_dtt(2:nx-1,2:ny-1) + calv_dtt(2:nx-1,2:ny-1) - ablbord_dtt(2:nx-1,2:ny-1)) / dtt_flux
147   
148! bilan d'eau sur la grille :
149  water_bilan=sum(tot_water(:,:))
150  diff_H = diff_H/dtt_flux
151
152!999 format(f0.2,1x,e15.8,1x,i10,8(1x,e15.8))
153!       write(6,999),time,sum_H,count(ice(:,:)==1),diff_H,water_bilan,sum(calv_dtt(:,:))/dtt,sum(ablbord_dtt(:,:))/dtt,sum(bmelt_dtt(:,:),mask=ice(:,:)==1)/dtt,sum(bm(:,:),mask=ice(:,:)==1),sum(Bm_dtt(:,:))/dtt,sum(bmelt_dtt(:,:))/dtt
154  diff_H_water_bilan(2:nx-1,2:ny-1)=tot_water(2:nx-1,2:ny-1)-diff_H_2D(2:nx-1,2:ny-1)
155
156endif
157sum_H_old=sum_H
158
159H_beau_old(:,:)=H(:,:)
160
161end subroutine bilan_eau
162
163end module bilan_eau_mod
Note: See TracBrowser for help on using the repository browser.