source: trunk/LMDZ.MARS/libf/phymars/iniwrite.F @ 1233

Last change on this file since 1233 was 1226, checked in by aslmd, 11 years ago

LMDZ.MARS : Replaced comcstfi and planete includes by modules.

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