source: trunk/SOURCES/ablation_ann_july_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: 3.9 KB
Line 
1!> \file ablation_ann_july_mod.f90
2!! Module pour l'initialisation et le calcul de l'ablation annuelle
3!<
4
5!> \namespace ablation_ann
6!! Initialisation et calcul de l'ablation annuelle
7!! \author ...
8!! \date ...
9!! @note Used module
10!! @note use module3D_phy
11!! @todo essayer de  mettre certaines autres variables en local.summ, tt,pdd,pds,simax,pdsi, .....
12!<
13     
14module ablation_ann             ! positive degree day  base annuelle
15
16                                ! remarque : on pourrait mettre certaines autres variables en local.
17                                ! summ, tt,pdd,pds,simax,pdsi, .....
18
19use module3d_phy
20implicit none
21
22! dimensionnements locaux
23real ::  dtp                   !< integrating step for positive degree days (degrees)
24real ::  pddct                 !< ct for PDD calculation     
25real ::  py                    !< ct for PDD calculation
26real ::  s22                   !< ct for PDD calculation
27real ::  nyear                 !< number of months in 1 year, st. dev. for temp *)
28
29
30! parametres qui seront lus dans le fichier param
31
32real ::  sigma                 !< variabilite Tday
33real ::  CSI                   !< proportion of melted water that can refreeze *)
34real ::  Cice                  !< melting factors for ice
35real ::  Csnow                 !< melting factors for snow
36
37contains
38!-------------------------------------------------------------------------------
39
40!> SUBROUTINE: init_ablation
41!! @note Attribution des parametres lies a la methode de calcul
42!>
43subroutine init_ablation
44
45namelist/ablation_ann/Cice,Csnow,csi,sigma
46
47! lecture des parametres
48
49rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
50read(num_param,ablation_ann)
51
52write(num_rep_42,*)'!___________________________________________________________' 
53write(num_rep_42,*)'!             module ablation_ann     PDD base Tann et Tjuly'
54write(num_rep_42,ablation_ann)                   ! pour ecrire les parametres lus
55write(num_rep_42,*)
56write(num_rep_42,*)
57write(num_rep_42,*)'! Cice and Csnow, melting factors for ice and snow '
58write(num_rep_42,*)'! sigma variabilite Tday'
59write(num_rep_42,*)'! csi proportion of melted water that can refreeze '
60
61
62! attribution des parametres lies a la methode de calcul
63dtp=2.0
64nyear=12
65s22=0.5/sigma/sigma
66py=2*pi/nyear
67pddct=dtp/sigma/sqrt(2.*pi)/nyear*365.
68
69
70end subroutine init_ablation
71
72!-----------------------------------------------------------------------------
73!> SUBROUTINE: ablation
74!!@note Calcul de l'ablation annuelle july
75!>
76subroutine ablation
77
78! positive degree day
79do j=1,ny
80   do i=1,nx
81        summ=0.0
82        do nday=1,nyear   
83           tt(nday)=tann(i,j)+(tjuly(i,j)-tann(i,j))*cos(py*nday) 
84           temp=0.0 
85           do k=1,50
86              if (temp.gt.tt(nday)+2.5*sigma) goto 205
87              summ=summ+temp*exp(-(temp-tt(nday))*(temp-tt(nday))*s22) 
88              temp=temp+dtp
89           end do
90205        continue
91        end do
92        pdd(i,j)=summ*pddct 
93     end do
94  end do
95
96do j=1,ny
97   do i=1,nx
98
99      pds=acc(i,j)/csnow        !  positive degrees required to melt the snow laye
100      simax=acc(i,j)*csi        !  maximum amount of super. ice that can be formed
101      pdsi=simax/cice           !  pos. degrees required to melt the superimposed ice
102
103      if (pdd(i,j).le.csi*pds) then
104         sif=pdd(i,j)*csnow
105      else
106         sif=simax 
107      endif
108
109!     impact of refreesing on surface temperature
110
111        ts(i,j)=(tann(i,j)+26.6*sif)
112        ts(i,j)=min(0.0,ts(i,j))
113        tshelf(i,j)=ts(i,j)
114
115!      mass balance
116        if (pdd(i,j).le.csi*pds)  bm(i,j)=acc(i,j) 
117
118        if ((csi*pds.lt.pdd(i,j)).and.(pdd(i,j).le.pds))       &
119             bm(i,j)=acc(i,j)+simax-pdd(i,j)*csnow 
120
121        if ((pds.lt.pdd(i,j)).and.(pdd(i,j).le.pds+pdsi))      &
122             bm(i,j)=simax-(pdd(i,j)-pds)*cice
123
124        if (pds+pdsi.le.pdd(i,j))   bm(i,j)=(pds+pdsi-pdd(i,j))*cice
125
126     end do
127  end do
128
129debug_3D(:,:,31)=bm(:,:)-bmelt(:,:)
130
131end subroutine ablation
132
133end module ablation_ann
Note: See TracBrowser for help on using the repository browser.