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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 109 - (hide annotations)
Wed Sep 17 10:08:00 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/bibio/Writefield/createnewfield.f
File size: 1297 byte(s)
Moved a call to writefield from guide to tau2alpha. (dxdys does not
change with itau.) So dxdys does not need to be a module variable any
longer. Other variables of modules tau2alpha_m downgraded to local
variables of tau2alpha, since they were not used elsewhere.

Procedures write_field[13]d and formcoord were never called. Could
then remove int2str.

Inline writefield_gen into writefield.

CreateNewField takes an integer array argument instead of 3 scalar
integers. CreateNewField now creates a number of dimensions adapted to
the rank of the output field, instead of always 4 dimensions.

Changed names of variables of module write_field: fieldid to
ncid, fieldindex to record, fieldvarid to varid.

In writefield_gen, if index == -1, no need to call GetFieldIndex
again, we know that the result is nbfield.

In guide, moved calls to writefield for some variables inside if
first_call: those variables do not change with time. Removed ztau:
computed only to be output, does not seem meaningful. Removed
writefield for aire: does not change with time and is already in
"grilles_gcm.nc".

1 guez 108 module createnewfield_m
2    
3     implicit none
4    
5     contains
6    
7 guez 109 subroutine CreateNewField(name, my_shape)
8 guez 108
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 guez 109 USE write_field, ONLY: ncid, record, fieldname, varid, nbfield
12 guez 108
13     character(len = *), intent(in):: name
14 guez 109 integer, intent(in):: my_shape(:)
15 guez 108
16     ! Local:
17 guez 109 integer Dimid(size(my_shape) + 1), n_dim
18 guez 108
19     !--------------------------------------------------------------------
20    
21     NbField = NbField + 1
22 guez 109 FieldName(NbField) = ADJUSTL(name)
23     Record(NbField) = 1
24     n_dim = size(my_shape)
25 guez 108
26 guez 109 call NF95_CREATE(TRIM(FieldName(NbField)) // '.nc', NF90_CLOBBER, &
27     Ncid(NbField))
28     call NF95_DEF_DIM(Ncid(NbField), 'X', my_shape(1), Dimid(1))
29     if (n_dim >= 2) then
30     call NF95_DEF_DIM(Ncid(NbField), 'Y', my_shape(2), Dimid(2))
31     if (n_dim >= 3) then
32     call NF95_DEF_DIM(Ncid(NbField), 'Z', my_shape(3), Dimid(3))
33     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 guez 108
40     end subroutine CreateNewField
41    
42     end module createnewfield_m

  ViewVC Help
Powered by ViewVC 1.1.21