LMDZ
wstats.F90
Go to the documentation of this file.
1 subroutine wstats(ngrid,nom,titre,unite,dim,px)
2 
3 implicit none
4 
5 #include "dimensions.h"
6 #include "statto.h"
7 #include "netcdf.inc"
8 
9 integer,intent(in) :: ngrid
10 character (len=*),intent(in) :: nom,titre,unite
11 integer,intent(in) :: dim
12 real, dimension(ngrid,llm),intent(in) :: px
13 integer,parameter :: iip1=iim+1
14 integer,parameter :: jjp1=jjm+1
15 real, dimension(iip1,jjp1,llm) :: mean3d,sd3d,dx3
16 real, dimension(iip1,jjp1) :: mean2d,sd2d,dx2
17 character (len=50) :: namebis
18 character (len=50), save :: firstvar
19 integer :: ierr,varid,nbdim,nid
20 integer :: meanid,sdid
21 integer, dimension(4) :: id,start,size
22 logical, save :: firstcall=.true.
23 integer :: l,i,j,ig0
24 integer,save :: index
25 
26 integer, save :: step=0
27 
28 
29 if (firstcall) then
30  firstcall=.false.
31  firstvar=trim((nom))
32  call inistats(ierr)
33 endif
34 
35 if (firstvar==nom) then ! If we're back to the first variable
36  step = step + 1
37 endif
38 
39 if (mod(step,istats).ne.0) then
40  RETURN
41 endif
42 
43 ierr = nf_open("stats.nc",nf_write,nid)
44 
45 namebis=trim(nom)
46 ierr= nf_inq_varid(nid,namebis,meanid)
47 
48 if (ierr.ne.nf_noerr) then
49 
50  if (firstvar==nom) then
51  index=1
52  count=0
53  endif
54 
55 
56 !declaration de la variable
57 
58 ! choix du nom des coordonnees
59  ierr= nf_inq_dimid(nid,"longitude",id(1))
60  ierr= nf_inq_dimid(nid,"latitude",id(2))
61  if (dim.eq.3) then
62  ierr= nf_inq_dimid(nid,"altitude",id(3))
63  ierr= nf_inq_dimid(nid,"Time",id(4))
64  nbdim=4
65  else if (dim==2) then
66  ierr= nf_inq_dimid(nid,"Time",id(3))
67  nbdim=3
68  endif
69 
70  write (*,*) "====================="
71  write (*,*) "STATS: creation de ",nom
72  namebis=trim(nom)
73  call def_var_stats(nid,namebis,titre,unite,nbdim,id,meanid,ierr)
74  call inivar(nid,meanid,ngrid,dim,index,px,ierr)
75  namebis=trim(nom)//"_sd"
76  call def_var_stats(nid,namebis,trim(titre)//" total standard deviation over the season",unite,nbdim,id,sdid,ierr)
77  call inivar(nid,sdid,ngrid,dim,index,px,ierr)
78 
79  ierr= nf_close(nid)
80  return
81 
82 else
83  namebis=trim(nom)//"_sd"
84  ierr= nf_inq_varid(nid,namebis,sdid)
85 
86 endif
87 
88 if (firstvar==nom) then
89  count(index)=count(int(index))+1
90  index=index+1
91  if (index>istime) then
92  index=1
93  endif
94 endif
95 
96 if (count(index)==0) then
97  if (dim.eq.3) then
98  start=(/1,1,1,index/)
99  size=(/iip1,jjp1,llm,1/)
100  mean3d=0
101  sd3d=0
102  else if (dim.eq.2) then
103  start=(/1,1,index,0/)
104  size=(/iip1,jjp1,1,0/)
105  mean2d=0
106  sd2d=0
107  endif
108 else
109  if (dim.eq.3) then
110  start=(/1,1,1,index/)
111  size=(/iip1,jjp1,llm,1/)
112 #ifdef NC_DOUBLE
113  ierr = nf_get_vara_double(nid,meanid,start,size,mean3d)
114  ierr = nf_get_vara_double(nid,sdid,start,size,sd3d)
115 #else
116  ierr = nf_get_vara_real(nid,meanid,start,size,mean3d)
117  ierr = nf_get_vara_real(nid,sdid,start,size,sd3d)
118 #endif
119  if (ierr.ne.nf_noerr) then
120  write (*,*) nf_strerror(ierr)
121  stop ""
122  endif
123 
124  else if (dim.eq.2) then
125  start=(/1,1,index,0/)
126  size=(/iip1,jjp1,1,0/)
127 #ifdef NC_DOUBLE
128  ierr = nf_get_vara_double(nid,meanid,start,size,mean2d)
129  ierr = nf_get_vara_double(nid,sdid,start,size,sd2d)
130 #else
131  ierr = nf_get_vara_real(nid,meanid,start,size,mean2d)
132  ierr = nf_get_vara_real(nid,sdid,start,size,sd2d)
133 #endif
134  if (ierr.ne.nf_noerr) then
135  write (*,*) nf_strerror(ierr)
136  stop ""
137  endif
138  endif
139 endif
140 
141 if (dim.eq.3) then
142 
143 ! Passage variable physique --> variable dynamique
144 
145  DO l=1,llm
146  DO i=1,iip1
147  dx3(i,1,l)=px(1,l)
148  dx3(i,jjp1,l)=px(ngrid,l)
149  ENDDO
150  DO j=2,jjm
151  ig0= 1+(j-2)*iim
152  DO i=1,iim
153  dx3(i,j,l)=px(ig0+i,l)
154  ENDDO
155  dx3(iip1,j,l)=dx3(1,j,l)
156  ENDDO
157  ENDDO
158 
159  mean3d(:,:,:)= mean3d(:,:,:)+dx3(:,:,:)
160  sd3d(:,:,:)= sd3d(:,:,:)+dx3(:,:,:)**2
161 
162 #ifdef NC_DOUBLE
163  ierr = nf_put_vara_double(nid,meanid,start,size,mean3d)
164  ierr = nf_put_vara_double(nid,sdid,start,size,sd3d)
165 #else
166  ierr = nf_put_vara_real(nid,meanid,start,size,mean3d)
167  ierr = nf_put_vara_real(nid,sdid,start,size,sd3d)
168 #endif
169  if (ierr.ne.nf_noerr) then
170  write (*,*) nf_strerror(ierr)
171  stop ""
172  endif
173 
174 else if (dim.eq.2) then
175 
176 ! Passage variable physique --> physique dynamique
177 
178  DO i=1,iip1
179  dx2(i,1)=px(1,1)
180  dx2(i,jjp1)=px(ngrid,1)
181  ENDDO
182  DO j=2,jjm
183  ig0= 1+(j-2)*iim
184  DO i=1,iim
185  dx2(i,j)=px(ig0+i,1)
186  ENDDO
187  dx2(iip1,j)=dx2(1,j)
188  ENDDO
189 
190  mean2d(:,:)= mean2d(:,:)+dx2(:,:)
191  sd2d(:,:)= sd2d(:,:)+dx2(:,:)**2
192 
193 #ifdef NC_DOUBLE
194  ierr = nf_put_vara_double(nid,meanid,start,size,mean2d)
195  ierr = nf_put_vara_double(nid,sdid,start,size,sd2d)
196 #else
197  ierr = nf_put_vara_real(nid,meanid,start,size,mean2d)
198  ierr = nf_put_vara_real(nid,sdid,start,size,sd2d)
199 #endif
200  if (ierr.ne.nf_noerr) then
201  write (*,*) nf_strerror(ierr)
202  stop ""
203  endif
204 
205 endif
206 
207 ierr= nf_close(nid)
208 
209 end subroutine wstats
210 
211 !======================================================
212 subroutine inivar(nid,varid,ngrid,dim,index,px,ierr)
214 implicit none
215 
216 include "dimensions.h"
217 !!!!!!include "dimphys.h"
218 include "netcdf.inc"
219 
220 integer, intent(in) :: nid,varid,dim,index,ngrid
221 real, dimension(ngrid,llm), intent(in) :: px
222 integer, intent(out) :: ierr
223 
224 integer,parameter :: iip1=iim+1
225 integer,parameter :: jjp1=jjm+1
226 
227 integer :: l,i,j,ig0
228 integer, dimension(4) :: start,size
229 real, dimension(iip1,jjp1,llm) :: dx3
230 real, dimension(iip1,jjp1) :: dx2
231 
232 if (dim.eq.3) then
233 
234  start=(/1,1,1,index/)
235  size=(/iip1,jjp1,llm,1/)
236 
237 ! Passage variable physique --> variable dynamique
238 
239  DO l=1,llm
240  DO i=1,iip1
241  dx3(i,1,l)=px(1,l)
242  dx3(i,jjp1,l)=px(ngrid,l)
243  ENDDO
244  DO j=2,jjm
245  ig0= 1+(j-2)*iim
246  DO i=1,iim
247  dx3(i,j,l)=px(ig0+i,l)
248  ENDDO
249  dx3(iip1,j,l)=dx3(1,j,l)
250  ENDDO
251  ENDDO
252 
253 #ifdef NC_DOUBLE
254  ierr = nf_put_vara_double(nid,varid,start,size,dx3)
255 #else
256  ierr = nf_put_vara_real(nid,varid,start,size,dx3)
257 #endif
258 
259 else if (dim.eq.2) then
260 
261  start=(/1,1,index,0/)
262  size=(/iip1,jjp1,1,0/)
263 
264 ! Passage variable physique --> physique dynamique
265 
266  DO i=1,iip1
267  dx2(i,1)=px(1,1)
268  dx2(i,jjp1)=px(ngrid,1)
269  ENDDO
270  DO j=2,jjm
271  ig0= 1+(j-2)*iim
272  DO i=1,iim
273  dx2(i,j)=px(ig0+i,1)
274  ENDDO
275  dx2(iip1,j)=dx2(1,j)
276  ENDDO
277 
278 #ifdef NC_DOUBLE
279  ierr = nf_put_vara_double(nid,varid,start,size,dx2)
280 #else
281  ierr = nf_put_vara_real(nid,varid,start,size,dx2)
282 #endif
283 
284 endif
285 
286 end subroutine inivar
287 
288 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
289 
290 subroutine def_var_stats(nid,name,title,units,nbdim,dimids,nvarid,ierr)
292 ! This subroutine defines variable 'name' in a (pre-existing and opened)
293 ! NetCDF file (known from its NetCDF ID 'nid').
294 ! The number of dimensions 'nbdim' of the variable, as well as the IDs of
295 ! corresponding dimensions must be set (in array 'dimids').
296 ! Upon successfull definition of the variable, 'nvarid' contains the
297 ! NetCDF ID of the variable.
298 ! The variables' attributes 'title' (Note that 'long_name' would be more
299 ! appropriate) and 'units' are also set.
300 
301 implicit none
302 
303 #include "netcdf.inc"
304 
305 integer,intent(in) :: nid ! NetCDF file ID
306 character(len=*),intent(in) :: name ! the variable's name
307 character(len=*),intent(in) :: title ! 'title' attribute of variable
308 character(len=*),intent(in) :: units ! 'units' attribute of variable
309 integer,intent(in) :: nbdim ! number of dimensions of the variable
310 integer,dimension(nbdim),intent(in) :: dimids ! NetCDF IDs of the dimensions
311  ! the variable is defined along
312 integer,intent(out) :: nvarid ! NetCDF ID of the variable
313 integer,intent(out) :: ierr ! returned NetCDF staus code
314 
315 ! 1. Switch to NetCDF define mode
316 ierr=nf_redef(nid)
317 
318 ! 2. Define the variable
319 #ifdef NC_DOUBLE
320 ierr = nf_def_var(nid,adjustl(name),nf_double,nbdim,dimids,nvarid)
321 #else
322 ierr = nf_def_var(nid,adjustl(name),nf_float,nbdim,dimids,nvarid)
323 #endif
324 if(ierr/=nf_noerr) then
325  write(*,*) "def_var_stats: Failed defining variable "//trim(name)
326  write(*,*) nf_strerror(ierr)
327  stop ""
328 endif
329 
330 ! 3. Write attributes
331 ierr=nf_put_att_text(nid,nvarid,"title",&
332  len_trim(adjustl(title)),adjustl(title))
333 if(ierr/=nf_noerr) then
334  write(*,*) "def_var_stats: Failed writing title attribute for "//trim(name)
335  write(*,*) nf_strerror(ierr)
336  stop ""
337 endif
338 
339 ierr=nf_put_att_text(nid,nvarid,"units",&
340  len_trim(adjustl(units)),adjustl(units))
341 if(ierr/=nf_noerr) then
342  write(*,*) "def_var_stats: Failed writing units attribute for "//trim(name)
343  write(*,*) nf_strerror(ierr)
344  stop ""
345 endif
346 
347 ! 4. Switch out of NetCDF define mode
348 ierr = nf_enddef(nid)
349 
350 end subroutine def_var_stats
351 
subroutine wstats(ngrid, nom, titre, unite, dim, px)
Definition: wstats.F90:2
!$Id Turb_fcg_gcssold get_uvd hqturb_gcssold endif!large scale llm day day1 day day1 *dt_toga endif!time annee_ref dt_toga u_toga vq_toga w_prof vq_prof llm day day1 day day1 *dt_dice endif!time annee_ref dt_dice swup_dice vg_dice omega_dice tg_prof vg_profd w_profd omega_profd!do llm!print llm l llm
!$Id itapm1 ENDIF!IM on interpole les champs sur les niveaux STD de pression!IM a chaque pas de temps de la physique c!positionnement de l argument logique a false c!pour ne pas recalculer deux fois la meme chose!c!a cet effet un appel a plevel_new a ete deplace c!a la fin de la serie d appels c!la boucle DO nlevSTD a ete internalisee c!dans d ou la creation de cette routine c c!CALL false
Definition: calcul_STDlev.h:26
subroutine inivar(nid, varid, ngrid, dim, index, px, ierr)
Definition: wstats.F90:213
!$Id itapm1 ENDIF!IM on interpole les champs sur les niveaux STD de pression!IM a chaque pas de temps de la physique c!positionnement de l argument logique a false c!pour ne pas recalculer deux fois la meme chose!c!a cet effet un appel a plevel_new a ete deplace c!a la fin de la serie d appels c!la boucle DO nlevSTD a ete internalisee c!dans d ou la creation de cette routine c c!CALL ulevSTD CALL &zphi philevSTD CALL &zx_rh rhlevSTD!DO klev DO klon klev DO klon klev DO klon klev DO klon klev DO klon klev DO klon klev DO klon klev DO klon klev DO klon klev DO klon du jour ou toutes les read_climoz CALL true
c c zjulian c cym CALL iim cym klev iim
Definition: ini_bilKP_ave.h:24
subroutine inistats(ierr)
Definition: inistats.F90:2
subroutine def_var_stats(nid, name, title, units, nbdim, dimids, nvarid, ierr)
Definition: wstats.F90:291