source: trunk/LMDZ.MARS/libf/phymars/inistats.F @ 2461

Last change on this file since 2461 was 2399, checked in by emillour, 5 years ago

Mars GCM:
More code cleanup: use "call abort_physic()" instead of "stop" or "call abort"
EM

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