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

Last change on this file since 374 was 135, checked in by aslmd, 14 years ago

CHANGEMENT ARBORESCENCE ETAPE 2 -- NON COMPLET

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