1 | !======================================================================= |
---|
2 | |
---|
3 | ! Auteur: F. Hourdin |
---|
4 | ! ------- |
---|
5 | |
---|
6 | ! Objet: |
---|
7 | ! ------ |
---|
8 | ! Light interface for netcdf outputs. can be used outside LMDZ |
---|
9 | |
---|
10 | !======================================================================= |
---|
11 | |
---|
12 | MODULE lmdz_iotd |
---|
13 | IMPLICIT NONE; PRIVATE |
---|
14 | PUBLIC iotd_fin, iotd_ecrit, iotd_ini, imax, jmax |
---|
15 | |
---|
16 | INTEGER imax, jmax, lmax, nid |
---|
17 | INTEGER dim_coord(4) |
---|
18 | REAL iotd_ts, iotd_t0 |
---|
19 | INTEGER :: n_names_iotd_def |
---|
20 | CHARACTER*20, DIMENSION(200) :: names_iotd_def |
---|
21 | CHARACTER*20 :: un_nom |
---|
22 | |
---|
23 | !$OMP THREADPRIVATE(imax, jmax, lmax, nid, dim_coord, iotd_t0, iotd_ts) |
---|
24 | !$OMP THREADPRIVATE(n_names_iotd_def, names_iotd_def) |
---|
25 | CONTAINS |
---|
26 | SUBROUTINE iotd_fin |
---|
27 | USE netcdf, ONLY: nf90_close |
---|
28 | IMPLICIT NONE |
---|
29 | INTEGER ierr |
---|
30 | |
---|
31 | ierr = nf90_close(nid) |
---|
32 | END SUBROUTINE iotd_fin |
---|
33 | |
---|
34 | SUBROUTINE iotd_ecrit(nom, llm, titre, unite, px) |
---|
35 | !----------------------------------------------------------------------- |
---|
36 | ! ---------- |
---|
37 | ! nom : nom de la variable a sortir (chaine de caracteres) |
---|
38 | ! llm : nombre de couches |
---|
39 | ! titre: titre de la variable (chaine de caracteres) |
---|
40 | ! unite : unite de la variable (chaine de caracteres) |
---|
41 | ! px : variable a sortir |
---|
42 | !================================================================= |
---|
43 | |
---|
44 | USE netcdf, ONLY: nf90_put_var, nf90_inq_varid, nf90_enddef, nf90_redef, nf90_sync, nf90_noerr, & |
---|
45 | nf90_float, nf90_def_var |
---|
46 | IMPLICIT NONE |
---|
47 | |
---|
48 | ! Arguments on input: |
---|
49 | INTEGER llm |
---|
50 | CHARACTER (LEN = *) :: nom, titre, unite |
---|
51 | INTEGER imjmax |
---|
52 | parameter (imjmax = 100000) |
---|
53 | REAL px(imjmax * llm) |
---|
54 | |
---|
55 | ! Local variables: |
---|
56 | |
---|
57 | real(kind = 4) date |
---|
58 | real(kind = 4) zx(imjmax * llm) |
---|
59 | |
---|
60 | INTEGER ierr, ndim, dim_cc(4) |
---|
61 | INTEGER iq |
---|
62 | INTEGER i, j, l |
---|
63 | |
---|
64 | INTEGER zitau |
---|
65 | CHARACTER firstnom*20 |
---|
66 | SAVE firstnom |
---|
67 | SAVE zitau |
---|
68 | SAVE date |
---|
69 | DATA firstnom /'1234567890'/ |
---|
70 | DATA zitau /0/ |
---|
71 | |
---|
72 | ! Ajouts |
---|
73 | INTEGER, save :: ntime = 0 |
---|
74 | INTEGER :: idim, varid |
---|
75 | CHARACTER (LEN = 50) :: fichnom |
---|
76 | INTEGER, DIMENSION(4) :: id |
---|
77 | INTEGER, DIMENSION(4) :: edges, corner |
---|
78 | |
---|
79 | IF (n_names_iotd_def>0 .and..not.any(names_iotd_def==nom)) RETURN |
---|
80 | !*************************************************************** |
---|
81 | ! Initialisation of 'firstnom' and create/open the "diagfi.nc" NetCDF file |
---|
82 | ! ------------------------------------------------------------------------ |
---|
83 | ! (Au tout premier appel de la SUBROUTINE durant le run.) |
---|
84 | |
---|
85 | |
---|
86 | !-------------------------------------------------------- |
---|
87 | ! Write the variables to output file if it's time to do so |
---|
88 | !-------------------------------------------------------- |
---|
89 | |
---|
90 | |
---|
91 | ! Compute/write/extend 'time' coordinate (date given in days) |
---|
92 | ! (done every "first call" (at given time level) to writediagfi) |
---|
93 | ! Note: date is incremented as 1 step ahead of physics time |
---|
94 | !-------------------------------------------------------- |
---|
95 | |
---|
96 | zx(1:imax * jmax * llm) = px(1:imax * jmax * llm) |
---|
97 | IF (firstnom =='1234567890') THEN |
---|
98 | firstnom = nom |
---|
99 | endif |
---|
100 | |
---|
101 | !PRINT*,'nom ',nom,firstnom |
---|
102 | |
---|
103 | !! Quand on tombe sur la premiere variable on ajoute un pas de temps |
---|
104 | IF (nom==firstnom) THEN |
---|
105 | ! We have identified a "first call" (at given date) |
---|
106 | |
---|
107 | ntime = ntime + 1 ! increment # of stored time steps |
---|
108 | |
---|
109 | !! PRINT*,'ntime ',ntime |
---|
110 | date = iotd_t0 + ntime * iotd_ts |
---|
111 | !PRINT*,'iotd_ecrit ',iotd_ts,ntime, date |
---|
112 | ! date= float (zitau +1)/float (day_step) |
---|
113 | |
---|
114 | ! compute corresponding date (in days and fractions thereof) |
---|
115 | ! Get NetCDF ID of 'time' variable |
---|
116 | |
---|
117 | ierr = nf90_sync(nid) |
---|
118 | |
---|
119 | ierr = nf90_inq_varid(nid, "time", varid) |
---|
120 | ! Write (append) the new date to the 'time' array |
---|
121 | |
---|
122 | ierr = nf90_put_var(nid, varid, date, [ntime]) |
---|
123 | |
---|
124 | ! PRINT*,'date ',date,ierr,nid |
---|
125 | ! PRINT*,'IOTD Date ,varid,nid,ntime,date',varid,nid,ntime,date |
---|
126 | |
---|
127 | IF (ierr/=nf90_noerr) THEN |
---|
128 | WRITE(*, *) "***** PUT_VAR matter in writediagfi_nc" |
---|
129 | WRITE(*, *) "***** with time" |
---|
130 | WRITE(*, *) 'ierr=', ierr |
---|
131 | endif |
---|
132 | |
---|
133 | ! WRITE(6,*)'WRITEDIAGFI: date= ', date |
---|
134 | end if ! of if (nom.EQ.firstnom) |
---|
135 | |
---|
136 | |
---|
137 | !Case of a 3D variable |
---|
138 | !--------------------- |
---|
139 | IF (llm==lmax) THEN |
---|
140 | ndim = 4 |
---|
141 | corner(1) = 1 |
---|
142 | corner(2) = 1 |
---|
143 | corner(3) = 1 |
---|
144 | corner(4) = ntime |
---|
145 | edges(1) = imax |
---|
146 | edges(2) = jmax |
---|
147 | edges(3) = llm |
---|
148 | edges(4) = 1 |
---|
149 | dim_cc = dim_coord |
---|
150 | |
---|
151 | |
---|
152 | !Case of a 2D variable |
---|
153 | !--------------------- |
---|
154 | |
---|
155 | ELSE IF (llm==1) THEN |
---|
156 | ndim = 3 |
---|
157 | corner(1) = 1 |
---|
158 | corner(2) = 1 |
---|
159 | corner(3) = ntime |
---|
160 | corner(4) = 1 |
---|
161 | edges(1) = imax |
---|
162 | edges(2) = jmax |
---|
163 | edges(3) = 1 |
---|
164 | edges(4) = 1 |
---|
165 | dim_cc(1:2) = dim_coord(1:2) |
---|
166 | dim_cc(3) = dim_coord(4) |
---|
167 | |
---|
168 | END IF ! of if llm=1 ou llm |
---|
169 | |
---|
170 | ! AU premier pas de temps, on crée les variables |
---|
171 | !----------------------------------------------- |
---|
172 | |
---|
173 | IF (ntime==1) THEN |
---|
174 | ierr = nf90_redef (nid) |
---|
175 | ierr = nf90_def_var(nid, nom, nf90_float, dim_cc, varid) |
---|
176 | !PRINT*,'DEF ',nom,nid,varid |
---|
177 | ierr = nf90_enddef(nid) |
---|
178 | ELSE |
---|
179 | ierr = nf90_inq_varid(nid, nom, varid) |
---|
180 | !PRINT*,'INQ ',nom,nid,varid |
---|
181 | ! Commandes pour recuperer automatiquement les coordonnees |
---|
182 | ! ierr= nf90_inq_dimid(nid,"longitude",id(1)) |
---|
183 | END IF |
---|
184 | |
---|
185 | ierr = nf90_put_var(nid, varid, zx, corner, edges) |
---|
186 | |
---|
187 | IF (ierr/=nf90_noerr) THEN |
---|
188 | WRITE(*, *) "***** PUT_VAR problem in writediagfi" |
---|
189 | WRITE(*, *) "***** with ", nom |
---|
190 | WRITE(*, *) 'ierr=', ierr |
---|
191 | endif |
---|
192 | |
---|
193 | END |
---|
194 | |
---|
195 | SUBROUTINE iotd_ini(fichnom, iim, jjm, llm, prlon, prlat, pcoordv, jour0, mois0, an0, t0, timestep, calendrier) |
---|
196 | USE netcdf, ONLY: nf90_enddef, nf90_put_att, nf90_float, nf90_def_var, nf90_redef, & |
---|
197 | nf90_global, nf90_def_dim, nf90_create, nf90_clobber, nf90_unlimited, nf90_put_var |
---|
198 | IMPLICIT NONE |
---|
199 | |
---|
200 | INTEGER iim, jjm, llm |
---|
201 | REAL prlon(iim), prlat(jjm), pcoordv(llm), timestep, t0 |
---|
202 | INTEGER id_FOCE |
---|
203 | INTEGER jour0, mois0, an0 |
---|
204 | CHARACTER*(*) calendrier |
---|
205 | |
---|
206 | INTEGER corner(4), edges(4), ndim |
---|
207 | real px(1000) |
---|
208 | CHARACTER (LEN = 10) :: nom |
---|
209 | real(kind = 4) rlon(iim), rlat(jjm), coordv(llm) |
---|
210 | |
---|
211 | ! Local: |
---|
212 | ! ------ |
---|
213 | CHARACTER*3, DIMENSION(12) :: cmois = (/'JAN', 'FEB', 'MAR', 'APR', 'MAY', 'JUN', 'JUL', 'AUG', 'SEP', 'OCT', 'NOV', 'DEC'/) |
---|
214 | CHARACTER*10 date0 |
---|
215 | CHARACTER*11 date0b |
---|
216 | |
---|
217 | INTEGER :: ierr |
---|
218 | |
---|
219 | INTEGER :: nvarid |
---|
220 | INTEGER, DIMENSION(2) :: id |
---|
221 | |
---|
222 | CHARACTER*(*) fichnom |
---|
223 | |
---|
224 | REAL pi |
---|
225 | |
---|
226 | iotd_ts = timestep |
---|
227 | iotd_t0 = t0 |
---|
228 | PRINT*, 'iotd_ini, ', timestep, iotd_ts |
---|
229 | imax = iim |
---|
230 | jmax = jjm |
---|
231 | lmax = llm |
---|
232 | ! Utile pour passer en real*4 pour les ecritures |
---|
233 | rlon = prlon |
---|
234 | rlat = prlat |
---|
235 | coordv = pcoordv |
---|
236 | |
---|
237 | |
---|
238 | !----------------------------------------------------------------------- |
---|
239 | ! Possibilité de spécifier une liste de variables à sortir |
---|
240 | ! dans iotd.def |
---|
241 | ! Si iotd.def existe et est non vide, |
---|
242 | ! seules les variables faisant à la fois l'objet d'un CALL iotd_ecrit |
---|
243 | ! et étant spécifiées dans iotd.def sont sorties. |
---|
244 | ! Sinon, toutes les variables faisant l'objet d'un CALL iotd_ecrit |
---|
245 | ! sont sorties |
---|
246 | !----------------------------------------------------------------------- |
---|
247 | n_names_iotd_def = 0 |
---|
248 | open(99, file = 'iotd.def', form = 'formatted', status = 'old', iostat = ierr) |
---|
249 | IF (ierr==0) THEN |
---|
250 | ierr = 0 |
---|
251 | do while (ierr==0) |
---|
252 | read(99, *, iostat = ierr) un_nom |
---|
253 | IF (ierr==0) THEN |
---|
254 | n_names_iotd_def = n_names_iotd_def + 1 |
---|
255 | names_iotd_def(n_names_iotd_def) = un_nom |
---|
256 | endif |
---|
257 | enddo |
---|
258 | endif |
---|
259 | PRINT*, n_names_iotd_def, names_iotd_def(1:n_names_iotd_def) |
---|
260 | close(99) |
---|
261 | |
---|
262 | pi = 2. * asin(1.) |
---|
263 | |
---|
264 | ! Define dimensions |
---|
265 | |
---|
266 | ! Create the NetCDF file |
---|
267 | ierr = nf90_create(fichnom, nf90_clobber, nid) |
---|
268 | ierr = nf90_def_dim(nid, "lon", iim, dim_coord(1)) |
---|
269 | ierr = nf90_def_dim(nid, "lat", jjm, dim_coord(2)) |
---|
270 | ierr = nf90_def_dim(nid, "lev", llm, dim_coord(3)) |
---|
271 | ierr = nf90_def_dim(nid, "time", nf90_unlimited, dim_coord(4)) |
---|
272 | ierr = nf90_put_att(nid, nf90_global, 'Conventions', "CF-1.1") |
---|
273 | !ierr = nf90_put_att(nid,nf90_global,'file_name',TRIM(fname)) |
---|
274 | ierr = nf90_enddef(nid) |
---|
275 | |
---|
276 | ! Switch out of NetCDF Define mode |
---|
277 | |
---|
278 | ierr = nf90_enddef(nid) |
---|
279 | |
---|
280 | ! Contol parameters for this run |
---|
281 | ! ---- longitude ----------- |
---|
282 | |
---|
283 | ierr = nf90_redef(nid) |
---|
284 | ierr = nf90_def_var(nid, "lon", nf90_float, dim_coord(1), nvarid) |
---|
285 | ierr = nf90_put_att(nid, nvarid, 'axis', 'X') |
---|
286 | ierr = nf90_put_att(nid, nvarid, 'units', "degrees_east") |
---|
287 | ierr = nf90_enddef(nid) |
---|
288 | ierr = nf90_put_var(nid, nvarid, rlon) |
---|
289 | PRINT*, ierr |
---|
290 | |
---|
291 | ! ---- latitude ------------ |
---|
292 | ierr = nf90_redef(nid) |
---|
293 | ierr = nf90_def_var(nid, "lat", nf90_float, dim_coord(2), nvarid) |
---|
294 | ierr = nf90_put_att(nid, nvarid, 'axis', 'Y') |
---|
295 | ierr = nf90_put_att(nid, nvarid, 'units', "degrees_north") |
---|
296 | ierr = nf90_enddef(nid) |
---|
297 | ierr = nf90_put_var(nid, nvarid, rlat) |
---|
298 | |
---|
299 | ! ---- vertical ------------ |
---|
300 | ierr = nf90_redef(nid) |
---|
301 | ierr = nf90_def_var(nid, "lev", nf90_float, dim_coord(3), nvarid) |
---|
302 | ierr = nf90_put_att(nid, nvarid, "long_name", "vert level") |
---|
303 | IF (coordv(2)>coordv(1)) THEN |
---|
304 | ierr = nf90_put_att(nid, nvarid, "long_name", "pseudo-alt") |
---|
305 | ierr = nf90_put_att(nid, nvarid, 'positive', "up") |
---|
306 | else |
---|
307 | ierr = nf90_put_att(nid, nvarid, "long_name", "pressure") |
---|
308 | ierr = nf90_put_att(nid, nvarid, 'positive', "down") |
---|
309 | endif |
---|
310 | ierr = nf90_enddef(nid) |
---|
311 | ierr = nf90_put_var(nid, nvarid, coordv) |
---|
312 | |
---|
313 | ! ---- time ---------------- |
---|
314 | ierr = nf90_redef(nid) |
---|
315 | ! Define the 'time' variable |
---|
316 | ierr = nf90_def_var(nid, "time", nf90_float, dim_coord(4), nvarid) |
---|
317 | ! ! Add attributes |
---|
318 | ierr = nf90_put_att(nid, nvarid, 'axis', 'T') |
---|
319 | ierr = nf90_put_att(nid, nvarid, 'standard_name', 'time') |
---|
320 | WRITE(date0, '(i4.4,"-",i2.2,"-",i2.2)') an0, mois0, jour0 |
---|
321 | ierr = nf90_put_att(nid, nvarid, 'units', & |
---|
322 | "seconds since " // date0 // " 00:00:00") |
---|
323 | ierr = nf90_put_att(nid, nvarid, 'calendar', calendrier) |
---|
324 | !ierr = nf90_put_att(nid,nvarid,'calendar','360d') |
---|
325 | ierr = nf90_put_att(nid, nvarid, 'title', 'Time') |
---|
326 | ierr = nf90_put_att(nid, nvarid, 'long_name', 'Time axis') |
---|
327 | WRITE(date0b, '(i4.4,"-",a3,"-",i2.2)') an0, cmois(mois0), jour0 |
---|
328 | ierr = nf90_put_att(nid, nvarid, 'time_origin', & |
---|
329 | date0b // ' 00:00:00') |
---|
330 | ierr = nf90_enddef(nid) |
---|
331 | |
---|
332 | END |
---|
333 | |
---|
334 | END MODULE lmdz_iotd |
---|