source: LMDZ5/branches/LMDZ5_SPLA/libf/phylmd/inistats.F90 @ 2265

Last change on this file since 2265 was 1992, checked in by lguez, 11 years ago

Converted to free source form files in libf/phylmd which were still in
fixed source form. The conversion was done using the polish mode of
the NAG Fortran Compiler.

In addition to converting to free source form, the processing of the
files also:

-- indented the code (including comments);

-- set Fortran keywords to uppercase, and set all other identifiers
to lower case;

-- added qualifiers to end statements (for example "end subroutine
conflx", instead of "end");

-- changed the terminating statements of all DO loops so that each
loop ends with an ENDDO statement (instead of a labeled continue).

-- replaced #include by include.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 3.8 KB
Line 
1SUBROUTINE 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)>1.E-8*daysec) &
29    STOP 'Dans Instat:  1jour .ne. n pas physiques'
30
31  IF (mod(nsteppd,istime)/=0) STOP &
32    '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  END DO
44
45  ierr = nf_create('stats.nc', nf_clobber, nid)
46  IF (ierr/=nf_noerr) THEN
47    WRITE (*, *) nf_strerror(ierr)
48    STOP ''
49  END IF
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', 'hours since 0000-00-0 00:00:00', &
59    1, idim_time, nvarid, ierr)
60  ! Time is initialised later by mkstats subroutine
61
62  CALL def_var_stats(nid, 'latitude', 'latitude', 'degrees_north', 1, &
63    idim_lat, nvarid, ierr)
64#ifdef NC_DOUBLE
65  ierr = nf_put_var_double(nid, nvarid, rlatu/pi*180)
66#else
67  ierr = nf_put_var_real(nid, nvarid, rlatu/pi*180)
68#endif
69  CALL def_var_stats(nid, 'longitude', 'East longitude', 'degrees_east', 1, &
70    idim_lon, nvarid, ierr)
71#ifdef NC_DOUBLE
72  ierr = nf_put_var_double(nid, nvarid, rlonv/pi*180)
73#else
74  ierr = nf_put_var_real(nid, nvarid, rlonv/pi*180)
75#endif
76
77  ! Niveaux verticaux, aps et bps
78  ierr = nf_redef(nid)
79  ! presnivs
80#ifdef NC_DOUBLE
81  ierr = nf_def_var(nid, 'presnivs', nf_double, 1, idim_llm, nvarid)
82#else
83  ierr = nf_def_var(nid, 'presnivs', nf_float, 1, idim_llm, nvarid)
84#endif
85  ierr = nf_put_att_text(nid, nvarid, 'long_name', 15, 'Vertical levels')
86  ierr = nf_put_att_text(nid, nvarid, 'units', 2, 'Pa')
87  ierr = nf_put_att_text(nid, nvarid, 'positive', 4, 'down')
88  ierr = nf_enddef(nid)
89#ifdef NC_DOUBLE
90  ierr = nf_put_var_double(nid, nvarid, presnivs(1:llm))
91#else
92  ierr = nf_put_var_real(nid, nvarid, presnivs(1:llm))
93#endif
94  ! Pseudo alts
95#ifdef NC_DOUBLE
96  ierr = nf_def_var(nid, 'altitude', nf_double, 1, idim_llm, nvarid)
97#else
98  ierr = nf_def_var(nid, 'altitude', nf_float, 1, idim_llm, nvarid)
99#endif
100  ierr = nf_put_att_text(nid, nvarid, 'long_name', 8, 'altitude')
101  ierr = nf_put_att_text(nid, nvarid, 'units', 2, 'km')
102  ierr = nf_put_att_text(nid, nvarid, 'positive', 2, 'up')
103  ierr = nf_enddef(nid)
104#ifdef NC_DOUBLE
105  ierr = nf_put_var_double(nid, nvarid, pseudoalt)
106#else
107  ierr = nf_put_var_real(nid, nvarid, pseudoalt)
108#endif
109  ! call def_var_stats(nid,"aps","hybrid pressure at midlayers"," ",
110  ! &            1,idim_llm,nvarid,ierr)
111  ! #ifdef NC_DOUBLE
112  ! ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,aps)
113  ! #else
114  ! ierr = NF_PUT_VAR_REAL (nid,nvarid,aps)
115  ! #endif
116
117  ! call def_var_stats(nid,"bps","hybrid sigma at midlayers"," ",
118  ! &            1,idim_llm,nvarid,ierr)
119  ! #ifdef NC_DOUBLE
120  ! ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bps)
121  ! #else
122  ! ierr = NF_PUT_VAR_REAL (nid,nvarid,bps)
123  ! #endif
124
125  ierr = nf_close(nid)
126
127END SUBROUTINE inistats
128
Note: See TracBrowser for help on using the repository browser.