New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
module_io.F90 in utils/tools/ABL_TOOLS – NEMO

source: utils/tools/ABL_TOOLS/module_io.F90 @ 11872

Last change on this file since 11872 was 11589, checked in by gsamson, 5 years ago

dev_r11265_ABL : see #2131

  • ABL_TOOLS first working version (README empty and arch files ignored for now)
File size: 35.2 KB
Line 
1MODULE module_io
2   !!======================================================================
3   !!                   ***  MODULE  module_io  ***
4   !! ABL preprocessing tool netcdf I/O subroutines
5   !!=====================================================================
6   !! History : 2016-10  (F. Lemarié)  Original code
7   !!----------------------------------------------------------------------
8
9   !!----------------------------------------------------------------------
10   !!
11   !!   FUNCTIONS   : Var_Existence   
12   !!   SUBROUTINES : Read_Ncdf_dim, Write_Ncdf_dim, Read_Ncdf_var1d_Real,
13   !!                 Read_Ncdf_var4d_Real_nt,       Read_Ncdf_var4d_Real_t,   
14   !!                 Write_Ncdf_var1d_Real       , Write_Ncdf_var2d_Real  ,
15   !!                 Write_Ncdf_var4d_Real_t     , Write_Ncdf_var1d_Real_t,
16   !!                 Duplicate_lon_lat_time
17   !!----------------------------------------------------------------------
18   USE netcdf
19   IMPLICIT NONE
20   
21   INTERFACE Read_Ncdf_var
22          MODULE PROCEDURE Read_Ncdf_var1d_Real   ,       &
23             &             Read_Ncdf_var2d_Real   ,       &
24             &             Read_Ncdf_var2d_Real_t ,       &
25             &             Read_Ncdf_var3d_Real_t ,       &
26             &             Read_Ncdf_var4d_Real_t ,       &
27             &             Read_Ncdf_var2d_Real_nt,       &
28             &             Read_Ncdf_var4d_Real_nt
29   END INTERFACE       
30!
31   INTERFACE Write_Ncdf_var
32          MODULE PROCEDURE Write_Ncdf_var1d_Real  ,       &
33             &             Write_Ncdf_var2d_Real  ,       &         
34             &             Write_Ncdf_var1d_Real_t,       &
35             &             Write_Ncdf_var2d_Real_t,       &
36             &             Write_Ncdf_var4d_Real_t
37   END INTERFACE     
38!     
39CONTAINS
40
41
42
43
44
45   LOGICAL FUNCTION Var_Existence( varname , filename )
46      !!---------------------------------------------------------------------
47      !!                    ***  FUNCTION Var_Existence  ***
48      !!                   
49      !! ** Purpose : check if variables varname exists in filename file 
50      !!
51      !!----------------------------------------------------------------------
52      CHARACTER(*),   intent(in) :: varname,filename
53      INTEGER                    :: status,ncid,varid       
54     
55      status = nf90_open(TRIM(filename),NF90_NOWRITE,ncid)
56     
57      IF ( status/=nf90_noerr ) THEN   
58            WRITE(*,*) "*** Var_Existence: unable to open netcdf file : ",TRIM(filename)
59            stop
60      END IF
61     
62      status = nf90_inq_varid( ncid, varname, varid )
63     
64      IF ( status/=nf90_noerr ) THEN
65         Var_Existence = .false.
66      ELSE
67         Var_Existence = .true.
68      END IF   
69
70   END FUNCTION Var_Existence
71   
72   
73   
74   SUBROUTINE Read_Ncdf_dim ( dimname, file, dimval )         
75      !!---------------------------------------------------------------------
76      !!                    ***  SUBROUTINE Read_Ncdf_dim  ***
77      !!                   
78      !! ** Purpose : read the integer dimension dimname in input file and   
79      !!              store it in dimval
80      !!
81      !!----------------------------------------------------------------------
82
83      CHARACTER(*),INTENT(in) :: dimname,file   
84      INTEGER                 :: dimval
85      INTEGER                 :: ncid,status,dimid
86      !     
87      status = nf90_open(file,NF90_NOWRITE,ncid)
88     
89      IF ( status/=nf90_noerr ) THEN   
90         WRITE(*,*)"*** Read_Ncdf_dim: unable to open netcdf file : ",trim(file)
91         STOP
92      END IF     
93      !     
94      status = nf90_inq_dimid        ( ncid, dimname,     dimid  ) 
95      status = nf90_inquire_dimension( ncid, dimid,   len=dimval )
96      status = nf90_close( ncid )
97
98   END SUBROUTINE Read_Ncdf_dim 
99
100
101
102
103
104
105
106   SUBROUTINE Write_Ncdf_dim( dimname, file, dimval )   
107      !!---------------------------------------------------------------------
108      !!                    ***  SUBROUTINE Write_Ncdf_dim  ***
109      !!                   
110      !! ** Purpose : write the dimension dimname in the output file   
111      !!
112      !!----------------------------------------------------------------------               
113      CHARACTER(*), INTENT(in) :: dimname,file   
114      INTEGER                  :: dimval
115      INTEGER                  :: ncid,status,dimid
116      !     
117      status = nf90_open(file,NF90_WRITE,ncid)
118      IF ( status/=nf90_noerr ) THEN   
119         WRITE(*,*)"*** Write_Ncdf_dim: to open netcdf file : ",file
120         STOP
121      END IF         
122     
123      status = nf90_redef(ncid)
124     
125      If( dimval.eq.0 ) THEN
126         status = nf90_def_dim(ncid,dimname,nf90_unlimited,dimid)       
127      ELSE
128         status = nf90_def_dim(ncid,dimname,dimval,dimid) 
129      END If   
130      !!
131      status = nf90_enddef(ncid)
132      status = nf90_close(ncid)
133      !
134   END SUBROUTINE Write_Ncdf_dim   
135
136
137
138
139
140
141   SUBROUTINE Read_Ncdf_var1d_Real( varname, file, tabvar )
142      !!---------------------------------------------------------------------
143      !!                    ***  SUBROUTINE Read_Ncdf_var1d_Real  ***
144      !!                   
145      !! ** Purpose : read the 1D variable varname in the input file   
146      !!              and store it in tabvar
147      !!
148      !!---------------------------------------------------------------------- 
149      CHARACTER(*),              INTENT(in)    :: varname,file
150      REAL(8), DIMENSION(:),ALLOCATABLE        :: tabvar
151      INTEGER, DIMENSION(1)                    :: dimID
152      INTEGER                                  :: dim1
153      INTEGER                                  :: status,ncid
154      INTEGER                                  :: varid             
155      !
156      status = nf90_open(file,NF90_NOWRITE,ncid)   
157      IF ( status/=nf90_noerr ) THEN   
158         WRITE(*,*) "*** Read_Ncdf_var1d_Real: unable to open netcdf file : ",file
159         STOP
160      END IF         
161      status = nf90_inq_varid(ncid,varname,varid)
162      status = nf90_inquire_variable(ncid,varid,dimids=dimID)
163      status = nf90_inquire_dimension(ncid,dimID(1),len=dim1)               
164      IF( .not. allocated(tabvar) ) THEN
165         ALLOCATE( tabvar ( dim1 ) ) 
166      ELSE
167         IF ( any(shape(tabvar)/=(/dim1/)) ) THEN
168            DEALLOCATE ( tabvar          )   
169            ALLOCATE   ( tabvar ( dim1 ) )   
170            WRITE(*,*) 'Warning change shape of array for ',trim(varname)   
171         END IF     
172      END IF       
173      status = nf90_get_var(ncid,varid,tabvar)     
174      status = nf90_close(ncid)
175      !
176   END SUBROUTINE Read_Ncdf_var1d_Real
177
178
179   SUBROUTINE Read_Ncdf_var2d_Real( varname, file, tabvar )   
180      !!---------------------------------------------------------------------
181      !!                    ***  SUBROUTINE Read_Ncdf_var2d_Real  ***
182      !!                   
183      !! ** Purpose : read the 2D variable varname in the input file   
184      !!              and store it in 2D tabvar
185      !!
186      !!----------------------------------------------------------------------         
187      CHARACTER(*),               INTENT(in)    :: varname, file
188      REAL(8), DIMENSION(:,:),    INTENT(inout) :: tabvar
189      INTEGER, DIMENSION(4)                     :: dimIDS
190      INTEGER                                   :: dim1,dim2
191      INTEGER                                   :: status,ncid
192      INTEGER                                   :: varid             
193      !
194      status = nf90_open(file,NF90_NOWRITE,ncid)     
195      IF (status/=nf90_noerr) THEN   
196         WRITE(*,*)"*** Read_Ncdf_var2d_Real: unable to open netcdf file : ",file
197         STOP
198      END IF         
199      status = nf90_inq_varid        (ncid , varname, varid)       
200      status = nf90_inquire_variable (ncid , varid, dimids=dimIDS)
201      status = nf90_inquire_dimension(ncid , dimIDS(1), len=dim1)
202      status = nf90_inquire_dimension(ncid , dimIDS(2), len=dim2)
203      status = nf90_get_var( ncid, varid, tabvar(:,:),      &
204        &                start = (/1,1/), count=(/dim1,dim2/))                                           
205      IF (status/=nf90_noerr) THEN   
206         WRITE(*,*)"unable to retrieve netcdf variable : ",trim(varname)
207         WRITE(*,*)"in file : ",trim(file)
208         WRITE(*,*) "error code: ", status
209         STOP
210      END IF         
211      status = nf90_close(ncid)
212      !   
213   END SUBROUTINE Read_Ncdf_var2d_Real
214
215   
216   SUBROUTINE Read_Ncdf_var2d_Real_nt( varname, file, tabvar, time, level )   
217      !!---------------------------------------------------------------------
218      !!                    ***  SUBROUTINE Read_Ncdf_var2d_Real_nt  ***
219      !!                   
220      !! ** Purpose : read the 4D variable varname in the input file   
221      !!              for a given time and vertical level, and store it in 2D tabvar
222      !!
223      !!----------------------------------------------------------------------         
224      CHARACTER(*),               INTENT(in)    :: varname, file
225      INTEGER     ,               INTENT(in)    :: time, level
226      !REAL(8), DIMENSION(:,:),    ALLOCATABLE   :: tabvar
227      REAL(8), DIMENSION(:,:),    INTENT(inout) :: tabvar
228      INTEGER, DIMENSION(4)                     :: dimIDS
229      INTEGER                                   :: dim1,dim2
230      INTEGER                                   :: status,ncid
231      INTEGER                                   :: varid             
232      !
233      status = nf90_open(file,NF90_NOWRITE,ncid)     
234      IF (status/=nf90_noerr) THEN   
235         WRITE(*,*)"*** Read_Ncdf_var2d_Real_nt: unable to open netcdf file : ",file
236         STOP
237      END IF         
238      status = nf90_inq_varid        (ncid , varname, varid)       
239      status = nf90_inquire_variable (ncid , varid, dimids=dimIDS)
240      status = nf90_inquire_dimension(ncid , dimIDS(1), len=dim1)
241      status = nf90_inquire_dimension(ncid , dimIDS(2), len=dim2)
242      !IF( .not. allocated( tabvar ) ) then
243      !   ALLOCATE ( tabvar( dim1, dim2 ) ) 
244      !ELSE
245      !   IF ( (size(tabvar,1) /= dim1) .OR. (size(tabvar,2) /= dim2) ) THEN
246      !      DEALLOCATE( tabvar )
247      !      ALLOCATE  ( tabvar (dim1, dim2 ) )
248      !   END IF     
249      !END IF       
250      status = nf90_get_var( ncid, varid, tabvar(:,:),      &
251        &                start = (/1,1,level,time/), count=(/dim1,dim2,1,1/))                                           
252      IF (status/=nf90_noerr) THEN   
253         WRITE(*,*)"unable to retrieve netcdf variable : ",trim(varname)
254         WRITE(*,*)"in file : ",trim(file)
255         WRITE(*,*) "error code: ", status
256         STOP
257      END IF         
258      status = nf90_close(ncid)
259      !   
260   END SUBROUTINE Read_Ncdf_var2d_Real_nt
261
262
263 
264 
265   SUBROUTINE Read_Ncdf_var4d_Real_nt( varname, file, tabvar, time, level )   
266      !!---------------------------------------------------------------------
267      !!                    ***  SUBROUTINE Read_Ncdf_var4d_Real_nt  ***
268      !!                   
269      !! ** Purpose : read the 4D variable varname in the input file   
270      !!              for a given time and vertical level, and store it in 4D tabvar
271      !!
272      !!----------------------------------------------------------------------         
273      CHARACTER(*),               INTENT(in)    :: varname,file
274      INTEGER     ,               INTENT(in)    :: time,level
275      !REAL(8), DIMENSION(:,:,:,:),ALLOCATABLE   :: tabvar
276      REAL(8), DIMENSION(:,:,:,:),INTENT(inout) :: tabvar
277      INTEGER, DIMENSION(4)                     :: dimIDS
278      INTEGER                                   :: dim1,dim2
279      INTEGER                                   :: status,ncid
280      INTEGER                                   :: varid             
281      !
282      status = nf90_open(file,NF90_NOWRITE,ncid)     
283      IF (status/=nf90_noerr) THEN   
284         WRITE(*,*)"*** Read_Ncdf_var4d_Real_nt: unable to open netcdf file : ",file
285         STOP
286      END IF         
287      status = nf90_inq_varid        (ncid , varname, varid)       
288      status = nf90_inquire_variable (ncid , varid, dimids=dimIDS)
289      status = nf90_inquire_dimension(ncid , dimIDS(1), len=dim1)
290      status = nf90_inquire_dimension(ncid , dimIDS(2), len=dim2)
291      !IF( .not. allocated( tabvar ) ) then
292      !   ALLOCATE ( tabvar( dim1, dim2, 1, 1 ) ) 
293      !ELSE
294      !   IF ( (size(tabvar,1) /= dim1) .OR. (size(tabvar,2) /= dim2) ) THEN
295      !      DEALLOCATE( tabvar )
296      !      ALLOCATE  ( tabvar (dim1, dim2, 1, 1 ) )
297      !   END IF     
298      !END IF       
299      status = nf90_get_var( ncid, varid, tabvar(:,:,:,:),      &
300        &                start = (/1,1,level,time/), count=(/dim1,dim2,1,1/))                                           
301      IF (status/=nf90_noerr) THEN   
302         WRITE(*,*)"unable to retrieve netcdf variable : ",trim(varname)
303         WRITE(*,*)"in file : ",trim(file)
304         WRITE(*,*) "error code: ", status
305         STOP
306      END IF         
307      status = nf90_close(ncid)
308      !   
309   END SUBROUTINE Read_Ncdf_var4d_Real_nt
310     
311 
312 
313 
314 
315 
316   SUBROUTINE Read_Ncdf_var3d_Real_t( varname, file, tabvar, time)   
317      !!---------------------------------------------------------------------
318      !!                    ***  SUBROUTINE Read_Ncdf_var3d_Real_t  ***
319      !!                   
320      !! ** Purpose : read the 3D variable varname in the input file   
321      !!              for a given time level, and store it in tabvar
322      !!
323      !!----------------------------------------------------------------------           
324      CHARACTER(*),               INTENT(in)    :: varname,file
325      INTEGER     ,               INTENT(in)    :: time
326      !REAL(8), DIMENSION(:,:,:),ALLOCATABLE     :: tabvar
327      REAL(8), DIMENSION(:,:,:), INTENT(inout)  :: tabvar
328      INTEGER, DIMENSION(3)                     :: dimIDS
329      INTEGER                                   :: dim1,dim2
330      INTEGER                                   :: status,ncid
331      INTEGER                                   :: varid             
332      !
333      status = nf90_open(file,NF90_NOWRITE,ncid)     
334      IF (status/=nf90_noerr) then   
335         WRITE(*,*)"*** Read_Ncdf_var3d_Real_t: unable to open netcdf file : ",file
336         STOP
337      END IF
338
339      status = nf90_inq_varid        (ncid, varname,varid)       
340      status = nf90_inquire_variable (ncid, varid,dimids=dimIDS)
341      status = nf90_inquire_dimension(ncid, dimIDS(1),len=dim1)
342      status = nf90_inquire_dimension(ncid, dimIDS(2),len=dim2)
343
344      !IF( .not. allocated(tabvar) ) ALLOCATE(tabvar(dim1,dim2,1)) 
345
346      status=nf90_get_var(ncid,varid,tabvar,start=(/1,1,time/))                                             
347      IF ( status/=nf90_noerr ) THEN   
348         WRITE(*,*)"unable to retrieve netcdf variable : ",trim(varname)
349         WRITE(*,*)"in file : ",trim(file)
350         WRITE(*,*) "error code: ", status
351         STOP
352      END IF           
353      status = nf90_close(ncid)
354      !
355   END SUBROUTINE Read_Ncdf_var3d_Real_t       
356 
357     
358   SUBROUTINE Read_Ncdf_var2d_Real_t(varname,file,tabvar,time)   
359      !!---------------------------------------------------------------------
360      !!                    ***  SUBROUTINE Read_Ncdf_var2d_Real_t  ***
361      !!                   
362      !! ** Purpose : read the 3D variable varname in the input file   
363      !!              for a given time level, and store it in tabvar
364      !!
365      !!----------------------------------------------------------------------           
366      CHARACTER(*),               INTENT(in)    :: varname,file
367      INTEGER     ,               INTENT(in)    :: time
368      !REAL(8), DIMENSION(:,:),ALLOCATABLE       :: tabvar
369      REAL(8), DIMENSION(:,:), INTENT(inout)    :: tabvar
370      INTEGER, DIMENSION(3)                     :: dimIDS
371      INTEGER                                   :: dim1,dim2
372      INTEGER                                   :: status,ncid
373      INTEGER                                   :: varid             
374      !
375      status = nf90_open(file,NF90_NOWRITE,ncid)     
376      IF (status/=nf90_noerr) then   
377         WRITE(*,*)"*** Read_Ncdf_var2d_Real_t: unable to open netcdf file : ",file
378         STOP
379      END IF
380
381      status = nf90_inq_varid        (ncid, varname,varid)       
382      status = nf90_inquire_variable (ncid, varid,dimids=dimIDS)
383      status = nf90_inquire_dimension(ncid, dimIDS(1),len=dim1)
384      status = nf90_inquire_dimension(ncid, dimIDS(2),len=dim2)
385
386      !IF( .not. allocated(tabvar) ) ALLOCATE(tabvar(dim1,dim2)) 
387
388      status=nf90_get_var(ncid,varid,tabvar,start=(/1,1,time/))                                             
389      IF ( status/=nf90_noerr ) THEN   
390         WRITE(*,*)"unable to retrieve netcdf variable : ",trim(varname)
391         WRITE(*,*)"in file : ",trim(file)
392         STOP
393      END IF           
394      status = nf90_close(ncid)
395      !
396   END SUBROUTINE Read_Ncdf_var2d_Real_t
397
398
399   SUBROUTINE Read_Ncdf_var4d_Real_t(varname,file,tabvar,time)   
400      !!---------------------------------------------------------------------
401      !!                    ***  SUBROUTINE Read_Ncdf_var4d_Real_t  ***
402      !!                   
403      !! ** Purpose : read the 4D variable varname in the input file   
404      !!              for a given time level, and store it in tabvar
405      !!
406      !!----------------------------------------------------------------------           
407      CHARACTER(*),               INTENT(in)    :: varname,file
408      INTEGER     ,               INTENT(in)    :: time
409      !REAL(8), DIMENSION(:,:,:,:),ALLOCATABLE   :: tabvar
410      REAL(8), DIMENSION(:,:,:,:), INTENT(inout):: tabvar
411      INTEGER, DIMENSION(4)                     :: dimIDS
412      INTEGER                                   :: dim1,dim2,dim3
413      INTEGER                                   :: status,ncid
414      INTEGER                                   :: varid             
415      !
416      status = nf90_open(file,NF90_NOWRITE,ncid)     
417      IF (status/=nf90_noerr) then   
418         WRITE(*,*)"*** Read_Ncdf_var4d_Real_t: unable to open netcdf file : ",file
419         STOP
420      END IF         
421
422      status = nf90_inq_varid        (ncid, varname,varid)       
423      status = nf90_inquire_variable (ncid, varid,dimids=dimIDS)
424      status = nf90_inquire_dimension(ncid, dimIDS(1),len=dim1)
425      status = nf90_inquire_dimension(ncid, dimIDS(2),len=dim2)
426      status = nf90_inquire_dimension(ncid, dimIDS(3),len=dim3)
427
428      !IF( .not. allocated(tabvar) ) ALLOCATE(tabvar(dim1,dim2,dim3,1))
429
430      status=nf90_get_var(ncid,varid,tabvar,start=(/1,1,1,time/))                                             
431      IF ( status/=nf90_noerr ) THEN   
432         WRITE(*,*) "unable to retrieve netcdf variable : ",trim(varname)
433         WRITE(*,*) "in file : ",trim(file)
434         WRITE(*,*) "error code: ", status
435         STOP
436      END IF           
437      status = nf90_close(ncid)
438      !
439   END SUBROUTINE Read_Ncdf_var4d_Real_t       
440   
441
442
443
444
445
446   SUBROUTINE Write_Ncdf_var1d_Real( varname, dimname, file, tabvar, typevar )
447      !!---------------------------------------------------------------------
448      !!                    ***  SUBROUTINE Write_Ncdf_var1d_Real  ***
449      !!                   
450      !! ** Purpose : write the 1D variable varname stored in tabvar
451      !!              in the output file using typevar type (float or double)             
452      !!
453      !!----------------------------------------------------------------------             
454      CHARACTER(*),               INTENT(in)    :: varname,file,dimname,typevar
455      REAL(8), DIMENSION(:), INTENT(in)         :: tabvar
456      INTEGER                                   :: dimid
457      INTEGER                                   :: status,ncid
458      INTEGER                                   :: varid             
459!
460      status = nf90_open(file,NF90_WRITE,ncid)       
461      IF (status/=nf90_noerr) then   
462         WRITE(*,*)"*** Write_Ncdf_var1d_Real: unable to open netcdf file : ",file
463         STOP
464      END IF         
465      status = nf90_inq_dimid(ncid,dimname, dimid)
466      status = nf90_redef(ncid)
467      SELECT CASE( TRIM(typevar) )
468      CASE('double')
469         status = nf90_def_var(ncid,varname,nf90_double,(/dimid/),varid)
470      CASE('float')
471         status = nf90_def_var(ncid,varname,nf90_float,(/dimid/),varid)     
472      END SELECT
473      status = nf90_enddef ( ncid )
474      status = nf90_put_var( ncid ,varid,tabvar)       
475      status = nf90_close  ( ncid )
476      !
477   END SUBROUTINE Write_Ncdf_var1d_Real 
478     
479     
480     
481
482
483   SUBROUTINE Write_Ncdf_var2d_Real( varname, dimname, file, tabvar, typevar )
484      !!---------------------------------------------------------------------
485      !!                    ***  SUBROUTINE Write_Ncdf_var2d_Real  ***
486      !!                   
487      !! ** Purpose : write the 2D variable varname stored in tabvar
488      !!              in the output file using typevar type (float or double)             
489      !!
490      !!----------------------------------------------------------------------         
491      CHARACTER(*),               INTENT(in)    :: varname, file, typevar
492      CHARACTER(*), DIMENSION(2), INTENT(in)    :: dimname     
493      REAL(8), DIMENSION(:,:), INTENT(in)       :: tabvar
494      INTEGER                                   :: dimid1, dimid2
495      INTEGER                                   :: status, ncid
496      INTEGER                                   :: varid             
497!
498      status = nf90_open(file,NF90_WRITE,ncid)       
499      IF (status/=nf90_noerr) then   
500         WRITE(*,*)"*** Write_Ncdf_var2d_Real: unable to open netcdf file : ",file
501         WRITE(*,*)"*** Write_Ncdf_var2d_Real: variable : ", varname
502         STOP
503      END IF
504         status = nf90_inq_dimid(ncid,dimname(1), dimid1)
505         status = nf90_inq_dimid(ncid,dimname(2), dimid2)
506         status = nf90_inq_varid(ncid,varname,varid)
507         IF (status /= nf90_noerr) THEN
508            status = nf90_redef(ncid)   
509            SELECT CASE( TRIM(typevar) )
510            CASE('double')
511               status = nf90_def_var(ncid,varname,nf90_double,     &
512               &                    (/dimid1,dimid2/),varid)
513            CASE('float')
514               status = nf90_def_var(ncid,varname,nf90_float,      &
515               &                    (/dimid1,dimid2/),varid)     
516            END SELECT
517            status = nf90_enddef(ncid)
518         END IF
519      status = nf90_put_var(ncid,varid,tabvar)     
520      status = nf90_close(ncid)
521      !
522   END SUBROUTINE Write_Ncdf_var2d_Real 
523
524
525
526
527
528
529   SUBROUTINE Write_Ncdf_var4d_Real_t( varname, dimname, file, tabvar, time, typevar )
530      !!---------------------------------------------------------------------
531      !!                    ***  SUBROUTINE Write_Ncdf_var4d_Real_t  ***
532      !!                   
533      !! ** Purpose : write the 4D variable varname stored in tabvar
534      !!              in the output file at time level time using typevar type (float or double)             
535      !!
536      !!---------------------------------------------------------------------- 
537      CHARACTER(*),                INTENT(in)   :: varname,file,typevar
538      CHARACTER(*), DIMENSION(4) , INTENT(in)   :: dimname
539      INTEGER, INTENT(in)                       :: time
540      REAL(8), DIMENSION(:,:,:,:), INTENT(in)   :: tabvar
541      INTEGER                                   :: dimid1,dimid2,dimid3,dimid4
542      INTEGER                                   :: status,ncid
543      INTEGER                                   :: varid             
544!
545      status = nf90_open(file,NF90_WRITE,ncid)       
546      IF ( status/=nf90_noerr ) THEN   
547         WRITE(*,*)"*** Write_Ncdf_var4d_Real_t: unable to open netcdf file : ",file
548         STOP
549      END IF           
550      IF ( time==1 .and. (TRIM(typevar)=='double'  .or. TRIM(typevar)=='float') ) THEN     
551         status = nf90_inq_dimid(ncid,dimname(1), dimid1)
552         status = nf90_inq_dimid(ncid,dimname(2), dimid2)
553         status = nf90_inq_dimid(ncid,dimname(3), dimid3)
554         status = nf90_inq_dimid(ncid,dimname(4), dimid4)
555         status = nf90_redef(ncid)     
556         SELECT CASE(TRIM(typevar))
557         CASE('double')
558            status = nf90_def_var(ncid,TRIM(varname),nf90_double,     &
559            &                     (/dimid1,dimid2,dimid3,dimid4/),varid) 
560         CASE('float')
561            status = nf90_def_var(ncid,TRIM(varname),nf90_float,     &
562            &                     (/dimid1,dimid2,dimid3,dimid4/),varid)   
563         END SELECT
564         status = nf90_enddef(ncid)
565      ELSE
566         status = nf90_inq_varid(ncid, varname, varid)
567      END IF   
568      status = nf90_put_var(ncid,varid,tabvar,start=(/1,1,1,time/))
569      IF (status/=nf90_noerr) THEN   
570         WRITE(*,*)"unable to store variable ",varname, &
571                                          " in file ",file
572         WRITE(*,*)"erorr code: ", status
573         STOP
574      END IF       
575      status = nf90_close(ncid)
576      !
577   END SUBROUTINE Write_Ncdf_var4d_Real_t           
578
579
580
581
582
583
584
585   
586   SUBROUTINE Write_Ncdf_var1d_Real_t(varname,dimname,file,tabvar,time,typevar)
587      !!---------------------------------------------------------------------
588      !!                    ***  SUBROUTINE Write_Ncdf_var1d_Real_t  ***
589      !!                   
590      !! ** Purpose : write the 1D variable varname stored in tabvar
591      !!              in the output file at time level time using typevar type (float or double)             
592      !!
593      !!---------------------------------------------------------------------- 
594      CHARACTER(*),               INTENT(in)    :: varname,file,typevar
595      CHARACTER(*), DIMENSION(1) ,INTENT(in)    :: dimname
596      INTEGER, INTENT(in)                       :: time
597      REAL(8), DIMENSION(:), INTENT(in)         :: tabvar
598      INTEGER                                   :: dimid1
599      INTEGER                                   :: status,ncid
600      INTEGER                                   :: varid             
601      !
602      status = nf90_open(file,NF90_WRITE,ncid)       
603      IF (status/=nf90_noerr) THEN   
604         WRITE(*,*) "*** Write_Ncdf_var1d_Real_t: unable to open netcdf file : ",file
605         STOP
606      END IF         
607           
608      IF( time==1 ) THEN
609         status = nf90_inq_dimid(ncid,dimname(1), dimid1) 
610         status = nf90_redef(ncid)
611         SELECT CASE ( TRIM(typevar) )
612         CASE ('double')
613            status = nf90_def_var(ncid,varname,nf90_double,     &
614             &                                (/dimid1/),varid) 
615         CASE ('float')
616            status = nf90_def_var(ncid,varname,nf90_float,     &
617             &                                (/dimid1/),varid)   
618         END SELECT     
619         status = nf90_enddef(ncid)
620      ELSE
621         status = nf90_inq_varid(ncid, varname, varid)
622      END IF 
623         
624      status = nf90_put_var(ncid,varid,tabvar,start=(/time/))
625     
626      IF (status/=nf90_noerr) THEN   
627         WRITE(*,*) "unable to store variable ",varname, &
628                                        " in file ",file
629         STOP
630      END IF       
631      status = nf90_close(ncid)
632     
633   END SUBROUTINE Write_Ncdf_var1d_Real_t 
634
635
636
637   SUBROUTINE Write_Ncdf_var2d_Real_t( varname, dimname, file, tabvar, time, typevar)
638      !!---------------------------------------------------------------------
639      !!                    ***  SUBROUTINE Write_Ncdf_var2d_Real_t  ***
640      !!                   
641      !! ** Purpose : write the 1D variable varname stored in tabvar
642      !!              in the output file at time level time using typevar type (float or double)             
643      !!
644      !!---------------------------------------------------------------------- 
645      CHARACTER(*),               INTENT(in)    :: varname, file, typevar
646      CHARACTER(*), DIMENSION(2), INTENT(in)    :: dimname
647      INTEGER, INTENT(in)                       :: time
648      REAL(8), DIMENSION(:,:), INTENT(in)       :: tabvar
649      INTEGER                                   :: dimid1, dimid2, dimidt
650      INTEGER                                   :: status,ncid
651      INTEGER                                   :: varid
652      INTEGER, DIMENSION(2)                     :: dim_size
653      !
654      status = nf90_open( file, NF90_WRITE, ncid)       
655      IF (status/=nf90_noerr) THEN   
656         WRITE(*,*) "*** Write_Ncdf_var2d_Real_t: unable to open netcdf file : ",file
657         STOP
658      END IF         
659
660      IF( time==1 ) THEN
661         status = nf90_inq_dimid( ncid, dimname(1), dimid1) 
662         status = nf90_inq_dimid( ncid, dimname(2), dimid2)
663         status = nf90_inq_dimid( ncid, "time"    , dimidt)
664         !status = nf90_def_dim( ncid, "time", NF90_UNLIMITED, dimidt)
665         status = nf90_redef(ncid)
666         SELECT CASE ( TRIM(typevar) )
667         CASE ('double')
668            status = nf90_def_var( ncid, varname, nf90_double,     &
669             &                                (/dimid1,dimid2,dimidt/), varid) 
670         CASE ('float')
671            status = nf90_def_var( ncid, varname, nf90_float,     &
672             &                                (/dimid1,dimid2,dimidt/), varid)   
673         END SELECT
674         status = nf90_enddef(ncid)
675      ELSE
676         status = nf90_inq_varid( ncid, varname, varid )
677      END IF 
678
679      dim_size = SHAPE(tabvar)
680      status = nf90_put_var( ncid, varid, tabvar(:,:), start=(/1,1,time/), count=(/dim_size(1),dim_size(2),1/) )
681     
682      IF (status/=nf90_noerr) THEN   
683         WRITE(*,*) "unable to store variable ",varname, &
684                                        " in file ",file
685         WRITE(*,*) "error code: ", status
686         STOP
687      END IF       
688      status = nf90_close(ncid)
689     
690   END SUBROUTINE Write_Ncdf_var2d_Real_t 
691
692
693
694   SUBROUTINE add_globatt_real( file, att_name, att_value )
695
696      INTEGER      :: ncid, status
697      CHARACTER(*) :: file, att_name
698      REAL(8)      :: att_value
699
700      status = nf90_open( file, nf90_write, ncid)
701
702      ! Enter define mode so we can add the attribute
703      status = nf90_redef( ncid )
704
705      ! ...  put the range attribute, setting it to eight byte reals...
706      status = nf90_put_att( ncid, NF90_GLOBAL, att_name, att_value )
707
708      ! Leave define mode
709      status = nf90_enddef(ncid)
710
711   END SUBROUTINE add_globatt_real
712
713
714
715   SUBROUTINE Duplicate_lon_lat_time( file_in, file_out )
716      !!---------------------------------------------------------------------
717      !!                    ***  SUBROUTINE Duplicate_lon_lat_time  ***
718      !!                   
719      !! ** Purpose : duplicate the attribute of lon, lat and time variables           
720      !!              from file_in in file_out
721      !!---------------------------------------------------------------------- 
722      CHARACTER(*),               INTENT(in)    ::  file_in, file_out
723      INTEGER                                   ::  status
724      INTEGER                                   ::  ncid_in, ncid_out   
725      INTEGER                                   :: varid_in,varid_out 
726     
727      status = nf90_open(file_in ,NF90_NOWRITE,ncid_in )
728      status = nf90_open(file_out,NF90_WRITE  ,ncid_out)
729       
730      status = nf90_inq_varid(ncid_in,'lon',varid_in)
731      status = nf90_inq_varid(ncid_out,'lon',varid_out)
732      status = nf90_redef(ncid_out)     
733      status = nf90_copy_att(ncid_in,varid_in,'standard_name',ncid_out,varid_out)
734      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)
735      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out) 
736      status = nf90_copy_att(ncid_in,varid_in,'axis',ncid_out,varid_out)     
737      status = nf90_enddef(ncid_out) 
738     
739      status = nf90_inq_varid(ncid_in,'lat',varid_in)
740      status = nf90_inq_varid(ncid_out,'lat',varid_out)
741      status = nf90_redef(ncid_out)     
742      status = nf90_copy_att(ncid_in,varid_in,'standard_name',ncid_out,varid_out)
743      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)
744      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out) 
745      status = nf90_copy_att(ncid_in,varid_in,'axis',ncid_out,varid_out)     
746      status = nf90_enddef(ncid_out)       
747     
748      status = nf90_inq_varid(ncid_in,'time',varid_in)
749      status = nf90_inq_varid(ncid_out,'time',varid_out)
750      status = nf90_redef(ncid_out)     
751      status = nf90_copy_att(ncid_in,varid_in,'standard_name',ncid_out,varid_out)
752      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out)
753      status = nf90_copy_att(ncid_in,varid_in,'calendar',ncid_out,varid_out)   
754      status = nf90_copy_att(ncid_in,varid_in,'axis',ncid_out,varid_out) 
755      status = nf90_enddef(ncid_out)       
756     
757      status = nf90_close(ncid_in)
758      status = nf90_close(ncid_out)
759      !   
760   END SUBROUTINE Duplicate_lon_lat_time 
761
762
763
764   SUBROUTINE Duplicate_lev_hyb( file_in, file_out )
765      !!---------------------------------------------------------------------
766      !!                    ***  SUBROUTINE Duplicate_lon_lat_time  ***
767      !!                   
768      !! ** Purpose : duplicate the attribute of lon, lat and time variables           
769      !!              from file_in in file_out
770      !!---------------------------------------------------------------------- 
771      CHARACTER(*),               INTENT(in)    ::  file_in, file_out
772      INTEGER                                   ::  status
773      INTEGER                                   ::  ncid_in, ncid_out   
774      INTEGER                                   :: varid_in,varid_out 
775     
776      status = nf90_open(file_in ,NF90_NOWRITE,ncid_in )
777      status = nf90_open(file_out,NF90_WRITE  ,ncid_out)
778       
779      status = nf90_inq_varid(ncid_in,'hyai',varid_in)
780      status = nf90_inq_varid(ncid_out,'hyai',varid_out)
781      status = nf90_redef(ncid_out)     
782      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)
783      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out)     
784      status = nf90_enddef(ncid_out) 
785     
786      status = nf90_inq_varid(ncid_in,'hybi',varid_in)
787      status = nf90_inq_varid(ncid_out,'hybi',varid_out)
788      status = nf90_redef(ncid_out)     
789      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)
790      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out)   
791      status = nf90_enddef(ncid_out)       
792
793      status = nf90_inq_varid(ncid_in,'hyam',varid_in)
794      status = nf90_inq_varid(ncid_out,'hyam',varid_out)
795      status = nf90_redef(ncid_out)     
796      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)
797      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out)     
798      status = nf90_enddef(ncid_out) 
799     
800      status = nf90_inq_varid(ncid_in,'hybm',varid_in)
801      status = nf90_inq_varid(ncid_out,'hybm',varid_out)
802      status = nf90_redef(ncid_out)     
803      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)
804      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out)   
805      status = nf90_enddef(ncid_out) 
806     
807      status = nf90_inq_varid(ncid_in,'lev',varid_in)
808      status = nf90_inq_varid(ncid_out,'lev',varid_out)
809      status = nf90_redef(ncid_out)     
810      status = nf90_copy_att(ncid_in,varid_in,'standard_name',ncid_out,varid_out)
811      status = nf90_copy_att(ncid_in,varid_in,'long_name',ncid_out,varid_out)     
812      status = nf90_copy_att(ncid_in,varid_in,'formula',ncid_out,varid_out)   
813      status = nf90_copy_att(ncid_in,varid_in,'formula_terms',ncid_out,varid_out)               
814      status = nf90_copy_att(ncid_in,varid_in,'units',ncid_out,varid_out)
815      status = nf90_copy_att(ncid_in,varid_in,'positive',ncid_out,varid_out)   
816      status = nf90_enddef(ncid_out)       
817     
818      status = nf90_close(ncid_in)
819      status = nf90_close(ncid_out)
820      !   
821   END SUBROUTINE Duplicate_lev_hyb
822
823END MODULE module_io
Note: See TracBrowser for help on using the repository browser.