Actual source code: vecio.c

petsc-3.4.5 2014-06-29
  2: /*
  3:    This file contains simple binary input routines for vectors.  The
  4:    analogous output routines are within each vector implementation's
  5:    VecView (with viewer types PETSCVIEWERBINARY)
  6:  */

  8: #include <petscsys.h>
  9: #include <petscvec.h>         /*I  "petscvec.h"  I*/
 10: #include <petsc-private/vecimpl.h>
 11: #include <petscmat.h> /* so that MAT_FILE_CLASSID is defined */
 12: #include <petscviewerhdf5.h>

 16: static PetscErrorCode PetscViewerBinaryReadVecHeader_Private(PetscViewer viewer,PetscInt *rows)
 17: {
 19:   MPI_Comm       comm;
 20:   PetscInt       tr[2],type;

 23:   PetscObjectGetComm((PetscObject)viewer,&comm);
 24:   /* Read vector header */
 25:   PetscViewerBinaryRead(viewer,tr,2,PETSC_INT);
 26:   type = tr[0];
 27:   if (type != VEC_FILE_CLASSID) {
 28:     PetscLogEventEnd(VEC_Load,viewer,0,0,0);
 29:     if (type == MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix is next in file, not a vector as you requested");
 30:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a vector next in file");
 31:   }
 32:   *rows = tr[1];
 33:   return(0);
 34: }

 36: #if defined(PETSC_HAVE_MPIIO)
 39: static PetscErrorCode VecLoad_Binary_MPIIO(Vec vec, PetscViewer viewer)
 40: {
 42:   PetscMPIInt    lsize;
 43:   PetscScalar    *avec;
 44:   MPI_File       mfdes;
 45:   MPI_Offset     off;

 48:   VecGetArray(vec,&avec);
 49:   PetscMPIIntCast(vec->map->n,&lsize);

 51:   PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
 52:   PetscViewerBinaryGetMPIIOOffset(viewer,&off);
 53:   off += vec->map->rstart*sizeof(PetscScalar);
 54:   MPI_File_set_view(mfdes,off,MPIU_SCALAR,MPIU_SCALAR,(char*)"native",MPI_INFO_NULL);
 55:   MPIU_File_read_all(mfdes,avec,lsize,MPIU_SCALAR,MPI_STATUS_IGNORE);
 56:   PetscViewerBinaryAddMPIIOOffset(viewer,vec->map->N*sizeof(PetscScalar));

 58:   VecRestoreArray(vec,&avec);
 59:   VecAssemblyBegin(vec);
 60:   VecAssemblyEnd(vec);
 61:   return(0);
 62: }
 63: #endif

 67: PetscErrorCode VecLoad_Binary(Vec vec, PetscViewer viewer)
 68: {
 69:   PetscMPIInt    size,rank,tag;
 70:   int            fd;
 71:   PetscInt       i,rows = 0,n,*range,N,bs;
 73:   PetscBool      flag;
 74:   PetscScalar    *avec,*avecwork;
 75:   MPI_Comm       comm;
 76:   MPI_Request    request;
 77:   MPI_Status     status;
 78: #if defined(PETSC_HAVE_MPIIO)
 79:   PetscBool      useMPIIO;
 80: #endif

 83:   PetscObjectGetComm((PetscObject)viewer,&comm);
 84:   MPI_Comm_rank(comm,&rank);
 85:   MPI_Comm_size(comm,&size);

 87:   PetscViewerBinaryGetDescriptor(viewer,&fd);
 88:   PetscViewerBinaryReadVecHeader_Private(viewer,&rows);
 89:   /* Set Vec sizes,blocksize,and type if not already set. Block size first so that local sizes will be compatible. */
 90:   PetscOptionsGetInt(((PetscObject)vec)->prefix, "-vecload_block_size", &bs, &flag);
 91:   if (flag) {
 92:     VecSetBlockSize(vec, bs);
 93:   }
 94:   if (vec->map->n < 0 && vec->map->N < 0) {
 95:     VecSetSizes(vec,PETSC_DECIDE,rows);
 96:   }

 98:   /* If sizes and type already set,check if the vector global size is correct */
 99:   VecGetSize(vec, &N);
100:   if (N != rows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Vector in file different length (%d) then input vector (%d)", rows, N);

102: #if defined(PETSC_HAVE_MPIIO)
103:   PetscViewerBinaryGetMPIIO(viewer,&useMPIIO);
104:   if (useMPIIO) {
105:     VecLoad_Binary_MPIIO(vec, viewer);
106:     return(0);
107:   }
108: #endif

110:   VecGetLocalSize(vec,&n);
111:   PetscObjectGetNewTag((PetscObject)viewer,&tag);
112:   VecGetArray(vec,&avec);
113:   if (!rank) {
114:     PetscBinaryRead(fd,avec,n,PETSC_SCALAR);

116:     if (size > 1) {
117:       /* read in other chuncks and send to other processors */
118:       /* determine maximum chunck owned by other */
119:       range = vec->map->range;
120:       n = 1;
121:       for (i=1; i<size; i++) n = PetscMax(n,range[i+1] - range[i]);

123:       PetscMalloc(n*sizeof(PetscScalar),&avecwork);
124:       for (i=1; i<size; i++) {
125:         n    = range[i+1] - range[i];
126:         PetscBinaryRead(fd,avecwork,n,PETSC_SCALAR);
127:         MPI_Isend(avecwork,n,MPIU_SCALAR,i,tag,comm,&request);
128:         MPI_Wait(&request,&status);
129:       }
130:       PetscFree(avecwork);
131:     }
132:   } else {
133:     MPI_Recv(avec,n,MPIU_SCALAR,0,tag,comm,&status);
134:   }

136:   VecRestoreArray(vec,&avec);
137:   VecAssemblyBegin(vec);
138:   VecAssemblyEnd(vec);
139:   return(0);
140: }

142: #if defined(PETSC_HAVE_HDF5)
145: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
146: {
147:   hid_t          file_id, group;
148:   const char     *groupName = NULL;

152:   PetscViewerHDF5GetFileId(viewer, &file_id);
153:   PetscViewerHDF5GetGroup(viewer, &groupName);
154:   /* Open group */
155:   if (groupName) {
156:     PetscBool root;

158:     PetscStrcmp(groupName, "/", &root);
159:     if (!root && !H5Lexists(file_id, groupName, H5P_DEFAULT)) {
160: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
161:       group = H5Gcreate2(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT);
162: #else /* deprecated HDF5 1.6 API */
163:       group = H5Gcreate(file_id, groupName, 0);
164: #endif
165:       if (group < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Could not create group %s", groupName);
166:       H5Gclose(group);
167:     }
168: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
169:     group = H5Gopen2(file_id, groupName, H5P_DEFAULT);
170: #else
171:     group = H5Gopen(file_id, groupName);
172: #endif
173:     if (group < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Could not open group %s", groupName);
174:   } else group = file_id;

176:   *fileId  = file_id;
177:   *groupId = group;
178:   return(0);
179: }

183: /*
184:      This should handle properly the cases where PetscInt is 32 or 64 and hsize_t is 32 or 64. These means properly casting with
185:    checks back and forth between the two types of variables.
186: */
187: PetscErrorCode VecLoad_HDF5(Vec xin, PetscViewer viewer)
188: {
189:   hid_t          file_id, group, dset_id, filespace, memspace, plist_id;
190:   hsize_t        rdim, dim;
191:   hsize_t        dims[4], count[4], offset[4];
192:   herr_t         status;
193:   PetscInt       n, N, bs = 1, bsInd, lenInd, low, timestep;
194:   PetscScalar    *x;
195:   const char     *vecname;

199:   PetscViewerHDF5OpenGroup(viewer, &file_id, &group);
200:   PetscViewerHDF5GetTimestep(viewer, &timestep);
201:   VecGetBlockSize(xin,&bs);
202:   /* Create the dataset with default properties and close filespace */
203:   PetscObjectGetName((PetscObject)xin,&vecname);
204: #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
205:   dset_id = H5Dopen2(group, vecname, H5P_DEFAULT);
206: #else
207:   dset_id = H5Dopen(group, vecname);
208: #endif
209:   if (dset_id == -1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Dopen() with Vec named %s",vecname);
210:   /* Retrieve the dataspace for the dataset */
211:   filespace = H5Dget_space(dset_id);
212:   if (filespace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Dget_space()");
213:   dim = 0;
214:   if (timestep >= 0) ++dim;
215:   ++dim;
216:   if (bs >= 1) ++dim;
217: #if defined(PETSC_USE_COMPLEX)
218:   ++dim;
219: #endif
220:   rdim = H5Sget_simple_extent_dims(filespace, dims, NULL);
221: #if defined(PETSC_USE_COMPLEX)
222:   bsInd = rdim-2;
223: #else
224:   bsInd = rdim-1;
225: #endif
226:   lenInd = timestep >= 0 ? 1 : 0;
227:   if (rdim != dim) {
228:     if (rdim == dim+1 && bs == -1) bs = dims[bsInd];
229:     else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected",rdim,dim);
230:   } else if (bs >= 1 && bs != (PetscInt) dims[bsInd]) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Block size %d specified for vector does not match blocksize in file %d",bs,dims[bsInd]);

232:   /* Set Vec sizes,blocksize,and type if not already set */
233:   if ((xin)->map->n < 0 && (xin)->map->N < 0) {
234:     VecSetSizes(xin, PETSC_DECIDE, dims[lenInd]*bs);
235:   }
236:   /* If sizes and type already set,check if the vector global size is correct */
237:   VecGetSize(xin, &N);
238:   if (N/bs != (PetscInt) dims[lenInd]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Vector in file different length (%d) then input vector (%d)", (PetscInt) dims[lenInd], N/bs);

240:   /* Each process defines a dataset and reads it from the hyperslab in the file */
241:   VecGetLocalSize(xin, &n);
242:   dim  = 0;
243:   if (timestep >= 0) {
244:     count[dim] = 1;
245:     ++dim;
246:   }
247:   PetscHDF5IntCast(n/bs,count + dim);
248:   ++dim;
249:   if (bs >= 1) {
250:     count[dim] = bs;
251:     ++dim;
252:   }
253: #if defined(PETSC_USE_COMPLEX)
254:   count[dim] = 2;
255:   ++dim;
256: #endif
257:   memspace = H5Screate_simple(dim, count, NULL);
258:   if (memspace == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Screate_simple()");

260:   /* Select hyperslab in the file */
261:   VecGetOwnershipRange(xin, &low, NULL);
262:   dim  = 0;
263:   if (timestep >= 0) {
264:     offset[dim] = timestep;
265:     ++dim;
266:   }
267:   PetscHDF5IntCast(low/bs,offset + dim);
268:   ++dim;
269:   if (bs >= 1) {
270:     offset[dim] = 0;
271:     ++dim;
272:   }
273: #if defined(PETSC_USE_COMPLEX)
274:   offset[dim] = 0;
275:   ++dim;
276: #endif
277:   status = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offset, NULL, count, NULL);CHKERRQ(status);

279:   /* Create property list for collective dataset read */
280:   plist_id = H5Pcreate(H5P_DATASET_XFER);
281:   if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Pcreate()");
282: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
283:   status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
284: #endif
285:   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */

287:   VecGetArray(xin, &x);
288:   status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
289:   VecRestoreArray(xin, &x);

291:   /* Close/release resources */
292:   if (group != file_id) {
293:     status = H5Gclose(group);CHKERRQ(status);
294:   }
295:   status = H5Pclose(plist_id);CHKERRQ(status);
296:   status = H5Sclose(filespace);CHKERRQ(status);
297:   status = H5Sclose(memspace);CHKERRQ(status);
298:   status = H5Dclose(dset_id);CHKERRQ(status);

300:   VecAssemblyBegin(xin);
301:   VecAssemblyEnd(xin);
302:   return(0);
303: }
304: #endif


309: PetscErrorCode  VecLoad_Default(Vec newvec, PetscViewer viewer)
310: {
312:   PetscBool      isbinary;
313: #if defined(PETSC_HAVE_HDF5)
314:   PetscBool      ishdf5;
315: #endif

318:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
319: #if defined(PETSC_HAVE_HDF5)
320:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
321: #endif

323: #if defined(PETSC_HAVE_HDF5)
324:   if (ishdf5) {
325:     if (!((PetscObject)newvec)->name) {
326:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Since HDF5 format gives ASCII name for each object in file; must use VecLoad() after setting name of Vec with PetscObjectSetName()");
327:       PetscLogEventEnd(VEC_Load,viewer,0,0,0);
328:     }
329:     VecLoad_HDF5(newvec, viewer);
330:   } else
331: #endif
332:   {
333:     VecLoad_Binary(newvec, viewer);
334:   }
335:   return(0);
336: }

340: /*@
341:   VecChop - Set all values in the vector less than the tolerance to zero

343:   Input Parameters:
344: + v   - The vector
345: - tol - The zero tolerance

347:   Output Parameters:
348: . v - The chopped vector

350:   Level: intermediate

352: .seealso: VecCreate(), VecSet()
353: @*/
354: PetscErrorCode VecChop(Vec v, PetscReal tol)
355: {
356:   PetscScalar    *a;
357:   PetscInt       n, i;

361:   VecGetLocalSize(v, &n);
362:   VecGetArray(v, &a);
363:   for (i = 0; i < n; ++i) {
364:     if (PetscAbsScalar(a[i]) < tol) a[i] = 0.0;
365:   }
366:   VecRestoreArray(v, &a);
367:   return(0);
368: }