source: trunk/LMDZ.GENERIC/libf/phystd/iniwrite.F

Last change on this file was 3773, checked in by afalco, 6 weeks ago

Generic: allow to write controle in XIOS output.
AF

File size: 11.7 KB
Line 
1      SUBROUTINE iniwrite(nid,idayref,phis,area,nbplon,nbplat)
2
3      use comsoil_h, only: mlayer, nsoilmx
4      USE comcstfi_mod, only: g, mugaz, omeg, rad, rcp, pi
5      USE vertical_layers_mod, ONLY: ap,bp,aps,bps,pseudoalt
6!      USE logic_mod, ONLY: fxyhypb,ysinus
7!      USE serre_mod, ONLY: clon,clat,grossismx,grossismy,dzoomx,dzoomy
8      USE time_phylmdz_mod, ONLY: daysec, dtphys
9!      USE ener_mod, ONLY: etot0,ptot0,ztot0,stot0,ang0
10      USE regular_lonlat_mod, ONLY: lon_reg, lat_reg
11      USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, nbp_lev
12      use tabfi_mod, only: tab_cntrl_mod
13      IMPLICIT NONE
14
15c=======================================================================
16c
17c   Auteur:  L. Fairhead  ,  P. Le Van, Y. Wanherdrick, F. Forget
18c   -------
19c
20c   Objet:
21c   ------
22c
23c   'Initialize' the diagfi.nc file: write down dimensions as well
24c   as time-independent fields (e.g: geopotential, mesh area, ...)
25c
26c=======================================================================
27c-----------------------------------------------------------------------
28c   Declarations:
29c   -------------
30
31      include "netcdf.inc"
32
33c   Arguments:
34c   ----------
35
36      integer,intent(in) :: nid        ! NetCDF file ID
37      INTEGER*4,intent(in) :: idayref  ! date (initial date for this run)
38      real,intent(in) :: phis(nbplon,nbp_lat) ! surface geopotential
39      real,intent(in) :: area(nbplon,nbp_lat) ! mesh area (m2)
40      integer,intent(in) :: nbplon,nbplat ! sizes of area and phis arrays
41
42c   Local:
43c   ------
44    !   INTEGER length,l
45    !   parameter (length = 100)
46    !   REAL tab_cntrl(length) ! run parameters are stored in this array
47      INTEGER ierr
48      REAl,ALLOCATABLE :: lon_reg_ext(:) ! extended longitudes
49
50      integer :: nvarid,idim_index,idim_rlonv
51      integer :: idim_rlatu,idim_llmp1,idim_llm
52      integer :: idim_nsoilmx ! "subsurface_layers" dimension ID #
53      integer, dimension(2) :: id
54c-----------------------------------------------------------------------
55
56      IF (nbp_lon*nbp_lat==1) THEN
57        ! 1D model
58        ALLOCATE(lon_reg_ext(1))
59      ELSE
60        ! 3D model
61        ALLOCATE(lon_reg_ext(nbp_lon+1))
62      ENDIF
63
64    !   DO l=1,length
65    !      tab_cntrl(l)=0.
66    !   ENDDO
67    !   tab_cntrl(1)  = real(nbp_lon)
68    !   tab_cntrl(2)  = real(nbp_lat-1)
69    !   tab_cntrl(3)  = real(nbp_lev)
70    !   tab_cntrl(4)  = real(idayref)
71    !   tab_cntrl(5)  = rad
72    !   tab_cntrl(6)  = omeg
73    !   tab_cntrl(7)  = g
74    !   tab_cntrl(8)  = mugaz
75    !   tab_cntrl(9)  = rcp
76    !   tab_cntrl(10) = daysec
77    !   tab_cntrl(11) = dtphys
78!      tab_cntrl(12) = etot0
79!      tab_cntrl(13) = ptot0
80!      tab_cntrl(14) = ztot0
81!      tab_cntrl(15) = stot0
82!      tab_cntrl(16) = ang0
83c
84c    ..........    P.Le Van  ( ajout le 8/04/96 )    .........
85c         .....        parametres  pour le zoom          ......
86!      tab_cntrl(17)  = clon
87!      tab_cntrl(18)  = clat
88!      tab_cntrl(19)  = grossismx
89!      tab_cntrl(20)  = grossismy
90c
91c     .....   ajout  le 6/05/97 et le 15/10/97  .......
92c
93!      IF ( fxyhypb )   THEN
94!        tab_cntrl(21) = 1.
95!        tab_cntrl(22) = dzoomx
96!        tab_cntrl(23) = dzoomy
97!      ELSE
98!        tab_cntrl(21) = 0.
99!        tab_cntrl(22) = dzoomx
100!        tab_cntrl(23) = dzoomy
101!        tab_cntrl(24) = 0.
102!        IF( ysinus )  tab_cntrl(24) = 1.
103!      ENDIF
104
105c    .........................................................
106
107! Define dimensions
108
109      ierr = NF_REDEF (nid)
110
111      ierr = NF_DEF_DIM (nid, "index", SIZE(tab_cntrl_mod), idim_index)
112!      ierr = NF_DEF_DIM (nid, "rlonu", iip1, idim_rlonu)
113      ierr = NF_DEF_DIM (nid, "latitude", nbp_lat, idim_rlatu)
114      IF (nbp_lon*nbp_lat==1) THEN
115        ierr = NF_DEF_DIM (nid, "longitude", 1, idim_rlonv)
116      ELSE
117        ierr = NF_DEF_DIM (nid, "longitude", nbp_lon+1, idim_rlonv)
118      ENDIF
119!      ierr = NF_DEF_DIM (nid, "rlatv", jjm, idim_rlatv)
120      ierr = NF_DEF_DIM (nid, "interlayer", (nbp_lev+1), idim_llmp1)
121      ierr = NF_DEF_DIM (nid, "altitude", nbp_lev, idim_llm)
122      ierr = NF_DEF_DIM (nid,"subsurface_layers",nsoilmx,idim_nsoilmx)
123c
124      ierr = NF_ENDDEF(nid)
125
126c  Contol parameters for this run
127      ierr = NF_REDEF (nid)
128#ifdef NC_DOUBLE
129      ierr = NF_DEF_VAR (nid, "controle", NF_DOUBLE, 1,
130     .       idim_index,nvarid)
131#else
132      ierr = NF_DEF_VAR (nid, "controle", NF_FLOAT, 1,
133     .       idim_index,nvarid)
134#endif
135      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 18,
136     .                       "Control parameters")
137      ierr = NF_ENDDEF(nid)
138#ifdef NC_DOUBLE
139      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,tab_cntrl_mod)
140#else
141      ierr = NF_PUT_VAR_REAL (nid,nvarid,tab_cntrl_mod)
142#endif
143
144c --------------------------
145c  longitudes and latitudes
146!
147!      ierr = NF_REDEF (nid)
148!#ifdef NC_DOUBLE
149!      ierr = NF_DEF_VAR (nid, "rlonu", NF_DOUBLE, 1, idim_rlonu,nvarid)
150!#else
151!      ierr = NF_DEF_VAR (nid, "rlonu", NF_FLOAT, 1, idim_rlonu,nvarid)
152!#endif
153!      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 21,
154!     .                       "Longitudes at u nodes")
155!      ierr = NF_ENDDEF(nid)
156!#ifdef NC_DOUBLE
157!      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlonu/pi*180)
158!#else
159!      ierr = NF_PUT_VAR_REAL (nid,nvarid,rlonu/pi*180)
160!#endif
161c
162c --------------------------
163      ierr = NF_REDEF (nid)
164#ifdef NC_DOUBLE
165      ierr =NF_DEF_VAR(nid, "latitude", NF_DOUBLE, 1, idim_rlatu,nvarid)
166#else
167      ierr =NF_DEF_VAR (nid, "latitude", NF_FLOAT, 1, idim_rlatu,nvarid)
168#endif
169      ierr =NF_PUT_ATT_TEXT(nid,nvarid,'units',13,"degrees_north")
170      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name", 14,
171     .      "North latitude")
172      ierr = NF_ENDDEF(nid)
173#ifdef NC_DOUBLE
174      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lat_reg/pi*180)
175#else
176      ierr = NF_PUT_VAR_REAL (nid,nvarid,lat_reg/pi*180)
177#endif
178c
179c --------------------------
180
181      lon_reg_ext(1:nbp_lon)=lon_reg(1:nbp_lon)
182      IF (nbp_lon*nbp_lat/=1) THEN
183        ! In 3D, add extra redundant point (180 degrees,
184        ! since lon_reg starts at -180)
185        lon_reg_ext(nbp_lon+1)=-lon_reg_ext(1)
186      ENDIF
187
188      ierr = NF_REDEF (nid)
189#ifdef NC_DOUBLE
190      ierr =NF_DEF_VAR(nid,"longitude", NF_DOUBLE, 1, idim_rlonv,nvarid)
191#else
192      ierr = NF_DEF_VAR(nid,"longitude", NF_FLOAT, 1, idim_rlonv,nvarid)
193#endif
194      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name", 14,
195     .      "East longitude")
196      ierr = NF_PUT_ATT_TEXT(nid,nvarid,'units',12,"degrees_east")
197      ierr = NF_ENDDEF(nid)
198#ifdef NC_DOUBLE
199      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,lon_reg_ext/pi*180)
200#else
201      ierr = NF_PUT_VAR_REAL (nid,nvarid,lon_reg_ext/pi*180)
202#endif
203c
204c --------------------------
205      ierr = NF_REDEF (nid)
206#ifdef NC_DOUBLE
207      ierr = NF_DEF_VAR (nid, "altitude", NF_DOUBLE, 1,
208     .       idim_llm,nvarid)
209#else
210      ierr = NF_DEF_VAR (nid, "altitude", NF_FLOAT, 1,
211     .       idim_llm,nvarid)
212#endif
213      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"long_name",10,"pseudo-alt")
214      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'units',2,"km")
215      ierr = NF_PUT_ATT_TEXT (nid,nvarid,'positive',2,"up")
216
217      ierr = NF_ENDDEF(nid)
218#ifdef NC_DOUBLE
219      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,pseudoalt)
220#else
221      ierr = NF_PUT_VAR_REAL (nid,nvarid,pseudoalt)
222#endif
223c
224c --------------------------
225!      ierr = NF_REDEF (nid)
226!#ifdef NC_DOUBLE
227!      ierr = NF_DEF_VAR (nid, "rlatv", NF_DOUBLE, 1, idim_rlatv,nvarid)
228!#else
229!      ierr = NF_DEF_VAR (nid, "rlatv", NF_FLOAT, 1, idim_rlatv,nvarid)
230!#endif
231!      ierr = NF_PUT_ATT_TEXT (nid,nvarid,"title", 20,
232!     .                       "Latitudes at v nodes")
233!      ierr = NF_ENDDEF(nid)
234!#ifdef NC_DOUBLE
235!      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,rlatv/pi*180)
236!#else
237!      ierr = NF_PUT_VAR_REAL (nid,nvarid,rlatv/pi*180)
238!#endif
239c
240c --------------------------
241c  Vertical levels
242      call def_var(nid,"aps","hybrid pressure at midlayers ","Pa",
243     .            1,idim_llm,nvarid,ierr)
244#ifdef NC_DOUBLE
245      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,aps)
246#else
247      ierr = NF_PUT_VAR_REAL (nid,nvarid,aps)
248#endif
249
250      call def_var(nid,"bps","hybrid sigma at midlayers"," ",
251     .            1,idim_llm,nvarid,ierr)
252#ifdef NC_DOUBLE
253      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bps)
254#else
255      ierr = NF_PUT_VAR_REAL (nid,nvarid,bps)
256#endif
257
258      call def_var(nid,"ap","hybrid pressure at interlayers","Pa",
259     .            1,idim_llmp1,nvarid,ierr)
260#ifdef NC_DOUBLE
261      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,ap)
262#else
263      ierr = NF_PUT_VAR_REAL (nid,nvarid,ap)
264#endif
265
266      call def_var(nid,"bp","hybrid sigma at interlayers"," ",
267     .            1,idim_llmp1,nvarid,ierr)
268#ifdef NC_DOUBLE
269      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,bp)
270#else
271      ierr = NF_PUT_VAR_REAL (nid,nvarid,bp)
272#endif
273
274!-------------------------------
275! (soil) depth variable mlayer() (known from comsoil_h)
276!-------------------------------
277      ierr=NF_REDEF (nid) ! Enter NetCDF (re-)define mode
278      ! define variable
279#ifdef NC_DOUBLE
280      ierr=NF_DEF_VAR(nid,"soildepth",NF_DOUBLE,1,idim_nsoilmx,nvarid)
281#else
282      ierr=NF_DEF_VAR(nid,"soildepth",NF_FLOAT,1,idim_nsoilmx,nvarid)
283#endif
284      ierr=NF_PUT_ATT_TEXT (nid,nvarid,"long_name", 20,
285     .                        "Soil mid-layer depth")
286      ierr=NF_PUT_ATT_TEXT (nid,nvarid,"units",1,"m")
287      ierr=NF_PUT_ATT_TEXT (nid,nvarid,"positive",4,"down")
288      ierr=NF_ENDDEF(nid) ! Leave NetCDF define mode
289      ! write variable
290#ifdef NC_DOUBLE
291      ierr=NF_PUT_VAR_DOUBLE (nid,nvarid,mlayer)
292#else
293      ierr=NF_PUT_VAR_REAL (nid,nvarid,mlayer)
294#endif
295
296c
297c --------------------------
298c  Mesh area and conversion coefficients cov. <-> contra. <--> natural
299
300!      id(1)=idim_rlonu
301!      id(2)=idim_rlatu
302c
303!      ierr = NF_REDEF (nid)
304!#ifdef NC_DOUBLE
305!      ierr = NF_DEF_VAR (nid, "cu", NF_DOUBLE, 2, id,nvarid)
306!#else
307!      ierr = NF_DEF_VAR (nid, "cu", NF_FLOAT, 2, id,nvarid)
308!#endif
309!      ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 40,
310!     .             "Conversion coefficients cov <--> natural")
311!      ierr = NF_ENDDEF(nid)
312!#ifdef NC_DOUBLE
313!      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,cu)
314!#else
315!      ierr = NF_PUT_VAR_REAL (nid,nvarid,cu)
316!#endif
317c
318!      id(1)=idim_rlonv
319!      id(2)=idim_rlatv
320c
321c --------------------------
322!      ierr = NF_REDEF (nid)
323!#ifdef NC_DOUBLE
324!      ierr = NF_DEF_VAR (nid, "cv", NF_DOUBLE, 2, id,nvarid)
325!#else
326!      ierr = NF_DEF_VAR (nid, "cv", NF_FLOAT, 2, id,nvarid)
327!#endif
328!      ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 40,
329!     .             "Conversion coefficients cov <--> natural")
330!      ierr = NF_ENDDEF(nid)
331!#ifdef NC_DOUBLE
332!      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,cv)
333!#else
334!      ierr = NF_PUT_VAR_REAL (nid,nvarid,cv)
335!#endif
336c
337      id(1)=idim_rlonv
338      id(2)=idim_rlatu
339c
340c --------------------------
341      ierr = NF_REDEF (nid)
342#ifdef NC_DOUBLE
343      ierr = NF_DEF_VAR (nid, "aire", NF_DOUBLE, 2, id,nvarid)
344#else
345      ierr = NF_DEF_VAR (nid, "aire", NF_FLOAT, 2, id,nvarid)
346#endif
347      ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 9,
348     .                       "Mesh area")
349      ierr = NF_ENDDEF(nid)
350#ifdef NC_DOUBLE
351      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,area)
352#else
353      ierr = NF_PUT_VAR_REAL (nid,nvarid,area)
354#endif
355c
356c  Surface geopotential
357      id(1)=idim_rlonv
358      id(2)=idim_rlatu
359c
360      ierr = NF_REDEF (nid)
361#ifdef NC_DOUBLE
362      ierr = NF_DEF_VAR (nid, "phisinit", NF_DOUBLE, 2, id,nvarid)
363#else
364      ierr = NF_DEF_VAR (nid, "phisinit", NF_FLOAT, 2, id,nvarid)
365#endif
366      ierr = NF_PUT_ATT_TEXT (nid, nvarid, "title", 27,
367     .                       "Geopotential at the surface")
368      ierr = NF_ENDDEF(nid)
369#ifdef NC_DOUBLE
370      ierr = NF_PUT_VAR_DOUBLE (nid,nvarid,phis)
371#else
372      ierr = NF_PUT_VAR_REAL (nid,nvarid,phis)
373#endif
374c
375
376      write(*,*)'iniwrite: nbp_lon,nbp_lat,nbp_lev,idayref',
377     & nbp_lon,nbp_lat,nbp_lev,idayref
378      write(*,*)'iniwrite: rad,omeg,g,mugaz,rcp',
379     & rad,omeg,g,mugaz,rcp
380      write(*,*)'iniwrite: daysec,dtphys',daysec,dtphys
381
382      END
Note: See TracBrowser for help on using the repository browser.