1 | !> \file output_anta_mod-0.4.f90 |
---|
2 | !! Module de |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace output_antarcti_mod |
---|
6 | !! Module de |
---|
7 | !! \author ... |
---|
8 | !! \date ... |
---|
9 | !! @note Used module |
---|
10 | !! @note - use module3D_phy |
---|
11 | !< |
---|
12 | |
---|
13 | module output_antarcti_mod |
---|
14 | |
---|
15 | USE module3D_phy |
---|
16 | use bilan_eau_mod |
---|
17 | use netcdf |
---|
18 | use io_netcdf_grisli |
---|
19 | use bilan_flux_mod |
---|
20 | |
---|
21 | implicit none |
---|
22 | |
---|
23 | !real :: vol ; !integer :: np |
---|
24 | integer :: nflot !< nbr de point flottant |
---|
25 | real :: bmean !< |
---|
26 | real :: accmean !< accumulation moyenne |
---|
27 | real :: ablmean !< ablation moyenne |
---|
28 | real :: calvmean !< moyenne calving |
---|
29 | real :: ablbordmean !< |
---|
30 | real :: bmeltmean !< moyenne fusion basale |
---|
31 | real :: tbmean !< temperature basale moyenne |
---|
32 | real :: tbdotmean !< moyenne variation / temps temperature basale |
---|
33 | real :: vsmean !< vitesse de surface moyenne |
---|
34 | !real :: vsdotmean !< moyenne variation / temps vitesse de surface |
---|
35 | !real :: uzsmean !!!! utilise ? !< vitesse verticale de surface moyenne |
---|
36 | real :: uzsdotmean !< moyenne variation / temps vitesse verticale de surface |
---|
37 | real :: uzkmean !< moyenne vitesse verticale de surface |
---|
38 | real :: hdotmean !< moyenne derivee / temps de h |
---|
39 | real :: bdotmean !< moyenne bedrock derive / temps |
---|
40 | real :: volf !< volume au dessus de la flottaison |
---|
41 | |
---|
42 | real :: lim !< Total ice mass |
---|
43 | real :: iareag !< surface posee |
---|
44 | real :: iareaf !< surface flottante |
---|
45 | real :: tendacabf !< Total SMB flux |
---|
46 | real :: tendlibmassbf !< Total Basal mass balance flux |
---|
47 | real :: tendlibmassbffl !< Total Basal mass balance flux on floating ice |
---|
48 | real :: tendlicalvf !< Total calving flux |
---|
49 | real :: tendlifmassbf !< Total calving and ice front melting flux |
---|
50 | real :: tendligroundf !< Total grounding line flux |
---|
51 | |
---|
52 | real,dimension(nx,ny) :: corrsurf !< facteur de correction de la surface |
---|
53 | !real,parameter :: ice_density=910. !< densite de la glace pour conversion en masse |
---|
54 | |
---|
55 | ! variables netcdf |
---|
56 | integer,parameter :: ncshortout=1 ! 1 sorties netcdf short initMIP |
---|
57 | integer,parameter :: nvar=10 ! nombre de variables dans le fichier de sortie temporel netcdf |
---|
58 | integer :: ncid |
---|
59 | integer :: status |
---|
60 | integer :: timeDimID |
---|
61 | integer :: timeVarID |
---|
62 | integer,dimension(nvar) :: varID |
---|
63 | integer :: nbtimeout ! index time output |
---|
64 | real,dimension(nvar) :: var_shortoutput |
---|
65 | |
---|
66 | ! variables netcdf icecore (ic) |
---|
67 | integer,parameter :: ncicecoreout=1 ! 1 sorties netcdf icecore / 0 no output |
---|
68 | integer :: nbic ! number of ice core site in output |
---|
69 | integer,parameter :: nvaric=10 ! nombre de variables dans le fichier ice core |
---|
70 | real,parameter :: dticecoreout=50. ! sorties ice core tous les dticecoreout |
---|
71 | integer :: ncidic |
---|
72 | integer :: timeDimIDic |
---|
73 | integer :: timeVarIDic |
---|
74 | integer :: stationDimIDic |
---|
75 | integer :: stationVarIDic |
---|
76 | integer :: station_strlenDimIDic |
---|
77 | integer :: lonVarIDic |
---|
78 | integer :: latVarIDic |
---|
79 | integer,dimension(nvaric) :: varIDic |
---|
80 | integer :: nbtimeoutic ! index time output |
---|
81 | real,allocatable,dimension(:,:) :: var_icoutput |
---|
82 | character(len=100), allocatable, dimension(:) :: site_nameic |
---|
83 | integer, allocatable, dimension(:) :: iic,jic ! coordonnees i,j des sites ice core |
---|
84 | |
---|
85 | CONTAINS |
---|
86 | |
---|
87 | subroutine init_outshort |
---|
88 | |
---|
89 | double precision,dimension(:,:),pointer :: tab !< tableau 2d real pointer |
---|
90 | character(len=100),dimension(nvar) :: namevar ! name, standard_name, long_name, units |
---|
91 | character(len=100),dimension(nvar) :: standard_name |
---|
92 | character(len=100),dimension(nvar) :: long_name |
---|
93 | character(len=100),dimension(nvar) :: units |
---|
94 | |
---|
95 | ! ice core variables |
---|
96 | character(len=100),allocatable,dimension(:) :: ic_name ! ice core name |
---|
97 | real,allocatable,dimension(:) :: ic_lon ! ice core lon |
---|
98 | real,allocatable,dimension(:) :: ic_lat ! ice core lat |
---|
99 | character(len=100),dimension(nvaric) :: namevaric ! name, standard_name, long_name, units |
---|
100 | character(len=100),dimension(nvaric) :: standard_nameic |
---|
101 | character(len=100),dimension(nvaric) :: long_nameic |
---|
102 | character(len=100),dimension(nvaric) :: unitsic |
---|
103 | integer :: err ! recuperation erreur |
---|
104 | integer :: ii,jj |
---|
105 | integer :: nloop ! nbr de tour dans la boucle while (evite boucle folle) |
---|
106 | real :: dist_lonlat, dist_lonlat_new, distance ! distance to lon lat ice core site |
---|
107 | |
---|
108 | |
---|
109 | |
---|
110 | !ndisp sorite courte tous les ndisp |
---|
111 | NDISP=100 |
---|
112 | |
---|
113 | |
---|
114 | if (ncshortout.eq.1) then ! ecriture netcdf |
---|
115 | |
---|
116 | ! lecture du fichier avec les corrections de surface |
---|
117 | if (geoplace.eq.'ant16km') then |
---|
118 | call Read_Ncdf_var('z',trim(DIRNAMEINP)//'/corrsurf-initMIP-16km.grd',tab) |
---|
119 | corrsurf(:,:)=tab(:,:) |
---|
120 | else |
---|
121 | corrsurf(:,:)= 1. |
---|
122 | end if |
---|
123 | |
---|
124 | open(568,file=trim(dirsource)//'/Fichiers-parametres/short-ISMIPnc.dat',status='old') |
---|
125 | ! lecture en-tete |
---|
126 | read(568,*) |
---|
127 | read(568,*) |
---|
128 | read(568,*) |
---|
129 | read(568,*) |
---|
130 | read(568,*) |
---|
131 | ! lecture des infos sur les variables : |
---|
132 | do k=1,nvar |
---|
133 | read(568,'(a100)') namevar(k) |
---|
134 | read(568,'(a100)') standard_name(k) |
---|
135 | read(568,'(a100)') long_name(k) |
---|
136 | read(568,'(a100)') units(k) |
---|
137 | read(568,*) |
---|
138 | enddo |
---|
139 | close(568) |
---|
140 | |
---|
141 | ! Fichier Netcdf initMIP |
---|
142 | ! creation du fichier Netcdf : |
---|
143 | status=nf90_create(path = trim(dirnameout)//'short'//runname//'.nc', cmode = nf90_clobber, ncid = ncid) |
---|
144 | if (status /= nf90_noerr) call handle_err(status) |
---|
145 | |
---|
146 | ! definition des dimension : |
---|
147 | status=nf90_def_dim(ncid, name="time", len=NF90_UNLIMITED, dimid=timeDimID) |
---|
148 | if (status /= nf90_noerr) call handle_err(status) |
---|
149 | status=nf90_def_var(ncid, name="time", xtype=nf90_float, dimids=(/ timeDimID/), varid=timeVarID) |
---|
150 | if (status /= nf90_noerr) call handle_err(status) |
---|
151 | status=nf90_put_att(ncid, timeVarID, "standard_name", "time") |
---|
152 | if (status /= nf90_noerr) call handle_err(status) |
---|
153 | status=nf90_put_att(ncid, timeVarID,"units", "years since 1995-01-01 00:00:00") |
---|
154 | if (status /= nf90_noerr) call handle_err(status) |
---|
155 | status=nf90_put_att(ncid, timeVarID,"calendar", "360_day") |
---|
156 | if (status /= nf90_noerr) call handle_err(status) |
---|
157 | |
---|
158 | ! definition des variables de sortie : |
---|
159 | do k=1,nvar ! boucle sur le nbr de variable a definir |
---|
160 | status=nf90_def_var(ncid, name=trim(namevar(k)), xtype=nf90_float, dimids= & |
---|
161 | (/ timeDimID /), varid=varID(k)) |
---|
162 | if (status /= nf90_noerr) call handle_err(status) |
---|
163 | status=nf90_put_att(ncid, varID(k), "standard_name", trim(standard_name(k))) |
---|
164 | if (status /= nf90_noerr) call handle_err(status) |
---|
165 | status=nf90_put_att(ncid, varID(k), "long_name", trim(long_name(k))) |
---|
166 | if (status /= nf90_noerr) call handle_err(status) |
---|
167 | status=nf90_put_att(ncid, varID(k), "units", trim(units(k))) |
---|
168 | if (status /= nf90_noerr) call handle_err(status) |
---|
169 | enddo |
---|
170 | |
---|
171 | ! fin de la definition du fichier : |
---|
172 | status=nf90_enddef(ncid) |
---|
173 | if (status /= nf90_noerr) call handle_err(status) |
---|
174 | nbtimeout = 0 ! initialisation compteur sorties axe time |
---|
175 | else ! pas de sortie netcdf et sans correction de surface |
---|
176 | corrsurf(:,:)=1. |
---|
177 | endif |
---|
178 | |
---|
179 | |
---|
180 | if (ncicecoreout.eq.1) then ! ecriture netcdf Ice Core |
---|
181 | ! Fichier Netcdf ice core Antarctique |
---|
182 | ! Lecture du fichier avec la liste des sites de forages |
---|
183 | ! infos sur les variables : |
---|
184 | open(568,file=trim(dirsource)//'/Fichiers-parametres/icecore-Ant.dat',status='old') |
---|
185 | ! lecture en-tete |
---|
186 | read(568,*) |
---|
187 | read(568,*) |
---|
188 | read(568,*) |
---|
189 | read(568,*) |
---|
190 | read(568,*) |
---|
191 | ! lectuire du nombre de sites de forages: |
---|
192 | read(568,*) nbic |
---|
193 | read(568,*) |
---|
194 | ! allocation des variables dependant du nombre de site de forage (nbic) |
---|
195 | if (.not.allocated(iic)) then |
---|
196 | allocate(iic(nbic),stat=err) |
---|
197 | if (err/=0) then |
---|
198 | print *,"erreur a l'allocation du tableau iic dans init_outshort",err |
---|
199 | stop 4 |
---|
200 | end if |
---|
201 | end if |
---|
202 | if (.not.allocated(jic)) then |
---|
203 | allocate(jic(nbic),stat=err) |
---|
204 | if (err/=0) then |
---|
205 | print *,"erreur a l'allocation du tableau jic dans init_outshort",err |
---|
206 | stop 4 |
---|
207 | end if |
---|
208 | end if |
---|
209 | if (.not.allocated(ic_name)) then |
---|
210 | allocate(ic_name(nbic),stat=err) |
---|
211 | if (err/=0) then |
---|
212 | print *,"erreur a l'allocation du tableau ic_name dans init_outshort",err |
---|
213 | stop 4 |
---|
214 | end if |
---|
215 | end if |
---|
216 | if (.not.allocated(ic_lon)) then |
---|
217 | allocate(ic_lon(nbic),stat=err) |
---|
218 | if (err/=0) then |
---|
219 | print *,"erreur a l'allocation du tableau ic_lon dans init_outshort",err |
---|
220 | stop 4 |
---|
221 | end if |
---|
222 | end if |
---|
223 | if (.not.allocated(ic_lat)) then |
---|
224 | allocate(ic_lat(nbic),stat=err) |
---|
225 | if (err/=0) then |
---|
226 | print *,"erreur a l'allocation du tableau ic_lat dans init_outshort",err |
---|
227 | stop 4 |
---|
228 | end if |
---|
229 | end if |
---|
230 | if (.not.allocated(var_icoutput)) then |
---|
231 | allocate(var_icoutput(nbic,nvaric),stat=err) |
---|
232 | if (err/=0) then |
---|
233 | print *,"erreur a l'allocation du tableau var_icoutput dans init_outshort",err |
---|
234 | stop 4 |
---|
235 | end if |
---|
236 | end if |
---|
237 | ! lecture des infos sur les variables : |
---|
238 | do k=1,nbic |
---|
239 | read(568,'(a100)') ic_name(k) |
---|
240 | !ic_name(k)=trim(ic_name(k)) |
---|
241 | read(568,*) ic_lon(k) |
---|
242 | read(568,*) ic_lat(k) |
---|
243 | read(568,*) |
---|
244 | enddo |
---|
245 | close(568) |
---|
246 | ! recherche du point i,j correspondant aux coordonnées lon,lat de chaque site de forage |
---|
247 | do k=1,nbic |
---|
248 | dist_lonlat = 1.e7 ! initialisation |
---|
249 | dist_lonlat_new = 1.e6 ! initialisation |
---|
250 | ! on part du centre de la grille |
---|
251 | ! iic(k)=int(nx/2) |
---|
252 | ! jic(k)=int(ny/2) |
---|
253 | iic(k)=6 |
---|
254 | jic(k)=6 |
---|
255 | nloop=0 |
---|
256 | do while (dist_lonlat_new .LT. dist_lonlat .AND. nloop .LT. 10000) |
---|
257 | nloop=nloop+1 |
---|
258 | dist_lonlat = dist_lonlat_new |
---|
259 | ii = iic(k) |
---|
260 | jj = jic(k) |
---|
261 | do j=max(jj-5,1),min(jj+5,ny) |
---|
262 | do i=max(ii-5,1),min(ii+5,nx) |
---|
263 | distance = abs(xlong(i,j) - ic_lon(k)) + abs(ylat(i,j) - ic_lat(k)) |
---|
264 | if (distance .LT. dist_lonlat) then |
---|
265 | dist_lonlat_new = distance |
---|
266 | iic(k) = i |
---|
267 | jic(k) = j |
---|
268 | endif |
---|
269 | enddo |
---|
270 | enddo |
---|
271 | enddo |
---|
272 | !print*,'init_outshort : ice core ', trim(ic_name(k)), nloop |
---|
273 | !print*,'init_outshort lon,lat', ic_lon(k),ic_lat(k) |
---|
274 | !print*,'init_outshort i,j',iic(k),jic(k) |
---|
275 | !print*,'init_outshort lon lat GRISLI',xlong(iic(k),jic(k)),ylat(iic(k),jic(k)) |
---|
276 | enddo |
---|
277 | |
---|
278 | |
---|
279 | |
---|
280 | ! creation du fichier Netcdf : |
---|
281 | status=nf90_create(path = trim(dirnameout)//'icecore'//runname//'.nc', cmode = nf90_clobber, ncid = ncidic) |
---|
282 | if (status /= nf90_noerr) call handle_err(status) |
---|
283 | |
---|
284 | ! definition des dimension : |
---|
285 | status=nf90_def_dim(ncidic, name="time", len=NF90_UNLIMITED, dimid=timeDimIDic) |
---|
286 | if (status /= nf90_noerr) call handle_err(status) |
---|
287 | status=nf90_def_var(ncidic, name="time", xtype=nf90_float, dimids=(/ timeDimIDic/), varid=timeVarIDic) |
---|
288 | if (status /= nf90_noerr) call handle_err(status) |
---|
289 | status=nf90_put_att(ncidic, timeVarIDic, "standard_name", "time") |
---|
290 | if (status /= nf90_noerr) call handle_err(status) |
---|
291 | status=nf90_put_att(ncidic, timeVarIDic,"units", "years since 0000-01-01 00:00:00") |
---|
292 | if (status /= nf90_noerr) call handle_err(status) |
---|
293 | status=nf90_put_att(ncidic, timeVarIDic,"calendar", "360_day") |
---|
294 | if (status /= nf90_noerr) call handle_err(status) |
---|
295 | |
---|
296 | ! dimension "station" Sites de forages |
---|
297 | status=nf90_def_dim(ncidic, name="station", len=nbic, dimid=stationDimIDic) |
---|
298 | if (status /= nf90_noerr) call handle_err(status) |
---|
299 | status=nf90_def_dim(ncidic, name="name_strlen", len=len(ic_name(1)), dimid=station_strlenDimIDic) |
---|
300 | if (status /= nf90_noerr) call handle_err(status) |
---|
301 | status=nf90_def_var(ncidic, name="station_name", xtype=nf90_char, dimids=(/ station_strlenDimIDic, stationDimIDic /), varid=stationVarIDic) |
---|
302 | if (status /= nf90_noerr) call handle_err(status) |
---|
303 | status=nf90_put_att(ncidic, stationVarIDic, "standard_name", "station_name") |
---|
304 | if (status /= nf90_noerr) call handle_err(status) |
---|
305 | status=nf90_put_att(ncidic, stationVarIDic, "long_name", "station name") |
---|
306 | if (status /= nf90_noerr) call handle_err(status) |
---|
307 | status=nf90_def_var(ncidic, name="lon", xtype=nf90_float, dimids=(/ stationDimIDic/), varid=lonVarIDic) |
---|
308 | if (status /= nf90_noerr) call handle_err(status) |
---|
309 | status=nf90_put_att(ncidic, lonVarIDic, "standard_name", "longitude") |
---|
310 | if (status /= nf90_noerr) call handle_err(status) |
---|
311 | status=nf90_put_att(ncidic, lonVarIDic, "long_name", "station longitude") |
---|
312 | if (status /= nf90_noerr) call handle_err(status) |
---|
313 | status=nf90_def_var(ncidic, name="lat", xtype=nf90_float, dimids=(/ stationDimIDic/), varid=latVarIDic) |
---|
314 | if (status /= nf90_noerr) call handle_err(status) |
---|
315 | status=nf90_put_att(ncidic, latVarIDic, "standard_name", "latitude") |
---|
316 | if (status /= nf90_noerr) call handle_err(status) |
---|
317 | status=nf90_put_att(ncidic, latVarIDic, "long_name", "station latitude") |
---|
318 | if (status /= nf90_noerr) call handle_err(status) |
---|
319 | |
---|
320 | |
---|
321 | ! infos sur les variables : |
---|
322 | open(568,file=trim(dirsource)//'/Fichiers-parametres/icecore-var.dat',status='old') |
---|
323 | ! lecture en-tete |
---|
324 | read(568,*) |
---|
325 | read(568,*) |
---|
326 | read(568,*) |
---|
327 | read(568,*) |
---|
328 | read(568,*) |
---|
329 | ! lecture des infos sur les variables : |
---|
330 | do l=1,nvaric |
---|
331 | read(568,'(a100)') namevaric(l) |
---|
332 | read(568,'(a100)') standard_nameic(l) |
---|
333 | read(568,'(a100)') long_nameic(l) |
---|
334 | read(568,'(a100)') unitsic(l) |
---|
335 | read(568,*) |
---|
336 | enddo |
---|
337 | close(568) |
---|
338 | |
---|
339 | ! definition des variables de sortie : |
---|
340 | do l=1,nvaric ! boucle sur le nbr de variable a definir |
---|
341 | status=nf90_def_var(ncidic, name=trim(namevaric(l)), xtype=nf90_float, dimids= & |
---|
342 | (/ stationDimIDic, timeDimIDic /), varid=varIDic(l)) |
---|
343 | if (status /= nf90_noerr) call handle_err(status) |
---|
344 | status=nf90_put_att(ncidic, varIDic(l), "standard_name", trim(standard_nameic(l))) |
---|
345 | if (status /= nf90_noerr) call handle_err(status) |
---|
346 | status=nf90_put_att(ncidic, varIDic(l), "long_name", trim(long_nameic(l))) |
---|
347 | if (status /= nf90_noerr) call handle_err(status) |
---|
348 | status=nf90_put_att(ncidic, varIDic(l), "units", trim(unitsic(l))) |
---|
349 | if (status /= nf90_noerr) call handle_err(status) |
---|
350 | status=nf90_put_att(ncidic, varIDic(l), "coordinates", "lon lat") |
---|
351 | if (status /= nf90_noerr) call handle_err(status) |
---|
352 | enddo |
---|
353 | |
---|
354 | ! fin de la definition du fichier : |
---|
355 | status=nf90_enddef(ncidic) |
---|
356 | if (status /= nf90_noerr) call handle_err(status) |
---|
357 | |
---|
358 | nbtimeoutic = 0 ! initialisation compteur sorties axe time |
---|
359 | ! ecriture des noms des sites de forage dans la variable site : |
---|
360 | ! do k=1,nbic |
---|
361 | ! status=nf90_put_var(ncidic, stationVarIDic, trim(ic_name(k)),start=(/k/),count=(/1/)) |
---|
362 | ! if (status /= nf90_noerr) call handle_err(status) |
---|
363 | ! enddo |
---|
364 | !print*,'ic_name',ic_name |
---|
365 | status=nf90_put_var(ncidic, stationVarIDic, ic_name) !,start=(/1/),count=(/nbic/)) |
---|
366 | if (status /= nf90_noerr) call handle_err(status) |
---|
367 | |
---|
368 | status=nf90_put_var(ncidic, lonVarIDic, ic_lon,start=(/1/),count=(/nbic/)) |
---|
369 | if (status /= nf90_noerr) call handle_err(status) |
---|
370 | |
---|
371 | status=nf90_put_var(ncidic, latVarIDic, ic_lat,start=(/1/),count=(/nbic/)) |
---|
372 | if (status /= nf90_noerr) call handle_err(status) |
---|
373 | ! print*, 'debug 2 init_outshort' |
---|
374 | endif |
---|
375 | |
---|
376 | |
---|
377 | |
---|
378 | end subroutine init_outshort |
---|
379 | |
---|
380 | |
---|
381 | |
---|
382 | !_________________________________________________________________________ |
---|
383 | subroutine shortoutput |
---|
384 | |
---|
385 | ! 1_initialization |
---|
386 | !------------------ |
---|
387 | real :: smax |
---|
388 | |
---|
389 | vol=0. |
---|
390 | np=0 |
---|
391 | nflot=0 |
---|
392 | hmax=0. |
---|
393 | smax=0. |
---|
394 | bmean=0. |
---|
395 | accmean=0. |
---|
396 | ablmean=0. |
---|
397 | calvmean=0. |
---|
398 | ablbordmean=0. |
---|
399 | bmeltmean=0. |
---|
400 | tbmean=0. |
---|
401 | tbdotmean=0. |
---|
402 | vsmean=0. |
---|
403 | ! vsdotmean=0. |
---|
404 | ! uzsmean=0. |
---|
405 | uzsdotmean=0. |
---|
406 | uzkmean=0. |
---|
407 | hdotmean=0. |
---|
408 | bdotmean=0. |
---|
409 | volf=0. |
---|
410 | lim=0. |
---|
411 | tendacabf=0. |
---|
412 | iareag=0. |
---|
413 | iareaf=0. |
---|
414 | tendlicalvf=0. |
---|
415 | tendligroundf=0. |
---|
416 | tendlibmassbf=0. |
---|
417 | tendlibmassbffl=0. |
---|
418 | tendlifmassbf=0. |
---|
419 | ! 2_preparing outputs |
---|
420 | !-------------------- |
---|
421 | do j=1,ny |
---|
422 | do i=1,nx |
---|
423 | if (ice(i,j).eq.1) then ! point englace |
---|
424 | if (.not.flot(i,j)) then ! point pose |
---|
425 | np=np+1 |
---|
426 | vol=vol+h(i,j) |
---|
427 | iareag=iareag+1.*corrsurf(i,j) ! surface englacee posee |
---|
428 | |
---|
429 | ! calcul de la hauteur au dessus de la flottaison |
---|
430 | if (sealevel_2d(i,j)-B(i,j).le.0.) then ! socle au dessus du niveau des mers |
---|
431 | volf=volf+h(i,j)*corrsurf(i,j) ! volume au-dessus de la flottaison |
---|
432 | else |
---|
433 | volf=volf+(h(i,j)-row/ro*(sealevel_2d(i,j)-b(i,j)))*corrsurf(i,j) ! volume au-dessus de la flottaison |
---|
434 | endif |
---|
435 | |
---|
436 | if (h(i,j).gt.hmax) hmax=h(i,j) |
---|
437 | if (s(i,j).gt.smax) smax=s(i,j) |
---|
438 | bmean=bm(i,j)+bmean |
---|
439 | accmean=acc(i,j)+accmean |
---|
440 | tbmean=tbmean+t(i,j,nz) |
---|
441 | tbdotmean=tbdotmean+tbdot(i,j) |
---|
442 | vsmean=vsmean+sqrt(ux(i,j,1)**2+uy(i,j,1)**2) |
---|
443 | ! vsdotmean=vsdotmean+vsdot(i,j) |
---|
444 | ! uzsmean=uzsmean+uz(i,j,1) |
---|
445 | uzsdotmean=uzsdotmean+uzsdot(i,j) |
---|
446 | uzkmean=uzkmean+uzk(i,j) |
---|
447 | hdotmean=hdotmean+abs(hdot(i,j)) |
---|
448 | bdotmean=bdotmean+abs(bdot(i,j)) |
---|
449 | bmeltmean=bmeltmean+bmelt(i,j) |
---|
450 | else ! point flottant |
---|
451 | iareaf=iareaf+1.*corrsurf(i,j) ! surface flottante |
---|
452 | tendlibmassbffl=tendlibmassbffl-bmelt_dtt(i,j)*corrsurf(i,j)/dtt_flux ! fonte basale sur zone flottante |
---|
453 | endif |
---|
454 | lim=lim+h(i,j)*corrsurf(i,j) ! volume total de glace |
---|
455 | tendacabf=tendacabf+bm_dtt(i,j)*corrsurf(i,j)/dtt_flux ! smb surface |
---|
456 | tendlibmassbf=tendlibmassbf-bmelt_dtt(i,j)*corrsurf(i,j)/dtt_flux ! fonte basale |
---|
457 | endif |
---|
458 | tendlicalvf=tendlicalvf+calv_dtt(i,j)*corrsurf(i,j)/dtt_flux ! calving |
---|
459 | tendlifmassbf=tendlifmassbf+(calv_dtt(i,j)-ablbord_dtt(i,j))*corrsurf(i,j)/dtt_flux ! calving + ablbord (ablbord est positif) |
---|
460 | tendligroundf=tendligroundf+grline_dtt(i,j)*corrsurf(i,j)/dtt_flux ! flux a la grounding line |
---|
461 | end do |
---|
462 | end do |
---|
463 | |
---|
464 | |
---|
465 | if (np.ne.0) then |
---|
466 | hmean=vol/np |
---|
467 | vol=vol*dx*dy |
---|
468 | volf=volf*dx*dy*ice_density |
---|
469 | bmean=bmean/np |
---|
470 | accmean=accmean/np |
---|
471 | ablmean=bmean-accmean |
---|
472 | calvmean=calvmean/np |
---|
473 | bmeltmean=bmeltmean/np |
---|
474 | ablbordmean=ablbordmean/np |
---|
475 | tbmean=tbmean/np |
---|
476 | tbdotmean=tbdotmean/np |
---|
477 | vsmean=vsmean/np |
---|
478 | ! vsdotmean=vsdotmean/np |
---|
479 | ! uzsmean=uzsmean/np |
---|
480 | uzsdotmean=uzsdotmean/np |
---|
481 | uzkmean=uzkmean/np |
---|
482 | hdotmean=hdotmean/np |
---|
483 | endif |
---|
484 | lim=lim*dx*dy*ice_density |
---|
485 | iareag=iareag*dx*dy |
---|
486 | iareaf=iareaf*dx*dy |
---|
487 | tendacabf=tendacabf*dx*dy*ice_density/secyear |
---|
488 | tendlibmassbf=tendlibmassbf*dx*dy*ice_density/secyear |
---|
489 | tendlibmassbffl=tendlibmassbffl*dx*dy*ice_density/secyear |
---|
490 | tendlicalvf=tendlicalvf*dx*dy*ice_density/secyear |
---|
491 | tendlifmassbf=tendlifmassbf*dx*dy*ice_density/secyear |
---|
492 | tendligroundf=tendligroundf*ice_density/secyear |
---|
493 | |
---|
494 | |
---|
495 | bdotmean=bdotmean/nx/ny |
---|
496 | |
---|
497 | |
---|
498 | ! 2_writing outputs |
---|
499 | !------------------ |
---|
500 | ! **** short display **** |
---|
501 | |
---|
502 | write(num_ritz,904) nt,time,tafor,sealevel,vol,volf,np, & |
---|
503 | nint(hmean),nint(smax), & |
---|
504 | bmean,tbmean,nint(vsmean), & |
---|
505 | tbdotmean,hdotmean,dt,accmean, & |
---|
506 | diff_H,water_bilan,sum(calv_dtt(:,:))/dtt_flux, & |
---|
507 | sum(ablbord_dtt(:,:))/dtt_flux, & |
---|
508 | sum(Bm_dtt(:,:))/dtt_flux,sum(bmelt_dtt(:,:))/dtt_flux |
---|
509 | |
---|
510 | |
---|
511 | 903 format(i8,1x,f0.2,1x,f0.4,1x,f0.2,1x,2(e10.5,1x),i6,1x,i4,1x,i5,1x, & |
---|
512 | f0.4,1x,f0.3,1x,i3,4(1x,e8.2),1x,f0.4,1x,f0.4) |
---|
513 | 904 format(i8,1x,f0.2,4(1x,e15.8),3(1x,i8),2(1x,e15.8),1x,i8,10(1x,e15.8)) |
---|
514 | !940 format('%%%% ',a,' time=',f8.0,' %%%%') |
---|
515 | |
---|
516 | |
---|
517 | if (ncshortout.eq.1) then ! ecriture netcdf |
---|
518 | ! Total ice mass |
---|
519 | var_shortoutput(1)=lim |
---|
520 | ! Mass above floatation |
---|
521 | var_shortoutput(2)=volf |
---|
522 | ! Grounded ice area |
---|
523 | var_shortoutput(3)=iareag |
---|
524 | ! Floating ice area |
---|
525 | var_shortoutput(4)=iareaf |
---|
526 | ! Total SMB flux |
---|
527 | var_shortoutput(5)=tendacabf |
---|
528 | ! Total Basal mass balance flux |
---|
529 | var_shortoutput(6)=tendlibmassbf |
---|
530 | ! Total Basal mass balance flux beneath floating ice |
---|
531 | var_shortoutput(7)=tendlibmassbffl |
---|
532 | ! Total calving flux |
---|
533 | var_shortoutput(8)=tendlicalvf |
---|
534 | ! Total calving and ice front melting flux |
---|
535 | var_shortoutput(9)=tendlifmassbf |
---|
536 | ! Total grounding line flux |
---|
537 | var_shortoutput(10)=tendligroundf |
---|
538 | |
---|
539 | nbtimeout=nbtimeout+1 |
---|
540 | |
---|
541 | status=nf90_put_var(ncid, timeVarID, time, start=(/nbtimeout/)) |
---|
542 | if (status /= nf90_noerr) call handle_err(status) |
---|
543 | |
---|
544 | do k=1,nvar ! boucle sur le nbr de variable a ecrire |
---|
545 | status=nf90_put_var(ncid, varID(k), var_shortoutput(k),start=(/nbtimeout/)) |
---|
546 | if (status /= nf90_noerr) call handle_err(status) |
---|
547 | enddo |
---|
548 | status=nf90_sync(ncid) |
---|
549 | if (status /= nf90_noerr) call handle_err(status) |
---|
550 | endif |
---|
551 | |
---|
552 | |
---|
553 | ! ecriture netcdf ice core output (tous les 50 ans) |
---|
554 | if (ncicecoreout.eq.1 .and. (mod(abs(time),dticecoreout).lt.dtmin)) then ! ecriture netcdf |
---|
555 | nbtimeoutic=nbtimeoutic+1 |
---|
556 | do k=1,nbic ! boucle sur les stations |
---|
557 | ! Surface elevation |
---|
558 | var_icoutput(k,1)=S(iic(k),jic(k)) |
---|
559 | ! Ice thickness |
---|
560 | var_icoutput(k,2)=H(iic(k),jic(k)) |
---|
561 | ! Bedrock elevation |
---|
562 | var_icoutput(k,3)=Bsoc(iic(k),jic(k)) |
---|
563 | ! Surface Mass Balance |
---|
564 | var_icoutput(k,4)=Bm(iic(k),jic(k)) |
---|
565 | ! Basal melt |
---|
566 | var_icoutput(k,5)=Bmelt(iic(k),jic(k)) |
---|
567 | ! Surface velocity in x |
---|
568 | var_icoutput(k,6)=Uxbar(iic(k),jic(k)) |
---|
569 | ! Surface velocity in y |
---|
570 | var_icoutput(k,7)=Uybar(iic(k),jic(k)) |
---|
571 | ! Velocity magnitude |
---|
572 | var_icoutput(k,8)=sqrt(Uxbar(iic(k),jic(k))**2 + Uybar(iic(k),jic(k))**2) |
---|
573 | ! Surface Temperature |
---|
574 | var_icoutput(k,9)=Ts(iic(k),jic(k)) |
---|
575 | ! Accumulation |
---|
576 | var_icoutput(k,10)=Acc(iic(k),jic(k)) |
---|
577 | enddo |
---|
578 | ! ecriture de time |
---|
579 | status=nf90_put_var(ncidic, timeVarIDic, time, start=(/nbtimeoutic/)) |
---|
580 | if (status /= nf90_noerr) call handle_err(status) |
---|
581 | |
---|
582 | do l=1,nvaric ! boucle sur le nbr de variables a ecrire |
---|
583 | status=nf90_put_var(ncidic, varIDic(l), var_icoutput(:,l),start=(/1,nbtimeoutic/)) |
---|
584 | if (status /= nf90_noerr) call handle_err(status) |
---|
585 | enddo |
---|
586 | status=nf90_sync(ncidic) |
---|
587 | if (status /= nf90_noerr) call handle_err(status) |
---|
588 | endif |
---|
589 | |
---|
590 | |
---|
591 | |
---|
592 | end subroutine shortoutput |
---|
593 | |
---|
594 | subroutine handle_err(status) |
---|
595 | integer, intent(in) :: status |
---|
596 | if (status /= nf90_noerr) then |
---|
597 | print*,trim(nf90_strerror(status)) |
---|
598 | stop "stopped" |
---|
599 | end if |
---|
600 | end subroutine handle_err |
---|
601 | |
---|
602 | end module output_antarcti_mod |
---|