source: trunk/SOURCES/massb_perturb_mois.f90 @ 25

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

initial import GRISLI trunk

File size: 2.8 KB
Line 
1!> \file massb_perturb_mois.f90
2!!Module pour le calcul du bilan de masse en mode perturbation
3!<
4
5!> SUBROUTINE: massb_perturb_mois()
6!! \author ...
7!! \date ...
8!! @note Cette routine permet de calculer le mass balance en mode perturbation
9!! @note Used modules:
10!! @note    - use climat_perturb_mois_mod
11!! @note    - use declare_month
12!<
13subroutine massb_perturb_mois              ! calcule le mass balance en mode perturbation
14                                         
15  use climat_perturb_mois_mod
16  use declare_month
17
18  implicit none
19
20
21! variable de travail local
22!  real,dimension(nx,ny,mois+1) :: Tm_time
23!  real,dimension(nx,ny,mois+1) :: Pm_time ! Pm_time n'est pas utilise
24
25
26  real zel     ! variable de travail
27  real T_surf_ref
28
29
30
31   do j=1,ny
32      do i=1,nx
33! Zs est l'altitude de la surface qu'elle soit mer, glace ou lac
34         ZS(i,j)=max(slv(i,j),S(I,J))
35
36         do mo=1,mois
37
38! on corrige la temperature par rapport a l'alti de ref
39            T_surf_ref = Tm_0(i,j,mo) + lapserate(i,j,mo) * (ZS(i,j)-S0CLIM(i,j))
40
41! on ajoute le terme de climat-perturb
42            T_surf_ref = T_surf_ref+Tafor
43
44! on corrige egalement les precipes si on considere les retroactions
45            if(retroac.eq.1) then      ! full retroaction acc.
46               Pm_surf(i,j,mo)= Pm_surf(i,j,mo)*exp(rappact*(T_surf_ref-Tm_surf(i,j,mo)))
47            end if
48            Tm_surf(i,j,mo)=T_surf_ref
49
50         end do
51      end do
52   end do
53
54! Pour les modulos
55   Tm_surf(:,:,mois+1)=Tm_surf(:,:,1)
56   Pm_surf(:,:,mois+1)=Pm_surf(:,:,1)
57
58
59   do j=1,ny
60      do i=1,nx
61         Tann(i,j)=(sum(Tm_surf(i,j,:))-Tm_surf(i,j,1))/12   ! janvier compte en double
62         !Tann(i,j)=sum(Tm_surf(i,j,:))/12
63         Tjuly(i,j)= Tm_surf(i,j,7)
64      end do
65   end do
66
67   debug_3D(:,:,59)=Tjuly(:,:)
68! aurel : certainement pas a mettre la mais pour le moment...
69   do j=1,ny
70      do i=1,nx
71         debug_3D(i,j,60)=sqrt(ux(i,j,1)**2+uy(i,j,1)**2)
72      end do
73   end do
74
75
76call accum_month()    ! on calcule l'accumulation deduite de ce champ (temp, precipes)
77
78
79
80!     ablation et bilan de masse sont maintenant appeles par forclim
81
82
83
84!--------------------------------------------------------------------------------------
85! etape 3 : calcul tann et tjuly // sortie de debug3D
86
87!!$do j=1,ny
88!!$   do i=1,nx
89!!$      Tann(i,j)=sum(Tm_surf(i,j,:))/12
90!!$      Tjuly(i,j)=(Tm_surf(i,j,6)+Tm_surf(i,j,7)+Tm_surf(i,j,8))/3
91!!$      zel=max(S(i,j),20.*(ylat(i,j)-65.))
92!!$      Tann_par(i,j)=49.13-0.007992*Zel-0.7576*ylat(i,j)
93!!$      Tjuly_par(i,j)= 30.78-0.006277*S(i,j)-0.3262*ylat(i,j)
94!!$   end do
95!!$end do
96
97!!$debug_3D(:,:,29)=Tann(:,:)-Ta0(:,:)
98!!$debug_3D(:,:,30)=Acc(:,:)-precip(:,:)
99!!$debug_3D(:,:,31)=Acc(:,:)-bmelt(:,:)
100!!$debug_3D(:,:,46)=Tann(:,:)-Tann_par(:,:)
101!!$debug_3D(:,:,47)=Tjuly(:,:)-Tjuly_par(:,:)
102!!$debug_3D(:,:,48)=Tm_surf(:,:,7)-Tjuly_par(:,:)
103
104return
105end subroutine massb_perturb_mois
106
Note: See TracBrowser for help on using the repository browser.