/[lmdze]/trunk/Sources/misc/createnewfield.f
ViewVC logotype

Diff of /trunk/Sources/misc/createnewfield.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 108 by guez, Tue Sep 16 14:00:41 2014 UTC revision 109 by guez, Wed Sep 17 10:08:00 2014 UTC
# Line 4  module createnewfield_m Line 4  module createnewfield_m
4    
5  contains  contains
6    
7    subroutine CreateNewField(name, dimx, dimy, dimz)    subroutine CreateNewField(name, my_shape)
8    
9      USE netcdf, ONLY: nf90_clobber, nf90_double, nf90_unlimited      USE netcdf, ONLY: nf90_clobber, nf90_double, nf90_unlimited
10      USE netcdf95, ONLY: nf95_create, nf95_def_dim, nf95_def_var, nf95_enddef      USE netcdf95, ONLY: nf95_create, nf95_def_dim, nf95_def_var, nf95_enddef
11      USE write_field, ONLY: fieldid, fieldindex, fieldname, fieldvarid, nbfield      USE write_field, ONLY: ncid, record, fieldname, varid, nbfield
12    
13      character(len = *), intent(in):: name      character(len = *), intent(in):: name
14      integer, intent(in):: dimx, dimy, dimz      integer, intent(in):: my_shape(:)
15    
16      ! Local:      ! Local:
17      integer TabDim(4)      integer Dimid(size(my_shape) + 1), n_dim
18    
19      !--------------------------------------------------------------------      !--------------------------------------------------------------------
20    
21      NbField = NbField + 1      NbField = NbField + 1
22      FieldName(NbField) = TRIM(ADJUSTL(name))      FieldName(NbField) = ADJUSTL(name)
23      FieldIndex(NbField) = 1      Record(NbField) = 1
24        n_dim = size(my_shape)
25      call NF95_CREATE(TRIM(ADJUSTL(name)) // '.nc', NF90_CLOBBER, &  
26           FieldId(NbField))      call NF95_CREATE(TRIM(FieldName(NbField)) // '.nc', NF90_CLOBBER, &
27      call NF95_DEF_DIM(FieldId(NbField), 'X', dimx, TabDim(1))           Ncid(NbField))
28      call NF95_DEF_DIM(FieldId(NbField), 'Y', dimy, TabDim(2))      call NF95_DEF_DIM(Ncid(NbField), 'X', my_shape(1), Dimid(1))
29      call NF95_DEF_DIM(FieldId(NbField), 'Z', dimz, TabDim(3))      if (n_dim >= 2) then
30      call NF95_DEF_DIM(FieldId(NbField), 'iter', NF90_UNLIMITED, TabDim(4))         call NF95_DEF_DIM(Ncid(NbField), 'Y', my_shape(2), Dimid(2))
31      call NF95_DEF_VAR(FieldId(NbField), FieldName(NbField), NF90_DOUBLE, &         if (n_dim >= 3) then
32           TabDim, FieldVarId(NbField))            call NF95_DEF_DIM(Ncid(NbField), 'Z', my_shape(3), Dimid(3))
33      call NF95_ENDDEF(FieldId(NbField))         end if
34        end if
35        call NF95_DEF_DIM(Ncid(NbField), 'iter', NF90_UNLIMITED, Dimid(n_dim + 1))
36        call NF95_DEF_VAR(Ncid(NbField), FieldName(NbField), NF90_DOUBLE, &
37             Dimid, Varid(NbField))
38        call NF95_ENDDEF(Ncid(NbField))
39    
40    end subroutine CreateNewField    end subroutine CreateNewField
41    

Legend:
Removed from v.108  
changed lines
  Added in v.109

  ViewVC Help
Powered by ViewVC 1.1.21