source: trunk/SOURCES/output_global_mod.f90 @ 4

Last change on this file since 4 was 4, checked in by dumas, 10 years ago

initial import GRISLI trunk

File size: 4.3 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!! @note   - use param_phy_mod
12!<
13
14
15module  output_antarcti_mod
16
17USE module3D_phy
18USE param_phy_mod
19
20implicit none
21     
22!real ::  vol ; !integer :: np
23real ::  bmean                        !<
24real ::  accmean                      !< accumulation moyenne
25real ::  ablmean                      !< ablation moyenne
26real ::  calvmean                     !< moyenne calving
27real ::  ablbordmean                  !<
28real ::  bmeltmean                    !< moyenne fusion basale
29real ::  tbmean                       !< temperature basale moyenne
30real ::  tbdotmean                    !< moyenne variation / temps temperature basale
31real ::  vsmean                       !< vitesse de surface moyenne
32!real ::  vsdotmean                   ! moyenne variation / temps vitesse de surface
33!real ::  uzsmean   !!!! utilise ?    ! vitesse verticale de surface moyenne
34real ::  uzsdotmean                   !< moyenne variation / temps vitesse verticale de surface
35real ::  uzkmean                      !< moyenne vitesse verticale de surface
36real ::  hdotmean                     !< moyenne derivee / temps de h
37real ::  bdotmean                     !< moyenne bedrock derive / temps
38real :: volf                          !< volume au dessus de la flottaison
39
40
41CONTAINS
42
43subroutine init_outshort
44
45!ndisp sorite courte tous les ndisp
46NDISP=100
47end subroutine init_outshort
48
49
50
51!_________________________________________________________________________
52subroutine shortoutput
53
54! 1_initialization
55!------------------
56real ::  smax
57      vol=0. 
58      np=0
59      hmax=0. 
60      smax=0.
61      bmean=0. 
62      accmean=0. 
63      ablmean=0. 
64      calvmean=0. 
65      ablbordmean=0.
66      bmeltmean=0.
67      tbmean=0.
68      tbdotmean=0.
69      vsmean=0.
70!      vsdotmean=0.
71!      uzsmean=0.
72      uzsdotmean=0.
73      uzkmean=0.
74      hdotmean=0.
75      bdotmean=0.
76      volf=0.
77 
78! 2_preparing outputs
79!--------------------     
80    do i=1,nx 
81      do j=1,ny
82        if (.not.flot(i,j)) then 
83!       if (h(i,j).gt.1.) then
84          np=np+1
85          vol=vol+h(i,j) 
86
87!        calcul de la hauteur au dessus de la flottaison
88         if (sealevel-B(i,j).le.0.) then    ! socle au dessus du niveau des mers
89               volf=volf+h(i,j)
90         else
91            volf=volf+h(i,j)-row/ro*(sealevel-b(i,j))
92         endif
93
94
95          if (h(i,j).gt.hmax) hmax=h(i,j) 
96          if (s(i,j).gt.smax) smax=s(i,j) 
97          bmean=bm(i,j)+bmean 
98          accmean=acc(i,j)+accmean
99          tbmean=tbmean+t(i,j,nz)
100          tbdotmean=tbdotmean+tbdot(i,j)
101          vsmean=vsmean+sqrt(ux(i,j,1)**2+uy(i,j,1)**2)
102!          vsdotmean=vsdotmean+vsdot(i,j)
103!          uzsmean=uzsmean+uz(i,j,1)
104          uzsdotmean=uzsdotmean+uzsdot(i,j)
105          uzkmean=uzkmean+uzk(i,j)
106          hdotmean=hdotmean+abs(hdot(i,j)) 
107          bdotmean=bdotmean+abs(bdot(i,j)) 
108          bmeltmean=bmeltmean+bmelt(i,j)
109        endif
110        calvmean=calvmean+calv(i,j) 
111        ablbordmean=ablbordmean+ablbord(i,j)
112      end do
113      end do
114
115
116      if (np.ne.0) then
117        hmean=vol/np 
118        vol=vol*dx*dy 
119        volf=volf*dx*dy
120        bmean=bmean/np 
121        accmean=accmean/np 
122        ablmean=bmean-accmean 
123        calvmean=calvmean/np 
124        bmeltmean=bmeltmean/np
125        ablbordmean=ablbordmean/np
126        tbmean=tbmean/np
127        tbdotmean=tbdotmean/np
128        vsmean=vsmean/np
129!        vsdotmean=vsdotmean/np
130!        uzsmean=uzsmean/np
131        uzsdotmean=uzsdotmean/np
132        uzkmean=uzkmean/np
133        hdotmean=hdotmean/np 
134      endif
135
136      bdotmean=bdotmean/nx/ny 
137
138
139! 2_writing outputs
140!------------------   
141!     **** short display ****
142
143        write(num_ritz,903) nt,time,tafor,sealevel,vol,volf,np, &
144          nint(hmean),nint(smax),                    &
145          bmean,tbmean,nint(vsmean),                 &
146!         tbdotmean,vsdotmean,hdotmean,bdotmean,    &
147!          tbdotmean,hdotmean,dt,bmeltmean,accmean
148
149          tbdotmean,hdotmean,dt,bmelt(3,3),accmean 
150
151
152903   format(i8,1x,f0.2,1x,f0.4,1x,f0.2,1x,2(e15.8,1x),i0,1x,i0,1x,i0,1x, &
153             f0.4,1x,f0.3,1x,i0,4(1x,e11.4),1x,f0.4,1x,f0.4) 
154!940   format('%%%% ',a,'   time=',f8.0,' %%%%')
155
156end subroutine shortoutput
157end module  output_antarcti_mod
Note: See TracBrowser for help on using the repository browser.