source: LMDZ5/branches/LF-private/libf/phylmd/inistats.F @ 5456

Last change on this file since 5456 was 1654, checked in by aslmd, 12 years ago

  • Changements dans calfis pour appeler la physique generique planetaire NB: dans le futur il faudrait un appel commun mais actuellement differences dans les arguments...
  • Changements dans dynetat0 pour s'accommoder des tableaux "controle" terrestre et planeto
  • Passage de la variable pseudoalt dans comvert.h (compatibilite planeto)
  • arch : ajout des fichiers pour faire fonctionner sur ciclad avec ifort

verifie : compatibilite bit-par-bit des resultats et initialisations, ainsi que 1+1=2 avec unpun.sh


  • Changes in calfis so that generic physics for planets could be interfaced NB: ideally a common call should be implemented but for now there are differences in arguments...
  • Changes in dynetat0 to cope with difference in "controle" between earth and planeto
  • Transfer of variable pseudoalt into comvert.h (this is only used by planeto and done for compatibility)
  • arch: added files to allow use of ifort compiler on ciclad cluster

checked : bit-by-bit compatibility of results and starts, as well as 1+1=2 with unpun.sh


Author : AS equipe planeto

File size: 4.1 KB
Line 
1      subroutine inistats(ierr)
2
3      implicit none
4
5#include "dimensions.h"
6#include "paramet.h"
7#include "comgeom.h"
8#include "comvert.h"
9#include "comconst.h"
10#include "statto.h"
11#include "netcdf.inc"
12
13      integer,intent(out) :: ierr
14      integer :: nid
15      integer :: l,nsteppd
16      real, dimension(llm) ::  sig_s
17      integer :: idim_lat,idim_lon,idim_llm,idim_llmp1,idim_time
18      real, dimension(istime) :: lt
19      integer :: nvarid
20
21      write (*,*)
22      write (*,*) '                        || STATS ||'
23      write (*,*)
24      write (*,*) 'daysec',daysec
25      write (*,*) 'dtphys',dtphys
26      nsteppd=nint(daysec/dtphys)
27      write (*,*) 'nsteppd=',nsteppd
28      if (abs(float(nsteppd)-daysec/dtphys).gt.1.e-8*daysec)
29     &   stop'Dans Instat:  1jour .ne. n pas physiques'
30
31      if(mod(nsteppd,istime).ne.0)
32     &   stop'Dans Instat:  1jour .ne. n*istime pas physiques'
33
34      istats=nsteppd/istime
35      write (*,*) 'istats=',istats
36      write (*,*) 'Storing ',istime,'times per day'
37      write (*,*) 'thus every ',istats,'physical timestep '
38      write (*,*)
39
40      do l= 1, llm
41         sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2.
42         pseudoalt(l)=log(preff/presnivs(l))*8. 
43      enddo
44
45      ierr = NF_CREATE("stats.nc",NF_CLOBBER,nid)
46      if (ierr.ne.NF_NOERR) then
47         write (*,*) NF_STRERROR(ierr)
48         stop ""
49      endif
50
51      ierr = NF_DEF_DIM (nid, "latitude", jjp1, idim_lat)
52      ierr = NF_DEF_DIM (nid, "longitude", iip1, idim_lon)
53      ierr = NF_DEF_DIM (nid, "altitude", llm, idim_llm)
54      ierr = NF_DEF_DIM (nid, "llmp1", llm+1, idim_llmp1)
55      ierr = NF_DEF_DIM (nid, "Time", NF_UNLIMITED, idim_time)
56
57      ierr = NF_ENDDEF(nid)
58      call def_var_stats(nid,"Time","Time",
59     &            "hours since 0000-00-0 00:00:00",1,
60     &            idim_time,nvarid,ierr)
61! Time is initialised later by mkstats subroutine
62
63      call def_var_stats(nid,"latitude","latitude",
64     &            "degrees_north",1,idim_lat,nvarid,ierr)
65#ifdef NC_DOUBLE
66      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlatu/pi*180)
67#else
68      ierr = NF_PUT_VAR_REAL (nid,nvarid,rlatu/pi*180)
69#endif
70      call def_var_stats(nid,"longitude","East longitude",
71     &            "degrees_east",1,idim_lon,nvarid,ierr)
72#ifdef NC_DOUBLE
73      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlonv/pi*180)
74#else
75      ierr = NF_PUT_VAR_REAL (nid,nvarid,rlonv/pi*180)
76#endif
77
78! Niveaux verticaux, aps et bps
79      ierr = NF_REDEF (nid)
80! presnivs
81#ifdef NC_DOUBLE
82      ierr = NF_DEF_VAR (nid,"presnivs", NF_DOUBLE, 1,idim_llm,nvarid)
83#else
84      ierr = NF_DEF_VAR (nid,"presnivs", NF_FLOAT, 1,idim_llm,nvarid)
85#endif
86      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name",15,
87     &                        "Vertical levels")
88      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'units',2,"Pa")
89      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'positive',4,"down")
90      ierr = NF_ENDDEF(nid)
91#ifdef NC_DOUBLE
92      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,
93     &                          presnivs(1:llm))
94#else
95      ierr = NF_PUT_VAR_REAL (nid,nvarid,
96     &                        presnivs(1:llm))
97#endif
98! Pseudo alts
99#ifdef NC_DOUBLE
100      ierr = NF_DEF_VAR (nid,"altitude", NF_DOUBLE, 1,idim_llm,nvarid)
101#else
102      ierr = NF_DEF_VAR (nid,"altitude", NF_FLOAT, 1,idim_llm,nvarid)
103#endif
104      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name",8,"altitude")
105      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'units',2,"km")
106      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'positive',2,"up")
107      ierr = NF_ENDDEF(nid)
108#ifdef NC_DOUBLE
109      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,pseudoalt)
110#else
111      ierr = NF_PUT_VAR_REAL (nid,nvarid,pseudoalt)
112#endif
113!      call def_var_stats(nid,"aps","hybrid pressure at midlayers"," ",
114!     &            1,idim_llm,nvarid,ierr)
115!#ifdef NC_DOUBLE
116!      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,aps)
117!#else
118!      ierr = NF_PUT_VAR_REAL (nid,nvarid,aps)
119!#endif
120
121!      call def_var_stats(nid,"bps","hybrid sigma at midlayers"," ",
122!     &            1,idim_llm,nvarid,ierr)
123!#ifdef NC_DOUBLE
124!      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bps)
125!#else
126!      ierr = NF_PUT_VAR_REAL (nid,nvarid,bps)
127!#endif
128
129      ierr=NF_CLOSE(nid)
130
131      end subroutine inistats
132
Note: See TracBrowser for help on using the repository browser.