1 | !> \file climat-perturb_mod-0.4.f90 |
---|
2 | !! Module pour les variations temporelles des variables climatiques |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace climat_perturb_mod |
---|
6 | !!Module pour les variations temporelles des variables climatiques |
---|
7 | !! \author C. Dumas |
---|
8 | !! \date 12/2015 |
---|
9 | !! @note Used modules: |
---|
10 | !! @note - use module3D_phy |
---|
11 | !< |
---|
12 | module climat_perturb_mod |
---|
13 | |
---|
14 | |
---|
15 | use module3d_phy,only: S,S0,Tann,Tjuly,acc,Ylat,num_forc,num_param,num_rep_42,tafor,time,sealevel,sealevel_2d,coefbmshelf |
---|
16 | use geography, only: nx,ny,dirforcage,dirnameinp |
---|
17 | use io_netcdf_grisli, only: read_ncdf_var |
---|
18 | |
---|
19 | implicit none |
---|
20 | |
---|
21 | integer nft !< nombre de lignes a lire dans le fichier forcage climatique |
---|
22 | real,dimension(:),allocatable :: tdate !< time for climate forcing |
---|
23 | real,dimension(:),allocatable :: tpert !< temperature for climate forcing |
---|
24 | real,dimension(:),allocatable :: spert !< sea surface perturbation |
---|
25 | |
---|
26 | real,dimension(nx,ny) :: ta0 !< initial air temperature at sea level annual |
---|
27 | real,dimension(nx,ny) :: precip |
---|
28 | |
---|
29 | real :: coefT !< pour modifier l'amplitude de la perturb. T |
---|
30 | real :: rappact !< pour le calcul du rapport 'accumulation |
---|
31 | integer :: retroac !< 1-> full retroactions accum |
---|
32 | real :: rapbmshelf !< pour calcul coefbmshelf |
---|
33 | real :: mincoefbmelt !< butoirs mini |
---|
34 | real :: maxcoefbmelt !< butoirs maxi de coefbmelt |
---|
35 | character(len=80) :: filforc !< nom du fichier forcage |
---|
36 | |
---|
37 | ! Pour l'instant tafor est global meme si inutile si on utilise |
---|
38 | ! un forcage (variation spatiales) |
---|
39 | |
---|
40 | contains |
---|
41 | |
---|
42 | !-------------------------------------------------------------------------------- |
---|
43 | !> SUBROUTINE: input_clim |
---|
44 | !! Routine qui permet d'initialiser les variations temporelles des variables climatiques |
---|
45 | !> |
---|
46 | |
---|
47 | subroutine input_clim |
---|
48 | |
---|
49 | implicit none |
---|
50 | character(len=100) :: precip_file ! fichier precipitations |
---|
51 | character(len=100) :: temp_annual_file ! fichier temperature annuelle |
---|
52 | real :: coef_dens ! pour corriger si donnees en eq. eau |
---|
53 | logical :: temp_param ! si utilisation de temperature parametree |
---|
54 | real*8, dimension(:,:), pointer :: data_2D => null() ! donnees lues dans le netcdf |
---|
55 | |
---|
56 | |
---|
57 | character(len=8) :: control !label to check clim. forc. file (filin) is usable |
---|
58 | character(len=80):: filin |
---|
59 | integer :: err !< pour l'allocation des tableaux |
---|
60 | integer :: i,j |
---|
61 | |
---|
62 | namelist/climat_acc_T_gen/precip_file,coef_dens,temp_annual_file |
---|
63 | |
---|
64 | 428 format(A) |
---|
65 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
66 | read(num_param,climat_acc_T_gen) |
---|
67 | |
---|
68 | write(num_rep_42,428)'!___________________________________________________________' |
---|
69 | write(num_rep_42,428)'! module lect_clim_acc_T_ant_gen ' |
---|
70 | write(num_rep_42,climat_acc_T_gen) |
---|
71 | write(num_rep_42,428)'!___________________________________________________________' |
---|
72 | |
---|
73 | |
---|
74 | ! precipitation |
---|
75 | precip_file = trim(dirnameinp)//trim(precip_file) |
---|
76 | |
---|
77 | call Read_ncdf_var('precip',trim(precip_file),data_2D) |
---|
78 | precip(:,:)=data_2D(:,:) |
---|
79 | !call lect_datfile(nx,ny,precip,1,precip_file) |
---|
80 | |
---|
81 | precip(:,:)=precip(:,:)*coef_dens |
---|
82 | acc(:,:)=precip(:,:) |
---|
83 | |
---|
84 | if ((trim(temp_annual_file).eq.'no').or.(trim(temp_annual_file).eq.'NO')) then |
---|
85 | temp_param=.true. |
---|
86 | else |
---|
87 | temp_param=.false. |
---|
88 | end if |
---|
89 | |
---|
90 | ! temperature en surface |
---|
91 | |
---|
92 | test_param: if (.not.temp_param) then |
---|
93 | temp_annual_file = trim(dirnameinp)//trim(temp_annual_file) |
---|
94 | |
---|
95 | |
---|
96 | call Read_ncdf_var('Tann',trim(temp_annual_file),data_2D) |
---|
97 | Tann(:,:)=data_2D(:,:) |
---|
98 | ! call lect_input(3,'Tann',1,Tann,temp_annual_file,trim(dirnameinp)//trim(runname)//'.nc') |
---|
99 | !call lect_datfile(nx,ny,Tann,1,temp_annual_file) ! temperature annuelle |
---|
100 | |
---|
101 | else ! parametrisation de Fortuin pour la temperature annuelle. |
---|
102 | |
---|
103 | do j=1,ny |
---|
104 | do i=1,nx |
---|
105 | |
---|
106 | 7 if (s0(i,j).le.200.) then ! shelfs |
---|
107 | tann(i,j)=49.642-0.943*abs(ylat(i,j)) |
---|
108 | else if ((s0(i,j).gt.200.).and.(s0(i,j).lt.1500.)) then ! pente |
---|
109 | tann(i,j)=36.689-0.005102*s0(i,j)-0.725*abs(ylat(i,j)) |
---|
110 | else if (s0(i,j).ge.1500.) then ! plateau |
---|
111 | tann(i,j)=7.405-0.014285*s0(i,j)-0.180*abs(ylat(i,j)) |
---|
112 | endif |
---|
113 | end do |
---|
114 | end do |
---|
115 | end if test_param |
---|
116 | |
---|
117 | ta0(:,:)=tann(:,:) |
---|
118 | |
---|
119 | |
---|
120 | ! pour la temperature d'ete, idem parametrisation huybrechts |
---|
121 | do j=1,ny |
---|
122 | do i=1,nx |
---|
123 | |
---|
124 | tjuly(i,j)=tann(i,j)-17.65+0.00222*s0(i,j)& |
---|
125 | +0.40802*abs(ylat(i,j)) |
---|
126 | end do |
---|
127 | end do |
---|
128 | |
---|
129 | |
---|
130 | |
---|
131 | |
---|
132 | |
---|
133 | ! Lecture du forcage |
---|
134 | !----------------------- |
---|
135 | ! Le fichier de forcage est lu dans le fichier entree parametres |
---|
136 | |
---|
137 | filin=trim(dirforcage)//trim(filforc) |
---|
138 | |
---|
139 | open(num_forc,file=filin,status='old') |
---|
140 | |
---|
141 | read(num_forc,*) control,nft |
---|
142 | |
---|
143 | ! Determination of file size (line nb), allocation of perturbation array |
---|
144 | |
---|
145 | if (control.ne.'nb_lines') then |
---|
146 | write(6,*) filin,'indiquer le nb de ligne en debut de fichier:' |
---|
147 | write(6,*) 'le nb de lignes et le label de control nb_lines' |
---|
148 | stop |
---|
149 | endif |
---|
150 | |
---|
151 | ! Dimensionnement des tableaux tdate, .... |
---|
152 | |
---|
153 | if (.not.allocated(tdate)) then |
---|
154 | allocate(tdate(nft),stat=err) |
---|
155 | if (err/=0) then |
---|
156 | print *,"erreur a l'allocation du tableau Tdate",err |
---|
157 | stop 4 |
---|
158 | end if |
---|
159 | end if |
---|
160 | |
---|
161 | if (.not.allocated(spert)) then |
---|
162 | allocate(spert(nft),stat=err) |
---|
163 | if (err/=0) then |
---|
164 | print *,"erreur a l'allocation du tableau Spert",err |
---|
165 | stop 4 |
---|
166 | end if |
---|
167 | end if |
---|
168 | |
---|
169 | |
---|
170 | if (.not.allocated(tpert)) then |
---|
171 | allocate(tpert(nft),stat=err) |
---|
172 | if (err/=0) then |
---|
173 | print *,"erreur a l'allocation du tableau Tpert",err |
---|
174 | stop 4 |
---|
175 | end if |
---|
176 | end if |
---|
177 | |
---|
178 | do i=1,nft |
---|
179 | read(num_forc,*) tdate(i),spert(i),tpert(i) |
---|
180 | end do |
---|
181 | close(num_forc) |
---|
182 | |
---|
183 | tpert(:)=tpert(:)*coefT |
---|
184 | |
---|
185 | |
---|
186 | end subroutine input_clim |
---|
187 | |
---|
188 | !-------------------------------------------------------------------------------- |
---|
189 | !> SUBROUTINE: init_forclim |
---|
190 | !! Routine qui permet d'initialiser les variables climatiques au cours du temps |
---|
191 | !> |
---|
192 | subroutine init_forclim |
---|
193 | |
---|
194 | |
---|
195 | character(len=80):: filin |
---|
196 | integer :: i,j |
---|
197 | integer :: ivo,jvo |
---|
198 | |
---|
199 | namelist/clim_pert/coefT,rappact,retroac,rapbmshelf,mincoefbmelt,maxcoefbmelt,filforc |
---|
200 | |
---|
201 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
202 | read(num_param,clim_pert) |
---|
203 | |
---|
204 | ! formats pour les ecritures dans 42 |
---|
205 | 428 format(A) |
---|
206 | |
---|
207 | write(num_rep_42,428)'!___________________________________________________________' |
---|
208 | write(num_rep_42,428) '&clim_pert ! module climat_perturb_mod ' |
---|
209 | write(num_rep_42,*) |
---|
210 | write(num_rep_42,*) 'coefT = ', coefT |
---|
211 | write(num_rep_42,*) 'rappact = ', rappact |
---|
212 | write(num_rep_42,*) 'rretroac = ', retroac |
---|
213 | write(num_rep_42,*) 'rapbmshelf = ', rapbmshelf |
---|
214 | write(num_rep_42,*) 'mincoefbmelt = ', mincoefbmelt |
---|
215 | write(num_rep_42,*) 'maxcoefbmelt = ', maxcoefbmelt |
---|
216 | write(num_rep_42,'(A,A)') ' filforc = ', filforc |
---|
217 | write(num_rep_42,*)'/' |
---|
218 | write(num_rep_42,*) |
---|
219 | |
---|
220 | |
---|
221 | |
---|
222 | !cdc Commente pour être compatible avec lecture fichiers Cat Schoofing |
---|
223 | |
---|
224 | !~ !!!!!!!! ATTENTION AJOUTE POUR TEST MAIS A REMETTRE AU PROPRE PLUS TARD C. DUMAS !!!!!!!!!! |
---|
225 | !~ !!!!!!!! ancien input_climat_ref de lect_clim_act_anteis |
---|
226 | !~ ! accumulation de Philippe |
---|
227 | !~ filin='accumHUY40km.dat' |
---|
228 | !~ call lect_eis(nx,ny,precip,filin,DIRNAMEINP) |
---|
229 | !~ !====================================== La reponse est 42 =========== |
---|
230 | !~ write(num_rep_42,*) 'fichier accum : ', filin |
---|
231 | |
---|
232 | !~ ! cas particulier de Vostok |
---|
233 | !~ ivo=101 |
---|
234 | !~ jvo=62 |
---|
235 | !~ do j=jvo-1,jvo+1 |
---|
236 | !~ do i=ivo-1,ivo+1 |
---|
237 | !~ precip(i,j)=0.02 ! valeur plus faible a Vostok. |
---|
238 | !~ end do |
---|
239 | !~ end do |
---|
240 | !~ acc(:,:)=precip(:,:) |
---|
241 | |
---|
242 | !~ ! temperature en surface : |
---|
243 | !~ ! parametrisation de Fortuin pour la temperature annuelle. |
---|
244 | !~ do j=1,ny |
---|
245 | !~ do i=1,nx |
---|
246 | |
---|
247 | !~ if (s0(i,j).le.200.) then ! shelfs |
---|
248 | !~ tann(i,j)=49.642-0.943*abs(ylat(i,j)) |
---|
249 | !~ else if ((s0(i,j).gt.200.).and.(s0(i,j).lt.1500.)) then ! pente |
---|
250 | !~ tann(i,j)=36.689-0.005102*s0(i,j)-0.725*abs(ylat(i,j)) |
---|
251 | !~ else if (s0(i,j).ge.1500.) then ! plateau |
---|
252 | !~ tann(i,j)=7.405-0.014285*s0(i,j)-0.180*abs(ylat(i,j)) |
---|
253 | !~ endif |
---|
254 | |
---|
255 | !~ ta0(i,j)=tann(i,j) |
---|
256 | !~ ! pour la temperature d'ete, idem parametrisation huybrechts |
---|
257 | !~ tjuly(i,j)=tann(i,j)-17.65+0.00222*s0(i,j)& |
---|
258 | !~ +0.40802*abs(ylat(i,j)) |
---|
259 | !~ end do |
---|
260 | !~ end do |
---|
261 | !~ !!!!!!!! FIN MODIF TEMPORAIRE !!!!!!!!!! |
---|
262 | !cdc fin Commente pour être compatible avec lecture fichiers Cat Schoofing |
---|
263 | |
---|
264 | return |
---|
265 | end subroutine init_forclim |
---|
266 | |
---|
267 | !-------------------------------------------------------------------------------- |
---|
268 | !> SUBROUTINE: forclim |
---|
269 | !! |
---|
270 | !! Routine qui permet le calcule climatiques au cours du temps |
---|
271 | !! @note Au temps considere (time) attribue les scalaires |
---|
272 | !! @note - tafor : forcage en temperature |
---|
273 | !! @note - sealevel : forcage niveau des mers |
---|
274 | !! @note - coefbmelt : forcage fusion basale ice shelves |
---|
275 | !> |
---|
276 | subroutine forclim ! au temps considere (time) attribue les scalaires |
---|
277 | ! tafor : forcage en temperature |
---|
278 | ! sealevel : forcage niveau des mers |
---|
279 | ! coefbmelt : forcage fusion basale ice shelves |
---|
280 | |
---|
281 | ! use module3d_phy |
---|
282 | implicit none |
---|
283 | integer :: i,ift |
---|
284 | |
---|
285 | ! time en dehors des limites du fichier forcage |
---|
286 | if(time.lt.tdate(1)) then |
---|
287 | tafor=tpert(1) |
---|
288 | sealevel=spert(1) |
---|
289 | ift=1 |
---|
290 | |
---|
291 | else if (time.ge.tdate(nft)) then |
---|
292 | tafor=tpert(nft) |
---|
293 | sealevel=spert(nft) |
---|
294 | ift=nft |
---|
295 | |
---|
296 | else |
---|
297 | do i=1,nft-1 |
---|
298 | if((time.ge.tdate(i)).and.(time.lt.tdate(i+1))) then ! entre i et i+1 : cas general |
---|
299 | tafor=tpert(i)+(tpert(i+1)-tpert(i))* & |
---|
300 | (time-tdate(i))/(tdate(i+1)-tdate(i)) |
---|
301 | sealevel=spert(i)+(spert(i+1)-spert(i))* & |
---|
302 | (time-tdate(i))/(tdate(i+1)-tdate(i)) |
---|
303 | ift=i |
---|
304 | goto 100 |
---|
305 | endif |
---|
306 | end do |
---|
307 | endif |
---|
308 | 100 continue |
---|
309 | |
---|
310 | sealevel_2d(:,:)=sealevel |
---|
311 | |
---|
312 | ! coefmshelf est un coefficient qui fait varier bmgrz et bmshelf |
---|
313 | ! en fonction de tafor |
---|
314 | |
---|
315 | ! coefbmshelf=(1.+tafor/10.) ! coefbmshelf=0 pour tafor=-10deg |
---|
316 | |
---|
317 | if ((tafor.le.0).and.(tafor.gt.-5.)) then |
---|
318 | coefbmshelf=(1.+tafor/rapbmshelf) ! coefbmshelf=0 pour tafor=-7 standard precedent |
---|
319 | |
---|
320 | else if (tafor.le.-5) then ! lineaire en 0 a -10 et la valeur a -5 |
---|
321 | coefbmshelf=(1.-5./rapbmshelf)/5.*(tafor+10) |
---|
322 | |
---|
323 | else |
---|
324 | coefbmshelf=(1.+5.*tafor/rapbmshelf) ! 5 fois plus efficace vers le chaud |
---|
325 | endif |
---|
326 | |
---|
327 | coefbmshelf=max(coefbmshelf,mincoefbmelt) |
---|
328 | coefbmshelf=min(coefbmshelf,maxcoefbmelt) |
---|
329 | |
---|
330 | call massb_perturb_Tparam |
---|
331 | ! massb_perturb_Tparam remplace desormais massb_anteis_perturb |
---|
332 | ! le code est le meme pour l'Antarctique (simplement passe en format libre) |
---|
333 | ! dans le fichier massb-ant_perturb_Tparam.f90 |
---|
334 | |
---|
335 | |
---|
336 | end subroutine forclim |
---|
337 | |
---|
338 | subroutine massb_perturb_Tparam ! calcule le mass balance en mode perturbation |
---|
339 | ! avec la temperature parametree |
---|
340 | ! version pour l'antarctique |
---|
341 | ! simple copie nettoyee de massb_anteis_perturb |
---|
342 | implicit none |
---|
343 | |
---|
344 | integer :: i,j |
---|
345 | |
---|
346 | ! surface temperature et accumulation |
---|
347 | |
---|
348 | do j=1,ny |
---|
349 | do i=1,nx |
---|
350 | |
---|
351 | if(retroac.eq.1) then |
---|
352 | |
---|
353 | tann(i,j)=ta0(i,j)-0.00914*(s(i,j)-s0(i,j))+tafor |
---|
354 | tjuly(i,j)=tann(i,j)-17.65+0.00222*s(i,j) & |
---|
355 | +0.40802*abs(ylat(i,j)) |
---|
356 | |
---|
357 | acc(i,j)=precip(i,j)*exp(rappact*(tann(i,j)-ta0(i,j))) |
---|
358 | |
---|
359 | else if(retroac.eq.0) then |
---|
360 | |
---|
361 | tann(i,j)=ta0(i,j) |
---|
362 | tjuly(i,j)=tann(i,j) |
---|
363 | acc(i,j)=precip(i,j)*exp(rappact*(tann(i,j)-ta0(i,j))) |
---|
364 | endif |
---|
365 | end do |
---|
366 | end do |
---|
367 | |
---|
368 | end subroutine massb_perturb_Tparam |
---|
369 | |
---|
370 | |
---|
371 | end module climat_perturb_mod |
---|