Changeset 5270 for LMDZ6/trunk/libf/misc
- Timestamp:
- Oct 24, 2024, 1:55:38 PM (2 months ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ6/trunk/libf/misc/write_field.f90
r5268 r5270 3 3 ! 4 4 module write_field 5 USE netcdf, ONLY: nf90_sync, nf90_put_var, nf90_enddef, nf90_def_dim, nf90_unlimited, & 6 nf90_clobber, nf90_create, nf90_def_var, nf90_double 5 7 USE strings_mod, ONLY: int2str 6 8 IMPLICIT NONE; PRIVATE … … 75 77 subroutine WriteField_gen(name,Field,dimx,dimy,dimz) 76 78 implicit none 77 include 'netcdf.inc'78 79 character(len=*) :: name 79 80 integer :: dimx,dimy,dimz … … 104 105 count(4)=1 105 106 106 status = NF_PUT_VARA_DOUBLE(FieldId(Index),FieldVarId(Index),start,count,Field)107 status = NF_SYNC(FieldId(Index))107 status = nf90_put_var(FieldId(Index),FieldVarId(Index),Field,start,count) 108 status = nf90_sync(FieldId(Index)) 108 109 109 110 end subroutine WriteField_gen … … 111 112 subroutine CreateNewField(name,dimx,dimy,dimz) 112 113 implicit none 113 include 'netcdf.inc'114 114 character(len=*) :: name 115 115 integer :: dimx,dimy,dimz … … 123 123 124 124 125 status = NF_CREATE(TRIM(ADJUSTL(name))//'.nc', NF_CLOBBER, FieldId(NbField))126 status = NF_DEF_DIM(FieldId(NbField),'X',dimx,TabDim(1))127 status = NF_DEF_DIM(FieldId(NbField),'Y',dimy,TabDim(2))128 status = NF_DEF_DIM(FieldId(NbField),'Z',dimz,TabDim(3))129 status = NF_DEF_DIM(FieldId(NbField),'iter',NF_UNLIMITED,TabDim(4))130 status = NF_DEF_VAR(FieldId(NbField),FieldName(NbField),NF_DOUBLE,4,TabDim,FieldVarId(NbField))131 status = NF_ENDDEF(FieldId(NbField))125 status = nf90_create(TRIM(ADJUSTL(name))//'.nc', nf90_clobber, FieldId(NbField)) 126 status = nf90_def_dim(FieldId(NbField),'X',dimx,TabDim(1)) 127 status = nf90_def_dim(FieldId(NbField),'Y',dimy,TabDim(2)) 128 status = nf90_def_dim(FieldId(NbField),'Z',dimz,TabDim(3)) 129 status = nf90_def_dim(FieldId(NbField),'iter',nf90_unlimited,TabDim(4)) 130 status = nf90_def_var(FieldId(NbField),FieldName(NbField),nf90_double,TabDim,FieldVarId(NbField)) 131 status = nf90_enddef(FieldId(NbField)) 132 132 133 133 end subroutine CreateNewField
Note: See TracChangeset
for help on using the changeset viewer.