source: branches/iLoveclim/SOURCES/bilan_flux_output_mod.f90 @ 146

Last change on this file since 146 was 146, checked in by aquiquet, 7 years ago

Grisli-iLoveclim branch: merged to trunk at revision 145

File size: 4.1 KB
Line 
1!> \file bilan_flux_output_mod.f90
2!! Clacul du bilan des flux : sorties initMIP
3!<
4     
5
6!> MODULE: bilan_flux
7!! Calcule le bilan des flux : sortie initMIP
8!! \author C. Dumas
9!! \date 10/10/2017
10!! @note Calcul des flux output initMIP
11!>
12
13module bilan_flux_mod
14
15
16use module3D_phy
17implicit none
18
19real :: dtoutflux                       ! pas de temps des sorties flux
20real :: dtsum                           ! somme des pas de temps (pour verification)
21!real :: alpha_flot                      ! ro/row
22
23! variables flux :
24real,dimension(nx,ny) :: bm_flux        ! surface smb on ice mask
25real,dimension(nx,ny) :: bmelt_flux     ! bmelt on ice mask
26real,dimension(nx,ny) :: dhdt_flux      ! dhdt on ice mask
27real,dimension(nx,ny) :: calving_flux   ! calving on ice mask
28real,dimension(nx,ny) :: grline_flux    ! grounding line flux
29
30real,dimension(nx,ny) :: vieuxH_flux    ! H at previous time step
31!real,dimension(nx,ny) :: archimtab      ! point pose si > 0
32!integer,dimension(nx,ny) :: gr_line        ! point grounding line (pose sous niveau de la mer avec point flottant a cote)
33
34real,parameter :: ice_density=910.    !< densite de la glace pour conversion en masse
35
36contains
37
38subroutine init_bilan_flux
39 ! initialisation des variables
40!       diff_H=0.
41!       sum_H_old = sum(H(2:nx-1,2:ny-1),mask=ice(2:nx-1,2:ny-1)==1)
42!       tot_water(:,:)=0.
43!       bm_dt(:,:)=0.
44!       bmelt_dt(:,:)=0.
45
46  dtoutflux=5. ! initMIP sorties flux tous les 5 ans
47  dtsum=0.
48  bm_flux(:,:)=0.
49  bmelt_flux(:,:)=0.
50  dhdt_flux(:,:)=0.
51  calving_flux(:,:)=0.
52  grline_flux(:,:)=0.
53  vieuxH_flux(:,:) = H(:,:)
54 
55end subroutine init_bilan_flux
56
57       
58
59
60subroutine bilan_flux_output
61
62dtsum = dtsum + dt
63
64
65where (ice(2:nx-1,2:ny-1).eq.1)
66   bm_flux(2:nx-1,2:ny-1) = bm_flux(2:nx-1,2:ny-1) + bm(2:nx-1,2:ny-1)*dt          ! somme Bm sur dt
67   bmelt_flux(2:nx-1,2:ny-1) = bmelt_flux(2:nx-1,2:ny-1) + bmelt(2:nx-1,2:ny-1)*dt + ablbord(2:nx-1,2:ny-1)  ! somme bmelt sur dt
68   dhdt_flux(2:nx-1,2:ny-1) = dhdt_flux(2:nx-1,2:ny-1) + (H(2:nx-1,2:ny-1)-vieuxH_flux(2:nx-1,2:ny-1))       ! somme dhdt  sur dt   
69endwhere
70
71calving_flux(2:nx-1,2:ny-1) = calving_flux(2:nx-1,2:ny-1) + calv(2:nx-1,2:ny-1)
72
73where (gr_line(:,:).eq.1)
74!~   grline_flux(:,:)= (((uxbar(:,:)+eoshift(uxbar(:,:),shift=1,boundary=0.,dim=1))**2+ &
75!~                     (uybar(:,:)+eoshift(uybar(:,:),shift=1,boundary=0.,dim=2))**2)**0.5)*0.5 &
76!~                     *H(:,:) + grline_flux(:,:)
77   grline_flux(:,:)= - sqrt(                                                                       &
78                            ( (uxbar(:,:)+eoshift(uxbar(:,:),shift=1,boundary=0.,dim=1))*dy/2. )**2 +  &
79                            ( (uybar(:,:)+eoshift(uybar(:,:),shift=1,boundary=0.,dim=2))*dx/2. )**2 )  &
80                           * H(:,:) * dt /(dx*dy) + grline_flux(:,:)
81endwhere
82
83
84! calcul de la moyenne des flux sur dtoutflux
85!if ((mod(abs(time),dtoutflux).lt.dtmin).and.(abs(dtsum-dtoutflux).lt.dtmin)) then
86if (mod(abs(time),dtoutflux).lt.dtmin) then
87  bm_flux(:,:)=(bm_flux(:,:)/dtsum)*ice_density/secyear
88  bmelt_flux(:,:)=(bmelt_flux(:,:)/dtsum)*ice_density/secyear
89  dhdt_flux(:,:)=(dhdt_flux(:,:)/dtsum)/secyear
90  calving_flux(:,:)=(calving_flux(:,:)/dtsum)*ice_density/secyear
91  grline_flux(:,:)= (grline_flux(:,:)/dtsum)*ice_density/secyear 
92
93! sorties dans debug_3d:
94  debug_3d(:,:,106) = bm_flux(:,:)       ! acabf
95  debug_3d(:,:,107) = -bmelt_flux(:,:)   ! libmassbf
96  debug_3d(:,:,108) = dhdt_flux(:,:)     ! dlithkdt
97  debug_3d(:,:,109) = calving_flux(:,:)  ! licalvf
98  debug_3d(:,:,110) = grline_flux(:,:)   ! ligroundf
99
100! remise a 0
101  dtsum = 0.
102  bm_flux(:,:)      = 0.
103  bmelt_flux(:,:)   = 0.
104  dhdt_flux(:,:)    = 0.
105  calving_flux(:,:) = 0.
106  grline_flux(:,:)  = 0.
107 
108!~ else if (abs(dtsum-dtoutflux).lt.dtmin) then
109!~   print*,'pb sorties flux dtoutflux : time=', time
110!~   print*,'dtsum=',dtsum,' dtoutflux', dtoutflux
111!~ else if (mod(abs(time),dtoutflux).lt.dtmin) then
112!~   print*,'pb sorties flux dtsum : time=', time
113!~   print*,'dtsum=',dtsum,' dtoutflux', dtoutflux
114endif
115
116! reinitialisation des variables
117vieuxH_flux(:,:)=H(:,:)
118ablbord(:,:)=0.
119
120end subroutine bilan_flux_output
121
122end module bilan_flux_mod
Note: See TracBrowser for help on using the repository browser.