1 | !> \file io_netcdf_GRISLI.f90 |
---|
2 | !! NetCDF Fortran 90 read/write interface |
---|
3 | !< |
---|
4 | |
---|
5 | !> \namespace io_netcdf_grisli |
---|
6 | !! NetCDF Fortran 90 read/write interface |
---|
7 | !! @note Use input/output functions provided by unidata |
---|
8 | !! |
---|
9 | !! @note Glossary of variables |
---|
10 | !! @note varname = name of variable to trait |
---|
11 | !! @note file,file_in,file_out = netcdf file name |
---|
12 | !! @note tabvar = array containing values of the required variable |
---|
13 | !! @note typevar = type of the variable to trait |
---|
14 | !! @note dimIDs = array containing dimensions of variables |
---|
15 | !! @note dim'i' = dimensions of variables |
---|
16 | !! @note varid,varid_in,varid_out = variables identificator |
---|
17 | !! @note ncid,ncid_in,ncid_out = netcdf file identificator |
---|
18 | !! @note time = time corresponding to the values to read/write |
---|
19 | !! @note level = level corresponding to the values to read/write |
---|
20 | !! @note dimval = value of the required dimension |
---|
21 | !! @note dimname = name of dimension to retrieve |
---|
22 | !! |
---|
23 | !! http://my.unidata.ucar.edu/content/software/netcdf/docs/netcdf-f90/index.html |
---|
24 | !! |
---|
25 | !< |
---|
26 | |
---|
27 | module io_netcdf_grisli |
---|
28 | ! |
---|
29 | use netcdf |
---|
30 | use runparam,only:dirsource |
---|
31 | !> \interface Read_Ncdf_var |
---|
32 | !! Interfaces of functions and subroutines to read netcdf variables |
---|
33 | !!\author ... |
---|
34 | !!\date 2010 |
---|
35 | !! |
---|
36 | !< |
---|
37 | integer :: ncdf_type = 0 ! pour changer eventuellement de netcdf |
---|
38 | |
---|
39 | interface Read_Ncdf_var |
---|
40 | module procedure Read_Ncdf_var1d_Real, & ! Reading a 1d real variable |
---|
41 | Read_Ncdf_var2d_Real, & ! Reading a 2d real variable |
---|
42 | Read_Ncdf_var2d_Real_bis, & ! Reading a 2d real var with selection of data |
---|
43 | Read_Ncdf_var3d_Real, & ! Reading a 3d real variable |
---|
44 | Read_Ncdf_var4d_Real, & ! Reading a 4d real variable |
---|
45 | Read_Ncdf_var1d_Int, & ! Reading a 1d integer variable |
---|
46 | Read_Ncdf_var2d_Int, & ! Reading a 2d integer variable |
---|
47 | Read_Ncdf_var3d_Int, & ! Reading a 3d integer variable |
---|
48 | Read_Ncdf_var4d_Int, & ! Reading a 4d integer variable |
---|
49 | Read_Ncdf_var3d_Real_t, & ! Reading a 2d variable depending on time |
---|
50 | Read_Ncdf_var4d_Real_t, & ! Reading a 3d variable depending on time |
---|
51 | Read_Ncdf_var4d_Real_nt ! Reading a 2d variable at a time t and a level n |
---|
52 | end interface |
---|
53 | !************************************************************************************************************ |
---|
54 | |
---|
55 | !> \interface Write_Ncdf_var |
---|
56 | !! Interfaces of functions and subroutine to write netcdf variables |
---|
57 | !!\author ... |
---|
58 | !!\date 2010 |
---|
59 | !! |
---|
60 | !< |
---|
61 | interface Write_Ncdf_var |
---|
62 | module procedure Write_Ncdf_var1d_Real, & ! Writing a 1d real (simple/double) variable |
---|
63 | Write_Ncdf_var2d_Real, & ! Writing a 2d real (simple/double) variable |
---|
64 | Write_Ncdf_var3d_Real, & ! Writing a 3d real (simple/double) variable |
---|
65 | Write_Ncdf_var4d_Real, & ! Writing a 4d real (simple/double) variable |
---|
66 | Write_Ncdf_var1d_Int, & ! Writing a 1d integer variable |
---|
67 | Write_Ncdf_var2d_Int, & ! Writing a 2d integer variable |
---|
68 | Write_Ncdf_var3d_Int, & ! Writing a 3d integer variable |
---|
69 | Write_Ncdf_var4d_Int, & ! Writing a 4d integer variable |
---|
70 | Write_Ncdf_var1d_Int_t, & ! Writing a 1d integer variable depending on time |
---|
71 | Write_Ncdf_var1d_Realbis_t, & ! Writing a 1d real (simple/double) variable depending on time |
---|
72 | Write_Ncdf_var1d_Real_t, & ! Writing a 1d real (simple/double) variable depending on time |
---|
73 | Write_Ncdf_var2d_Int_t, & ! Writing a 2d integer variable depending on time |
---|
74 | Write_Ncdf_var2d_Real_t, & ! Writing a 2d real (simple/double) variable depending on time |
---|
75 | Write_Ncdf_var3d_Real_t, & ! Writing a 3d real (simple/double) variable depending on time |
---|
76 | Write_Ncdf_var4d_Real_t, & ! Writing a 4d real (simple/double) variable depending on time |
---|
77 | Write_Ncdf_var4d_Real_nt ! Writing a 2d variable at a time t and a level n |
---|
78 | |
---|
79 | end interface |
---|
80 | !************************************************************************************************************* |
---|
81 | |
---|
82 | !> \interface Copy_Ncdf_att |
---|
83 | !! Interfaces of functions and subroutines to copy netcdf variables |
---|
84 | !!\author ... |
---|
85 | !!\date 2010 |
---|
86 | !! |
---|
87 | !< |
---|
88 | interface Copy_Ncdf_att |
---|
89 | module procedure Copy_Ncdf_att_latlon,Copy_Ncdf_att_var |
---|
90 | end interface |
---|
91 | |
---|
92 | contains |
---|
93 | |
---|
94 | ! subroutine pour lire netcdf_type |
---|
95 | subroutine lect_netcdf_type |
---|
96 | |
---|
97 | implicit none |
---|
98 | ! open(22,file='../SOURCES/Fichiers-parametres/netcdf_type.dat') |
---|
99 | !open(22,file=trim(dirsource)//'/Fichiers-parametres/netcdf_type.dat') |
---|
100 | open(22,file=trim(dirsource)//'/netcdf_type.dat') |
---|
101 | |
---|
102 | read(22,'(i3)') ncdf_type |
---|
103 | close(22) |
---|
104 | write(6,*) ' lecture de ncdft_type',ncdf_type |
---|
105 | return |
---|
106 | end subroutine lect_netcdf_type |
---|
107 | |
---|
108 | !> Subroutine to retrieve values of a given 1D real variable |
---|
109 | !! @param[in] varname : name of variable to retrieve |
---|
110 | !! @param[in] file : netcdf file name |
---|
111 | !! @param[out] tabvar : array containing values of the required variable |
---|
112 | !! @return tabvar |
---|
113 | |
---|
114 | subroutine Read_Ncdf_var1d_Real(varname,file,tabvar) |
---|
115 | ! |
---|
116 | implicit none |
---|
117 | ! |
---|
118 | Character(*),intent(in) :: varname,file |
---|
119 | Real*8, dimension(:), pointer :: tabvar |
---|
120 | ! |
---|
121 | !local variables |
---|
122 | ! |
---|
123 | Integer, dimension(1) :: dimID |
---|
124 | Integer :: dim1 |
---|
125 | Integer :: status,ncid |
---|
126 | Integer :: varid |
---|
127 | ! |
---|
128 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
129 | |
---|
130 | if (status/=nf90_noerr) then |
---|
131 | write(*,*)"unable to open netcdf file : ",file |
---|
132 | stop |
---|
133 | endif |
---|
134 | ! |
---|
135 | status = nf90_inq_varid(ncid,varname,varid) |
---|
136 | |
---|
137 | if (status/=nf90_noerr) then |
---|
138 | write(*,*)"pb de nom de variable :", varname, file |
---|
139 | stop |
---|
140 | endif |
---|
141 | |
---|
142 | ! recupere l'identite des dimensions de la variable |
---|
143 | status = nf90_inquire_variable(ncid,varid,dimids=dimID) |
---|
144 | |
---|
145 | ! recupere la taille de la dimension de la variable |
---|
146 | status = nf90_inquire_dimension(ncid,dimID(1),len=dim1) |
---|
147 | ! |
---|
148 | if(.not. associated(tabvar)) then |
---|
149 | Allocate(tabvar(dim1)) |
---|
150 | else |
---|
151 | if( any(shape(tabvar)/=(/dim1/)) ) then |
---|
152 | deallocate(tabvar) |
---|
153 | Allocate(tabvar(dim1)) |
---|
154 | endif |
---|
155 | endif |
---|
156 | |
---|
157 | ! lit la variable |
---|
158 | status=nf90_get_var(ncid,varid,tabvar) |
---|
159 | |
---|
160 | ! referme le fichier |
---|
161 | status = nf90_close(ncid) |
---|
162 | ! |
---|
163 | end subroutine Read_Ncdf_var1d_Real |
---|
164 | |
---|
165 | |
---|
166 | |
---|
167 | !> Subroutine to retrieve values of a given 2D real variable |
---|
168 | !! @param[in] varname : name of variable to retrieve |
---|
169 | !! @param[in] file : netcdf file name |
---|
170 | !! @param[out] tabvar : array containing values of the required variable |
---|
171 | !! @return tabvar |
---|
172 | |
---|
173 | subroutine Read_Ncdf_var2d_Real(varname,file,tabvar) |
---|
174 | ! |
---|
175 | implicit none |
---|
176 | ! |
---|
177 | Character(*),intent(in) :: varname,file |
---|
178 | Real*8, dimension(:,:), pointer :: tabvar |
---|
179 | !local variables |
---|
180 | Integer, dimension(10) :: dimIDS |
---|
181 | Integer :: dim1,dim2 |
---|
182 | Integer :: status,ncid |
---|
183 | Integer :: varid |
---|
184 | ! |
---|
185 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
186 | ! |
---|
187 | if (status/=nf90_noerr) then |
---|
188 | write(*,*)"unable to open netcdf file : ",file |
---|
189 | stop |
---|
190 | endif |
---|
191 | ! |
---|
192 | status = nf90_inq_varid(ncid,varname,varid) |
---|
193 | if (status/=nf90_noerr) then |
---|
194 | write(*,*)"pb de nom de variable :", varname, file |
---|
195 | stop |
---|
196 | endif |
---|
197 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
198 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
199 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
200 | ! |
---|
201 | if(.not. associated(tabvar)) then |
---|
202 | Allocate(tabvar(dim1,dim2)) |
---|
203 | else |
---|
204 | if( any(shape(tabvar)/=(/dim1,dim2/)) ) then |
---|
205 | tabvar => null() |
---|
206 | Allocate(tabvar(dim1,dim2)) |
---|
207 | endif |
---|
208 | endif |
---|
209 | ! |
---|
210 | status=nf90_get_var(ncid,varid,tabvar) |
---|
211 | ! |
---|
212 | status = nf90_close(ncid) |
---|
213 | ! |
---|
214 | end subroutine Read_Ncdf_var2d_Real |
---|
215 | |
---|
216 | |
---|
217 | |
---|
218 | !> Subroutine to retrieve values of a given 2D real variable |
---|
219 | !! @param[in] varname : name of variable to retrieve |
---|
220 | !! @param[in] file : netcdf file name |
---|
221 | !! @param[out] tabvar : array containing values of the required variable |
---|
222 | !! @param[in] strt : array of integers specifiying the index in the |
---|
223 | !! variable from which the first (or only) of the |
---|
224 | !! data values will be read. The elements of start |
---|
225 | !! correspond, in order, to the variable's dimensions |
---|
226 | !! @param[in] cnt : A vector of integers specifying the number of |
---|
227 | !! indices selected along each dimension. |
---|
228 | !! @return tabvar |
---|
229 | |
---|
230 | |
---|
231 | subroutine Read_Ncdf_var2d_Real_bis(varname,file,tabvar,strt,cnt) |
---|
232 | ! |
---|
233 | implicit none |
---|
234 | ! |
---|
235 | Character(*),intent(in) :: varname,file |
---|
236 | Real*8, dimension(:,:), pointer :: tabvar |
---|
237 | !local variables |
---|
238 | Integer, dimension(10) :: dimIDS |
---|
239 | Integer, dimension(2) :: strt,cnt |
---|
240 | Integer :: dim1,dim2 |
---|
241 | Integer :: status,ncid |
---|
242 | Integer :: varid |
---|
243 | ! |
---|
244 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
245 | ! |
---|
246 | if (status/=nf90_noerr) then |
---|
247 | write(*,*)"unable to open netcdf file : ",file |
---|
248 | stop |
---|
249 | endif |
---|
250 | ! |
---|
251 | dim1 = cnt(1) |
---|
252 | dim2 = cnt(2) |
---|
253 | ! |
---|
254 | status = nf90_inq_varid(ncid,varname,varid) |
---|
255 | if (status/=nf90_noerr) then |
---|
256 | write(*,*)"pb de nom de variable :", varname, file |
---|
257 | stop |
---|
258 | endif |
---|
259 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
260 | ! |
---|
261 | if(.not. associated(tabvar)) then |
---|
262 | Allocate(tabvar(dim1,dim2)) |
---|
263 | else |
---|
264 | if( any(shape(tabvar)/=(/dim1,dim2/)) ) then |
---|
265 | deallocate(tabvar) |
---|
266 | Allocate(tabvar(dim1,dim2)) |
---|
267 | endif |
---|
268 | endif |
---|
269 | ! |
---|
270 | status=nf90_get_var(ncid,varid,tabvar,start = strt,count = cnt) |
---|
271 | ! |
---|
272 | status = nf90_close(ncid) |
---|
273 | ! |
---|
274 | end subroutine Read_Ncdf_var2d_Real_bis |
---|
275 | |
---|
276 | |
---|
277 | !> Subroutine to retrieve values of a given 3D real variable |
---|
278 | !! @param[in] varname : name of variable to retrieve |
---|
279 | !! @param[in] file : netcdf file name |
---|
280 | !! @param[out] tabvar : array containing values of the required variable |
---|
281 | !! @return tabvar |
---|
282 | subroutine Read_Ncdf_var3d_Real(varname,file,tabvar) |
---|
283 | ! |
---|
284 | implicit none |
---|
285 | ! |
---|
286 | Character(*),intent(in) :: varname,file |
---|
287 | Real*8, dimension(:,:,:), pointer :: tabvar |
---|
288 | ! |
---|
289 | !local variables |
---|
290 | ! |
---|
291 | Integer, dimension(10) :: dimIDS |
---|
292 | Integer :: dim1,dim2,dim3 |
---|
293 | Integer :: status,ncid |
---|
294 | Integer :: varid |
---|
295 | ! |
---|
296 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
297 | ! |
---|
298 | if (status/=nf90_noerr) then |
---|
299 | write(*,*)"unable to open netcdf file : ",file |
---|
300 | stop |
---|
301 | endif |
---|
302 | ! |
---|
303 | status = nf90_inq_varid(ncid,varname,varid) |
---|
304 | if (status/=nf90_noerr) then |
---|
305 | write(*,*)"pb de nom de variable :", varname, file |
---|
306 | stop |
---|
307 | endif |
---|
308 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
309 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
310 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
311 | status=nf90_inquire_dimension(ncid,dimIDS(3),len=dim3) |
---|
312 | ! |
---|
313 | if(.not. associated(tabvar)) then |
---|
314 | Allocate(tabvar(dim1,dim2,dim3)) |
---|
315 | else |
---|
316 | if( any(shape(tabvar) /= (/dim1,dim2,dim3/)) ) then |
---|
317 | deallocate(tabvar) |
---|
318 | Allocate(tabvar(dim1,dim2,dim3)) |
---|
319 | endif |
---|
320 | endif |
---|
321 | ! |
---|
322 | status=nf90_get_var(ncid,varid,tabvar) |
---|
323 | if (status/=nf90_noerr) then |
---|
324 | write(*,*)"unable to retrieve netcdf variable : ",trim(varname) |
---|
325 | stop |
---|
326 | endif |
---|
327 | ! |
---|
328 | status = nf90_close(ncid) |
---|
329 | ! |
---|
330 | end subroutine Read_Ncdf_var3d_Real |
---|
331 | |
---|
332 | |
---|
333 | !> Subroutine to retrieve values of a given 4D real variable |
---|
334 | !! @param[in] varname : name of variable to retrieve |
---|
335 | !! @param[in] file : netcdf file name |
---|
336 | !! @param[out] tabvar : array containing values of the required variable |
---|
337 | !! @return tabvar |
---|
338 | subroutine Read_Ncdf_var4d_Real(varname,file,tabvar) |
---|
339 | ! |
---|
340 | implicit none |
---|
341 | ! |
---|
342 | Character(*),intent(in) :: varname,file |
---|
343 | Real*8, dimension(:,:,:,:), pointer :: tabvar |
---|
344 | ! |
---|
345 | !local variables |
---|
346 | ! |
---|
347 | Integer, dimension(10) :: dimIDS |
---|
348 | Integer :: dim1,dim2,dim3,dim4 |
---|
349 | Integer :: status,ncid |
---|
350 | Integer :: varid |
---|
351 | ! |
---|
352 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
353 | ! |
---|
354 | if (status/=nf90_noerr) then |
---|
355 | write(*,*)"unable to open netcdf file : ",file |
---|
356 | stop |
---|
357 | endif |
---|
358 | ! |
---|
359 | status = nf90_inq_varid(ncid,varname,varid) |
---|
360 | if (status/=nf90_noerr) then |
---|
361 | write(*,*)"pb de nom de variable :", varname, file |
---|
362 | stop |
---|
363 | endif |
---|
364 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
365 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
366 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
367 | status=nf90_inquire_dimension(ncid,dimIDS(3),len=dim3) |
---|
368 | status=nf90_inquire_dimension(ncid,dimIDS(4),len=dim4) |
---|
369 | ! |
---|
370 | if(.not. associated(tabvar)) then |
---|
371 | Allocate(tabvar(dim1,dim2,dim3,dim4)) |
---|
372 | else |
---|
373 | if( any(shape(tabvar) /= (/dim1,dim2,dim3,dim4/)) ) then |
---|
374 | deallocate(tabvar) |
---|
375 | Allocate(tabvar(dim1,dim2,dim3,dim4)) |
---|
376 | endif |
---|
377 | endif |
---|
378 | ! |
---|
379 | status=nf90_get_var(ncid,varid,tabvar) |
---|
380 | ! |
---|
381 | status = nf90_close(ncid) |
---|
382 | ! |
---|
383 | end subroutine Read_Ncdf_var4d_Real |
---|
384 | |
---|
385 | !> Subroutine to retrieve values of a given 1D integer variable |
---|
386 | !! @param[in] varname : name of variable to retrieve |
---|
387 | !! @param[in] file : netcdf file name |
---|
388 | !! @param[out] tabvar : array containing values of the required variable |
---|
389 | !! @return tabvar |
---|
390 | subroutine Read_Ncdf_var1d_Int(varname,file,tabvar) |
---|
391 | ! |
---|
392 | implicit none |
---|
393 | ! |
---|
394 | Character(*),intent(in) :: varname,file |
---|
395 | Integer, dimension(:), pointer :: tabvar |
---|
396 | ! |
---|
397 | !local variables |
---|
398 | ! |
---|
399 | Integer,dimension(10) :: dimID |
---|
400 | Integer :: dim1 |
---|
401 | Integer :: status,ncid |
---|
402 | Integer :: varid |
---|
403 | ! |
---|
404 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
405 | ! |
---|
406 | if (status/=nf90_noerr) then |
---|
407 | write(*,*)"unable to open netcdf file : ",file |
---|
408 | stop |
---|
409 | endif |
---|
410 | ! |
---|
411 | status = nf90_inq_varid(ncid,varname,varid) |
---|
412 | if (status/=nf90_noerr) then |
---|
413 | write(*,*)"pb de nom de variable :", varname, file |
---|
414 | stop |
---|
415 | endif |
---|
416 | status=nf90_inquire_variable(ncid,varid,dimids=dimID) |
---|
417 | status=nf90_inquire_dimension(ncid,dimID(1),len=dim1) |
---|
418 | ! |
---|
419 | if(.not. associated(tabvar)) then |
---|
420 | Allocate(tabvar(dim1)) |
---|
421 | else |
---|
422 | if( any(shape(tabvar) /= (/dim1/)) ) then |
---|
423 | deallocate(tabvar) |
---|
424 | Allocate(tabvar(dim1)) |
---|
425 | endif |
---|
426 | endif |
---|
427 | ! |
---|
428 | status=nf90_get_var(ncid,varid,tabvar) |
---|
429 | ! |
---|
430 | status = nf90_close(ncid) |
---|
431 | ! |
---|
432 | end subroutine Read_Ncdf_var1d_Int |
---|
433 | |
---|
434 | |
---|
435 | !> Subroutine to retrieve values of a given 2D integer variable |
---|
436 | !! @param[in] varname : name of variable to retrieve |
---|
437 | !! @param[in] file : netcdf file name |
---|
438 | !! @param[out] tabvar : array containing values of the required variable |
---|
439 | !! @return tabvar |
---|
440 | subroutine Read_Ncdf_var2d_Int(varname,file,tabvar) |
---|
441 | ! |
---|
442 | implicit none |
---|
443 | ! |
---|
444 | Character(*),intent(in) :: varname,file |
---|
445 | Integer, dimension(:,:), pointer :: tabvar |
---|
446 | !local variables |
---|
447 | Integer, dimension(10) :: dimIDS |
---|
448 | Integer :: dim1,dim2 |
---|
449 | Integer :: status,ncid |
---|
450 | Integer :: varid |
---|
451 | ! |
---|
452 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
453 | ! |
---|
454 | if (status/=nf90_noerr) then |
---|
455 | write(*,*)"unable to open netcdf file : ",file |
---|
456 | stop |
---|
457 | endif |
---|
458 | ! |
---|
459 | status = nf90_inq_varid(ncid,varname,varid) |
---|
460 | if (status/=nf90_noerr) then |
---|
461 | write(*,*)"pb de nom de variable :", varname, file |
---|
462 | stop |
---|
463 | endif |
---|
464 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
465 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
466 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
467 | ! |
---|
468 | if(.not. associated(tabvar)) then |
---|
469 | Allocate(tabvar(dim1,dim2)) |
---|
470 | else |
---|
471 | if( any(shape(tabvar) /= (/dim1,dim2/)) ) then |
---|
472 | deallocate(tabvar) |
---|
473 | Allocate(tabvar(dim1,dim2)) |
---|
474 | endif |
---|
475 | endif |
---|
476 | ! |
---|
477 | status=nf90_get_var(ncid,varid,tabvar) |
---|
478 | ! |
---|
479 | status = nf90_close(ncid) |
---|
480 | ! |
---|
481 | end subroutine Read_Ncdf_var2d_Int |
---|
482 | |
---|
483 | !> Subroutine to retrieve values of a given 3D integer variable |
---|
484 | !! @param[in] varname : name of variable to retrieve |
---|
485 | !! @param[in] file : netcdf file name |
---|
486 | !! @param[out] tabvar : array containing values of the required variable |
---|
487 | !! @return tabvar |
---|
488 | subroutine Read_Ncdf_var3d_Int(varname,file,tabvar) |
---|
489 | ! |
---|
490 | implicit none |
---|
491 | ! |
---|
492 | Character(*),intent(in) :: varname,file |
---|
493 | Integer, dimension(:,:,:), pointer :: tabvar |
---|
494 | ! |
---|
495 | !local variables |
---|
496 | ! |
---|
497 | Integer, dimension(10) :: dimIDS |
---|
498 | Integer :: dim1,dim2,dim3 |
---|
499 | Integer :: status,ncid |
---|
500 | Integer :: varid |
---|
501 | ! |
---|
502 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
503 | ! |
---|
504 | if (status/=nf90_noerr) then |
---|
505 | write(*,*)"unable to open netcdf file : ",file |
---|
506 | stop |
---|
507 | endif |
---|
508 | ! |
---|
509 | status = nf90_inq_varid(ncid,varname,varid) |
---|
510 | if (status/=nf90_noerr) then |
---|
511 | write(*,*)"pb de nom de variable :", varname, file |
---|
512 | stop |
---|
513 | endif |
---|
514 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
515 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
516 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
517 | status=nf90_inquire_dimension(ncid,dimIDS(3),len=dim3) |
---|
518 | ! |
---|
519 | if(.not. associated(tabvar)) then |
---|
520 | Allocate(tabvar(dim1,dim2,dim3)) |
---|
521 | else |
---|
522 | if( any(shape(tabvar) /= (/dim1,dim2,dim3/)) ) then |
---|
523 | deallocate(tabvar) |
---|
524 | Allocate(tabvar(dim1,dim2,dim3)) |
---|
525 | endif |
---|
526 | endif |
---|
527 | ! |
---|
528 | status=nf90_get_var(ncid,varid,tabvar) |
---|
529 | ! |
---|
530 | status = nf90_close(ncid) |
---|
531 | ! |
---|
532 | end subroutine Read_Ncdf_var3d_Int |
---|
533 | |
---|
534 | !> Subroutine to retrieve values of a given 4D integer variable |
---|
535 | !! @param[in] varname : name of variable to retrieve |
---|
536 | !! @param[in] file : netcdf file name |
---|
537 | !! @param[out] tabvar : array containing values of the required variable |
---|
538 | !! @return tabvar |
---|
539 | |
---|
540 | subroutine Read_Ncdf_var4d_Int(varname,file,tabvar) |
---|
541 | ! |
---|
542 | implicit none |
---|
543 | ! |
---|
544 | Character(*),intent(in) :: varname,file |
---|
545 | Integer, dimension(:,:,:,:), pointer :: tabvar |
---|
546 | ! |
---|
547 | !local variables |
---|
548 | ! |
---|
549 | Integer, dimension(10) :: dimIDS |
---|
550 | Integer :: dim1,dim2,dim3,dim4 |
---|
551 | Integer :: status,ncid |
---|
552 | Integer :: varid |
---|
553 | ! |
---|
554 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
555 | ! |
---|
556 | if (status/=nf90_noerr) then |
---|
557 | write(*,*)"unable to open netcdf file : ",file |
---|
558 | stop |
---|
559 | endif |
---|
560 | ! |
---|
561 | status = nf90_inq_varid(ncid,varname,varid) |
---|
562 | if (status/=nf90_noerr) then |
---|
563 | write(*,*)"pb de nom de variable :", varname, file |
---|
564 | stop |
---|
565 | endif |
---|
566 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
567 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
568 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
569 | status=nf90_inquire_dimension(ncid,dimIDS(3),len=dim3) |
---|
570 | status=nf90_inquire_dimension(ncid,dimIDS(4),len=dim4) |
---|
571 | ! |
---|
572 | if(.not. associated(tabvar)) then |
---|
573 | Allocate(tabvar(dim1,dim2,dim3,dim4)) |
---|
574 | else |
---|
575 | if( any(shape(tabvar) /= (/dim1,dim2,dim3,dim4/)) ) then |
---|
576 | deallocate(tabvar) |
---|
577 | Allocate(tabvar(dim1,dim2,dim3,dim4)) |
---|
578 | endif |
---|
579 | endif |
---|
580 | ! |
---|
581 | status=nf90_get_var(ncid,varid,tabvar) |
---|
582 | ! |
---|
583 | status = nf90_close(ncid) |
---|
584 | ! |
---|
585 | end subroutine Read_Ncdf_var4d_Int |
---|
586 | |
---|
587 | |
---|
588 | !> Subroutine to retrieve values of a given 4D real variable |
---|
589 | !! @param[in] varname : name of variable to retrieve |
---|
590 | !! @param[in] file : netcdf file name |
---|
591 | !! @param[out] tabvar : array containing values of the required variable |
---|
592 | !! @param[in] time : time corresponding to the values to read |
---|
593 | !! @param[in] level : level corresponding to the values to read |
---|
594 | !! @return tabvar |
---|
595 | |
---|
596 | subroutine Read_Ncdf_var4d_Real_nt(varname,file,tabvar,time,level) |
---|
597 | implicit none |
---|
598 | ! |
---|
599 | Character(*),intent(in) :: varname,file |
---|
600 | Integer,intent(in) :: time,level |
---|
601 | Real*8, dimension(:,:,:,:), pointer :: tabvar |
---|
602 | ! |
---|
603 | !local variables |
---|
604 | ! |
---|
605 | Integer, dimension(4) :: dimIDS |
---|
606 | Integer :: dim1,dim2 |
---|
607 | Integer :: status,ncid |
---|
608 | Integer :: varid |
---|
609 | ! |
---|
610 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
611 | ! |
---|
612 | if (status/=nf90_noerr) then |
---|
613 | write(*,*)"unable to open netcdf file : ",file |
---|
614 | stop |
---|
615 | endif |
---|
616 | ! |
---|
617 | status = nf90_inq_varid(ncid,varname,varid) |
---|
618 | if (status/=nf90_noerr) then |
---|
619 | write(*,*)"pb de nom de variable :", varname, file |
---|
620 | stop |
---|
621 | endif |
---|
622 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
623 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
624 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
625 | ! |
---|
626 | if(.not. associated(tabvar)) then |
---|
627 | Allocate(tabvar(dim1,dim2,1,1)) |
---|
628 | else |
---|
629 | if( any(shape(tabvar) /= (/dim1,dim2,1,1/)) ) then |
---|
630 | deallocate(tabvar) |
---|
631 | Allocate(tabvar(dim1,dim2,1,1)) |
---|
632 | endif |
---|
633 | endif |
---|
634 | ! |
---|
635 | status=nf90_get_var(ncid,varid,tabvar,start=(/1,1,level,time/)) |
---|
636 | ! |
---|
637 | if (status/=nf90_noerr) then |
---|
638 | write(*,*)"unable to retrieve netcdf variable : ",trim(varname) |
---|
639 | stop |
---|
640 | endif |
---|
641 | ! |
---|
642 | status = nf90_close(ncid) |
---|
643 | ! |
---|
644 | end subroutine Read_Ncdf_var4d_Real_nt |
---|
645 | |
---|
646 | |
---|
647 | !> Subroutine to retrieve values of a given 3D real variable withe time dependent |
---|
648 | !! @param[in] varname : name of variable to retrieve |
---|
649 | !! @param[in] file : netcdf file name |
---|
650 | !! @param[out] tabvar : array containing values of the required variable |
---|
651 | !! @param[in] time : time corresponding to the values to read |
---|
652 | !! @return tabvar |
---|
653 | subroutine Read_Ncdf_var4d_Real_t(varname,file,tabvar,time) |
---|
654 | implicit none |
---|
655 | ! |
---|
656 | Character(*),intent(in) :: varname,file |
---|
657 | Integer,intent(in) :: time |
---|
658 | Real*8, dimension(:,:,:,:), pointer :: tabvar |
---|
659 | ! |
---|
660 | !local variables |
---|
661 | ! |
---|
662 | Integer, dimension(4) :: dimIDS |
---|
663 | Integer :: dim1,dim2,dim3 |
---|
664 | Integer :: status,ncid |
---|
665 | Integer :: varid |
---|
666 | ! |
---|
667 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
668 | ! |
---|
669 | if (status/=nf90_noerr) then |
---|
670 | write(*,*)"unable to open netcdf file : ",file |
---|
671 | stop |
---|
672 | endif |
---|
673 | ! |
---|
674 | status = nf90_inq_varid(ncid,varname,varid) |
---|
675 | if (status/=nf90_noerr) then |
---|
676 | write(*,*)"pb de nom de variable :", varname, file |
---|
677 | stop |
---|
678 | endif |
---|
679 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
680 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
681 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
682 | status=nf90_inquire_dimension(ncid,dimIDS(3),len=dim3) |
---|
683 | ! |
---|
684 | if(.not. associated(tabvar)) Allocate(tabvar(dim1,dim2,dim3,1)) |
---|
685 | status=nf90_get_var(ncid,varid,tabvar,start=(/1,1,1,time/)) |
---|
686 | |
---|
687 | if (status/=nf90_noerr) then |
---|
688 | write(*,*)"unable to retrieve netcdf variable : ",trim(varname) |
---|
689 | stop |
---|
690 | endif |
---|
691 | ! |
---|
692 | status = nf90_close(ncid) |
---|
693 | ! |
---|
694 | end subroutine Read_Ncdf_var4d_Real_t |
---|
695 | |
---|
696 | !> Subroutine to read a 3D real variable in a given file for time t |
---|
697 | !! @param[in] varname : name of variable to retrieve |
---|
698 | !! @param[in] file : netcdf file name |
---|
699 | !! @param[out] tabvar : array containing values of the required variable |
---|
700 | !! @param[in] time : time corresponding to the values to read |
---|
701 | !! @return tabvar |
---|
702 | |
---|
703 | |
---|
704 | subroutine Read_Ncdf_var3d_Real_t(varname,file,tabvar,time) |
---|
705 | implicit none |
---|
706 | ! |
---|
707 | Character(*),intent(in) :: varname,file |
---|
708 | Integer,intent(in) :: time |
---|
709 | Real*8, dimension(:,:,:), pointer :: tabvar |
---|
710 | ! |
---|
711 | !local variables |
---|
712 | ! |
---|
713 | Integer, dimension(3) :: dimIDS |
---|
714 | Integer :: dim1,dim2 |
---|
715 | Integer :: status,ncid |
---|
716 | Integer :: varid |
---|
717 | ! |
---|
718 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
719 | ! |
---|
720 | if (status/=nf90_noerr) then |
---|
721 | write(*,*)"unable to open netcdf file : ",file |
---|
722 | stop |
---|
723 | endif |
---|
724 | ! |
---|
725 | status = nf90_inq_varid(ncid,varname,varid) |
---|
726 | ! |
---|
727 | status=nf90_inquire_variable(ncid,varid,dimids=dimIDS) |
---|
728 | status=nf90_inquire_dimension(ncid,dimIDS(1),len=dim1) |
---|
729 | status=nf90_inquire_dimension(ncid,dimIDS(2),len=dim2) |
---|
730 | ! |
---|
731 | if(.not. associated(tabvar)) then |
---|
732 | Allocate(tabvar(dim1,dim2,1)) |
---|
733 | else |
---|
734 | if( any(shape(tabvar) /= (/dim1,dim2,1/)) ) then |
---|
735 | deallocate(tabvar) |
---|
736 | Allocate(tabvar(dim1,dim2,1)) |
---|
737 | endif |
---|
738 | endif |
---|
739 | |
---|
740 | status=nf90_get_var(ncid,varid,tabvar,start=(/1,1,time/)) |
---|
741 | |
---|
742 | if (status/=nf90_noerr) then |
---|
743 | write(*,*)"unable to retrieve netcdf variable : ",trim(varname) |
---|
744 | stop |
---|
745 | endif |
---|
746 | ! |
---|
747 | status = nf90_close(ncid) |
---|
748 | ! |
---|
749 | end subroutine Read_Ncdf_var3d_Real_t |
---|
750 | |
---|
751 | !> Subroutine to write a 1D real variable in a given file |
---|
752 | !!@param[in] varname : name of variable to store |
---|
753 | !!@param[in] dimname : name of dimensions of the given variable |
---|
754 | !!@param[in] file : netcdf file name |
---|
755 | !!@param[in] tabvar : values of the variable to write |
---|
756 | !!@param[in] typevar : type of the variable to write |
---|
757 | |
---|
758 | subroutine Write_Ncdf_var1d_Real(varname,dimname,file,tabvar,typevar) |
---|
759 | ! |
---|
760 | implicit none |
---|
761 | ! |
---|
762 | Character(*),intent(in) :: varname,file,dimname,typevar |
---|
763 | Real*8, dimension(:), pointer :: tabvar |
---|
764 | ! |
---|
765 | ! local variables |
---|
766 | ! |
---|
767 | Integer :: dimid |
---|
768 | Integer :: status,ncid |
---|
769 | Integer :: varid |
---|
770 | ! |
---|
771 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
772 | if (status/=nf90_noerr) then |
---|
773 | write(*,*)"unable to open netcdf file : ",file |
---|
774 | stop |
---|
775 | endif |
---|
776 | |
---|
777 | |
---|
778 | status = nf90_inq_dimid(ncid,dimname, dimid) |
---|
779 | |
---|
780 | status = nf90_inq_varid(ncid,varname,varid) |
---|
781 | ! write(6,*) 'varid status',status,'no error :',nf90_NoErr,' varid :',varid,' varname: ',varname |
---|
782 | ! write(6,*) 'erreur',trim(nf90_strerror(status)) |
---|
783 | |
---|
784 | status = nf90_redef(ncid) |
---|
785 | |
---|
786 | select case(TRIM(typevar)) |
---|
787 | case('double') |
---|
788 | status = nf90_def_var(ncid,varname,nf90_double,(/dimid/),varid) |
---|
789 | |
---|
790 | case('float') |
---|
791 | status = nf90_def_var(ncid,varname,nf90_float,(/dimid/),varid) |
---|
792 | end select |
---|
793 | |
---|
794 | status = nf90_enddef(ncid) |
---|
795 | status = nf90_put_var(ncid,varid,tabvar) |
---|
796 | |
---|
797 | ! |
---|
798 | status = nf90_close(ncid) |
---|
799 | ! |
---|
800 | end subroutine Write_Ncdf_var1d_Real |
---|
801 | |
---|
802 | !> Subroutine to write a 2D real variable in a given file |
---|
803 | !!@param[in] varname : name of variable to store |
---|
804 | !!@param[in] dimname : name of dimensions of the given variable |
---|
805 | !!@param[in] file : netcdf file name |
---|
806 | !!@param[in] tabvar : values of the variable to write |
---|
807 | !!@param[in] typevar : type of the variable to write |
---|
808 | |
---|
809 | |
---|
810 | subroutine Write_Ncdf_var2d_Real(varname,dimname,file,tabvar,typevar) |
---|
811 | ! |
---|
812 | ! implicit none |
---|
813 | ! |
---|
814 | Character(*),intent(in) :: varname,file,typevar |
---|
815 | Character(*), dimension(2) :: dimname |
---|
816 | Real*8, dimension(:,:), pointer :: tabvar |
---|
817 | ! |
---|
818 | ! local variables |
---|
819 | ! |
---|
820 | Integer :: dimid1,dimid2 |
---|
821 | Integer :: status,ncid |
---|
822 | Integer :: varid |
---|
823 | ! |
---|
824 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
825 | if (status/=nf90_noerr) then |
---|
826 | write(*,*)"unable to open netcdf file : ",file |
---|
827 | stop |
---|
828 | endif |
---|
829 | ! |
---|
830 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
831 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
832 | |
---|
833 | status = nf90_inq_varid(ncid,varname,varid) |
---|
834 | status = nf90_redef(ncid) |
---|
835 | |
---|
836 | select case(TRIM(typevar)) |
---|
837 | case('double') |
---|
838 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
839 | (/dimid1,dimid2/),varid) |
---|
840 | case('float') |
---|
841 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
842 | (/dimid1,dimid2/),varid) |
---|
843 | end select |
---|
844 | ! |
---|
845 | status = nf90_enddef(ncid) |
---|
846 | status = nf90_put_var(ncid,varid,tabvar) |
---|
847 | ! |
---|
848 | status = nf90_close(ncid) |
---|
849 | ! |
---|
850 | end subroutine Write_Ncdf_var2d_Real |
---|
851 | |
---|
852 | !> Subroutine to write a 3D real variable in a given file |
---|
853 | !!@param[in] varname : name of variable to store |
---|
854 | !!@param[in] dimname : name of dimensions of the given variable |
---|
855 | !!@param[in] file : netcdf file name |
---|
856 | !!@param[in] tabvar : values of the variable to write |
---|
857 | !!@param[in] typevar : type of the variable to write |
---|
858 | |
---|
859 | subroutine Write_Ncdf_var3d_Real(varname,dimname,file,tabvar,typevar) |
---|
860 | ! |
---|
861 | implicit none |
---|
862 | ! |
---|
863 | Character(*),intent(in) :: varname,file,typevar |
---|
864 | Character(*),dimension(3),intent(in) :: dimname |
---|
865 | Real*8, dimension(:,:,:), pointer :: tabvar |
---|
866 | ! |
---|
867 | ! local variables |
---|
868 | ! |
---|
869 | Integer :: dimid1,dimid2,dimid3 |
---|
870 | Integer :: status,ncid |
---|
871 | Integer :: varid |
---|
872 | ! |
---|
873 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
874 | if (status/=nf90_noerr) then |
---|
875 | write(*,*)"unable to open netcdf file : ",file |
---|
876 | stop |
---|
877 | endif |
---|
878 | ! |
---|
879 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
880 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
881 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
882 | status = nf90_inq_varid(ncid,varname,varid) |
---|
883 | status = nf90_redef(ncid) |
---|
884 | |
---|
885 | select case(TRIM(typevar)) |
---|
886 | case('double') |
---|
887 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
888 | (/dimid1,dimid2,dimid3/),varid) |
---|
889 | case('float') |
---|
890 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
891 | (/dimid1,dimid2,dimid3/),varid) |
---|
892 | end select |
---|
893 | ! |
---|
894 | status = nf90_enddef(ncid) |
---|
895 | status = nf90_put_var(ncid,varid,tabvar) |
---|
896 | ! |
---|
897 | status = nf90_close(ncid) |
---|
898 | ! |
---|
899 | end subroutine Write_Ncdf_var3d_Real |
---|
900 | |
---|
901 | !> Subroutine to write a 4D real variable in a given file |
---|
902 | !!@param[in] varname : name of variable to store |
---|
903 | !!@param[in] dimname : name of dimensions of the given variable |
---|
904 | !!@param[in] file : netcdf file name |
---|
905 | !!@param[in] tabvar : values of the variable to write |
---|
906 | !!@param[in] typevar : type of the variable to write |
---|
907 | |
---|
908 | subroutine Write_Ncdf_var4d_Real(varname,dimname,file,tabvar,typevar) |
---|
909 | ! |
---|
910 | implicit none |
---|
911 | ! |
---|
912 | Character(*),intent(in) :: varname,file,typevar |
---|
913 | Character(*),dimension(4),intent(in) :: dimname |
---|
914 | Real*8, dimension(:,:,:,:), pointer :: tabvar |
---|
915 | ! |
---|
916 | ! local variables |
---|
917 | ! |
---|
918 | Integer :: dimid1,dimid2,dimid3,dimid4 |
---|
919 | Integer :: status,ncid |
---|
920 | Integer :: varid |
---|
921 | ! |
---|
922 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
923 | if (status/=nf90_noerr) then |
---|
924 | write(*,*)"unable to open netcdf file : ",file |
---|
925 | stop |
---|
926 | endif |
---|
927 | ! |
---|
928 | status = nf90_inq_varid(ncid,varname,varid) |
---|
929 | ! |
---|
930 | if(status/=nf90_noerr) then |
---|
931 | ! |
---|
932 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
933 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
934 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
935 | status = nf90_inq_dimid(ncid,dimname(4), dimid4) |
---|
936 | status = nf90_redef(ncid) |
---|
937 | ! |
---|
938 | select case(TRIM(typevar)) |
---|
939 | case('double') |
---|
940 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
941 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
942 | case('float') |
---|
943 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
944 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
945 | end select |
---|
946 | ! |
---|
947 | status = nf90_enddef(ncid) |
---|
948 | endif |
---|
949 | ! |
---|
950 | status = nf90_put_var(ncid,varid,tabvar) |
---|
951 | ! |
---|
952 | status = nf90_close(ncid) |
---|
953 | ! |
---|
954 | end subroutine Write_Ncdf_var4d_Real |
---|
955 | |
---|
956 | |
---|
957 | !> Subroutine to write a 1D integer variable in a given file |
---|
958 | !!@param[in] varname : name of variable to store |
---|
959 | !!@param[in] dimname : name of dimensions of the given variable |
---|
960 | !!@param[in] file : netcdf file name |
---|
961 | !!@param[in] tabvar : values of the variable to write |
---|
962 | |
---|
963 | |
---|
964 | subroutine Write_Ncdf_var1d_Int(varname,dimname,file,tabvar) |
---|
965 | ! |
---|
966 | implicit none |
---|
967 | ! |
---|
968 | Character(*),intent(in) :: varname,file,dimname |
---|
969 | Integer, dimension(:), pointer :: tabvar |
---|
970 | ! |
---|
971 | ! local variables |
---|
972 | ! |
---|
973 | Integer :: dimid |
---|
974 | Integer :: status,ncid |
---|
975 | Integer :: varid |
---|
976 | ! |
---|
977 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
978 | if (status/=nf90_noerr) then |
---|
979 | write(*,*)"unable to open netcdf file : ",file |
---|
980 | stop |
---|
981 | endif |
---|
982 | ! |
---|
983 | status = nf90_inq_dimid(ncid,dimname, dimid) |
---|
984 | status = nf90_inq_varid(ncid,varname,varid) |
---|
985 | status = nf90_redef(ncid) |
---|
986 | status = nf90_def_var(ncid,varname,nf90_int,(/dimid/),varid) |
---|
987 | status = nf90_enddef(ncid) |
---|
988 | status = nf90_put_var(ncid,varid,tabvar) |
---|
989 | ! |
---|
990 | status = nf90_close(ncid) |
---|
991 | ! |
---|
992 | end subroutine Write_Ncdf_var1d_Int |
---|
993 | |
---|
994 | !> Subroutine to write a 2D integer variable in a given file |
---|
995 | !!@param[in] varname : name of variable to store |
---|
996 | !!@param[in] dimname : name of dimensions of the given variable |
---|
997 | !!@param[in] file : netcdf file name |
---|
998 | !!@param[in] tabvar : values of the variable to write |
---|
999 | |
---|
1000 | subroutine Write_Ncdf_var2d_Int(varname,dimname,file,tabvar) |
---|
1001 | ! |
---|
1002 | implicit none |
---|
1003 | ! |
---|
1004 | Character(*),intent(in) :: varname,file |
---|
1005 | Character(*), dimension(2) :: dimname |
---|
1006 | Integer, dimension(:,:), pointer :: tabvar |
---|
1007 | ! |
---|
1008 | ! local variables |
---|
1009 | ! |
---|
1010 | Integer :: dimid1,dimid2 |
---|
1011 | Integer :: status,ncid |
---|
1012 | Integer :: varid |
---|
1013 | ! |
---|
1014 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1015 | if (status/=nf90_noerr) then |
---|
1016 | write(*,*)"unable to open netcdf file : ",file |
---|
1017 | stop |
---|
1018 | endif |
---|
1019 | ! |
---|
1020 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1021 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1022 | status = nf90_inq_varid(ncid,varname,varid) |
---|
1023 | status = nf90_redef(ncid) |
---|
1024 | status = nf90_def_var(ncid,varname,nf90_int, & |
---|
1025 | (/dimid1,dimid2/),varid) |
---|
1026 | status = nf90_enddef(ncid) |
---|
1027 | status = nf90_put_var(ncid,varid,tabvar) |
---|
1028 | ! |
---|
1029 | status = nf90_close(ncid) |
---|
1030 | ! |
---|
1031 | end subroutine Write_Ncdf_var2d_Int |
---|
1032 | |
---|
1033 | !> Subroutine to write a 3D integer variable in a given file |
---|
1034 | !!@param[in] varname : name of variable to store |
---|
1035 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1036 | !!@param[in] file : netcdf file name |
---|
1037 | !!@param[in] tabvar : values of the variable to write |
---|
1038 | |
---|
1039 | subroutine Write_Ncdf_var3d_Int(varname,dimname,file,tabvar) |
---|
1040 | ! |
---|
1041 | implicit none |
---|
1042 | ! |
---|
1043 | Character(*),intent(in) :: varname,file |
---|
1044 | Character(*),dimension(3),intent(in) :: dimname |
---|
1045 | Integer, dimension(:,:,:), pointer :: tabvar |
---|
1046 | ! |
---|
1047 | ! local variables |
---|
1048 | ! |
---|
1049 | Integer :: dimid1,dimid2,dimid3 |
---|
1050 | Integer :: status,ncid |
---|
1051 | Integer :: varid |
---|
1052 | ! |
---|
1053 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1054 | if (status/=nf90_noerr) then |
---|
1055 | write(*,*)"unable to open netcdf file : ",file |
---|
1056 | stop |
---|
1057 | endif |
---|
1058 | ! |
---|
1059 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1060 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1061 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1062 | status = nf90_inq_varid(ncid,varname,varid) |
---|
1063 | status = nf90_redef(ncid) |
---|
1064 | status = nf90_def_var(ncid,varname,nf90_int, & |
---|
1065 | (/dimid1,dimid2,dimid3/),varid) |
---|
1066 | status = nf90_enddef(ncid) |
---|
1067 | status = nf90_put_var(ncid,varid,tabvar) |
---|
1068 | ! |
---|
1069 | status = nf90_close(ncid) |
---|
1070 | ! |
---|
1071 | end subroutine Write_Ncdf_var3d_Int |
---|
1072 | |
---|
1073 | !> Subroutine to write a 4D integer variable in a given file |
---|
1074 | !!@param[in] varname : name of variable to store |
---|
1075 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1076 | !!@param[in] file : netcdf file name |
---|
1077 | !!@param[in] tabvar : values of the variable to write |
---|
1078 | |
---|
1079 | |
---|
1080 | subroutine Write_Ncdf_var4d_Int(varname,dimname,file,tabvar) |
---|
1081 | ! |
---|
1082 | implicit none |
---|
1083 | ! |
---|
1084 | Character(*),intent(in) :: varname,file |
---|
1085 | Character(*),dimension(4),intent(in) :: dimname |
---|
1086 | Integer, dimension(:,:,:,:), pointer :: tabvar |
---|
1087 | ! |
---|
1088 | ! local variables |
---|
1089 | ! |
---|
1090 | Integer :: dimid1,dimid2,dimid3,dimid4 |
---|
1091 | Integer :: status,ncid |
---|
1092 | Integer :: varid |
---|
1093 | ! |
---|
1094 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1095 | if (status/=nf90_noerr) then |
---|
1096 | write(*,*)"unable to open netcdf file : ",file |
---|
1097 | stop |
---|
1098 | endif |
---|
1099 | ! |
---|
1100 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1101 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1102 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1103 | status = nf90_inq_dimid(ncid,dimname(4), dimid4) |
---|
1104 | status = nf90_inq_varid(ncid,varname,varid) |
---|
1105 | status = nf90_redef(ncid) |
---|
1106 | status = nf90_def_var(ncid,varname,nf90_int, & |
---|
1107 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1108 | status = nf90_enddef(ncid) |
---|
1109 | status = nf90_put_var(ncid,varid,tabvar) |
---|
1110 | ! |
---|
1111 | status = nf90_close(ncid) |
---|
1112 | ! |
---|
1113 | end subroutine Write_Ncdf_var4d_Int |
---|
1114 | |
---|
1115 | !> Subroutine to write a 2D real variable with time dependent in a given file |
---|
1116 | !!@param[in] varname : name of variable to store |
---|
1117 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1118 | !!@param[in] file : netcdf file name |
---|
1119 | !!@param[in] tabvar : values of the variable to write |
---|
1120 | !!@param[in] time : time corresponding to the values to write |
---|
1121 | !!@param[inout] idef : |
---|
1122 | !!@param[in] typevar : type of the variable to write |
---|
1123 | !!@return idef |
---|
1124 | subroutine Write_Ncdf_var2d_Real_t(varname,dimname,file,tabvar,time,idef,typevar) |
---|
1125 | ! |
---|
1126 | implicit none |
---|
1127 | ! |
---|
1128 | Character(*),intent(in) :: varname,file,typevar |
---|
1129 | Character(*),dimension(3),intent(in) :: dimname |
---|
1130 | Integer :: time |
---|
1131 | Integer,intent(inout) :: idef |
---|
1132 | Real*8, dimension(:,:), pointer :: tabvar |
---|
1133 | ! |
---|
1134 | ! local variables |
---|
1135 | ! |
---|
1136 | Integer :: dimid1,dimid2,dimid3 |
---|
1137 | Integer :: status,ncid |
---|
1138 | Integer :: varid |
---|
1139 | ! |
---|
1140 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1141 | if (status/=nf90_noerr) then |
---|
1142 | write(*,*)"unable to open netcdf file : ",file |
---|
1143 | stop |
---|
1144 | endif |
---|
1145 | ! |
---|
1146 | if(time==1 .OR. idef==0) then |
---|
1147 | ! |
---|
1148 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1149 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1150 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1151 | status = nf90_redef(ncid) |
---|
1152 | ! |
---|
1153 | select case(TRIM(typevar)) |
---|
1154 | case('double') |
---|
1155 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
1156 | (/dimid1,dimid2,dimid3/),varid) |
---|
1157 | case('float') |
---|
1158 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
1159 | (/dimid1,dimid2,dimid3/),varid) |
---|
1160 | end select |
---|
1161 | ! |
---|
1162 | status = nf90_enddef(ncid) |
---|
1163 | idef = 1 |
---|
1164 | else |
---|
1165 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1166 | endif |
---|
1167 | ! |
---|
1168 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,1,time/)) |
---|
1169 | if (status/=nf90_noerr) then |
---|
1170 | write(*,*)"unable to store variable ",varname, & |
---|
1171 | " in file ",file |
---|
1172 | stop |
---|
1173 | endif |
---|
1174 | ! |
---|
1175 | status = nf90_close(ncid) |
---|
1176 | ! |
---|
1177 | end subroutine Write_Ncdf_var2d_Real_t |
---|
1178 | |
---|
1179 | !> Subroutine to write a 2D integer variable with time dependent in a given file |
---|
1180 | !!@param[in] varname : name of variable to store |
---|
1181 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1182 | !!@param[in] file : netcdf file name |
---|
1183 | !!@param[in] tabvar : values of the variable to write |
---|
1184 | !!@param[in] time : time corresponding to the values to write |
---|
1185 | !!@param[inout] idef : |
---|
1186 | !!@return idef |
---|
1187 | |
---|
1188 | subroutine Write_Ncdf_var2d_Int_t(varname,dimname,file,tabvar,time,idef) |
---|
1189 | ! |
---|
1190 | implicit none |
---|
1191 | ! |
---|
1192 | Character(*),intent(in) :: varname,file |
---|
1193 | Character(*),dimension(3),intent(in) :: dimname |
---|
1194 | Integer :: time |
---|
1195 | Integer,intent(inout) :: idef |
---|
1196 | Integer, dimension(:,:), pointer :: tabvar |
---|
1197 | ! |
---|
1198 | ! local variables |
---|
1199 | ! |
---|
1200 | Integer :: dimid1,dimid2,dimid3 |
---|
1201 | Integer :: status,ncid |
---|
1202 | Integer :: varid |
---|
1203 | ! |
---|
1204 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1205 | if (status/=nf90_noerr) then |
---|
1206 | write(*,*)"unable to open netcdf file : ",file |
---|
1207 | stop |
---|
1208 | endif |
---|
1209 | ! |
---|
1210 | if(time==1 .OR. idef==0 ) then |
---|
1211 | ! |
---|
1212 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1213 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1214 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1215 | status = nf90_redef(ncid) |
---|
1216 | ! |
---|
1217 | status = nf90_def_var(ncid,varname,nf90_int, & |
---|
1218 | (/dimid1,dimid2,dimid3/),varid) |
---|
1219 | ! |
---|
1220 | status = nf90_enddef(ncid) |
---|
1221 | idef=1 |
---|
1222 | else |
---|
1223 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1224 | endif |
---|
1225 | ! |
---|
1226 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,1,time/)) |
---|
1227 | if (status/=nf90_noerr) then |
---|
1228 | write(*,*)"unable to store variable ",varname, & |
---|
1229 | " in file ",file |
---|
1230 | stop |
---|
1231 | endif |
---|
1232 | ! |
---|
1233 | status = nf90_close(ncid) |
---|
1234 | ! |
---|
1235 | end subroutine Write_Ncdf_var2d_Int_t |
---|
1236 | |
---|
1237 | |
---|
1238 | !> Subroutine to write a 1D real variable with time dependent in a given file |
---|
1239 | !!@param[in] varname : name of variable to store |
---|
1240 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1241 | !!@param[in] file : netcdf file name |
---|
1242 | !!@param[in] tabvar : values of the variable to write |
---|
1243 | !!@param[in] time : time corresponding to the values to write |
---|
1244 | !!@param[inou] idef : |
---|
1245 | !!@param[in] typevar : type of the variable to write |
---|
1246 | !!@return idef |
---|
1247 | |
---|
1248 | |
---|
1249 | subroutine Write_Ncdf_var1d_Real_t(varname,dimname,file,tabvar,time,idef,typevar) |
---|
1250 | ! |
---|
1251 | implicit none |
---|
1252 | ! |
---|
1253 | Character(*),intent(in) :: varname,file,typevar |
---|
1254 | Character(*),dimension(2),intent(in) :: dimname |
---|
1255 | Integer :: time |
---|
1256 | Integer,intent(inout) :: idef |
---|
1257 | Real*8, dimension(:) :: tabvar |
---|
1258 | ! |
---|
1259 | ! local variables |
---|
1260 | ! |
---|
1261 | Integer :: dimid1,dimid2 |
---|
1262 | Integer :: status,ncid |
---|
1263 | Integer :: varid |
---|
1264 | ! |
---|
1265 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1266 | if (status/=nf90_noerr) then |
---|
1267 | write(*,*)"unable to open netcdf file : ",file |
---|
1268 | stop |
---|
1269 | endif |
---|
1270 | ! |
---|
1271 | if(time==1 .OR. idef==0) then |
---|
1272 | ! |
---|
1273 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1274 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1275 | status = nf90_redef(ncid) |
---|
1276 | ! |
---|
1277 | select case(TRIM(typevar)) |
---|
1278 | case('double') |
---|
1279 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
1280 | (/dimid1,dimid2/),varid) |
---|
1281 | case('float') |
---|
1282 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
1283 | (/dimid1,dimid2/),varid) |
---|
1284 | end select |
---|
1285 | ! |
---|
1286 | status = nf90_enddef(ncid) |
---|
1287 | idef=1 |
---|
1288 | else |
---|
1289 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1290 | endif |
---|
1291 | ! |
---|
1292 | |
---|
1293 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,time/)) |
---|
1294 | if (status/=nf90_noerr) then |
---|
1295 | write(*,*)"unable to store variable ",varname, & |
---|
1296 | " in file ",file |
---|
1297 | stop |
---|
1298 | endif |
---|
1299 | ! |
---|
1300 | status = nf90_close(ncid) |
---|
1301 | ! |
---|
1302 | end subroutine Write_Ncdf_var1d_Real_t |
---|
1303 | |
---|
1304 | !> Subroutine to write a 1D real variable with time dependent in a given file |
---|
1305 | !!@param[in] varname : name of variable to store |
---|
1306 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1307 | !!@param[in] file : netcdf file name |
---|
1308 | !!@param[in] tabvar : values of the variable to write |
---|
1309 | !!@param[in] time : time corresponding to the values to write |
---|
1310 | !!@param[inout] idef : |
---|
1311 | !!@param[in] typevar : type of the variable to write |
---|
1312 | !!@return idef |
---|
1313 | |
---|
1314 | subroutine Write_Ncdf_var1d_Realbis_t(varname,dimname,file,tabvar,time,idef,typevar) |
---|
1315 | ! |
---|
1316 | implicit none |
---|
1317 | ! |
---|
1318 | Character(*),intent(in) :: varname,file,typevar |
---|
1319 | Character(*),intent(in) :: dimname |
---|
1320 | Integer :: time |
---|
1321 | Integer,intent(inout) :: idef |
---|
1322 | Real*8, dimension(:), pointer :: tabvar |
---|
1323 | ! |
---|
1324 | ! local variables |
---|
1325 | ! |
---|
1326 | Integer :: dimid1 |
---|
1327 | Integer :: status,ncid |
---|
1328 | Integer :: varid |
---|
1329 | ! |
---|
1330 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1331 | if (status/=nf90_noerr) then |
---|
1332 | write(*,*)"unable to open netcdf file : ",file |
---|
1333 | stop |
---|
1334 | endif |
---|
1335 | ! |
---|
1336 | if(time==1 .OR. idef==0) then |
---|
1337 | ! |
---|
1338 | status = nf90_inq_dimid(ncid,dimname, dimid1) |
---|
1339 | status = nf90_redef(ncid) |
---|
1340 | ! |
---|
1341 | select case(TRIM(typevar)) |
---|
1342 | case('double') |
---|
1343 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
1344 | (/dimid1/),varid) |
---|
1345 | case('float') |
---|
1346 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
1347 | (/dimid1/),varid) |
---|
1348 | end select |
---|
1349 | ! |
---|
1350 | status = nf90_enddef(ncid) |
---|
1351 | idef=1 |
---|
1352 | else |
---|
1353 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1354 | endif |
---|
1355 | ! |
---|
1356 | |
---|
1357 | status = nf90_put_var(ncid,varid,tabvar,start=(/time/)) |
---|
1358 | if (status/=nf90_noerr) then |
---|
1359 | write(*,*)"unable to store variable ",varname, & |
---|
1360 | " in file ",file |
---|
1361 | stop |
---|
1362 | endif |
---|
1363 | ! |
---|
1364 | status = nf90_close(ncid) |
---|
1365 | ! |
---|
1366 | end subroutine Write_Ncdf_var1d_Realbis_t |
---|
1367 | |
---|
1368 | |
---|
1369 | !> Subroutine to write a 1D integer variable with time dependent in a given file |
---|
1370 | !!@param[in] varname : name of variable to store |
---|
1371 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1372 | !!@param[in] file : netcdf file name |
---|
1373 | !!@param[in] tabvar : values of the variable to write |
---|
1374 | !!@param[in] time : time corresponding to the values to write |
---|
1375 | !!@param[inout] idef : |
---|
1376 | !!@param[in] typevar : type of the variable to write |
---|
1377 | !!@return idef |
---|
1378 | subroutine Write_Ncdf_var1d_Int_t(varname,dimname,file,tabvar,time,idef,typevar) |
---|
1379 | ! |
---|
1380 | implicit none |
---|
1381 | ! |
---|
1382 | Character(*),intent(in) :: varname,file,typevar |
---|
1383 | Character(*),dimension(2),intent(in) :: dimname |
---|
1384 | Integer :: time |
---|
1385 | Integer,intent(inout) :: idef |
---|
1386 | !Integer, dimension(:), pointer :: tabvar |
---|
1387 | Integer, dimension(:) :: tabvar |
---|
1388 | ! |
---|
1389 | ! local variables |
---|
1390 | ! |
---|
1391 | Integer :: dimid1,dimid2 |
---|
1392 | Integer :: status,ncid |
---|
1393 | Integer :: varid |
---|
1394 | ! |
---|
1395 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1396 | if (status/=nf90_noerr) then |
---|
1397 | write(*,*)"unable to open netcdf file : ",file |
---|
1398 | stop |
---|
1399 | endif |
---|
1400 | ! |
---|
1401 | if(time==1 .OR. idef==0) then |
---|
1402 | ! |
---|
1403 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1404 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1405 | status = nf90_redef(ncid) |
---|
1406 | |
---|
1407 | ! |
---|
1408 | |
---|
1409 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
1410 | (/dimid1,dimid2/),varid) |
---|
1411 | |
---|
1412 | ! |
---|
1413 | status = nf90_enddef(ncid) |
---|
1414 | idef=1 |
---|
1415 | else |
---|
1416 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1417 | endif |
---|
1418 | ! |
---|
1419 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,time/)) |
---|
1420 | if (status/=nf90_noerr) then |
---|
1421 | write(*,*)"unable to store variable ",varname, & |
---|
1422 | " in file ",file |
---|
1423 | stop |
---|
1424 | endif |
---|
1425 | ! |
---|
1426 | status = nf90_close(ncid) |
---|
1427 | ! |
---|
1428 | end subroutine Write_Ncdf_var1d_Int_t |
---|
1429 | |
---|
1430 | !> Subroutine to write a 2D real variable with time dependent in a given file |
---|
1431 | !!@param[in] varname : name of variable to store |
---|
1432 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1433 | !!@param[in] file : netcdf file name |
---|
1434 | !!@param[in] tabvar : values of the variable to write |
---|
1435 | !!@param[in] time : time corresponding to the values to write |
---|
1436 | !!@param[inout] idef : |
---|
1437 | !!@param[in] typevar : type of the variable to write |
---|
1438 | !!@return idef |
---|
1439 | subroutine Write_Ncdf_var3d_Real_t(varname,dimname,file,tabvar,time,idef,typevar) |
---|
1440 | ! |
---|
1441 | implicit none |
---|
1442 | ! |
---|
1443 | Character(*),intent(in) :: varname,file,typevar |
---|
1444 | Character(*),dimension(4),intent(in) :: dimname |
---|
1445 | Integer :: time |
---|
1446 | Integer,intent(inout) :: idef |
---|
1447 | Real*8, dimension(:,:,:), pointer :: tabvar |
---|
1448 | ! |
---|
1449 | ! local variables |
---|
1450 | ! |
---|
1451 | Integer :: dimid1,dimid2,dimid3,dimid4 |
---|
1452 | Integer :: status,ncid |
---|
1453 | Integer :: varid |
---|
1454 | ! |
---|
1455 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1456 | if (status/=nf90_noerr) then |
---|
1457 | write(*,*)"unable to open netcdf file : ",file |
---|
1458 | stop |
---|
1459 | endif |
---|
1460 | ! |
---|
1461 | if(time==1 .OR. idef==0) then |
---|
1462 | ! |
---|
1463 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1464 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1465 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1466 | status = nf90_inq_dimid(ncid,dimname(4), dimid4) |
---|
1467 | status = nf90_redef(ncid) |
---|
1468 | ! |
---|
1469 | select case(TRIM(typevar)) |
---|
1470 | case('double') |
---|
1471 | status = nf90_def_var(ncid,varname,nf90_double, & |
---|
1472 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1473 | case('float') |
---|
1474 | status = nf90_def_var(ncid,varname,nf90_float, & |
---|
1475 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1476 | end select |
---|
1477 | ! |
---|
1478 | status = nf90_enddef(ncid) |
---|
1479 | idef=1 |
---|
1480 | else |
---|
1481 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1482 | endif |
---|
1483 | ! |
---|
1484 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,1,1,time/)) |
---|
1485 | if (status/=nf90_noerr) then |
---|
1486 | write(*,*)"unable to store variable ",varname, & |
---|
1487 | " in file ",file |
---|
1488 | stop |
---|
1489 | endif |
---|
1490 | ! |
---|
1491 | status = nf90_close(ncid) |
---|
1492 | ! |
---|
1493 | end subroutine Write_Ncdf_var3d_Real_t |
---|
1494 | |
---|
1495 | !> Subroutine to write a 3D real variable with time dependent in a given file |
---|
1496 | !!@param[in] varname : name of variable to store |
---|
1497 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1498 | !!@param[in] file : netcdf file name |
---|
1499 | !!@param[in] tabvar : values of the variable to write |
---|
1500 | !!@param[in] time : time corresponding to the values to write |
---|
1501 | !!@param[inout] idef : |
---|
1502 | !!@param[in] typevar : type of the variable to write |
---|
1503 | !!@return idef |
---|
1504 | subroutine Write_Ncdf_var4d_Real_t(varname,dimname,file,tabvar,time,idef,typevar) |
---|
1505 | ! |
---|
1506 | implicit none |
---|
1507 | ! |
---|
1508 | Character(*),intent(in) :: varname,file,typevar |
---|
1509 | Character(*),dimension(4),intent(in) :: dimname |
---|
1510 | Integer :: time |
---|
1511 | Integer,intent(inout) :: idef |
---|
1512 | Real*8, dimension(:,:,:,:), pointer :: tabvar |
---|
1513 | ! |
---|
1514 | ! local variables |
---|
1515 | ! |
---|
1516 | Integer :: dimid1,dimid2,dimid3,dimid4 |
---|
1517 | Integer :: status,ncid |
---|
1518 | Integer :: varid |
---|
1519 | ! |
---|
1520 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1521 | if (status/=nf90_noerr) then |
---|
1522 | write(*,*)"unable to open netcdf file : ",file |
---|
1523 | stop |
---|
1524 | endif |
---|
1525 | ! |
---|
1526 | if(time==1 .OR. idef==0) then |
---|
1527 | ! |
---|
1528 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1529 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1530 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1531 | status = nf90_inq_dimid(ncid,dimname(4), dimid4) |
---|
1532 | status = nf90_redef(ncid) |
---|
1533 | ! |
---|
1534 | select case(TRIM(typevar)) |
---|
1535 | case('double') |
---|
1536 | status = nf90_def_var(ncid,TRIM(varname),nf90_double, & |
---|
1537 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1538 | case('float') |
---|
1539 | status = nf90_def_var(ncid,TRIM(varname),nf90_float, & |
---|
1540 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1541 | end select |
---|
1542 | ! |
---|
1543 | status = nf90_enddef(ncid) |
---|
1544 | idef=1 |
---|
1545 | else |
---|
1546 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1547 | endif |
---|
1548 | ! |
---|
1549 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,1,1,time/)) |
---|
1550 | if (status/=nf90_noerr) then |
---|
1551 | write(*,*)"unable to store variable ",varname, & |
---|
1552 | " in file ",file |
---|
1553 | stop |
---|
1554 | endif |
---|
1555 | ! |
---|
1556 | status = nf90_close(ncid) |
---|
1557 | ! |
---|
1558 | end subroutine Write_Ncdf_var4d_Real_t |
---|
1559 | |
---|
1560 | |
---|
1561 | !> Subroutine to write a 3D real variable with time dependent in a given file |
---|
1562 | !!@param[in] varname : name of variable to store |
---|
1563 | !!@param[in] dimname : name of dimensions of the given variable |
---|
1564 | !!@param[in] file : netcdf file name |
---|
1565 | !!@param[in] tabvar : values of the variable to write |
---|
1566 | !!@param[in] time : time corresponding to the values to write |
---|
1567 | !!@param[in] level : level corresponding to the values to read |
---|
1568 | !!@param[inout] idef : |
---|
1569 | !!@param[in] typevar : type of the variable to write |
---|
1570 | !!@return idef |
---|
1571 | |
---|
1572 | subroutine Write_Ncdf_var4d_Real_nt(varname,dimname,file,tabvar,time,level,idef,typevar) |
---|
1573 | ! |
---|
1574 | implicit none |
---|
1575 | ! |
---|
1576 | Character(*),intent(in) :: varname,file,typevar |
---|
1577 | Character(*),dimension(4),intent(in) :: dimname |
---|
1578 | Integer :: time,level |
---|
1579 | Integer,intent(inout) :: idef |
---|
1580 | Real*8, dimension(:,:,:,:), pointer :: tabvar |
---|
1581 | ! |
---|
1582 | ! local variables |
---|
1583 | ! |
---|
1584 | Integer :: dimid1,dimid2,dimid3,dimid4 |
---|
1585 | Integer :: status,ncid |
---|
1586 | Integer :: varid |
---|
1587 | ! |
---|
1588 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1589 | if (status/=nf90_noerr) then |
---|
1590 | write(*,*)"unable to open netcdf file : ",file |
---|
1591 | stop |
---|
1592 | endif |
---|
1593 | ! |
---|
1594 | if((time==1.and.level==1) .OR. idef==0 ) then |
---|
1595 | ! |
---|
1596 | status = nf90_inq_dimid(ncid,dimname(1), dimid1) |
---|
1597 | status = nf90_inq_dimid(ncid,dimname(2), dimid2) |
---|
1598 | status = nf90_inq_dimid(ncid,dimname(3), dimid3) |
---|
1599 | status = nf90_inq_dimid(ncid,dimname(4), dimid4) |
---|
1600 | status = nf90_redef(ncid) |
---|
1601 | ! |
---|
1602 | select case(TRIM(typevar)) |
---|
1603 | case('double') |
---|
1604 | status = nf90_def_var(ncid,TRIM(varname),nf90_double, & |
---|
1605 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1606 | case('float') |
---|
1607 | status = nf90_def_var(ncid,TRIM(varname),nf90_float, & |
---|
1608 | (/dimid1,dimid2,dimid3,dimid4/),varid) |
---|
1609 | end select |
---|
1610 | ! |
---|
1611 | status = nf90_enddef(ncid) |
---|
1612 | idef=1 |
---|
1613 | else |
---|
1614 | status = nf90_inq_varid(ncid, varname, varid) |
---|
1615 | endif |
---|
1616 | ! |
---|
1617 | status = nf90_put_var(ncid,varid,tabvar,start=(/1,1,level,time/)) |
---|
1618 | if (status/=nf90_noerr) then |
---|
1619 | write(*,*)"unable to store variable ",varname, & |
---|
1620 | " in file ",file |
---|
1621 | stop |
---|
1622 | endif |
---|
1623 | ! |
---|
1624 | status = nf90_close(ncid) |
---|
1625 | ! |
---|
1626 | end subroutine Write_Ncdf_var4d_Real_nt |
---|
1627 | |
---|
1628 | !> Subroutine to retrieve names of all variables included in a given file |
---|
1629 | !!@param[in] filename : netcdf file name |
---|
1630 | !!@param[out] tabvarname : array containing various variables names |
---|
1631 | !!@return tabvarname |
---|
1632 | |
---|
1633 | subroutine Read_Ncdf_VarName(filename,tabvarname) |
---|
1634 | ! |
---|
1635 | Character(*),intent(in) :: filename |
---|
1636 | Character*20,dimension(:),pointer :: tabvarname |
---|
1637 | Integer :: nDimensions,nVariables |
---|
1638 | Integer :: nAttributes,unlimitedDimId,i |
---|
1639 | Integer :: ncid,status |
---|
1640 | ! |
---|
1641 | status = nf90_open(filename,NF90_NOWRITE,ncid) |
---|
1642 | if (status/=nf90_noerr) then |
---|
1643 | write(*,*)"unable to open netcdf file : ",filename |
---|
1644 | stop |
---|
1645 | endif |
---|
1646 | ! |
---|
1647 | status = nf90_inquire(ncid,nDimensions,nVariables,nAttributes, & |
---|
1648 | unlimitedDimId) |
---|
1649 | ! |
---|
1650 | Allocate(tabvarname(nVariables)) |
---|
1651 | ! |
---|
1652 | Do i=1,nVariables |
---|
1653 | status = nf90_inquire_variable(ncid,i,tabvarname(i)) |
---|
1654 | End do |
---|
1655 | |
---|
1656 | end subroutine Read_Ncdf_Varname |
---|
1657 | |
---|
1658 | |
---|
1659 | !> Subroutine copy a netcdf varaiable included in a given file to an other file |
---|
1660 | !!@param[in] varname : array containing various variables names |
---|
1661 | !!@param[in] filein : netcdf entry file name |
---|
1662 | !!@param[in] fileout : netcdf result file name |
---|
1663 | !!@return |
---|
1664 | |
---|
1665 | subroutine Copy_Ncdf_att_var(varname,filein,fileout) |
---|
1666 | ! |
---|
1667 | Character(*),intent(in) :: filein,fileout |
---|
1668 | Character(*),intent(in) :: varname |
---|
1669 | Integer :: ncid_in,ncid_out,status,varid_in,varid_out |
---|
1670 | ! |
---|
1671 | ! print *,'filein = ',filein,fileout |
---|
1672 | status = nf90_open(filein,NF90_NOWRITE,ncid_in) |
---|
1673 | if (status/=nf90_noerr) then |
---|
1674 | write(*,*)"unable to open input netcdf file : ",filein |
---|
1675 | stop |
---|
1676 | endif |
---|
1677 | ! |
---|
1678 | status = nf90_open(fileout,NF90_WRITE,ncid_out) |
---|
1679 | if (status/=nf90_noerr) then |
---|
1680 | write(*,*)"unable to open output netcdf file : ",fileout |
---|
1681 | stop |
---|
1682 | endif |
---|
1683 | ! |
---|
1684 | ! print *,'ici1' |
---|
1685 | status = nf90_inq_varid(ncid_in,varname,varid_in) |
---|
1686 | status = nf90_inq_varid(ncid_out,varname,varid_out) |
---|
1687 | ! |
---|
1688 | status = nf90_redef(ncid_out) |
---|
1689 | ! |
---|
1690 | status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out) |
---|
1691 | status = nf90_copy_att(ncid_in,varid_in,'valid_min',ncid_out,varid_out) |
---|
1692 | status = nf90_copy_att(ncid_in,varid_in,'valid_max',ncid_out,varid_out) |
---|
1693 | status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out) |
---|
1694 | status = nf90_copy_att(ncid_in,varid_in,'calendar',ncid_out,varid_out) |
---|
1695 | status = nf90_copy_att(ncid_in,varid_in,'title',ncid_out,varid_out) |
---|
1696 | status = nf90_copy_att(ncid_in,varid_in,'time_origin',ncid_out,varid_out) |
---|
1697 | status = nf90_copy_att(ncid_in,varid_in,'positive',ncid_out,varid_out) |
---|
1698 | status = nf90_copy_att(ncid_in,varid_in,'tstep_sec',ncid_out,varid_out) |
---|
1699 | status = nf90_copy_att(ncid_in,varid_in,'nav_model',ncid_out,varid_out) |
---|
1700 | status = nf90_copy_att(ncid_in,varid_in,'Minvalue=',ncid_out,varid_out) |
---|
1701 | status = nf90_copy_att(ncid_in,varid_in,'Maxvalue=',ncid_out,varid_out) |
---|
1702 | status = nf90_copy_att(ncid_in,varid_in,'short_name',ncid_out,varid_out) |
---|
1703 | status = nf90_copy_att(ncid_in,varid_in,'online_operation',ncid_out,varid_out) |
---|
1704 | status = nf90_copy_att(ncid_in,varid_in,'axis',ncid_out,varid_out) |
---|
1705 | status = nf90_copy_att(ncid_in,varid_in,'interval_operation',ncid_out,varid_out) |
---|
1706 | status = nf90_copy_att(ncid_in,varid_in,'interval_write',ncid_out,varid_out) |
---|
1707 | status = nf90_copy_att(ncid_in,varid_in,'associate',ncid_out,varid_out) |
---|
1708 | status = nf90_copy_att(ncid_in,varid_in,'actual_range',ncid_out,varid_out) |
---|
1709 | status = nf90_copy_att(ncid_in,varid_in,'longitude',ncid_out,varid_out) |
---|
1710 | status = nf90_copy_att(ncid_in,varid_in,'latitude',ncid_out,varid_out) |
---|
1711 | status = nf90_copy_att(ncid_in,varid_in,'scale_factor',ncid_out,varid_out) |
---|
1712 | status = nf90_copy_att(ncid_in,varid_in,'add_offset',ncid_out,varid_out) |
---|
1713 | status = nf90_copy_att(ncid_in,varid_in,'missing_value',ncid_out,varid_out) |
---|
1714 | ! |
---|
1715 | status = nf90_enddef(ncid_out) |
---|
1716 | ! |
---|
1717 | status = nf90_close(ncid_in) |
---|
1718 | status = nf90_close(ncid_out) |
---|
1719 | ! print *,'ici2' |
---|
1720 | ! |
---|
1721 | end subroutine Copy_Ncdf_att_var |
---|
1722 | |
---|
1723 | !> Subroutine copy a netcdf varaiable included in a given file to an other file at latlon |
---|
1724 | !!@param[in] varname : array containing various variables names |
---|
1725 | !!@param[in] filein : netcdf entry file name |
---|
1726 | !!@param[in] fileout : netcdf result file name |
---|
1727 | !!@param[in] min : |
---|
1728 | !!@param[in] max : |
---|
1729 | !!@return |
---|
1730 | |
---|
1731 | |
---|
1732 | subroutine Copy_Ncdf_att_latlon(varname,filein,fileout,min,max) |
---|
1733 | ! |
---|
1734 | Character(*),intent(in) :: filein,fileout |
---|
1735 | Character(*),intent(in) :: varname |
---|
1736 | real*8 :: min,max |
---|
1737 | Integer :: ncid_in,ncid_out,status,varid_in,varid_out |
---|
1738 | ! |
---|
1739 | status = nf90_open(filein,NF90_NOWRITE,ncid_in) |
---|
1740 | if (status/=nf90_noerr) then |
---|
1741 | write(*,*)"unable to open netcdf file : ",filein |
---|
1742 | stop |
---|
1743 | endif |
---|
1744 | ! |
---|
1745 | status = nf90_open(fileout,NF90_WRITE,ncid_out) |
---|
1746 | if (status/=nf90_noerr) then |
---|
1747 | write(*,*)"unable to open netcdf file : ",fileout |
---|
1748 | stop |
---|
1749 | endif |
---|
1750 | ! |
---|
1751 | status = nf90_inq_varid(ncid_in,varname,varid_in) |
---|
1752 | status = nf90_inq_varid(ncid_out,varname,varid_out) |
---|
1753 | ! |
---|
1754 | status = nf90_redef(ncid_out) |
---|
1755 | ! |
---|
1756 | select case (varname) |
---|
1757 | ! |
---|
1758 | case('nav_lon') |
---|
1759 | status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out) |
---|
1760 | status = nf90_put_att(ncid_out,varid_out,'valid_min',real(min,4)) |
---|
1761 | status = nf90_put_att(ncid_out,varid_out,'valid_max',real(max,4)) |
---|
1762 | status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out) |
---|
1763 | status = nf90_copy_att(ncid_in,varid_in,'nav_model',ncid_out,varid_out) |
---|
1764 | status = nf90_copy_att(ncid_in,varid_in,'title',ncid_out,varid_out) |
---|
1765 | ! |
---|
1766 | case('nav_lat') |
---|
1767 | status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out) |
---|
1768 | status = nf90_put_att(ncid_out,varid_out,'valid_min',real(min,4)) |
---|
1769 | status = nf90_put_att(ncid_out,varid_out,'valid_max',real(max,4)) |
---|
1770 | status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out) |
---|
1771 | status = nf90_copy_att(ncid_in,varid_in,'nav_model',ncid_out,varid_out) |
---|
1772 | status = nf90_copy_att(ncid_in,varid_in,'title',ncid_out,varid_out) |
---|
1773 | ! |
---|
1774 | end select |
---|
1775 | ! |
---|
1776 | status = nf90_enddef(ncid_out) |
---|
1777 | ! |
---|
1778 | status = nf90_close(ncid_in) |
---|
1779 | status = nf90_close(ncid_out) |
---|
1780 | end subroutine Copy_Ncdf_att_latlon |
---|
1781 | |
---|
1782 | !> Subroutine to retrieve value of a given dimension |
---|
1783 | !!@param[in] dimname : name of dimension to retrieve |
---|
1784 | !!@param[in] file : netcdf file name |
---|
1785 | !!@param[out] dimval : value of the required dimension |
---|
1786 | !!@return dimval |
---|
1787 | |
---|
1788 | subroutine Read_Ncdf_dim(dimname,file,dimval) |
---|
1789 | ! |
---|
1790 | implicit none |
---|
1791 | ! |
---|
1792 | Character(*),intent(in) :: dimname,file |
---|
1793 | Integer :: dimval |
---|
1794 | ! |
---|
1795 | ! local variables |
---|
1796 | ! |
---|
1797 | integer ncid,status,dimid |
---|
1798 | ! |
---|
1799 | status = nf90_open(file,NF90_NOWRITE,ncid) |
---|
1800 | if (status/=nf90_noerr) then |
---|
1801 | write(*,*)"unable to open netcdf file : ",file |
---|
1802 | stop |
---|
1803 | endif |
---|
1804 | ! |
---|
1805 | status = nf90_inq_dimid(ncid,dimname,dimid) |
---|
1806 | status = nf90_inquire_dimension(ncid,dimid,len=dimval) |
---|
1807 | ! |
---|
1808 | status = nf90_close(ncid) |
---|
1809 | ! |
---|
1810 | end subroutine Read_Ncdf_dim |
---|
1811 | |
---|
1812 | !> Subroutine to write a dimension in a given file |
---|
1813 | !!@param[in] dimname : name of dimension to retrieve |
---|
1814 | !!@param[in] file : netcdf file name |
---|
1815 | !!@param[out] dimval : value of the required dimension |
---|
1816 | !!@return dimval |
---|
1817 | |
---|
1818 | subroutine Write_Ncdf_dim(dimname,file,dimval) |
---|
1819 | ! |
---|
1820 | implicit none |
---|
1821 | ! |
---|
1822 | Character(*),intent(in) :: dimname,file |
---|
1823 | Integer :: dimval |
---|
1824 | ! |
---|
1825 | ! local variables |
---|
1826 | ! |
---|
1827 | integer ncid,status,dimid |
---|
1828 | ! |
---|
1829 | status = nf90_open(file,NF90_WRITE,ncid) |
---|
1830 | if (status/=nf90_noerr) then |
---|
1831 | write(*,*)"unable to open netcdf file : ",file |
---|
1832 | stop |
---|
1833 | endif |
---|
1834 | ! |
---|
1835 | status = nf90_redef(ncid) |
---|
1836 | If(dimval.eq.0) then |
---|
1837 | status = nf90_def_dim(ncid,dimname,nf90_unlimited,dimid) |
---|
1838 | Else |
---|
1839 | status = nf90_def_dim(ncid,dimname,dimval,dimid) |
---|
1840 | End If |
---|
1841 | status = nf90_enddef(ncid) |
---|
1842 | ! |
---|
1843 | status = nf90_close(ncid) |
---|
1844 | ! |
---|
1845 | end subroutine Write_Ncdf_dim |
---|
1846 | |
---|
1847 | !> Function to get the dimension of a given vriable in a given file |
---|
1848 | !!@param[in] varname : name of variable to retrieve dimension |
---|
1849 | !!@param[in] filename : netcdf file name |
---|
1850 | !!@return |
---|
1851 | |
---|
1852 | Integer function Get_NbDims( varname , filename ) |
---|
1853 | ! |
---|
1854 | Character(*),intent(in) :: varname,filename |
---|
1855 | Integer :: status,ncid,varid |
---|
1856 | ! |
---|
1857 | status = nf90_open(TRIM(filename),NF90_NOWRITE,ncid) |
---|
1858 | if (status/=nf90_noerr) then |
---|
1859 | write(*,*)"unable to open netcdf file : ",TRIM(filename) |
---|
1860 | stop |
---|
1861 | endif |
---|
1862 | status = nf90_inq_varid(ncid,TRIM(varname),varid) |
---|
1863 | status = nf90_inquire_variable(ncid, varid , ndims = Get_NbDims) |
---|
1864 | ! |
---|
1865 | return |
---|
1866 | ! |
---|
1867 | end function Get_NbDims |
---|
1868 | |
---|
1869 | |
---|
1870 | !> Function to test exitence dimension in a given file |
---|
1871 | !!@param[in] varname : name of variable to retrieve dimension |
---|
1872 | !!@param[in] filename : netcdf file name |
---|
1873 | !!@return |
---|
1874 | |
---|
1875 | ! function Get_NbDims_Existence |
---|
1876 | !************************************************************** |
---|
1877 | ! |
---|
1878 | Logical function Dims_Existence( dimname , filename ) |
---|
1879 | ! |
---|
1880 | Character(*),intent(in) :: dimname,filename |
---|
1881 | Integer :: status,ncid,dimid |
---|
1882 | ! |
---|
1883 | status = nf90_open(TRIM(filename),NF90_NOWRITE,ncid) |
---|
1884 | if (status/=nf90_noerr) then |
---|
1885 | write(*,*)"unable to open netcdf file : ",TRIM(filename) |
---|
1886 | stop |
---|
1887 | endif |
---|
1888 | status = nf90_inq_dimid(ncid,dimname,dimid) |
---|
1889 | ! |
---|
1890 | if (status/=nf90_noerr) then |
---|
1891 | Dims_Existence = .false. |
---|
1892 | else |
---|
1893 | Dims_Existence = .true. |
---|
1894 | endif |
---|
1895 | ! |
---|
1896 | return |
---|
1897 | ! |
---|
1898 | end function Dims_Existence |
---|
1899 | ! |
---|
1900 | !************************************************************** |
---|
1901 | !************************************************************** |
---|
1902 | |
---|
1903 | end module io_netcdf_grisli |
---|
1904 | |
---|
1905 | |
---|
1906 | |
---|
1907 | |
---|
1908 | |
---|
1909 | |
---|