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

Last change on this file since 1477 was 1422, checked in by milmd, 10 years ago

In GENERIC, MARS and COMMON models replace some include files by modules (usefull for decoupling physics with dynamics).

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