Actual source code: vecio.c

petsc-dev 2014-04-22
Report Typos and Errors
  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:       PetscMalloc1(n,&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:   htri_t         found;
149:   const char     *groupName = NULL;

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

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

178:   *fileId  = file_id;
179:   *groupId = group;
180:   return(0);
181: }

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

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

238:   /* Set Vec sizes,blocksize,and type if not already set */
239:   if ((xin)->map->n < 0 && (xin)->map->N < 0) {
240:     VecSetSizes(xin, PETSC_DECIDE, dims[lenInd]*bs);
241:   }
242:   /* If sizes and type already set,check if the vector global size is correct */
243:   VecGetSize(xin, &N);
244:   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);

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

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

285:   /* Create property list for collective dataset read */
286:   plist_id = H5Pcreate(H5P_DATASET_XFER);
287:   if (plist_id == -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Could not H5Pcreate()");
288: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
289:   status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);CHKERRQ(status);
290: #endif
291:   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */

293:   VecGetArray(xin, &x);
294:   status = H5Dread(dset_id, H5T_NATIVE_DOUBLE, memspace, filespace, plist_id, x);CHKERRQ(status);
295:   VecRestoreArray(xin, &x);

297:   /* Close/release resources */
298:   if (group != file_id) {
299:     status = H5Gclose(group);CHKERRQ(status);
300:   }
301:   status = H5Pclose(plist_id);CHKERRQ(status);
302:   status = H5Sclose(filespace);CHKERRQ(status);
303:   status = H5Sclose(memspace);CHKERRQ(status);
304:   status = H5Dclose(dset_id);CHKERRQ(status);

306:   VecAssemblyBegin(xin);
307:   VecAssemblyEnd(xin);
308:   return(0);
309: }
310: #endif


315: PetscErrorCode  VecLoad_Default(Vec newvec, PetscViewer viewer)
316: {
318:   PetscBool      isbinary;
319: #if defined(PETSC_HAVE_HDF5)
320:   PetscBool      ishdf5;
321: #endif

324:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
325: #if defined(PETSC_HAVE_HDF5)
326:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
327: #endif

329: #if defined(PETSC_HAVE_HDF5)
330:   if (ishdf5) {
331:     if (!((PetscObject)newvec)->name) {
332:       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()");
333:       PetscLogEventEnd(VEC_Load,viewer,0,0,0);
334:     }
335:     VecLoad_HDF5(newvec, viewer);
336:   } else
337: #endif
338:   {
339:     VecLoad_Binary(newvec, viewer);
340:   }
341:   return(0);
342: }

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

349:   Input Parameters:
350: + v   - The vector
351: - tol - The zero tolerance

353:   Output Parameters:
354: . v - The chopped vector

356:   Level: intermediate

358: .seealso: VecCreate(), VecSet()
359: @*/
360: PetscErrorCode VecChop(Vec v, PetscReal tol)
361: {
362:   PetscScalar    *a;
363:   PetscInt       n, i;

367:   VecGetLocalSize(v, &n);
368:   VecGetArray(v, &a);
369:   for (i = 0; i < n; ++i) {
370:     if (PetscAbsScalar(a[i]) < tol) a[i] = 0.0;
371:   }
372:   VecRestoreArray(v, &a);
373:   return(0);
374: }