Actual source code: vecio.c

petsc-3.11.0 2019-03-29
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>
 10:  #include <petsc/private/vecimpl.h>
 11:  #include <petsc/private/viewerimpl.h>

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

 20:   PetscObjectGetComm((PetscObject)viewer,&comm);
 21:   /* Read vector header */
 22:   PetscViewerBinaryRead(viewer,tr,2,NULL,PETSC_INT);
 23:   type = tr[0];
 24:   if (type != VEC_FILE_CLASSID) {
 25:     PetscLogEventEnd(VEC_Load,viewer,0,0,0);
 26:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a vector next in file");
 27:   }
 28:   *rows = tr[1];
 29:   return(0);
 30: }

 32: #if defined(PETSC_HAVE_MPIIO)
 33: static PetscErrorCode VecLoad_Binary_MPIIO(Vec vec, PetscViewer viewer)
 34: {
 36:   PetscMPIInt    lsize;
 37:   PetscScalar    *avec;
 38:   MPI_File       mfdes;
 39:   MPI_Offset     off;

 42:   VecGetArray(vec,&avec);
 43:   PetscMPIIntCast(vec->map->n,&lsize);

 45:   PetscViewerBinaryGetMPIIODescriptor(viewer,&mfdes);
 46:   PetscViewerBinaryGetMPIIOOffset(viewer,&off);
 47:   off += vec->map->rstart*sizeof(PetscScalar);
 48:   MPI_File_set_view(mfdes,off,MPIU_SCALAR,MPIU_SCALAR,(char*)"native",MPI_INFO_NULL);
 49:   MPIU_File_read_all(mfdes,avec,lsize,MPIU_SCALAR,MPI_STATUS_IGNORE);
 50:   PetscViewerBinaryAddMPIIOOffset(viewer,vec->map->N*sizeof(PetscScalar));

 52:   VecRestoreArray(vec,&avec);
 53:   VecAssemblyBegin(vec);
 54:   VecAssemblyEnd(vec);
 55:   return(0);
 56: }
 57: #endif

 59: PetscErrorCode VecLoad_Binary(Vec vec, PetscViewer viewer)
 60: {
 61:   PetscMPIInt    size,rank,tag;
 62:   int            fd;
 63:   PetscInt       i,rows = 0,n,*range,N,bs;
 65:   PetscBool      flag,skipheader;
 66:   PetscScalar    *avec,*avecwork;
 67:   MPI_Comm       comm;
 68:   MPI_Request    request;
 69:   MPI_Status     status;
 70: #if defined(PETSC_HAVE_MPIIO)
 71:   PetscBool      useMPIIO;
 72: #endif

 75:   /* force binary viewer to load .info file if it has not yet done so */
 76:   PetscViewerSetUp(viewer);
 77:   PetscObjectGetComm((PetscObject)viewer,&comm);
 78:   MPI_Comm_rank(comm,&rank);
 79:   MPI_Comm_size(comm,&size);

 81:   PetscViewerBinaryGetDescriptor(viewer,&fd);
 82:   PetscViewerBinaryGetSkipHeader(viewer,&skipheader);
 83:   if (!skipheader) {
 84:     PetscViewerBinaryReadVecHeader_Private(viewer,&rows);
 85:   } else {
 86:     VecType vtype;
 87:     VecGetType(vec,&vtype);
 88:     if (!vtype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER, "Vector binary file header was skipped, thus the user must specify the type and length of input vector");
 89:     VecGetSize(vec,&N);
 90:     if (!N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER, "Vector binary file header was skipped, thus the user must specify the length of input vector");
 91:     rows = N;
 92:   }
 93:   /* Set Vec sizes,blocksize,and type if not already set. Block size first so that local sizes will be compatible. */
 94:   PetscOptionsGetInt(NULL,((PetscObject)vec)->prefix, "-vecload_block_size", &bs, &flag);
 95:   if (flag) {
 96:     VecSetBlockSize(vec, bs);
 97:   }
 98:   if (vec->map->n < 0 && vec->map->N < 0) {
 99:     VecSetSizes(vec,PETSC_DECIDE,rows);
100:   }

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

106: #if defined(PETSC_HAVE_MPIIO)
107:   PetscViewerBinaryGetUseMPIIO(viewer,&useMPIIO);
108:   if (useMPIIO) {
109:     VecLoad_Binary_MPIIO(vec, viewer);
110:     return(0);
111:   }
112: #endif

114:   VecGetLocalSize(vec,&n);
115:   PetscObjectGetNewTag((PetscObject)viewer,&tag);
116:   VecGetArray(vec,&avec);
117:   if (!rank) {
118:     PetscBinaryRead(fd,avec,n,PETSC_SCALAR);

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

127:       PetscMalloc1(n,&avecwork);
128:       for (i=1; i<size; i++) {
129:         n    = range[i+1] - range[i];
130:         PetscBinaryRead(fd,avecwork,n,PETSC_SCALAR);
131:         MPI_Isend(avecwork,n,MPIU_SCALAR,i,tag,comm,&request);
132:         MPI_Wait(&request,&status);
133:       }
134:       PetscFree(avecwork);
135:     }
136:   } else {
137:     MPI_Recv(avec,n,MPIU_SCALAR,0,tag,comm,&status);
138:   }

140:   VecRestoreArray(vec,&avec);
141:   VecAssemblyBegin(vec);
142:   VecAssemblyEnd(vec);
143:   return(0);
144: }

146: #if defined(PETSC_HAVE_HDF5)
147: /*
148:      This should handle properly the cases where PetscInt is 32 or 64 and hsize_t is 32 or 64. These means properly casting with
149:    checks back and forth between the two types of variables.
150: */
151: PetscErrorCode VecLoad_HDF5(Vec xin, PetscViewer viewer)
152: {
153:   hid_t          scalartype; /* scalar type (H5T_NATIVE_FLOAT or H5T_NATIVE_DOUBLE) */
154:   PetscScalar    *x;
155:   const char     *vecname;

159:   if (!((PetscObject)xin)->name) SETERRQ(PetscObjectComm((PetscObject)xin), PETSC_ERR_SUP, "Vec name must be set with PetscObjectSetName() before VecLoad()");
160: #if defined(PETSC_USE_REAL_SINGLE)
161:   scalartype = H5T_NATIVE_FLOAT;
162: #elif defined(PETSC_USE_REAL___FLOAT128)
163: #error "HDF5 output with 128 bit floats not supported."
164: #elif defined(PETSC_USE_REAL___FP16)
165: #error "HDF5 output with 16 bit floats not supported."
166: #else
167:   scalartype = H5T_NATIVE_DOUBLE;
168: #endif
169:   PetscObjectGetName((PetscObject)xin, &vecname);
170:   PetscViewerHDF5Load(viewer, vecname, xin->map, scalartype, (void**)&x);
171:   VecSetUp(xin);
172:   VecReplaceArray(xin, x);
173:   return(0);
174: }
175: #endif

177: #if defined(PETSC_HAVE_ADIOS)
178: #include <adios.h>
179: #include <adios_read.h>
180:  #include <petsc/private/vieweradiosimpl.h>
181:  #include <petsc/private/viewerimpl.h>

183: PetscErrorCode VecLoad_ADIOS(Vec xin, PetscViewer viewer)
184: {
185:   PetscViewer_ADIOS *adios = (PetscViewer_ADIOS*)viewer->data;
186:   PetscErrorCode    ierr;
187:   PetscScalar       *x;
188:   PetscInt          Nfile,N,rstart,n;
189:   uint64_t          N_t,rstart_t;
190:   const char        *vecname;
191:   ADIOS_VARINFO     *v;
192:   ADIOS_SELECTION   *sel;

195:   PetscObjectGetName((PetscObject) xin, &vecname);

197:   v    = adios_inq_var(adios->adios_fp, vecname);
198:   if (v->ndim != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Vector in file is not of dimension 1 (%D)", v->ndim);
199:   Nfile = (PetscInt) v->dims[0];

201:   /* Set Vec sizes,blocksize,and type if not already set */
202:   if ((xin)->map->n < 0 && (xin)->map->N < 0) {
203:     VecSetSizes(xin, PETSC_DECIDE, Nfile);
204:   }
205:   /* If sizes and type already set,check if the vector global size is correct */
206:   VecGetSize(xin, &N);
207:   VecGetLocalSize(xin, &n);
208:   if (N != Nfile) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Vector in file different length (%D) then input vector (%D)", Nfile, N);
209: 
210:   VecGetOwnershipRange(xin,&rstart,NULL);
211:   rstart_t = rstart;
212:   N_t  = n;
213:   sel  = adios_selection_boundingbox (v->ndim, &rstart_t, &N_t);
214:   VecGetArray(xin,&x);
215:   adios_schedule_read(adios->adios_fp, sel, vecname, 0, 1, x);
216:   adios_perform_reads (adios->adios_fp, 1);
217:   VecRestoreArray(xin,&x);
218:   adios_selection_delete(sel);

220:   return(0);
221: }
222: #endif

224: #if defined(PETSC_HAVE_ADIOS2)
225: #include <adios2_c.h>
226: #include <petsc/private/vieweradios2impl.h>
227:  #include <petsc/private/viewerimpl.h>

229: PetscErrorCode VecLoad_ADIOS2(Vec xin, PetscViewer viewer)
230: {
231:   PetscViewer_ADIOS2 *adios2 = (PetscViewer_ADIOS2*)viewer->data;
232:   PetscErrorCode     ierr;
233:   PetscScalar        *x;
234:   PetscInt           Nfile,N,rstart,n;
235:   const char         *vecname;

238:   PetscObjectGetName((PetscObject) xin, &vecname);

240:   /* Set Vec sizes,blocksize,and type if not already set */
241:   if ((xin)->map->n < 0 && (xin)->map->N < 0) {
242:     VecSetSizes(xin, PETSC_DECIDE, Nfile);
243:   }
244:   /* If sizes and type already set,check if the vector global size is correct */
245:   VecGetSize(xin, &N);
246:   VecGetLocalSize(xin, &n);
247:   if (N != Nfile) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Vector in file different length (%D) then input vector (%D)", Nfile, N);

249:   VecGetOwnershipRange(xin,&rstart,NULL);
250:   VecGetArray(xin,&x);
251:   VecRestoreArray(xin,&x);
252:   return(0);
253: }
254: #endif

256: PetscErrorCode  VecLoad_Default(Vec newvec, PetscViewer viewer)
257: {
259:   PetscBool      isbinary;
260: #if defined(PETSC_HAVE_HDF5)
261:   PetscBool      ishdf5;
262: #endif
263: #if defined(PETSC_HAVE_ADIOS)
264:   PetscBool      isadios;
265: #endif
266: #if defined(PETSC_HAVE_ADIOS2)
267:   PetscBool      isadios2;
268: #endif

271:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
272: #if defined(PETSC_HAVE_HDF5)
273:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
274: #endif
275: #if defined(PETSC_HAVE_ADIOS)
276:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERADIOS,&isadios);
277: #endif
278: #if defined(PETSC_HAVE_ADIOS2)
279:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERADIOS2,&isadios2);
280: #endif

282: #if defined(PETSC_HAVE_HDF5)
283:   if (ishdf5) {
284:     if (!((PetscObject)newvec)->name) {
285:       PetscLogEventEnd(VEC_Load,viewer,0,0,0);
286:       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()");
287:     }
288:     VecLoad_HDF5(newvec, viewer);
289:   } else
290: #endif
291: #if defined(PETSC_HAVE_ADIOS)
292:   if (isadios) {
293:     if (!((PetscObject)newvec)->name) {
294:       PetscLogEventEnd(VEC_Load,viewer,0,0,0);
295:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Since ADIOS format gives ASCII name for each object in file; must use VecLoad() after setting name of Vec with PetscObjectSetName()");
296:     }
297:     VecLoad_ADIOS(newvec, viewer);
298:   } else
299: #endif
300: #if defined(PETSC_HAVE_ADIOS2)
301:   if (isadios2) {
302:     if (!((PetscObject)newvec)->name) {
303:       PetscLogEventEnd(VEC_Load,viewer,0,0,0);
304:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Since ADIOS2 format gives ASCII name for each object in file; must use VecLoad() after setting name of Vec with PetscObjectSetName()");
305:     }
306:     VecLoad_ADIOS2(newvec, viewer);
307:   } else
308: #endif
309:   {
310:     VecLoad_Binary(newvec, viewer);
311:   }
312:   return(0);
313: }

315: /*@
316:   VecChop - Set all values in the vector with an absolute value less than the tolerance to zero

318:   Input Parameters:
319: + v   - The vector
320: - tol - The zero tolerance

322:   Output Parameters:
323: . v - The chopped vector

325:   Level: intermediate

327: .seealso: VecCreate(), VecSet()
328: @*/
329: PetscErrorCode VecChop(Vec v, PetscReal tol)
330: {
331:   PetscScalar    *a;
332:   PetscInt       n, i;

336:   VecGetLocalSize(v, &n);
337:   VecGetArray(v, &a);
338:   for (i = 0; i < n; ++i) {
339:     if (PetscAbsScalar(a[i]) < tol) a[i] = 0.0;
340:   }
341:   VecRestoreArray(v, &a);
342:   return(0);
343: }