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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 108 - (hide annotations)
Tue Sep 16 14:00:41 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/bibio/Writefield/createnewfield.f
File size: 1192 byte(s)
Imported writefield from LMDZ. Close at the end of gcm the files which
were created by writefiled (not done in LMDZ).

Removed procedures for the output of Grads files. Removed calls to
dump2d. In guide, replaced calls to wrgrads by calls to writefield.

In vlspltqs, removed redundant programming of saturation
pressure. Call foeew from module FCTTRE instead.

Bug fix in interpre: size of w exceeding size of correponding actual
argument wg in advtrac.

In leapfrog, call guide until the end of the run, instead of six hours
before the end.

Bug fix in readsulfate_preind: type of arguments.

1 guez 108 module createnewfield_m
2    
3     implicit none
4    
5     contains
6    
7     subroutine CreateNewField(name, dimx, dimy, dimz)
8    
9     USE netcdf, ONLY: nf90_clobber, nf90_double, nf90_unlimited
10     USE netcdf95, ONLY: nf95_create, nf95_def_dim, nf95_def_var, nf95_enddef
11     USE write_field, ONLY: fieldid, fieldindex, fieldname, fieldvarid, nbfield
12    
13     character(len = *), intent(in):: name
14     integer, intent(in):: dimx, dimy, dimz
15    
16     ! Local:
17     integer TabDim(4)
18    
19     !--------------------------------------------------------------------
20    
21     NbField = NbField + 1
22     FieldName(NbField) = TRIM(ADJUSTL(name))
23     FieldIndex(NbField) = 1
24    
25     call NF95_CREATE(TRIM(ADJUSTL(name)) // '.nc', NF90_CLOBBER, &
26     FieldId(NbField))
27     call NF95_DEF_DIM(FieldId(NbField), 'X', dimx, TabDim(1))
28     call NF95_DEF_DIM(FieldId(NbField), 'Y', dimy, TabDim(2))
29     call NF95_DEF_DIM(FieldId(NbField), 'Z', dimz, TabDim(3))
30     call NF95_DEF_DIM(FieldId(NbField), 'iter', NF90_UNLIMITED, TabDim(4))
31     call NF95_DEF_VAR(FieldId(NbField), FieldName(NbField), NF90_DOUBLE, &
32     TabDim, FieldVarId(NbField))
33     call NF95_ENDDEF(FieldId(NbField))
34    
35     end subroutine CreateNewField
36    
37     end module createnewfield_m

  ViewVC Help
Powered by ViewVC 1.1.21