source: branches/GRISLIv3/SOURCES/Alps_files/output_alps_mod.f90 @ 465

Last change on this file since 465 was 465, checked in by aquiquet, 4 months ago

Cleaning branch: start cleaning module3D

File size: 11.2 KB
Line 
1!> \file output_alps.f90
2!! Module de
3!<
4
5!> \namespace output_alps_mod
6!! Module de
7!! \author ...
8!! \date ...
9!! @note Used module
10!! @note   - use module3D_phy
11!<
12
13module  output_alps_mod
14
15USE module3D_phy
16use bilan_eau_mod
17use netcdf
18use io_netcdf_grisli
19use bilan_flux_mod
20
21implicit none
22     
23!real ::  vol ; !integer :: np
24integer :: nflot                      !< nbr de point flottant
25real ::  bmean                        !<
26real ::  accmean                      !< accumulation moyenne
27real ::  ablmean                      !< ablation moyenne
28real ::  calvmean                     !< moyenne calving
29real ::  ablbordmean                  !<
30real ::  bmeltmean                    !< moyenne fusion basale
31real ::  tbmean                       !< temperature basale moyenne
32real ::  tbdotmean                    !< moyenne variation / temps temperature basale
33real ::  vsmean                       !< vitesse de surface moyenne
34!real ::  vsdotmean                    !< moyenne variation / temps vitesse de surface
35!real ::  uzsmean   !!!! utilise ?     !< vitesse verticale de surface moyenne
36real ::  uzsdotmean                   !< moyenne variation / temps vitesse verticale de surface
37real ::  uzkmean                      !< moyenne vitesse verticale de surface
38real ::  hdotmean                     !< moyenne derivee / temps de h
39real ::  bdotmean                     !< moyenne bedrock derive / temps
40real ::  volf                         !< volume au dessus de la flottaison
41
42real :: lim                           !< Total ice mass
43real :: iareag                        !< surface posee
44real :: iareaf                        !< surface flottante
45real :: tendacabf                     !< Total SMB flux
46real :: tendlibmassbf                 !< Total Basal mass balance flux
47real :: tendlibmassbffl               !< Total Basal mass balance flux on floating ice
48real :: tendlicalvf                   !< Total calving flux
49real :: tendlifmassbf                 !< Total calving and ice front melting flux
50real :: tendligroundf                 !< Total grounding line flux
51
52real,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
56integer,parameter :: ncshortout=1 ! 1 sorties netcdf short initMIP
57integer,parameter :: nvar=10 ! nombre de variables dans le fichier de sortie temporel netcdf
58integer :: ncid
59integer :: status
60integer :: timeDimID
61integer :: timeVarID
62integer,dimension(nvar) :: varID
63integer :: nbtimeout  ! index time output
64
65real,dimension(nvar) :: var_shortoutput
66
67
68CONTAINS
69
70subroutine init_outshort
71
72double precision,dimension(:,:),pointer      :: tab               !< tableau 2d real pointer
73character(len=100),dimension(nvar) :: namevar ! name, standard_name, long_name, units
74character(len=100),dimension(nvar) :: standard_name
75character(len=100),dimension(nvar) :: long_name
76character(len=100),dimension(nvar) :: units
77
78
79
80if (ncshortout.eq.1) then  ! ecriture netcdf
81
82! lecture du fichier avec les corrections de surface
83  if (geoplace.eq.'ant16km') then
84      call Read_Ncdf_var('z',trim(DIRNAMEINP)//'/corrsurf-initMIP-16km.grd',tab)
85      corrsurf(:,:)=tab(:,:)
86  else
87      corrsurf(:,:)= 1.
88  end if
89
90  open(568,file=trim(dirsource)//'/Fichiers-parametres/short-ISMIPnc.dat',status='old')
91! lecture en-tete
92  read(568,*)
93  read(568,*)
94  read(568,*)
95  read(568,*)
96  read(568,*)
97! lecture des infos sur les variables :
98  do k=1,nvar
99    read(568,'(a100)') namevar(k)
100    read(568,'(a100)') standard_name(k)
101    read(568,'(a100)') long_name(k)
102    read(568,'(a100)') units(k)
103    read(568,*)
104  enddo
105  close(568)
106! Fichier Netcdf initMIP
107! creation du fichier Netcdf :
108  status=nf90_create(path = trim(dirnameout)//'short'//runname//'.nc', cmode = nf90_clobber, ncid = ncid)
109  if (status /= nf90_noerr) call handle_err(status)
110
111! definition des dimension :
112  status=nf90_def_dim(ncid, name="time", len=NF90_UNLIMITED, dimid=timeDimID)
113  if (status /= nf90_noerr) call handle_err(status)
114  status=nf90_def_var(ncid, name="time", xtype=nf90_float, dimids=(/ timeDimID/), varid=timeVarID) 
115  if (status /= nf90_noerr) call handle_err(status)
116  status=nf90_put_att(ncid, timeVarID, "standard_name", "time")
117  if (status /= nf90_noerr) call handle_err(status)
118  status=nf90_put_att(ncid, timeVarID,"units", "years since 1995-01-01 00:00:00")
119  if (status /= nf90_noerr) call handle_err(status)
120  status=nf90_put_att(ncid, timeVarID,"calendar", "360_day")
121  if (status /= nf90_noerr) call handle_err(status)
122 
123! definition des variables de sortie :
124  do k=1,nvar  ! boucle sur le nbr de variable a definir
125    status=nf90_def_var(ncid, name=trim(namevar(k)), xtype=nf90_float, dimids= &
126          (/ timeDimID /), varid=varID(k))
127    if (status /= nf90_noerr) call handle_err(status)
128    status=nf90_put_att(ncid, varID(k), "standard_name", trim(standard_name(k)))
129    if (status /= nf90_noerr) call handle_err(status)
130    status=nf90_put_att(ncid, varID(k), "long_name", trim(long_name(k)))
131    if (status /= nf90_noerr) call handle_err(status)
132    status=nf90_put_att(ncid, varID(k), "units", trim(units(k)))
133    if (status /= nf90_noerr) call handle_err(status)
134  enddo
135
136! fin de la definition du fchier :
137  status=nf90_enddef(ncid)
138  if (status /= nf90_noerr) call handle_err(status)
139  nbtimeout = 0 ! initialisation compteur sorties axe time
140else ! pas de sortie netcdf et sans correction de surface
141  corrsurf(:,:)=1.
142endif
143
144end subroutine init_outshort
145
146
147
148!_________________________________________________________________________
149subroutine shortoutput
150
151! 1_initialization
152!------------------
153real ::  smax
154integer :: np
155
156      vol=0. 
157      np=0
158      nflot=0
159      hmax=0. 
160      smax=0.
161      bmean=0. 
162      accmean=0. 
163      ablmean=0. 
164      calvmean=0. 
165      ablbordmean=0.
166      bmeltmean=0.
167      tbmean=0.
168      tbdotmean=0.
169      vsmean=0.
170!      vsdotmean=0.
171!      uzsmean=0.
172      uzsdotmean=0.
173      uzkmean=0.
174      hdotmean=0.
175      bdotmean=0.
176      volf=0.
177      lim=0.
178      tendacabf=0.
179      iareag=0.
180      iareaf=0.
181      tendlicalvf=0.
182      tendligroundf=0.
183      tendlibmassbf=0.
184      tendlibmassbffl=0.
185      tendlifmassbf=0.
186! 2_preparing outputs
187!--------------------     
188    do j=1,ny
189      do i=1,nx
190        if (ice(i,j).eq.1) then ! point englace
191          if (.not.flot(i,j)) then ! point pose
192            np=np+1
193            vol=vol+h(i,j)
194            iareag=iareag+1.*corrsurf(i,j)                            ! surface englacee posee
195
196!         calcul de la hauteur au dessus de la flottaison
197            if (sealevel_2d(i,j)-B(i,j).le.0.) then    ! socle au dessus du niveau des mers
198              volf=volf+h(i,j)*corrsurf(i,j)                          ! volume au-dessus de la flottaison
199            else
200              volf=volf+(h(i,j)-row/ro*(sealevel_2d(i,j)-b(i,j)))*corrsurf(i,j) ! volume au-dessus de la flottaison
201            endif
202
203            if (h(i,j).gt.hmax) hmax=h(i,j) 
204            if (s(i,j).gt.smax) smax=s(i,j) 
205            bmean=bm(i,j)+bmean 
206            accmean=acc(i,j)+accmean
207            tbmean=tbmean+t(i,j,nz)
208            tbdotmean=tbdotmean+tbdot(i,j)
209            vsmean=vsmean+sqrt(ux(i,j,1)**2+uy(i,j,1)**2)
210!          vsdotmean=vsdotmean+vsdot(i,j)
211!          uzsmean=uzsmean+uz(i,j,1)
212            uzsdotmean=uzsdotmean+uzsdot(i,j)
213            uzkmean=uzkmean+uzk(i,j)
214            hdotmean=hdotmean+abs(hdot(i,j)) 
215            bdotmean=bdotmean+abs(bdot(i,j)) 
216            bmeltmean=bmeltmean+bmelt(i,j)
217          else ! point flottant
218            iareaf=iareaf+1.*corrsurf(i,j)                               ! surface flottante
219            tendlibmassbffl=tendlibmassbffl-bmelt_dtt(i,j)*corrsurf(i,j)/dtt_flux   ! fonte basale sur zone flottante
220          endif
221          lim=lim+h(i,j)*corrsurf(i,j)                                   ! volume total de glace
222          tendacabf=tendacabf+bm_dtt(i,j)*corrsurf(i,j)/dtt_flux              ! smb surface
223          tendlibmassbf=tendlibmassbf-bmelt_dtt(i,j)*corrsurf(i,j)/dtt_flux   ! fonte basale     
224        endif
225        tendlicalvf=tendlicalvf+calv_dtt(i,j)*corrsurf(i,j)/dtt_flux        ! calving
226        tendlifmassbf=tendlifmassbf+(calv_dtt(i,j)-ablbord_dtt(i,j))*corrsurf(i,j)/dtt_flux ! calving + ablbord (ablbord est positif)
227        tendligroundf=tendligroundf+grline_dtt(i,j)*corrsurf(i,j)/dtt_flux  ! flux a la grounding line
228      end do
229    end do
230
231
232      if (np.ne.0) then
233        hmean=vol/np 
234        vol=vol*dx*dy 
235        volf=volf*dx*dy*ice_density
236        bmean=bmean/np 
237        accmean=accmean/np 
238        ablmean=bmean-accmean 
239        calvmean=calvmean/np 
240        bmeltmean=bmeltmean/np
241        ablbordmean=ablbordmean/np
242        tbmean=tbmean/np
243        tbdotmean=tbdotmean/np
244        vsmean=vsmean/np
245!        vsdotmean=vsdotmean/np
246!        uzsmean=uzsmean/np
247        uzsdotmean=uzsdotmean/np
248        uzkmean=uzkmean/np
249        hdotmean=hdotmean/np
250      endif
251      lim=lim*dx*dy*ice_density
252      iareag=iareag*dx*dy 
253      iareaf=iareaf*dx*dy
254      tendacabf=tendacabf*dx*dy*ice_density/secyear
255      tendlibmassbf=tendlibmassbf*dx*dy*ice_density/secyear
256      tendlibmassbffl=tendlibmassbffl*dx*dy*ice_density/secyear
257      tendlicalvf=tendlicalvf*dx*dy*ice_density/secyear
258      tendlifmassbf=tendlifmassbf*dx*dy*ice_density/secyear
259      tendligroundf=tendligroundf*ice_density/secyear
260     
261     
262      bdotmean=bdotmean/nx/ny 
263
264
265! 2_writing outputs
266!------------------   
267!     **** short display ****
268
269        write(num_ritz,904) nt,time,tafor,sealevel,vol,volf,np, &
270          nint(hmean),nint(smax),                    &
271          bmean,tbmean,nint(vsmean),                 &
272                                  tbdotmean,hdotmean,dt,accmean, &
273                                  diff_H,water_bilan,sum(calv_dtt(:,:))/dtt_flux, &
274                                  sum(ablbord_dtt(:,:))/dtt_flux, &
275                                  sum(Bm_dtt(:,:))/dtt_flux,sum(bmelt_dtt(:,:))/dtt_flux
276
277
278903   format(i8,1x,f0.2,1x,f0.4,1x,f0.2,1x,2(e10.5,1x),i6,1x,i4,1x,i5,1x, &
279             f0.4,1x,f0.3,1x,i3,4(1x,e8.2),1x,f0.4,1x,f0.4) 
280904   format(i8,1x,f0.2,4(1x,e15.8),3(1x,i8),2(1x,e15.8),1x,i8,10(1x,e15.8))
281!940   format('%%%% ',a,'   time=',f8.0,' %%%%')
282
283
284if (ncshortout.eq.1) then  ! ecriture netcdf
285  ! Total ice mass
286  var_shortoutput(1)=lim
287  ! Mass above floatation
288  var_shortoutput(2)=volf
289  ! Grounded ice area
290  var_shortoutput(3)=iareag
291  ! Floating ice area
292  var_shortoutput(4)=iareaf
293  ! Total SMB flux
294  var_shortoutput(5)=tendacabf
295  ! Total Basal mass balance flux
296  var_shortoutput(6)=tendlibmassbf
297  ! Total Basal mass balance flux beneath floating ice
298  var_shortoutput(7)=tendlibmassbffl
299  ! Total calving flux
300  var_shortoutput(8)=tendlicalvf
301  ! Total calving and ice front melting flux
302  var_shortoutput(9)=tendlifmassbf
303  ! Total grounding line flux
304  var_shortoutput(10)=tendligroundf
305
306  nbtimeout=nbtimeout+1
307 
308  status=nf90_put_var(ncid, timeVarID, time, start=(/nbtimeout/))
309  if (status /= nf90_noerr) call handle_err(status)
310
311  do k=1,nvar  ! boucle sur le nbr de variable a ecrire
312    status=nf90_put_var(ncid, varID(k), var_shortoutput(k),start=(/nbtimeout/))
313    if (status /= nf90_noerr) call handle_err(status)
314  enddo
315  status=nf90_sync(ncid)
316  if (status /= nf90_noerr) call handle_err(status)
317endif
318
319
320end subroutine shortoutput
321
322subroutine handle_err(status)
323  integer, intent(in) :: status
324  if (status /= nf90_noerr) then
325     print*,trim(nf90_strerror(status))
326     stop "stopped"
327  end if
328end subroutine handle_err
329
330end module  output_alps_mod
Note: See TracBrowser for help on using the repository browser.