1 | subroutine mkstats(ierr) |
---|
2 | |
---|
3 | |
---|
4 | ! |
---|
5 | ! This program writes a stats.nc file from sums and sums of squares |
---|
6 | ! to means and standard deviations and also writes netcdf style |
---|
7 | ! file so that the data can be viewed easily. The data file is |
---|
8 | ! overwritten in place. |
---|
9 | ! SRL 21 May 1996 |
---|
10 | ! Yann W. july 2003 |
---|
11 | |
---|
12 | use statto_mod, only: istime,count |
---|
13 | use mod_phys_lmdz_para, only : is_master |
---|
14 | use mod_grid_phy_lmdz, only : nbp_lon, nbp_lat, nbp_lev, klon_glo |
---|
15 | |
---|
16 | implicit none |
---|
17 | |
---|
18 | include "netcdf.inc" |
---|
19 | |
---|
20 | integer :: ierr,nid,nbvar,i,ndims,lt,nvarid |
---|
21 | integer, dimension(4) :: id,varid,start,size |
---|
22 | integer, dimension(5) :: dimids |
---|
23 | character (len=50) :: name,nameout,units,title |
---|
24 | real,allocatable :: sum3d(:,:,:),square3d(:,:,:),mean3d(:,:,:),sd3d(:,:,:) |
---|
25 | real,allocatable :: sum2d(:,:),square2d(:,:),mean2d(:,:),sd2d(:,:) |
---|
26 | real, dimension(istime) :: time |
---|
27 | real, dimension(nbp_lat) :: lat |
---|
28 | real,allocatable :: lon(:) |
---|
29 | real, dimension(nbp_lev) :: alt |
---|
30 | logical :: lcopy=.true. |
---|
31 | !integer :: latid,lonid,altid,timeid |
---|
32 | integer :: meanid,sdid |
---|
33 | !integer, dimension(4) :: dimout |
---|
34 | |
---|
35 | ! Incrementation of count for the last step, which is not done in wstats |
---|
36 | count(istime)=count(istime)+1 |
---|
37 | |
---|
38 | if (is_master) then |
---|
39 | ! only the master needs do this |
---|
40 | if (klon_glo>1) then |
---|
41 | allocate(lon(nbp_lon+1)) |
---|
42 | allocate(sum3d(nbp_lon+1,nbp_lat,nbp_lev)) |
---|
43 | allocate(square3d(nbp_lon+1,nbp_lat,nbp_lev)) |
---|
44 | allocate(mean3d(nbp_lon+1,nbp_lat,nbp_lev)) |
---|
45 | allocate(sd3d(nbp_lon+1,nbp_lat,nbp_lev)) |
---|
46 | allocate(sum2d(nbp_lon+1,nbp_lat)) |
---|
47 | allocate(square2d(nbp_lon+1,nbp_lat)) |
---|
48 | allocate(mean2d(nbp_lon+1,nbp_lat)) |
---|
49 | allocate(sd2d(nbp_lon+1,nbp_lat)) |
---|
50 | else ! 1D model case |
---|
51 | allocate(lon(1)) |
---|
52 | allocate(sum3d(1,1,nbp_lev)) |
---|
53 | allocate(square3d(1,1,nbp_lev)) |
---|
54 | allocate(mean3d(1,1,nbp_lev)) |
---|
55 | allocate(sd3d(1,1,nbp_lev)) |
---|
56 | allocate(sum2d(1,1)) |
---|
57 | allocate(square2d(1,1)) |
---|
58 | allocate(mean2d(1,1)) |
---|
59 | allocate(sd2d(1,1)) |
---|
60 | endif |
---|
61 | |
---|
62 | ierr = NF_OPEN("stats.nc",NF_WRITE,nid) |
---|
63 | |
---|
64 | ! We catch the id of dimensions of the stats file |
---|
65 | |
---|
66 | ierr= NF_INQ_DIMID(nid,"latitude",id(1)) |
---|
67 | ierr= NF_INQ_DIMID(nid,"longitude",id(2)) |
---|
68 | ierr= NF_INQ_DIMID(nid,"altitude",id(3)) |
---|
69 | ierr= NF_INQ_DIMID(nid,"Time",id(4)) |
---|
70 | |
---|
71 | ierr= NF_INQ_VARID(nid,"latitude",varid(1)) |
---|
72 | ierr= NF_INQ_VARID(nid,"longitude",varid(2)) |
---|
73 | ierr= NF_INQ_VARID(nid,"altitude",varid(3)) |
---|
74 | ierr= NF_INQ_VARID(nid,"Time",varid(4)) |
---|
75 | |
---|
76 | ! Time initialisation |
---|
77 | |
---|
78 | do i=1,istime |
---|
79 | time(i)=i*24./istime |
---|
80 | #ifdef NC_DOUBLE |
---|
81 | ierr= NF_PUT_VARA_DOUBLE(nid,varid(4),i,1,time(i)) |
---|
82 | #else |
---|
83 | ierr= NF_PUT_VARA_REAL(nid,varid(4),i,1,time(i)) |
---|
84 | #endif |
---|
85 | enddo |
---|
86 | |
---|
87 | ! We catche the values of the variables |
---|
88 | |
---|
89 | #ifdef NC_DOUBLE |
---|
90 | ierr = NF_GET_VAR_DOUBLE(nid,varid(1),lat) |
---|
91 | ierr = NF_GET_VAR_DOUBLE(nid,varid(2),lon) |
---|
92 | ierr = NF_GET_VAR_DOUBLE(nid,varid(3),alt) |
---|
93 | #else |
---|
94 | ierr = NF_GET_VAR_REAL(nid,varid(1),lat) |
---|
95 | ierr = NF_GET_VAR_REAL(nid,varid(2),lon) |
---|
96 | ierr = NF_GET_VAR_REAL(nid,varid(3),alt) |
---|
97 | #endif |
---|
98 | |
---|
99 | ! We catch the number of variables in the stats file |
---|
100 | ierr = NF_INQ_NVARS(nid,nbvar) |
---|
101 | |
---|
102 | ! to catche the "real" number of variables (without the "additionnal variables") |
---|
103 | nbvar=(nbvar-4)/2 |
---|
104 | |
---|
105 | do i=1,nbvar |
---|
106 | varid=(i-1)*2+5 |
---|
107 | |
---|
108 | ! What's the variable's name? |
---|
109 | ierr=NF_INQ_VARNAME(nid,varid,name) |
---|
110 | write(*,*) "OK variable ",name |
---|
111 | ! Its units? |
---|
112 | units=" " |
---|
113 | ierr=NF_GET_ATT_TEXT(nid,varid,"units",units) |
---|
114 | ! Its title? |
---|
115 | title=" " |
---|
116 | ierr=NF_GET_ATT_TEXT(nid,varid,"title",title) |
---|
117 | ! Its number of dimensions? |
---|
118 | ierr=NF_INQ_VARNDIMS(nid,varid,ndims) |
---|
119 | ! Its values? |
---|
120 | |
---|
121 | if(ndims==4) then ! lat, lon, alt & time |
---|
122 | |
---|
123 | ! dimout(1)=lonid |
---|
124 | ! dimout(2)=latid |
---|
125 | ! dimout(3)=altid |
---|
126 | ! dimout(4)=timeid |
---|
127 | |
---|
128 | if (klon_glo>1) then ! general case |
---|
129 | size=(/nbp_lon+1,nbp_lat,nbp_lev,1/) |
---|
130 | else ! 1D model |
---|
131 | size=(/1,1,nbp_lev,1/) |
---|
132 | endif |
---|
133 | do lt=1,istime |
---|
134 | start=(/1,1,1,lt/) |
---|
135 | ! Extraction of the "source" variables |
---|
136 | #ifdef NC_DOUBLE |
---|
137 | ierr = NF_GET_VARA_DOUBLE(nid,varid,start,size,sum3d) |
---|
138 | ierr = NF_GET_VARA_DOUBLE(nid,varid+1,start,size,square3d) |
---|
139 | #else |
---|
140 | ierr = NF_GET_VARA_REAL(nid,varid,start,size,sum3d) |
---|
141 | ierr = NF_GET_VARA_REAL(nid,varid+1,start,size,square3d) |
---|
142 | #endif |
---|
143 | ! Calculation of these variables |
---|
144 | mean3d=sum3d/count(lt) |
---|
145 | sd3d=sqrt(max(0.,square3d/count(lt)-mean3d**2)) |
---|
146 | ! Writing of the variables |
---|
147 | #ifdef NC_DOUBLE |
---|
148 | ierr = NF_PUT_VARA_DOUBLE(nid,varid,start,size,mean3d) |
---|
149 | ierr = NF_PUT_VARA_DOUBLE(nid,varid+1,start,size,sd3d) |
---|
150 | #else |
---|
151 | ierr = NF_PUT_VARA_REAL(nid,varid,start,size,mean3d) |
---|
152 | ierr = NF_PUT_VARA_REAL(nid,varid+1,start,size,sd3d) |
---|
153 | #endif |
---|
154 | enddo |
---|
155 | |
---|
156 | else if (ndims.eq.3) then |
---|
157 | |
---|
158 | ! dimout(1)=lonid |
---|
159 | ! dimout(2)=latid |
---|
160 | ! dimout(3)=timeid |
---|
161 | |
---|
162 | if (klon_glo>1) then ! general case |
---|
163 | size=(/nbp_lon+1,nbp_lat,1,0/) |
---|
164 | else |
---|
165 | size=(/1,1,1,0/) |
---|
166 | endif |
---|
167 | do lt=1,istime |
---|
168 | start=(/1,1,lt,0/) |
---|
169 | ! Extraction of the "source" variables |
---|
170 | #ifdef NC_DOUBLE |
---|
171 | ierr = NF_GET_VARA_DOUBLE(nid,varid,start,size,sum2d) |
---|
172 | ierr = NF_GET_VARA_DOUBLE(nid,varid+1,start,size,square2d) |
---|
173 | #else |
---|
174 | ierr = NF_GET_VARA_REAL(nid,varid,start,size,sum2d) |
---|
175 | ierr = NF_GET_VARA_REAL(nid,varid+1,start,size,square2d) |
---|
176 | #endif |
---|
177 | ! Calculation of these variables |
---|
178 | mean2d=sum2d/count(lt) |
---|
179 | sd2d=sqrt(max(0.,square2d/count(lt)-mean2d**2)) |
---|
180 | ! Writing of the variables |
---|
181 | #ifdef NC_DOUBLE |
---|
182 | ierr = NF_PUT_VARA_DOUBLE(nid,varid,start,size,mean2d) |
---|
183 | ierr = NF_PUT_VARA_DOUBLE(nid,varid+1,start,size,sd2d) |
---|
184 | #else |
---|
185 | ierr = NF_PUT_VARA_REAL(nid,varid,start,size,mean2d) |
---|
186 | ierr = NF_PUT_VARA_REAL(nid,varid+1,start,size,sd2d) |
---|
187 | #endif |
---|
188 | enddo |
---|
189 | |
---|
190 | endif |
---|
191 | enddo |
---|
192 | |
---|
193 | ierr= NF_CLOSE(nid) |
---|
194 | |
---|
195 | endif ! of if (is_master) |
---|
196 | |
---|
197 | end |
---|