source: trunk/SOURCES/Ant40_files/output_anta40_mod-0.4.f90 @ 139

Last change on this file since 139 was 139, checked in by dumas, 7 years ago

small corrections for tendlibmassbf and tendlicalvf unity

File size: 10.9 KB
Line 
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
13module  output_antarcti_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 :: tendlicalvf                   !< Total calving flux
48real :: tendligroundf                 !< Total grounding line flux
49
50real,dimension(nx,ny) :: corrsurf     !< facteur de correction de la surface
51real,parameter :: ice_density=910.    !< densite de la glace pour conversion en masse
52
53! variables netcdf
54integer,parameter :: ncshortout=1 ! 1 sorties netcdf short initMIP
55integer,parameter :: nvar=8 ! nombre de variables dans le fichier de sortie temporel netcdf
56integer :: ncid
57integer :: status
58integer :: timeDimID
59integer :: timeVarID
60integer,dimension(nvar) :: varID
61integer :: nbtimeout  ! index time output
62
63real,dimension(nvar) :: var_shortoutput
64
65
66CONTAINS
67
68subroutine init_outshort
69
70double precision,dimension(:,:),pointer      :: tab               !< tableau 2d real pointer
71character(len=100),dimension(nvar) :: namevar ! name, standard_name, long_name, units
72character(len=100),dimension(nvar) :: standard_name
73character(len=100),dimension(nvar) :: long_name
74character(len=100),dimension(nvar) :: units
75
76!ndisp sorite courte tous les ndisp
77NDISP=100
78
79
80
81if (ncshortout.eq.1) then  ! ecriture netcdf
82
83! lecture du fichier avec les corrections de surface
84  call Read_Ncdf_var('z',trim(DIRNAMEINP)//'/corrsurf-initMIP-16km.grd',tab)
85  corrsurf(:,:)=tab(:,:)
86
87  open(568,file=trim(dirsource)//'/Fichiers-parametres/short-initMIPnc.dat',status='old')
88! lecture en-tete
89  read(568,*)
90  read(568,*)
91  read(568,*)
92  read(568,*)
93  read(568,*)
94! lecture des infos sur les variables :
95  do k=1,nvar
96    read(568,'(a100)') namevar(k)
97    print*,trim(namevar(k))
98    read(568,'(a100)') standard_name(k)
99    print*,trim(standard_name(k))
100    read(568,'(a100)') long_name(k)
101    print*,trim(long_name(k))
102    read(568,'(a100)') units(k)
103    print*,trim(units(k))
104    read(568,*)
105  enddo
106  close(568)
107  print*,'apres lecture fichier descriptions variables'
108! Fichier Netcdf initMIP
109! creation du fichier Netcdf :
110  status=nf90_create(path = 'short'//runname//'.nc', cmode = nf90_clobber, ncid = ncid)
111  if (status /= nf90_noerr) call handle_err(status)
112  print*,'apres nf90_create'
113
114! definition des dimension :
115  status=nf90_def_dim(ncid, name="time", len=NF90_UNLIMITED, dimid=timeDimID)
116  if (status /= nf90_noerr) call handle_err(status)
117  print*,'apres 2'
118  status=nf90_def_var(ncid, name="time", xtype=nf90_float, dimids=(/ timeDimID/), varid=timeVarID) 
119  if (status /= nf90_noerr) call handle_err(status)
120  status=nf90_put_att(ncid, timeVarID, "standard_name", "time")
121  if (status /= nf90_noerr) call handle_err(status)
122  print*,'apres 3'
123  status=nf90_put_att(ncid, timeVarID,"units", "years since 2007-01-01 00:00:00")
124  if (status /= nf90_noerr) call handle_err(status)
125  print*,'apres 4'
126
127! definition des variables de sortie :
128  do k=1,nvar  ! boucle sur le nbr de variable a definir
129    status=nf90_def_var(ncid, name=trim(namevar(k)), xtype=nf90_float, dimids= &
130          (/ timeDimID /), varid=varID(k))
131  print*,'apres 5'
132    if (status /= nf90_noerr) call handle_err(status)
133    status=nf90_put_att(ncid, varID(k), "standard_name", trim(standard_name(k)))
134    print*,'apres 6',trim(standard_name(k))
135    if (status /= nf90_noerr) call handle_err(status)
136    status=nf90_put_att(ncid, varID(k), "long_name", trim(long_name(k)))
137    print*,'apres 7',trim(long_name(k))
138    if (status /= nf90_noerr) call handle_err(status)
139    status=nf90_put_att(ncid, varID(k), "units", trim(units(k)))
140    print*,'apres 8', trim(units(k))
141    if (status /= nf90_noerr) call handle_err(status)
142  enddo
143  print*,'apres boucle'
144
145! fin de la definition du fchier :
146  status=nf90_enddef(ncid)
147  if (status /= nf90_noerr) call handle_err(status)
148  nbtimeout = 0 ! initialisation compteur sorties axe time
149else ! pas de sortie netcdf et sans correction de surface
150  corrsurf(:,:)=1.
151endif
152
153end subroutine init_outshort
154
155
156
157!_________________________________________________________________________
158subroutine shortoutput
159
160! 1_initialization
161!------------------
162real ::  smax
163
164      vol=0. 
165      np=0
166      nflot=0
167      hmax=0. 
168      smax=0.
169      bmean=0. 
170      accmean=0. 
171      ablmean=0. 
172      calvmean=0. 
173      ablbordmean=0.
174      bmeltmean=0.
175      tbmean=0.
176      tbdotmean=0.
177      vsmean=0.
178!      vsdotmean=0.
179!      uzsmean=0.
180      uzsdotmean=0.
181      uzkmean=0.
182      hdotmean=0.
183      bdotmean=0.
184      volf=0.
185      lim=0.
186      tendacabf=0.
187      iareag=0.
188      iareaf=0.
189      tendlicalvf=0.
190      tendligroundf=0.
191      tendlibmassbf=0.
192! 2_preparing outputs
193!--------------------     
194    do j=1,ny
195      do i=1,nx
196        if (ice(i,j).eq.1) then ! point englace
197          if (.not.flot(i,j)) then ! point pose
198            np=np+1
199            vol=vol+h(i,j)
200            iareag=iareag+1.*corrsurf(i,j)                            ! surface englacee posee
201
202!         calcul de la hauteur au dessus de la flottaison
203            if (sealevel-B(i,j).le.0.) then    ! socle au dessus du niveau des mers
204              volf=volf+h(i,j)*corrsurf(i,j)                          ! volume au-dessus de la flottaison
205            else
206              volf=volf+(h(i,j)-row/ro*(sealevel-b(i,j)))*corrsurf(i,j) ! volume au-dessus de la flottaison
207            endif
208
209            if (h(i,j).gt.hmax) hmax=h(i,j) 
210            if (s(i,j).gt.smax) smax=s(i,j) 
211            bmean=bm(i,j)+bmean 
212            accmean=acc(i,j)+accmean
213            tbmean=tbmean+t(i,j,nz)
214            tbdotmean=tbdotmean+tbdot(i,j)
215            vsmean=vsmean+sqrt(ux(i,j,1)**2+uy(i,j,1)**2)
216!          vsdotmean=vsdotmean+vsdot(i,j)
217!          uzsmean=uzsmean+uz(i,j,1)
218            uzsdotmean=uzsdotmean+uzsdot(i,j)
219            uzkmean=uzkmean+uzk(i,j)
220            hdotmean=hdotmean+abs(hdot(i,j)) 
221            bdotmean=bdotmean+abs(bdot(i,j)) 
222            bmeltmean=bmeltmean+bmelt(i,j)
223          else ! point flottant
224            iareaf=iareaf+1.*corrsurf(i,j)                               ! surface flottante
225          endif
226          lim=lim+h(i,j)*corrsurf(i,j)                                   ! volume total de glace
227          tendacabf=tendacabf+bm_dtt(i,j)*corrsurf(i,j)/dtt              ! smb surface
228          tendlibmassbf=tendlibmassbf-bmelt_dtt(i,j)*corrsurf(i,j)/dtt   ! fonte basale     
229        endif
230        tendlicalvf=tendlicalvf+calv_dtt(i,j)*corrsurf(i,j)/dtt        ! calving
231        tendlibmassbf=tendlibmassbf-ablbord_dtt(i,j)*corrsurf(i,j)/dtt ! partie ablbord de la fonte basale
232        tendligroundf=tendligroundf+grline_dtt(i,j)*corrsurf(i,j)/dtt  ! flux a la grounding line
233      end do
234    end do
235
236
237      if (np.ne.0) then
238        hmean=vol/np 
239        vol=vol*dx*dy 
240        volf=volf*dx*dy*ice_density
241        bmean=bmean/np 
242        accmean=accmean/np 
243        ablmean=bmean-accmean 
244        calvmean=calvmean/np 
245        bmeltmean=bmeltmean/np
246        ablbordmean=ablbordmean/np
247        tbmean=tbmean/np
248        tbdotmean=tbdotmean/np
249        vsmean=vsmean/np
250!        vsdotmean=vsdotmean/np
251!        uzsmean=uzsmean/np
252        uzsdotmean=uzsdotmean/np
253        uzkmean=uzkmean/np
254        hdotmean=hdotmean/np
255      endif
256      lim=lim*dx*dy*ice_density
257      iareag=iareag*dx*dy 
258      iareaf=iareaf*dx*dy
259      tendacabf=tendacabf*dx*dy*ice_density/secyear
260      tendlibmassbf=tendlibmassbf*dx*dy*ice_density/secyear
261      tendlicalvf=tendlicalvf*dx*dy*ice_density/secyear
262      tendligroundf=tendligroundf*ice_density/secyear
263     
264     
265      bdotmean=bdotmean/nx/ny 
266
267
268! 2_writing outputs
269!------------------   
270!     **** short display ****
271
272        write(num_ritz,904) nt,time,tafor,sealevel,vol,volf,np, &
273          nint(hmean),nint(smax),                    &
274          bmean,tbmean,nint(vsmean),                 &
275                                  tbdotmean,hdotmean,dt,accmean, &
276                                  diff_H,water_bilan,sum(calv_dtt(:,:))/dtt, &
277                                  sum(ablbord_dtt(:,:))/dtt, &
278                                  sum(Bm_dtt(:,:))/dtt,sum(bmelt_dtt(:,:))/dtt
279
280
281903   format(i8,1x,f0.2,1x,f0.4,1x,f0.2,1x,2(e10.5,1x),i6,1x,i4,1x,i5,1x, &
282             f0.4,1x,f0.3,1x,i3,4(1x,e8.2),1x,f0.4,1x,f0.4) 
283904   format(i8,1x,f0.2,4(1x,e15.8),3(1x,i8),2(1x,e15.8),1x,i8,10(1x,e15.8))
284!940   format('%%%% ',a,'   time=',f8.0,' %%%%')
285
286print*, 'dans shortoutput avant netcdf ncshortout', ncshortout
287
288if (ncshortout.eq.1) then  ! ecriture netcdf
289  ! Total ice mass
290  var_shortoutput(1)=lim
291  ! Mass above floatation
292  var_shortoutput(2)=volf
293  ! Grounded ice area
294  var_shortoutput(3)=iareag
295  ! Floating ice area
296  var_shortoutput(4)=iareaf
297  ! Total SMB flux
298  var_shortoutput(5)=tendacabf
299  ! Total Basal mass balance flux
300  var_shortoutput(6)=tendlibmassbf
301  ! Total calving flux
302  var_shortoutput(7)=tendlicalvf
303  ! Total grounding line flux
304  var_shortoutput(8)=tendligroundf
305
306  nbtimeout=nbtimeout+1
307  print*,'ecriture shortnc', nbtimeout
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_antarcti_mod
Note: See TracBrowser for help on using the repository browser.