source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/iniwrite.F @ 264

Last change on this file since 264 was 222, checked in by ymipsl, 10 years ago

Creating temporary dynamico/lmdz/saturn branche

YM

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