source: branches/GRISLIv3/SOURCES/ablation_mod.f90 @ 446

Last change on this file since 446 was 446, checked in by aquiquet, 8 months ago

Cleaning branch: use only statement in module3D

File size: 11.0 KB
Line 
1!********************************************************************
2! Module pour le calcul de l'ablation
3! lecture dans le fichier param de la methode de caclul du BM a utiliser
4! pdd_type : defini le type de pdd utilise
5! pdd_type=0  ! pdd reeh standard
6! pdd_type=1  ! pdd fausto
7! pdd_type=2  ! pdd Tarasov
8module ablation_mod
9
10  use geography, only : nx,ny
11
12  implicit none
13  integer :: pdd_type  ! type de calcul du PDD utilise
14  logical :: annual    ! T = annuel, F = mensuel
15
16  real :: sigma_ice                           ! variabilite Tday
17  real :: Cice                                ! melting factors for ice
18  real :: Csnow                               ! melting factors for snow
19  real :: csi                                 ! proportion of melted water that can refreeze *)
20
21  real, dimension(nx,ny) :: Cice_2D           ! melting factors for ice 2D
22  real, dimension(nx,ny) :: Csnow_2D          ! melting factors for snow 2D
23  real, dimension(nx,ny) :: csi_2D            ! proportion of melted water that can refreeze 2D
24  REAL, dimension(nx,ny) :: sigma_ice_2D      ! sigma fn altitude 2D
25
26  REAL, PARAMETER :: PI_L= 3.1415926
27  REAL, PARAMETER :: DTP=2.0                  ! integrating step for positive degree days (degrees)
28
29  REAL, dimension(nx,ny) :: S22
30  INTEGER, PARAMETER :: NYEAR = 12, NYEAR2 = 360 ! number of months/days in 1 year, st. dev. for temp *)
31  real, parameter :: PYG=2*PI_L/NYEAR            ! ct for PDD calculation (PYG remplace PY qui existe dans CLIMBER
32  REAL, dimension(nx,ny) :: PDDCT                ! ct for PDD calculation
33  REAL, dimension(nx,ny) :: PDDCT2               ! ct for PDD calculation
34  integer :: methode_abl=0                       ! selection methode pdd (0) ou van den Berg 2008 (1)
35
36contains
37
38
39subroutine init_ablation
40 
41  use module3d_phy,only:num_param,num_rep_42
42
43  real :: Cice, Csnow, csi
44 
45! fichiers snapshots
46  namelist/ablation/pdd_type,annual,Cice,Csnow,csi,sigma_ice
47
48! formats pour les ecritures dans 42
49428 format(A)
50  rewind(num_param)        ! pour revenir au debut du fichier param_list.dat
51  read(num_param,ablation)
52  write(num_rep_42,428)'!___________________________________________________________' 
53  write(num_rep_42,428) '&ablation                             ! module ablation_mod'
54  write(num_rep_42,'(A,I2)')       'pdd_type = ',pdd_type
55  write(num_rep_42,'(A,L7)')       'annual = ',annual
56  write(num_rep_42,'(A,F12.6)')    'Cice = ',Cice
57  write(num_rep_42,'(A,F12.6)')    'Csnow = ',Csnow
58  write(num_rep_42,'(A,F12.6)')    'csi = ',csi
59  write(num_rep_42,'(A,F12.6)')    'sigma_ice = ',sigma_ice
60  write(num_rep_42,*)'/'                     
61  write(num_rep_42,428) '! pdd_type : 0 reeh, 1 Fausto, 2 Tarasov'
62  write(num_rep_42,428) '! annual : T = annuel, F = mensuel'
63  write(num_rep_42,428)'! Cice and Csnow, melting factors for ice and snow '
64  write(num_rep_42,428)'! sigma variabilite Tday'
65  write(num_rep_42,428)'! csi proportion of melted water that can refreeze '
66  write(num_rep_42,*)
67
68! initialisation des variables 2D :
69  Cice_2D(:,:)=Cice
70  Csnow_2D(:,:)=Csnow
71  csi_2D(:,:)=csi
72  sigma_ice_2D(:,:)=sigma_ice
73  S22=0.5/sigma_ice/sigma_ice
74  PDDCT=DTP/sigma_ice/SQRT(2.*PI_L)/NYEAR*365.
75  PDDCT2=DTP/sigma_ice/SQRT(2.*PI_L)/NYEAR2*365.
76end subroutine init_ablation
77
78subroutine ablation
79! calcul de l'ablation avec methode annuelle ou mensuelle (flag dans fichier param)
80! remplace l'ancienne subroutine ablation-0.2.f
81!     *********************************************************
82!     (******** ABLATION ********)
83!     *********************************************************
84
85
86  use module3d_phy,only:Tjuly,Tann,Tmois,acc,pdd,TS,Tshelf,precip,BM,Abl,S,cl
87  use param_phy_mod,only:dice
88
89  IMPLICIT NONE
90
91  real, dimension(nx,ny) :: pds, simax, pdsi, sif
92  integer :: i,j,k,mo,nday
93  real :: summ
94  real :: temp
95  real,dimension(365) :: TT             !< air temperature yearly cycle, for PDD
96!- SC- Uniquement pour pdd Tarasov
97  REAL, DIMENSION(nx,ny) :: pr_ice_eq, totmelt, snowmelt, cpsurf
98  REAL, DIMENSION(nx,ny) :: refr2, refreezed_ice
99
100
101! pour pdd fausto calcul Cice_2D, csnow_2D, csi_2D,sigma_ice_2D
102  IF (pdd_type.eq.1) THEN
103      WHERE (Tjuly(:,:).GE.10.)
104          Cice_2D(:,:)=7.
105          Csnow_2D(:,:)=3.
106      ELSEWHERE ((-1.LE.Tjuly(:,:)).AND.(Tjuly(:,:).LT.10.))
107          Cice_2D(:,:)=7.+((15.-7.)/((10.+1)**3))*((10.-Tjuly(:,:))**3)
108          Csnow_2D(:,:)=3.
109      ELSEWHERE
110          Cice_2D(:,:)=15.
111          Csnow_2D(:,:)=3.
112      ENDWHERE
113! pour etre en m/an :
114      Cice_2D(:,:)=Cice_2D(:,:)/1000.
115      Csnow_2D(:,:)=Csnow_2D(:,:)/1000.
116!cdc version std 
117      sigma_ice_2D(:,:)=1.574+0.0012224*S(:,:)
118! calcul de S22, PDDCT et PDDCT2
119      S22(:,:)=0.5/SIGMA_ICE_2D(:,:)/SIGMA_ICE_2D(:,:)
120      PDDCT(:,:)=DTP/SIGMA_ICE_2D(:,:)/SQRT(2.*PI_L)/NYEAR*365.
121      PDDCT2(:,:)=DTP/SIGMA_ICE_2D(:,:)/SQRT(2.*PI_L)/NYEAR2*365.
122
123! calcul du tx de regel
124      WHERE (S(:,:).LE.800)
125          CSI_2D(:,:)=0.
126      ELSEWHERE ((800.LT.S(:,:)).AND.(S(:,:).LT.2000.))
127          CSI_2D(:,:)=(S(:,:)-800.)*0.000833
128      ELSEWHERE
129          CSI_2D(:,:)=1.
130      ENDWHERE
131
132  ELSEIF (pdd_type.EQ.2) THEN ! pdd Tarasov
133
134      where (TJULY(:,:).gt.10.)
135        Cice_2D(:,:) = 8.3*1e-3
136        Csnow_2D(:,:) = 4.3*1e-3
137      endwhere
138      where(TJULY(:,:).gt.-1..and.TJULY(:,:).lt.10.)
139        Cice_2D(:,:) = 1e-3*(8.3+0.0067*(10.-TJULY(:,:))**3)
140        Csnow_2D(:,:) = 1e-3*(2.8+0.15*TJULY(:,:))
141      endwhere
142      where(TJULY(:,:).le.-1.)
143        Cice_2D(:,:) = 17.22*1e-3
144        Csnow_2D(:,:) = 2.65*1e-3
145      endwhere
146
147      sigma_ice_2D(:,:)=sigma_ice
148  ENDIF
149
150
151
152  IF (annual) THEN
153!     (*** positive degree days Tann et Tjuly***)
154      DO i=1,nx
155        DO j=1,ny
156          summ=0.0
157          month_reconstr: DO nday=1,nyear   ! reconstitution du cycle annuel
158            tt(nday)=tann(i,j)+(tjuly(i,j)-tann(i,j))*COS(pyg*nday) 
159            temp=0.0
160            k=1
161            average_day_reconstr: DO WHILE ((temp.LE.tt(nday)+2.5*sigma_ice_2D(i,j)).AND.k.LE.50) ! pdd d'un jour par mois   
162              summ=summ+temp*EXP(-(temp-tt(nday))*(temp-tt(nday))*s22(i,j)) 
163              temp=temp+dtp
164              k=k+1
165            END DO average_day_reconstr
166          END DO month_reconstr
167          pdd(i,j)=summ*pddct(i,j) 
168        END DO
169      END DO
170  ELSE            ! pdd mensuel
171      DO i=1,nx
172        DO j=1,ny
173          summ=0.0
174! On a deja calcule Tmois(i,j,m) : temperature moyenne de chaque mois
175          month: DO mo=1,12 ! boucle sur les mois
176            temp=0.0          ! variable d'integration
177            k=1
178            average_day: DO WHILE ((temp.LE.Tmois(i,j,mo)+2.5*sigma_ice_2D(i,j)).AND.k.LE.50) ! pdd d'un jour par mois   
179              summ=summ+temp*EXP(-((temp-Tmois(i,j,mo))*(temp-Tmois(i,j,mo)))*s22(i,j))
180              temp=temp+dtp     ! pdtp pas d'integration
181              k= k+1
182            END DO average_day
183          END DO month          !  boucle sur le mois
184          pdd(i,j)=summ*pddct(i,j)   ! pdd pour toute l'annee
185        END DO
186      END DO
187  ENDIF
188
189
190!     calcul du Bilan de masse
191  IF (methode_abl.EQ.0) THEN
192     IF (pdd_type.EQ.1) THEN
193        PDS(:,:)=ACC(:,:)/Csnow_2D(:,:)
194        SIMAX(:,:)=ACC(:,:)*CSI_2D(:,:)
195        PDSI(:,:)=SIMAX(:,:)/Cice_2D(:,:)
196! avec regel de 60% puis fonte (2 premiers where) :
197!      WHERE (PDD(:,:).LE.CSI_2D*PDS(:,:))
198!          BM(:,:)=ACC(:,:)
199!          SIF(:,:)=PDD(:,:)*Csnow_2D(:,:)
200!      endwhere
201!      WHERE ((CSI_2D*PDS(:,:).LT.PDD(:,:)).AND.(PDD(:,:).LE.PDS(:,:)))
202!          BM(:,:)=ACC(:,:)+SIMAX(:,:)-PDD(:,:)*Csnow_2D
203!          SIF(:,:)=SIMAX(:,:)
204!      endwhere
205        WHERE (PDD(:,:).LE.PDS(:,:))  ! test avec regel de 60% progressif
206           BM(:,:)=ACC(:,:)-PDD(:,:)*Csnow_2D*(1-CSI_2D(:,:))
207           SIF(:,:)=PDD(:,:)*Csnow_2D(:,:)*(1-CSI_2D(:,:))
208        endwhere
209        WHERE ((PDS(:,:).LT.PDD(:,:)).AND.(PDD(:,:).LE.PDS(:,:)+PDSI(:,:)))
210           BM(:,:)=SIMAX(:,:)-(PDD(:,:)-PDS(:,:))*Cice_2D
211           SIF(:,:)=SIMAX(:,:)
212        endwhere
213        WHERE (PDS(:,:)+PDSI(:,:).LE.PDD(:,:))
214           BM(:,:)=(PDS(:,:)+PDSI(:,:)-PDD(:,:))*Cice_2D
215           SIF(:,:)=SIMAX(:,:)
216        endwhere
217     ELSEIF (PDD_type.EQ.2) THEN  ! PDD Tarasov
218        PDS(:,:)=ACC(:,:)/Csnow_2D(:,:)
219        pr_ice_eq(:,:) = amax1(0.,((PRECIP(:,:)/DICE)-ACC(:,:)))  ! precipe liquide (ice equivalent)
220
221        WHERE (PDD(:,:).LE.PDS(:,:))
222           totmelt(:,:) = Csnow_2D(:,:)*PDD(:,:)
223        ELSEWHERE
224           totmelt(:,:) = Csnow_2D(:,:)*PDS(:,:) + Cice_2D(:,:)*(PDD(:,:)-PDS(:,:))
225        ENDWHERE
226
227        snowmelt(:,:) = amin1(totmelt(:,:),ACC(:,:))
228
229! Deux formules possibles pour la capacité calorifique (en J/kg.K):
230! Formule 1 correspond à celle figurant dans prop-thermiques_mod.f90
231! Formule 2 : celle de Tarasov
232        cpsurf(:,:) = 2115.3+7.79293*TANN(:,:) ! Formule Catherine
233!      cpsurf(:,:) = 152.5 + 7.122*TANN(:,:) ! Formule Tarasov
234
235!        where(snowmelt(:,:).lt.ACC(:,:))
236        refr2(:,:) = 2.2*(ACC(:,:)-snowmelt(:,:))-(cpsurf(:,:)/CL)*amin1(TANN(:,:),0.)
237        refreezed_ice(:,:) = amin1(pr_ice_eq(:,:)+snowmelt(:,:),refr2(:,:))
238!        elsewhere
239!           refr2(:,:) = -(cpsurf(:,:)/CL)*amin1(TANN(:,:),0.)
240!           refreezed_ice(:,:) = amin1(pr_ice_eq(:,:)+snowmelt(:,:),refr2(:,:))
241!        endwhere
242
243        bm(:,:) = ACC(:,:)+refreezed_ice(:,:)-totmelt(:,:)
244        SIF(:,:)=refreezed_ice(:,:)
245
246     ELSE ! pdd standard reeh
247!       (* Positive degrees required to melt the snow layer *)
248        PDS(:,:)=ACC(:,:)/Csnow_2D(:,:)
249!       (* Maximum amount of super. ice that can be formed *)
250        SIMAX(:,:)=ACC(:,:)*CSI_2D(:,:)
251!       (* Pos. degrees required to melt the superimposed ice *)
252        PDSI(:,:)=SIMAX(:,:)/Cice_2D(:,:)
253! avec regel de 60% puis fonte (2 premiers where) :
254        WHERE (PDD(:,:).LE.CSI_2D*PDS(:,:))
255           BM(:,:)=ACC(:,:)
256           SIF(:,:)=PDD(:,:)*Csnow_2D(:,:)
257        endwhere
258        WHERE ((CSI_2D*PDS(:,:).LT.PDD(:,:)).AND.(PDD(:,:).LE.PDS(:,:)))
259           BM(:,:)=ACC(:,:)+SIMAX(:,:)-PDD(:,:)*Csnow_2D
260           SIF(:,:)=SIMAX(:,:)
261        endwhere
262!       WHERE (PDD(:,:).LE.PDS(:,:))                  ! test avec regel de 60% progressif|
263!           BM(:,:)=ACC(:,:)-PDD(:,:)*Csnow_2D*(1-CSI_2D)   ! remplace les 2 premiers where    |
264!           SIF(:,:)=PDD(:,:)*Csnow_2D(:,:)*(1-CSI_2D(:,:)) !----------------------------------|
265!      endwhere
266        WHERE ((PDS(:,:).LT.PDD(:,:)).AND.(PDD(:,:).LE.PDS(:,:)+PDSI(:,:)))
267           BM(:,:)=SIMAX(:,:)-(PDD(:,:)-PDS(:,:))*Cice_2D
268           SIF(:,:)=SIMAX(:,:)
269        endwhere
270        WHERE (PDS(:,:)+PDSI(:,:).LE.PDD(:,:))
271           BM(:,:)=(PDS(:,:)+PDSI(:,:)-PDD(:,:))*Cice_2D
272           SIF(:,:)=SIMAX(:,:)
273        endwhere
274     ENDIF
275  ELSEIF ( methode_abl.EQ.1 ) THEN ! Insolation Temperature Melt equation, van den Berg, 2008
276     SIMAX(:,:)=ACC(:,:)*CSI
277     SIF(:,:)=SIMAX(:,:)
278  ELSE
279     print*,'ablation.f90 : pb methode calcul BM'
280     print*,'ATTENTION methode_abl > 1 NON COMPATIBLE AVEC iLOVECLIM' 
281     stop
282  ENDIF
283
284! calcul de la temperature de surface (utilisee dans icetemp) :
285  TS(:,:)=(TANN(:,:)+26.6*SIF(:,:))
286  TS(:,:)=min(0.0,TS(:,:))
287  tshelf(:,:)=TS(:,:)
288  Abl(:,:)=BM(:,:)-Acc(:,:)
289
290END SUBROUTINE ABLATION
291
292end module ablation_mod
Note: See TracBrowser for help on using the repository browser.