1 | subroutine inistats(ierr) |
---|
2 | |
---|
3 | use mod_phys_lmdz_para, only : is_master |
---|
4 | USE comvert_mod, ONLY: ap,bp,aps,bps,preff,pseudoalt,presnivs |
---|
5 | USE comcstfi_h, ONLY: pi |
---|
6 | USE time_phylmdz_mod, ONLY: daysec,dtphys |
---|
7 | USE regular_lonlat_mod, ONLY: lon_reg, lat_reg |
---|
8 | USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, nbp_lev |
---|
9 | implicit none |
---|
10 | |
---|
11 | include "statto.h" |
---|
12 | include "netcdf.inc" |
---|
13 | |
---|
14 | integer,intent(out) :: ierr |
---|
15 | integer :: nid |
---|
16 | integer :: l,nsteppd |
---|
17 | real, dimension(nbp_lev) :: sig_s |
---|
18 | real :: lon_reg_ext(nbp_lon+1) ! extended longitudes |
---|
19 | integer :: idim_lat,idim_lon,idim_llm,idim_llmp1,idim_time |
---|
20 | real, dimension(istime) :: lt |
---|
21 | integer :: nvarid |
---|
22 | |
---|
23 | write (*,*) |
---|
24 | write (*,*) ' || STATS ||' |
---|
25 | write (*,*) |
---|
26 | write (*,*) 'daysec',daysec |
---|
27 | write (*,*) 'dtphys',dtphys |
---|
28 | nsteppd=nint(daysec/dtphys) |
---|
29 | write (*,*) 'nsteppd=',nsteppd |
---|
30 | if (abs(float(nsteppd)-daysec/dtphys).gt.1.e-8*daysec) |
---|
31 | & stop'Dans Instat: 1jour .ne. n pas physiques' |
---|
32 | |
---|
33 | if(mod(nsteppd,istime).ne.0) |
---|
34 | & stop'Dans Instat: 1jour .ne. n*istime pas physiques' |
---|
35 | |
---|
36 | istats=nsteppd/istime |
---|
37 | write (*,*) 'istats=',istats |
---|
38 | write (*,*) 'Storing ',istime,'times per day' |
---|
39 | write (*,*) 'thus every ',istats,'physical timestep ' |
---|
40 | write (*,*) |
---|
41 | |
---|
42 | do l= 1, nbp_lev |
---|
43 | sig_s(l)=((ap(l)+ap(l+1))/preff+bp(l)+bp(l+1))/2. |
---|
44 | pseudoalt(l)=-10.*log(presnivs(l)/preff) |
---|
45 | enddo |
---|
46 | |
---|
47 | lon_reg_ext(1:nbp_lon)=lon_reg(1:nbp_lon) |
---|
48 | !add extra redundant point (180 degrees, since lon_reg starts at -180 |
---|
49 | lon_reg_ext(nbp_lon+1)=-lon_reg_ext(1) |
---|
50 | |
---|
51 | if (is_master) then |
---|
52 | ! only the master needs do this |
---|
53 | |
---|
54 | ierr = NF_CREATE("stats.nc",IOR(NF_CLOBBER,NF_64BIT_OFFSET),nid) |
---|
55 | if (ierr.ne.NF_NOERR) then |
---|
56 | write (*,*) NF_STRERROR(ierr) |
---|
57 | stop "" |
---|
58 | endif |
---|
59 | |
---|
60 | ierr = NF_DEF_DIM (nid, "latitude", nbp_lat, idim_lat) |
---|
61 | ierr = NF_DEF_DIM (nid, "longitude", nbp_lon+1, idim_lon) |
---|
62 | ierr = NF_DEF_DIM (nid, "altitude", nbp_lev, idim_llm) |
---|
63 | ierr = NF_DEF_DIM (nid, "llmp1", nbp_lev+1, idim_llmp1) |
---|
64 | ierr = NF_DEF_DIM (nid, "Time", NF_UNLIMITED, idim_time) |
---|
65 | |
---|
66 | ierr = NF_ENDDEF(nid) |
---|
67 | call def_var_stats(nid,"Time","Time", |
---|
68 | & "days since 0000-00-0 00:00:00",1, |
---|
69 | & idim_time,nvarid,ierr) |
---|
70 | ! Time is initialised later by mkstats subroutine |
---|
71 | |
---|
72 | call def_var_stats(nid,"latitude","latitude", |
---|
73 | & "degrees_north",1,idim_lat,nvarid,ierr) |
---|
74 | #ifdef NC_DOUBLE |
---|
75 | ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lat_reg/pi*180) |
---|
76 | #else |
---|
77 | ierr = NF_PUT_VAR_REAL (nid,nvarid,lat_reg/pi*180) |
---|
78 | #endif |
---|
79 | call def_var_stats(nid,"longitude","East longitude", |
---|
80 | & "degrees_east",1,idim_lon,nvarid,ierr) |
---|
81 | #ifdef NC_DOUBLE |
---|
82 | ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lon_reg_ext/pi*180) |
---|
83 | #else |
---|
84 | ierr = NF_PUT_VAR_REAL (nid,nvarid,lon_reg_ext/pi*180) |
---|
85 | #endif |
---|
86 | |
---|
87 | ! Niveaux verticaux, aps et bps |
---|
88 | ierr = NF_REDEF (nid) |
---|
89 | #ifdef NC_DOUBLE |
---|
90 | ierr = NF_DEF_VAR (nid,"altitude", NF_DOUBLE, 1,idim_llm,nvarid) |
---|
91 | #else |
---|
92 | ierr = NF_DEF_VAR (nid,"altitude", NF_FLOAT, 1,idim_llm,nvarid) |
---|
93 | #endif |
---|
94 | ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name",8,"altitude") |
---|
95 | ierr = NF_PUT_ATT_TEXT (nid,nvarid,'units',2,"km") |
---|
96 | ierr = NF_PUT_ATT_TEXT (nid,nvarid,'positive',2,"up") |
---|
97 | ierr = NF_ENDDEF(nid) |
---|
98 | #ifdef NC_DOUBLE |
---|
99 | ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,pseudoalt) |
---|
100 | #else |
---|
101 | ierr = NF_PUT_VAR_REAL (nid,nvarid,pseudoalt) |
---|
102 | #endif |
---|
103 | call def_var_stats(nid,"aps","hybrid pressure at midlayers" |
---|
104 | & ," ",1,idim_llm,nvarid,ierr) |
---|
105 | #ifdef NC_DOUBLE |
---|
106 | ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,aps) |
---|
107 | #else |
---|
108 | ierr = NF_PUT_VAR_REAL (nid,nvarid,aps) |
---|
109 | #endif |
---|
110 | |
---|
111 | call def_var_stats(nid,"bps","hybrid sigma at midlayers" |
---|
112 | & ," ",1,idim_llm,nvarid,ierr) |
---|
113 | #ifdef NC_DOUBLE |
---|
114 | ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bps) |
---|
115 | #else |
---|
116 | ierr = NF_PUT_VAR_REAL (nid,nvarid,bps) |
---|
117 | #endif |
---|
118 | |
---|
119 | ierr=NF_CLOSE(nid) |
---|
120 | |
---|
121 | endif ! of if (is_master) |
---|
122 | end subroutine inistats |
---|
123 | |
---|