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 mod_phys_lmdz_para, only : is_master |
---|
13 | |
---|
14 | implicit none |
---|
15 | |
---|
16 | #include "dimensions.h" |
---|
17 | #include "statto.h" |
---|
18 | #include "netcdf.inc" |
---|
19 | |
---|
20 | integer,parameter :: iip1=iim+1 |
---|
21 | integer,parameter :: jjp1=jjm+1 |
---|
22 | integer :: ierr,nid,nbvar,i,ndims,lt,nvarid |
---|
23 | integer, dimension(4) :: id,varid,start,size |
---|
24 | integer, dimension(5) :: dimids |
---|
25 | character (len=50) :: name,nameout,units,title |
---|
26 | real, dimension(iip1,jjp1,llm) :: sum3d,square3d,mean3d,sd3d |
---|
27 | real, dimension(iip1,jjp1) :: sum2d,square2d,mean2d,sd2d |
---|
28 | real, dimension(istime) :: time |
---|
29 | real, dimension(jjp1) :: lat |
---|
30 | real, dimension(iip1) :: lon |
---|
31 | real, dimension(llm) :: alt |
---|
32 | logical :: lcopy=.true. |
---|
33 | !integer :: latid,lonid,altid,timeid |
---|
34 | integer :: meanid,sdid |
---|
35 | !integer, dimension(4) :: dimout |
---|
36 | |
---|
37 | ! Incrementation of count for the last step, which is not done in wstats |
---|
38 | count(istime)=count(istime)+1 |
---|
39 | |
---|
40 | if (is_master) then |
---|
41 | ! only the master needs do this |
---|
42 | |
---|
43 | ierr = NF_OPEN("stats.nc",NF_WRITE,nid) |
---|
44 | |
---|
45 | ! We catch the id of dimensions of the stats file |
---|
46 | |
---|
47 | ierr= NF_INQ_DIMID(nid,"latitude",id(1)) |
---|
48 | ierr= NF_INQ_DIMID(nid,"longitude",id(2)) |
---|
49 | ierr= NF_INQ_DIMID(nid,"altitude",id(3)) |
---|
50 | ierr= NF_INQ_DIMID(nid,"Time",id(4)) |
---|
51 | |
---|
52 | ierr= NF_INQ_VARID(nid,"latitude",varid(1)) |
---|
53 | ierr= NF_INQ_VARID(nid,"longitude",varid(2)) |
---|
54 | ierr= NF_INQ_VARID(nid,"altitude",varid(3)) |
---|
55 | ierr= NF_INQ_VARID(nid,"Time",varid(4)) |
---|
56 | |
---|
57 | ! Time initialisation |
---|
58 | |
---|
59 | do i=1,istime |
---|
60 | time(i)=i*24./istime |
---|
61 | #ifdef NC_DOUBLE |
---|
62 | ierr= NF_PUT_VARA_DOUBLE(nid,varid(4),i,1,time(i)) |
---|
63 | #else |
---|
64 | ierr= NF_PUT_VARA_REAL(nid,varid(4),i,1,time(i)) |
---|
65 | #endif |
---|
66 | enddo |
---|
67 | |
---|
68 | ! We catche the values of the variables |
---|
69 | |
---|
70 | #ifdef NC_DOUBLE |
---|
71 | ierr = NF_GET_VAR_DOUBLE(nid,varid(1),lat) |
---|
72 | ierr = NF_GET_VAR_DOUBLE(nid,varid(2),lon) |
---|
73 | ierr = NF_GET_VAR_DOUBLE(nid,varid(3),alt) |
---|
74 | #else |
---|
75 | ierr = NF_GET_VAR_REAL(nid,varid(1),lat) |
---|
76 | ierr = NF_GET_VAR_REAL(nid,varid(2),lon) |
---|
77 | ierr = NF_GET_VAR_REAL(nid,varid(3),alt) |
---|
78 | #endif |
---|
79 | |
---|
80 | ! We catch the number of variables in the stats file |
---|
81 | ierr = NF_INQ_NVARS(nid,nbvar) |
---|
82 | |
---|
83 | ! to catche the "real" number of variables (without the "additionnal variables") |
---|
84 | nbvar=(nbvar-4)/2 |
---|
85 | |
---|
86 | do i=1,nbvar |
---|
87 | varid=(i-1)*2+5 |
---|
88 | |
---|
89 | ! What's the variable's name? |
---|
90 | ierr=NF_INQ_VARNAME(nid,varid,name) |
---|
91 | write(*,*) "OK variable ",name |
---|
92 | ! Its units? |
---|
93 | units=" " |
---|
94 | ierr=NF_GET_ATT_TEXT(nid,varid,"units",units) |
---|
95 | ! Its title? |
---|
96 | title=" " |
---|
97 | ierr=NF_GET_ATT_TEXT(nid,varid,"title",title) |
---|
98 | ! Its number of dimensions? |
---|
99 | ierr=NF_INQ_VARNDIMS(nid,varid,ndims) |
---|
100 | ! Its values? |
---|
101 | |
---|
102 | if(ndims==4) then ! lat, lon, alt & time |
---|
103 | |
---|
104 | ! dimout(1)=lonid |
---|
105 | ! dimout(2)=latid |
---|
106 | ! dimout(3)=altid |
---|
107 | ! dimout(4)=timeid |
---|
108 | |
---|
109 | size=(/iip1,jjp1,llm,1/) |
---|
110 | do lt=1,istime |
---|
111 | start=(/1,1,1,lt/) |
---|
112 | ! Extraction of the "source" variables |
---|
113 | #ifdef NC_DOUBLE |
---|
114 | ierr = NF_GET_VARA_DOUBLE(nid,varid,start,size,sum3d) |
---|
115 | ierr = NF_GET_VARA_DOUBLE(nid,varid+1,start,size,square3d) |
---|
116 | #else |
---|
117 | ierr = NF_GET_VARA_REAL(nid,varid,start,size,sum3d) |
---|
118 | ierr = NF_GET_VARA_REAL(nid,varid+1,start,size,square3d) |
---|
119 | #endif |
---|
120 | ! Calculation of these variables |
---|
121 | mean3d=sum3d/count(lt) |
---|
122 | sd3d=sqrt(max(0.,square3d/count(lt)-mean3d**2)) |
---|
123 | ! Writing of the variables |
---|
124 | #ifdef NC_DOUBLE |
---|
125 | ierr = NF_PUT_VARA_DOUBLE(nid,varid,start,size,mean3d) |
---|
126 | ierr = NF_PUT_VARA_DOUBLE(nid,varid+1,start,size,sd3d) |
---|
127 | #else |
---|
128 | ierr = NF_PUT_VARA_REAL(nid,varid,start,size,mean3d) |
---|
129 | ierr = NF_PUT_VARA_REAL(nid,varid+1,start,size,sd3d) |
---|
130 | #endif |
---|
131 | enddo |
---|
132 | |
---|
133 | else if (ndims.eq.3) then |
---|
134 | |
---|
135 | ! dimout(1)=lonid |
---|
136 | ! dimout(2)=latid |
---|
137 | ! dimout(3)=timeid |
---|
138 | |
---|
139 | size=(/iip1,jjp1,1,0/) |
---|
140 | do lt=1,istime |
---|
141 | start=(/1,1,lt,0/) |
---|
142 | ! Extraction of the "source" variables |
---|
143 | #ifdef NC_DOUBLE |
---|
144 | ierr = NF_GET_VARA_DOUBLE(nid,varid,start,size,sum2d) |
---|
145 | ierr = NF_GET_VARA_DOUBLE(nid,varid+1,start,size,square2d) |
---|
146 | #else |
---|
147 | ierr = NF_GET_VARA_REAL(nid,varid,start,size,sum2d) |
---|
148 | ierr = NF_GET_VARA_REAL(nid,varid+1,start,size,square2d) |
---|
149 | #endif |
---|
150 | ! Calculation of these variables |
---|
151 | mean2d=sum2d/count(lt) |
---|
152 | sd2d=sqrt(max(0.,square2d/count(lt)-mean2d**2)) |
---|
153 | ! Writing of the variables |
---|
154 | #ifdef NC_DOUBLE |
---|
155 | ierr = NF_PUT_VARA_DOUBLE(nid,varid,start,size,mean2d) |
---|
156 | ierr = NF_PUT_VARA_DOUBLE(nid,varid+1,start,size,sd2d) |
---|
157 | #else |
---|
158 | ierr = NF_PUT_VARA_REAL(nid,varid,start,size,mean2d) |
---|
159 | ierr = NF_PUT_VARA_REAL(nid,varid+1,start,size,sd2d) |
---|
160 | #endif |
---|
161 | enddo |
---|
162 | |
---|
163 | endif |
---|
164 | enddo |
---|
165 | |
---|
166 | ierr= NF_CLOSE(nid) |
---|
167 | |
---|
168 | endif ! of if (is_master) |
---|
169 | |
---|
170 | end |
---|