Actual source code: bvec2.c

  1: /*
  2:    Implements the sequential vectors.
  3: */

  5: #include <../src/vec/vec/impls/dvecimpl.h>
  6: #include <../src/vec/vec/impls/mpi/pvecimpl.h>
  7: #include <petsc/private/glvisviewerimpl.h>
  8: #include <petsc/private/glvisvecimpl.h>
  9: #include <petscblaslapack.h>

 11: static PetscErrorCode VecPointwiseApply_Seq(Vec win, Vec xin, Vec yin, PetscScalar (*const func)(PetscScalar, PetscScalar))
 12: {
 13:   const PetscInt n = win->map->n;
 14:   PetscScalar   *ww, *xx, *yy; /* cannot make xx or yy const since might be ww */

 16:   PetscFunctionBegin;
 17:   PetscCall(VecGetArrayRead(xin, (const PetscScalar **)&xx));
 18:   PetscCall(VecGetArrayRead(yin, (const PetscScalar **)&yy));
 19:   PetscCall(VecGetArray(win, &ww));
 20:   for (PetscInt i = 0; i < n; ++i) ww[i] = func(xx[i], yy[i]);
 21:   PetscCall(VecRestoreArrayRead(xin, (const PetscScalar **)&xx));
 22:   PetscCall(VecRestoreArrayRead(yin, (const PetscScalar **)&yy));
 23:   PetscCall(VecRestoreArray(win, &ww));
 24:   PetscCall(PetscLogFlops(n));
 25:   PetscFunctionReturn(PETSC_SUCCESS);
 26: }

 28: static PetscScalar MaxRealPart(PetscScalar x, PetscScalar y)
 29: {
 30:   // use temporaries to avoid reevaluating side-effects
 31:   const PetscReal rx = PetscRealPart(x), ry = PetscRealPart(y);

 33:   return PetscMax(rx, ry);
 34: }

 36: PetscErrorCode VecPointwiseMax_Seq(Vec win, Vec xin, Vec yin)
 37: {
 38:   PetscFunctionBegin;
 39:   PetscCall(VecPointwiseApply_Seq(win, xin, yin, MaxRealPart));
 40:   PetscFunctionReturn(PETSC_SUCCESS);
 41: }

 43: static PetscScalar MinRealPart(PetscScalar x, PetscScalar y)
 44: {
 45:   // use temporaries to avoid reevaluating side-effects
 46:   const PetscReal rx = PetscRealPart(x), ry = PetscRealPart(y);

 48:   return PetscMin(rx, ry);
 49: }

 51: PetscErrorCode VecPointwiseMin_Seq(Vec win, Vec xin, Vec yin)
 52: {
 53:   PetscFunctionBegin;
 54:   PetscCall(VecPointwiseApply_Seq(win, xin, yin, MinRealPart));
 55:   PetscFunctionReturn(PETSC_SUCCESS);
 56: }

 58: static PetscScalar MaxAbs(PetscScalar x, PetscScalar y)
 59: {
 60:   return PetscMax(PetscAbsScalar(x), PetscAbsScalar(y));
 61: }

 63: PetscErrorCode VecPointwiseMaxAbs_Seq(Vec win, Vec xin, Vec yin)
 64: {
 65:   PetscFunctionBegin;
 66:   PetscCall(VecPointwiseApply_Seq(win, xin, yin, MaxAbs));
 67:   PetscFunctionReturn(PETSC_SUCCESS);
 68: }

 70: #include <../src/vec/vec/impls/seq/ftn-kernels/fxtimesy.h>

 72: PetscErrorCode VecPointwiseMult_Seq(Vec win, Vec xin, Vec yin)
 73: {
 74:   PetscInt     n = win->map->n, i;
 75:   PetscScalar *ww, *xx, *yy; /* cannot make xx or yy const since might be ww */

 77:   PetscFunctionBegin;
 78:   PetscCall(VecGetArrayRead(xin, (const PetscScalar **)&xx));
 79:   PetscCall(VecGetArrayRead(yin, (const PetscScalar **)&yy));
 80:   PetscCall(VecGetArray(win, &ww));
 81:   if (ww == xx) {
 82:     for (i = 0; i < n; i++) ww[i] *= yy[i];
 83:   } else if (ww == yy) {
 84:     for (i = 0; i < n; i++) ww[i] *= xx[i];
 85:   } else {
 86: #if defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
 87:     fortranxtimesy_(xx, yy, ww, &n);
 88: #else
 89:     for (i = 0; i < n; i++) ww[i] = xx[i] * yy[i];
 90: #endif
 91:   }
 92:   PetscCall(VecRestoreArrayRead(xin, (const PetscScalar **)&xx));
 93:   PetscCall(VecRestoreArrayRead(yin, (const PetscScalar **)&yy));
 94:   PetscCall(VecRestoreArray(win, &ww));
 95:   PetscCall(PetscLogFlops(n));
 96:   PetscFunctionReturn(PETSC_SUCCESS);
 97: }

 99: static PetscScalar ScalDiv(PetscScalar x, PetscScalar y)
100: {
101:   return y == 0.0 ? 0.0 : x / y;
102: }

104: PetscErrorCode VecPointwiseDivide_Seq(Vec win, Vec xin, Vec yin)
105: {
106:   PetscFunctionBegin;
107:   PetscCall(VecPointwiseApply_Seq(win, xin, yin, ScalDiv));
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: PetscErrorCode VecSetRandom_Seq(Vec xin, PetscRandom r)
112: {
113:   PetscScalar *xx;

115:   PetscFunctionBegin;
116:   PetscCall(VecGetArrayWrite(xin, &xx));
117:   PetscCall(PetscRandomGetValues(r, xin->map->n, xx));
118:   PetscCall(VecRestoreArrayWrite(xin, &xx));
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

122: PetscErrorCode VecGetSize_Seq(Vec vin, PetscInt *size)
123: {
124:   PetscFunctionBegin;
125:   *size = vin->map->n;
126:   PetscFunctionReturn(PETSC_SUCCESS);
127: }

129: PetscErrorCode VecConjugate_Seq(Vec xin)
130: {
131:   PetscFunctionBegin;
132:   if (PetscDefined(USE_COMPLEX)) {
133:     const PetscInt n = xin->map->n;
134:     PetscScalar   *x;

136:     PetscCall(VecGetArray(xin, &x));
137:     for (PetscInt i = 0; i < n; ++i) x[i] = PetscConj(x[i]);
138:     PetscCall(VecRestoreArray(xin, &x));
139:   }
140:   PetscFunctionReturn(PETSC_SUCCESS);
141: }

143: PetscErrorCode VecResetArray_Seq(Vec vin)
144: {
145:   Vec_Seq *v = (Vec_Seq *)vin->data;

147:   PetscFunctionBegin;
148:   v->array         = v->unplacedarray;
149:   v->unplacedarray = NULL;
150:   PetscFunctionReturn(PETSC_SUCCESS);
151: }

153: PetscErrorCode VecCopy_Seq(Vec xin, Vec yin)
154: {
155:   PetscFunctionBegin;
156:   if (xin != yin) {
157:     const PetscScalar *xa;
158:     PetscScalar       *ya;

160:     PetscCall(VecGetArrayRead(xin, &xa));
161:     PetscCall(VecGetArrayWrite(yin, &ya));
162:     PetscCall(PetscArraycpy(ya, xa, xin->map->n));
163:     PetscCall(VecRestoreArrayRead(xin, &xa));
164:     PetscCall(VecRestoreArrayWrite(yin, &ya));
165:   }
166:   PetscFunctionReturn(PETSC_SUCCESS);
167: }

169: PetscErrorCode VecSwap_Seq(Vec xin, Vec yin)
170: {
171:   PetscFunctionBegin;
172:   if (xin != yin) {
173:     const PetscBLASInt one = 1;
174:     PetscScalar       *ya, *xa;
175:     PetscBLASInt       bn;

177:     PetscCall(PetscBLASIntCast(xin->map->n, &bn));
178:     PetscCall(VecGetArray(xin, &xa));
179:     PetscCall(VecGetArray(yin, &ya));
180:     PetscCallBLAS("BLASswap", BLASswap_(&bn, xa, &one, ya, &one));
181:     PetscCall(VecRestoreArray(xin, &xa));
182:     PetscCall(VecRestoreArray(yin, &ya));
183:   }
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: PetscErrorCode VecNorm_Seq(Vec xin, NormType type, PetscReal *z)
188: {
189:   // use a local variable to ensure compiler doesn't think z aliases any of the other arrays
190:   PetscReal      ztmp[] = {0.0, 0.0};
191:   const PetscInt n      = xin->map->n;

193:   PetscFunctionBegin;
194:   if (n) {
195:     const PetscScalar *xx;
196:     const PetscBLASInt one = 1;
197:     PetscBLASInt       bn  = 0;

199:     PetscCall(PetscBLASIntCast(n, &bn));
200:     PetscCall(VecGetArrayRead(xin, &xx));
201:     if (type == NORM_2 || type == NORM_FROBENIUS) {
202:     NORM_1_AND_2_DOING_NORM_2:
203:       if (PetscDefined(USE_REAL___FP16)) {
204:         PetscCallBLAS("BLASnrm2", ztmp[type == NORM_1_AND_2] = BLASnrm2_(&bn, xx, &one));
205:       } else {
206:         PetscCallBLAS("BLASdot", ztmp[type == NORM_1_AND_2] = PetscSqrtReal(PetscRealPart(BLASdot_(&bn, xx, &one, xx, &one))));
207:       }
208:       PetscCall(PetscLogFlops(2.0 * n - 1));
209:     } else if (type == NORM_INFINITY) {
210:       for (PetscInt i = 0; i < n; ++i) {
211:         const PetscReal tmp = PetscAbsScalar(xx[i]);

213:         /* check special case of tmp == NaN */
214:         if ((tmp > ztmp[0]) || (tmp != tmp)) {
215:           ztmp[0] = tmp;
216:           if (tmp != tmp) break;
217:         }
218:       }
219:     } else if (type == NORM_1 || type == NORM_1_AND_2) {
220:       if (PetscDefined(USE_COMPLEX)) {
221:         // BLASasum() returns the nonstandard 1 norm of the 1 norm of the complex entries so we
222:         // provide a custom loop instead
223:         for (PetscInt i = 0; i < n; ++i) ztmp[0] += PetscAbsScalar(xx[i]);
224:       } else {
225:         PetscCallBLAS("BLASasum", ztmp[0] = BLASasum_(&bn, xx, &one));
226:       }
227:       PetscCall(PetscLogFlops(n - 1.0));
228:       /* slight reshuffle so we can skip getting the array again (but still log the flops) if we
229:          do norm2 after this */
230:       if (type == NORM_1_AND_2) goto NORM_1_AND_2_DOING_NORM_2;
231:     }
232:     PetscCall(VecRestoreArrayRead(xin, &xx));
233:   }
234:   z[0] = ztmp[0];
235:   if (type == NORM_1_AND_2) z[1] = ztmp[1];
236:   PetscFunctionReturn(PETSC_SUCCESS);
237: }

239: static PetscErrorCode VecView_Seq_ASCII(Vec xin, PetscViewer viewer)
240: {
241:   PetscInt           i, n = xin->map->n;
242:   const char        *name;
243:   PetscViewerFormat  format;
244:   const PetscScalar *xv;

246:   PetscFunctionBegin;
247:   PetscCall(VecGetArrayRead(xin, &xv));
248:   PetscCall(PetscViewerGetFormat(viewer, &format));
249:   if (format == PETSC_VIEWER_ASCII_MATLAB) {
250:     PetscCall(PetscObjectGetName((PetscObject)xin, &name));
251:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s = [\n", name));
252:     for (i = 0; i < n; i++) {
253: #if defined(PETSC_USE_COMPLEX)
254:       if (PetscImaginaryPart(xv[i]) > 0.0) {
255:         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e + %18.16ei\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
256:       } else if (PetscImaginaryPart(xv[i]) < 0.0) {
257:         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e - %18.16ei\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i])));
258:       } else {
259:         PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e\n", (double)PetscRealPart(xv[i])));
260:       }
261: #else
262:       PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e\n", (double)xv[i]));
263: #endif
264:     }
265:     PetscCall(PetscViewerASCIIPrintf(viewer, "];\n"));
266:   } else if (format == PETSC_VIEWER_ASCII_SYMMODU) {
267:     for (i = 0; i < n; i++) {
268: #if defined(PETSC_USE_COMPLEX)
269:       PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e %18.16e\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
270: #else
271:       PetscCall(PetscViewerASCIIPrintf(viewer, "%18.16e\n", (double)xv[i]));
272: #endif
273:     }
274:   } else if (format == PETSC_VIEWER_ASCII_PCICE) {
275:     PetscInt bs, b;

277:     PetscCall(VecGetBlockSize(xin, &bs));
278:     PetscCheck(bs >= 1 && bs <= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCICE can only handle up to 3D objects, but vector dimension is %" PetscInt_FMT, bs);
279:     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", xin->map->N / bs));
280:     for (i = 0; i < n / bs; i++) {
281:       PetscCall(PetscViewerASCIIPrintf(viewer, "%7" PetscInt_FMT "   ", i + 1));
282:       for (b = 0; b < bs; b++) {
283:         if (b > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " "));
284: #if !defined(PETSC_USE_COMPLEX)
285:         PetscCall(PetscViewerASCIIPrintf(viewer, "% 12.5E", (double)xv[i * bs + b]));
286: #endif
287:       }
288:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
289:     }
290:   } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
291:     /* GLVis ASCII visualization/dump: this function mimics mfem::GridFunction::Save() */
292:     const PetscScalar      *array;
293:     PetscInt                i, n, vdim, ordering = 1; /* mfem::FiniteElementSpace::Ordering::byVDIM */
294:     PetscContainer          glvis_container;
295:     PetscViewerGLVisVecInfo glvis_vec_info;
296:     PetscViewerGLVisInfo    glvis_info;

298:     /* mfem::FiniteElementSpace::Save() */
299:     PetscCall(VecGetBlockSize(xin, &vdim));
300:     PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
301:     PetscCall(PetscObjectQuery((PetscObject)xin, "_glvis_info_container", (PetscObject *)&glvis_container));
302:     PetscCheck(glvis_container, PetscObjectComm((PetscObject)xin), PETSC_ERR_PLIB, "Missing GLVis container");
303:     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_vec_info));
304:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", glvis_vec_info->fec_type));
305:     PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", vdim));
306:     PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: %" PetscInt_FMT "\n", ordering));
307:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
308:     /* mfem::Vector::Print() */
309:     PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
310:     PetscCheck(glvis_container, PetscObjectComm((PetscObject)viewer), PETSC_ERR_PLIB, "Missing GLVis container");
311:     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info));
312:     if (glvis_info->enabled) {
313:       PetscCall(VecGetLocalSize(xin, &n));
314:       PetscCall(VecGetArrayRead(xin, &array));
315:       for (i = 0; i < n; i++) {
316:         PetscCall(PetscViewerASCIIPrintf(viewer, glvis_info->fmt, (double)PetscRealPart(array[i])));
317:         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
318:       }
319:       PetscCall(VecRestoreArrayRead(xin, &array));
320:     }
321:   } else if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
322:     /* No info */
323:   } else {
324:     for (i = 0; i < n; i++) {
325:       if (format == PETSC_VIEWER_ASCII_INDEX) PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ": ", i));
326: #if defined(PETSC_USE_COMPLEX)
327:       if (PetscImaginaryPart(xv[i]) > 0.0) {
328:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g + %g i\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
329:       } else if (PetscImaginaryPart(xv[i]) < 0.0) {
330:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g - %g i\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i])));
331:       } else {
332:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g\n", (double)PetscRealPart(xv[i])));
333:       }
334: #else
335:       PetscCall(PetscViewerASCIIPrintf(viewer, "%g\n", (double)xv[i]));
336: #endif
337:     }
338:   }
339:   PetscCall(PetscViewerFlush(viewer));
340:   PetscCall(VecRestoreArrayRead(xin, &xv));
341:   PetscFunctionReturn(PETSC_SUCCESS);
342: }

344: #include <petscdraw.h>
345: static PetscErrorCode VecView_Seq_Draw_LG(Vec xin, PetscViewer v)
346: {
347:   PetscDraw          draw;
348:   PetscBool          isnull;
349:   PetscDrawLG        lg;
350:   PetscInt           i, c, bs = xin->map->bs, n = xin->map->n / bs;
351:   const PetscScalar *xv;
352:   PetscReal         *xx, *yy, xmin, xmax, h;
353:   int                colors[] = {PETSC_DRAW_RED};
354:   PetscViewerFormat  format;
355:   PetscDrawAxis      axis;
356:   const char        *name;

358:   PetscFunctionBegin;
359:   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
360:   PetscCall(PetscDrawIsNull(draw, &isnull));
361:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

363:   PetscCall(PetscObjectGetName((PetscObject)xin, &name));
364:   PetscCall(PetscDrawSetTitle(draw, name));
365:   PetscCall(PetscViewerGetFormat(v, &format));
366:   PetscCall(PetscMalloc2(n, &xx, n, &yy));
367:   PetscCall(VecGetArrayRead(xin, &xv));
368:   for (c = 0; c < bs; c++) {
369:     PetscCall(PetscViewerDrawGetDrawLG(v, c, &lg));
370:     PetscCall(PetscDrawLGReset(lg));
371:     PetscCall(PetscDrawLGSetDimension(lg, 1));
372:     PetscCall(PetscDrawLGSetColors(lg, colors));
373:     if (format == PETSC_VIEWER_DRAW_LG_XRANGE) {
374:       PetscCall(PetscDrawLGGetAxis(lg, &axis));
375:       PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, NULL, NULL));
376:       h = (xmax - xmin) / n;
377:       for (i = 0; i < n; i++) xx[i] = i * h + 0.5 * h; /* cell center */
378:     } else {
379:       for (i = 0; i < n; i++) xx[i] = (PetscReal)i;
380:     }
381:     for (i = 0; i < n; i++) yy[i] = PetscRealPart(xv[c + i * bs]);

383:     PetscCall(PetscDrawLGAddPoints(lg, n, &xx, &yy));
384:     PetscCall(PetscDrawLGDraw(lg));
385:     PetscCall(PetscDrawLGSave(lg));
386:   }
387:   PetscCall(VecRestoreArrayRead(xin, &xv));
388:   PetscCall(PetscFree2(xx, yy));
389:   PetscFunctionReturn(PETSC_SUCCESS);
390: }

392: static PetscErrorCode VecView_Seq_Draw(Vec xin, PetscViewer v)
393: {
394:   PetscDraw draw;
395:   PetscBool isnull;

397:   PetscFunctionBegin;
398:   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
399:   PetscCall(PetscDrawIsNull(draw, &isnull));
400:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

402:   PetscCall(VecView_Seq_Draw_LG(xin, v));
403:   PetscFunctionReturn(PETSC_SUCCESS);
404: }

406: static PetscErrorCode VecView_Seq_Binary(Vec xin, PetscViewer viewer)
407: {
408:   return VecView_Binary(xin, viewer);
409: }

411: #if defined(PETSC_HAVE_MATLAB)
412: #include <petscmatlab.h>
413:   #include <mat.h> /* MATLAB include file */
414: PetscErrorCode VecView_Seq_Matlab(Vec vec, PetscViewer viewer)
415: {
416:   PetscInt           n;
417:   const PetscScalar *array;

419:   PetscFunctionBegin;
420:   PetscCall(VecGetLocalSize(vec, &n));
421:   PetscCall(PetscObjectName((PetscObject)vec));
422:   PetscCall(VecGetArrayRead(vec, &array));
423:   PetscCall(PetscViewerMatlabPutArray(viewer, n, 1, array, ((PetscObject)vec)->name));
424:   PetscCall(VecRestoreArrayRead(vec, &array));
425:   PetscFunctionReturn(PETSC_SUCCESS);
426: }
427: #endif

429: PetscErrorCode VecView_Seq(Vec xin, PetscViewer viewer)
430: {
431:   PetscBool isdraw, iascii, issocket, isbinary;
432: #if defined(PETSC_HAVE_MATHEMATICA)
433:   PetscBool ismathematica;
434: #endif
435: #if defined(PETSC_HAVE_MATLAB)
436:   PetscBool ismatlab;
437: #endif
438: #if defined(PETSC_HAVE_HDF5)
439:   PetscBool ishdf5;
440: #endif
441:   PetscBool isglvis;
442: #if defined(PETSC_HAVE_ADIOS)
443:   PetscBool isadios;
444: #endif

446:   PetscFunctionBegin;
447:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
448:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
449:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
450:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
451: #if defined(PETSC_HAVE_MATHEMATICA)
452:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATHEMATICA, &ismathematica));
453: #endif
454: #if defined(PETSC_HAVE_HDF5)
455:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
456: #endif
457: #if defined(PETSC_HAVE_MATLAB)
458:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATLAB, &ismatlab));
459: #endif
460:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
461: #if defined(PETSC_HAVE_ADIOS)
462:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERADIOS, &isadios));
463: #endif

465:   if (isdraw) {
466:     PetscCall(VecView_Seq_Draw(xin, viewer));
467:   } else if (iascii) {
468:     PetscCall(VecView_Seq_ASCII(xin, viewer));
469:   } else if (isbinary) {
470:     PetscCall(VecView_Seq_Binary(xin, viewer));
471: #if defined(PETSC_HAVE_MATHEMATICA)
472:   } else if (ismathematica) {
473:     PetscCall(PetscViewerMathematicaPutVector(viewer, xin));
474: #endif
475: #if defined(PETSC_HAVE_HDF5)
476:   } else if (ishdf5) {
477:     PetscCall(VecView_MPI_HDF5(xin, viewer)); /* Reusing VecView_MPI_HDF5 ... don't want code duplication*/
478: #endif
479: #if defined(PETSC_HAVE_ADIOS)
480:   } else if (isadios) {
481:     PetscCall(VecView_MPI_ADIOS(xin, viewer)); /* Reusing VecView_MPI_ADIOS ... don't want code duplication*/
482: #endif
483: #if defined(PETSC_HAVE_MATLAB)
484:   } else if (ismatlab) {
485:     PetscCall(VecView_Seq_Matlab(xin, viewer));
486: #endif
487:   } else if (isglvis) PetscCall(VecView_GLVis(xin, viewer));
488:   PetscFunctionReturn(PETSC_SUCCESS);
489: }

491: PetscErrorCode VecGetValues_Seq(Vec xin, PetscInt ni, const PetscInt ix[], PetscScalar y[])
492: {
493:   const PetscBool    ignorenegidx = xin->stash.ignorenegidx;
494:   const PetscScalar *xx;

496:   PetscFunctionBegin;
497:   PetscCall(VecGetArrayRead(xin, &xx));
498:   for (PetscInt i = 0; i < ni; ++i) {
499:     if (ignorenegidx && (ix[i] < 0)) continue;
500:     if (PetscDefined(USE_DEBUG)) {
501:       PetscCheck(ix[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " cannot be negative", ix[i]);
502:       PetscCheck(ix[i] < xin->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT ", should be less than %" PetscInt_FMT, ix[i], xin->map->n);
503:     }
504:     y[i] = xx[ix[i]];
505:   }
506:   PetscCall(VecRestoreArrayRead(xin, &xx));
507:   PetscFunctionReturn(PETSC_SUCCESS);
508: }

510: PetscErrorCode VecSetValues_Seq(Vec xin, PetscInt ni, const PetscInt ix[], const PetscScalar y[], InsertMode m)
511: {
512:   const PetscBool ignorenegidx = xin->stash.ignorenegidx;
513:   PetscScalar    *xx;

515:   PetscFunctionBegin;
516:   // call to getarray (not e.g. getarraywrite() if m is INSERT_VALUES) is deliberate! If this
517:   // is secretly a VECSEQCUDA it may have values currently on the device, in which case --
518:   // unless we are replacing the entire array -- we need to copy them up
519:   PetscCall(VecGetArray(xin, &xx));
520:   for (PetscInt i = 0; i < ni; i++) {
521:     if (ignorenegidx && (ix[i] < 0)) continue;
522:     if (PetscDefined(USE_DEBUG)) {
523:       PetscCheck(ix[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " cannot be negative", ix[i]);
524:       PetscCheck(ix[i] < xin->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT ", should be less than %" PetscInt_FMT, ix[i], xin->map->n);
525:     }
526:     if (m == INSERT_VALUES) {
527:       xx[ix[i]] = y[i];
528:     } else {
529:       xx[ix[i]] += y[i];
530:     }
531:   }
532:   PetscCall(VecRestoreArray(xin, &xx));
533:   PetscFunctionReturn(PETSC_SUCCESS);
534: }

536: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin, PetscInt ni, const PetscInt ix[], const PetscScalar yin[], InsertMode m)
537: {
538:   PetscScalar *xx;
539:   PetscInt     bs;

541:   /* For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling */
542:   PetscFunctionBegin;
543:   PetscCall(VecGetBlockSize(xin, &bs));
544:   PetscCall(VecGetArray(xin, &xx));
545:   for (PetscInt i = 0; i < ni; ++i, yin += bs) {
546:     const PetscInt start = bs * ix[i];

548:     if (start < 0) continue;
549:     PetscCheck(start < xin->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT ", should be less than %" PetscInt_FMT, start, xin->map->n);
550:     for (PetscInt j = 0; j < bs; j++) {
551:       if (m == INSERT_VALUES) {
552:         xx[start + j] = yin[j];
553:       } else {
554:         xx[start + j] += yin[j];
555:       }
556:     }
557:   }
558:   PetscCall(VecRestoreArray(xin, &xx));
559:   PetscFunctionReturn(PETSC_SUCCESS);
560: }

562: static PetscErrorCode VecResetPreallocationCOO_Seq(Vec x)
563: {
564:   Vec_Seq *vs = (Vec_Seq *)x->data;

566:   PetscFunctionBegin;
567:   if (vs) {
568:     PetscCall(PetscFree(vs->jmap1)); /* Destroy old stuff */
569:     PetscCall(PetscFree(vs->perm1));
570:   }
571:   PetscFunctionReturn(PETSC_SUCCESS);
572: }

574: PetscErrorCode VecSetPreallocationCOO_Seq(Vec x, PetscCount coo_n, const PetscInt coo_i[])
575: {
576:   PetscInt    m, *i;
577:   PetscCount  k, nneg;
578:   PetscCount *perm1, *jmap1;
579:   Vec_Seq    *vs = (Vec_Seq *)x->data;

581:   PetscFunctionBegin;
582:   PetscCall(VecResetPreallocationCOO_Seq(x)); /* Destroy old stuff */
583:   PetscCall(PetscMalloc1(coo_n, &i));
584:   PetscCall(PetscArraycpy(i, coo_i, coo_n)); /* Make a copy since we'll modify it */
585:   PetscCall(PetscMalloc1(coo_n, &perm1));
586:   for (k = 0; k < coo_n; k++) perm1[k] = k;
587:   PetscCall(PetscSortIntWithCountArray(coo_n, i, perm1));
588:   for (k = 0; k < coo_n; k++) {
589:     if (i[k] >= 0) break;
590:   } /* Advance k to the first entry with a non-negative index */
591:   nneg = k;

593:   PetscCall(VecGetLocalSize(x, &m));
594:   PetscCheck(!nneg || x->stash.ignorenegidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Found a negative index in VecSetPreallocateCOO() but VEC_IGNORE_NEGATIVE_INDICES was not set");
595:   PetscCheck(!coo_n || i[coo_n - 1] < m, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Found index (%" PetscInt_FMT ") greater than the size of the vector (%" PetscInt_FMT ") in VecSetPreallocateCOO()", i[coo_n - 1], m);

597:   PetscCall(PetscCalloc1(m + 1, &jmap1));
598:   for (; k < coo_n; k++) jmap1[i[k] + 1]++;         /* Count repeats of each entry */
599:   for (k = 0; k < m; k++) jmap1[k + 1] += jmap1[k]; /* Transform jmap[] to CSR-like data structure */
600:   PetscCall(PetscFree(i));

602:   if (nneg) { /* Discard leading negative indices */
603:     PetscCount *perm1_new;
604:     PetscCall(PetscMalloc1(coo_n - nneg, &perm1_new));
605:     PetscCall(PetscArraycpy(perm1_new, perm1 + nneg, coo_n - nneg));
606:     PetscCall(PetscFree(perm1));
607:     perm1 = perm1_new;
608:   }

610:   /* Record COO fields */
611:   vs->coo_n = coo_n;
612:   vs->tot1  = coo_n - nneg;
613:   vs->jmap1 = jmap1; /* [m+1] */
614:   vs->perm1 = perm1; /* [tot] */
615:   PetscFunctionReturn(PETSC_SUCCESS);
616: }

618: PetscErrorCode VecSetValuesCOO_Seq(Vec x, const PetscScalar coo_v[], InsertMode imode)
619: {
620:   Vec_Seq          *vs    = (Vec_Seq *)x->data;
621:   const PetscCount *perm1 = vs->perm1, *jmap1 = vs->jmap1;
622:   PetscScalar      *xv;
623:   PetscInt          m;

625:   PetscFunctionBegin;
626:   PetscCall(VecGetLocalSize(x, &m));
627:   PetscCall(VecGetArray(x, &xv));
628:   for (PetscInt i = 0; i < m; i++) {
629:     PetscScalar sum = 0.0;
630:     for (PetscCount j = jmap1[i]; j < jmap1[i + 1]; j++) sum += coo_v[perm1[j]];
631:     xv[i] = (imode == INSERT_VALUES ? 0.0 : xv[i]) + sum;
632:   }
633:   PetscCall(VecRestoreArray(x, &xv));
634:   PetscFunctionReturn(PETSC_SUCCESS);
635: }

637: PetscErrorCode VecDestroy_Seq(Vec v)
638: {
639:   Vec_Seq *vs = (Vec_Seq *)v->data;

641:   PetscFunctionBegin;
642:   PetscCall(PetscLogObjectState((PetscObject)v, "Length=%" PetscInt_FMT, v->map->n));
643:   if (vs) PetscCall(PetscShmgetDeallocateArray((void **)&vs->array_allocated));
644:   PetscCall(VecResetPreallocationCOO_Seq(v));
645:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEnginePut_C", NULL));
646:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEngineGet_C", NULL));
647:   PetscCall(PetscFree(v->data));
648:   PetscFunctionReturn(PETSC_SUCCESS);
649: }

651: PetscErrorCode VecSetOption_Seq(Vec v, VecOption op, PetscBool flag)
652: {
653:   PetscFunctionBegin;
654:   if (op == VEC_IGNORE_NEGATIVE_INDICES) v->stash.ignorenegidx = flag;
655:   PetscFunctionReturn(PETSC_SUCCESS);
656: }

658: // duplicate w to v. v is half-baked, potentially already with arrays allocated.
659: static PetscErrorCode VecDuplicate_Seq_Private(Vec w, Vec v)
660: {
661:   PetscFunctionBegin;
662:   PetscCall(VecSetType(v, ((PetscObject)w)->type_name));
663:   PetscCall(PetscObjectListDuplicate(((PetscObject)w)->olist, &((PetscObject)v)->olist));
664:   PetscCall(PetscFunctionListDuplicate(((PetscObject)w)->qlist, &((PetscObject)v)->qlist));

666:   // Vec ops are not necessarily fully set by VecSetType(), e.g., see DMCreateGlobalVector_DA, so we copy w's to v
667:   v->ops[0] = w->ops[0];
668: #if defined(PETSC_HAVE_DEVICE)
669:   v->boundtocpu        = w->boundtocpu;
670:   v->bindingpropagates = w->bindingpropagates;
671: #endif
672:   PetscFunctionReturn(PETSC_SUCCESS);
673: }

675: PetscErrorCode VecDuplicate_Seq(Vec win, Vec *V)
676: {
677:   PetscFunctionBegin;
678:   PetscCall(VecCreateWithLayout_Private(win->map, V));
679:   PetscCall(VecDuplicate_Seq_Private(win, *V));
680:   PetscFunctionReturn(PETSC_SUCCESS);
681: }

683: PetscErrorCode VecReplaceArray_Default_GEMV_Error(Vec v, const PetscScalar *a)
684: {
685:   PetscFunctionBegin;
686:   PetscCheck(PETSC_FALSE, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "VecReplaceArray() is not supported on the first Vec obtained from VecDuplicateVecs(). \
687: You could either 1) use -vec_mdot_use_gemv 0 -vec_maxpy_use_gemv 0 to turn off an optimization to allow your current code to work or 2) use VecDuplicate() to duplicate the vector.");
688:   (void)a;
689:   PetscFunctionReturn(PETSC_SUCCESS);
690: }

692: static PetscErrorCode VecDuplicateVecs_Seq_GEMV(Vec w, PetscInt m, Vec *V[])
693: {
694:   PetscScalar *array;
695:   PetscInt64   lda; // use 64-bit as we will do "m * lda"

697:   PetscFunctionBegin;
698:   PetscCall(PetscMalloc1(m, V));
699:   VecGetLocalSizeAligned(w, 64, &lda); // get in lda the 64-bytes aligned local size
700:   PetscCall(PetscCalloc1(m * lda, &array));
701:   for (PetscInt i = 0; i < m; i++) {
702:     Vec v;
703:     PetscCall(VecCreateSeqWithLayoutAndArray_Private(w->map, PetscSafePointerPlusOffset(array, i * lda), &v));
704:     PetscCall(VecDuplicate_Seq_Private(w, v));
705:     (*V)[i] = v;
706:   }
707:   // so when the first vector is destroyed it will destroy the array
708:   if (m) ((Vec_Seq *)(*V)[0]->data)->array_allocated = array;
709:   // disable replacearray of the first vector, as freeing its memory also frees others in the group.
710:   // But replacearray of others is ok, as they don't own their array.
711:   if (m > 1) (*V)[0]->ops->replacearray = VecReplaceArray_Default_GEMV_Error;
712:   PetscFunctionReturn(PETSC_SUCCESS);
713: }

715: static struct _VecOps DvOps = {
716:   PetscDesignatedInitializer(duplicate, VecDuplicate_Seq), /* 1 */
717:   PetscDesignatedInitializer(duplicatevecs, VecDuplicateVecs_Default),
718:   PetscDesignatedInitializer(destroyvecs, VecDestroyVecs_Default),
719:   PetscDesignatedInitializer(dot, VecDot_Seq),
720:   PetscDesignatedInitializer(mdot, VecMDot_Seq),
721:   PetscDesignatedInitializer(norm, VecNorm_Seq),
722:   PetscDesignatedInitializer(tdot, VecTDot_Seq),
723:   PetscDesignatedInitializer(mtdot, VecMTDot_Seq),
724:   PetscDesignatedInitializer(scale, VecScale_Seq),
725:   PetscDesignatedInitializer(copy, VecCopy_Seq), /* 10 */
726:   PetscDesignatedInitializer(set, VecSet_Seq),
727:   PetscDesignatedInitializer(swap, VecSwap_Seq),
728:   PetscDesignatedInitializer(axpy, VecAXPY_Seq),
729:   PetscDesignatedInitializer(axpby, VecAXPBY_Seq),
730:   PetscDesignatedInitializer(maxpy, VecMAXPY_Seq),
731:   PetscDesignatedInitializer(aypx, VecAYPX_Seq),
732:   PetscDesignatedInitializer(waxpy, VecWAXPY_Seq),
733:   PetscDesignatedInitializer(axpbypcz, VecAXPBYPCZ_Seq),
734:   PetscDesignatedInitializer(pointwisemult, VecPointwiseMult_Seq),
735:   PetscDesignatedInitializer(pointwisedivide, VecPointwiseDivide_Seq),
736:   PetscDesignatedInitializer(setvalues, VecSetValues_Seq), /* 20 */
737:   PetscDesignatedInitializer(assemblybegin, NULL),
738:   PetscDesignatedInitializer(assemblyend, NULL),
739:   PetscDesignatedInitializer(getarray, NULL),
740:   PetscDesignatedInitializer(getsize, VecGetSize_Seq),
741:   PetscDesignatedInitializer(getlocalsize, VecGetSize_Seq),
742:   PetscDesignatedInitializer(restorearray, NULL),
743:   PetscDesignatedInitializer(max, VecMax_Seq),
744:   PetscDesignatedInitializer(min, VecMin_Seq),
745:   PetscDesignatedInitializer(setrandom, VecSetRandom_Seq),
746:   PetscDesignatedInitializer(setoption, VecSetOption_Seq), /* 30 */
747:   PetscDesignatedInitializer(setvaluesblocked, VecSetValuesBlocked_Seq),
748:   PetscDesignatedInitializer(destroy, VecDestroy_Seq),
749:   PetscDesignatedInitializer(view, VecView_Seq),
750:   PetscDesignatedInitializer(placearray, VecPlaceArray_Seq),
751:   PetscDesignatedInitializer(replacearray, VecReplaceArray_Seq),
752:   PetscDesignatedInitializer(dot_local, VecDot_Seq),
753:   PetscDesignatedInitializer(tdot_local, VecTDot_Seq),
754:   PetscDesignatedInitializer(norm_local, VecNorm_Seq),
755:   PetscDesignatedInitializer(mdot_local, VecMDot_Seq),
756:   PetscDesignatedInitializer(mtdot_local, VecMTDot_Seq), /* 40 */
757:   PetscDesignatedInitializer(load, VecLoad_Default),
758:   PetscDesignatedInitializer(reciprocal, VecReciprocal_Default),
759:   PetscDesignatedInitializer(conjugate, VecConjugate_Seq),
760:   PetscDesignatedInitializer(setlocaltoglobalmapping, NULL),
761:   PetscDesignatedInitializer(getlocaltoglobalmapping, NULL),
762:   PetscDesignatedInitializer(resetarray, VecResetArray_Seq),
763:   PetscDesignatedInitializer(setfromoptions, NULL),
764:   PetscDesignatedInitializer(maxpointwisedivide, VecMaxPointwiseDivide_Seq),
765:   PetscDesignatedInitializer(pointwisemax, VecPointwiseMax_Seq),
766:   PetscDesignatedInitializer(pointwisemaxabs, VecPointwiseMaxAbs_Seq),
767:   PetscDesignatedInitializer(pointwisemin, VecPointwiseMin_Seq),
768:   PetscDesignatedInitializer(getvalues, VecGetValues_Seq),
769:   PetscDesignatedInitializer(sqrt, NULL),
770:   PetscDesignatedInitializer(abs, NULL),
771:   PetscDesignatedInitializer(exp, NULL),
772:   PetscDesignatedInitializer(log, NULL),
773:   PetscDesignatedInitializer(shift, NULL),
774:   PetscDesignatedInitializer(create, NULL),
775:   PetscDesignatedInitializer(stridegather, VecStrideGather_Default),
776:   PetscDesignatedInitializer(stridescatter, VecStrideScatter_Default),
777:   PetscDesignatedInitializer(dotnorm2, NULL),
778:   PetscDesignatedInitializer(getsubvector, NULL),
779:   PetscDesignatedInitializer(restoresubvector, NULL),
780:   PetscDesignatedInitializer(getarrayread, NULL),
781:   PetscDesignatedInitializer(restorearrayread, NULL),
782:   PetscDesignatedInitializer(stridesubsetgather, VecStrideSubSetGather_Default),
783:   PetscDesignatedInitializer(stridesubsetscatter, VecStrideSubSetScatter_Default),
784:   PetscDesignatedInitializer(viewnative, VecView_Seq),
785:   PetscDesignatedInitializer(loadnative, NULL),
786:   PetscDesignatedInitializer(createlocalvector, NULL),
787:   PetscDesignatedInitializer(getlocalvector, NULL),
788:   PetscDesignatedInitializer(restorelocalvector, NULL),
789:   PetscDesignatedInitializer(getlocalvectorread, NULL),
790:   PetscDesignatedInitializer(restorelocalvectorread, NULL),
791:   PetscDesignatedInitializer(bindtocpu, NULL),
792:   PetscDesignatedInitializer(getarraywrite, NULL),
793:   PetscDesignatedInitializer(restorearraywrite, NULL),
794:   PetscDesignatedInitializer(getarrayandmemtype, NULL),
795:   PetscDesignatedInitializer(restorearrayandmemtype, NULL),
796:   PetscDesignatedInitializer(getarrayreadandmemtype, NULL),
797:   PetscDesignatedInitializer(restorearrayreadandmemtype, NULL),
798:   PetscDesignatedInitializer(getarraywriteandmemtype, NULL),
799:   PetscDesignatedInitializer(restorearraywriteandmemtype, NULL),
800:   PetscDesignatedInitializer(concatenate, NULL),
801:   PetscDesignatedInitializer(sum, NULL),
802:   PetscDesignatedInitializer(setpreallocationcoo, VecSetPreallocationCOO_Seq),
803:   PetscDesignatedInitializer(setvaluescoo, VecSetValuesCOO_Seq),
804:   PetscDesignatedInitializer(errorwnorm, NULL),
805:   PetscDesignatedInitializer(maxpby, NULL),
806: };

808: /*
809:   Create a VECSEQ with the given layout and array

811:   Input Parameter:
812: + map   - the layout
813: - array - the array on host

815:   Output Parameter:
816: . V  - The vector object
817: */
818: PetscErrorCode VecCreateSeqWithLayoutAndArray_Private(PetscLayout map, const PetscScalar array[], Vec *V)
819: {
820:   PetscMPIInt size;

822:   PetscFunctionBegin;
823:   PetscCall(VecCreateWithLayout_Private(map, V));
824:   PetscCallMPI(MPI_Comm_size(map->comm, &size));
825:   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQ on more than one process");
826:   PetscCall(VecCreate_Seq_Private(*V, array));
827:   PetscFunctionReturn(PETSC_SUCCESS);
828: }

830: /*
831:       This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
832: */
833: PetscErrorCode VecCreate_Seq_Private(Vec v, const PetscScalar array[])
834: {
835:   Vec_Seq  *s;
836:   PetscBool mdot_use_gemv  = PETSC_TRUE;
837:   PetscBool maxpy_use_gemv = PETSC_FALSE; // default is false as we saw bad performance with vendors' GEMV with tall skinny matrices.

839:   PetscFunctionBegin;
840:   PetscCall(PetscNew(&s));
841:   v->ops[0] = DvOps;

843:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_mdot_use_gemv", &mdot_use_gemv, NULL));
844:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_maxpy_use_gemv", &maxpy_use_gemv, NULL));

846:   // allocate multiple vectors together
847:   if (mdot_use_gemv || maxpy_use_gemv) v->ops[0].duplicatevecs = VecDuplicateVecs_Seq_GEMV;

849:   if (mdot_use_gemv) {
850:     v->ops[0].mdot        = VecMDot_Seq_GEMV;
851:     v->ops[0].mdot_local  = VecMDot_Seq_GEMV;
852:     v->ops[0].mtdot       = VecMTDot_Seq_GEMV;
853:     v->ops[0].mtdot_local = VecMTDot_Seq_GEMV;
854:   }
855:   if (maxpy_use_gemv) v->ops[0].maxpy = VecMAXPY_Seq_GEMV;

857:   v->data            = (void *)s;
858:   v->petscnative     = PETSC_TRUE;
859:   s->array           = (PetscScalar *)array;
860:   s->array_allocated = NULL;
861:   if (array) v->offloadmask = PETSC_OFFLOAD_CPU;

863:   PetscCall(PetscLayoutSetUp(v->map));
864:   PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQ));
865: #if defined(PETSC_HAVE_MATLAB)
866:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEnginePut_C", VecMatlabEnginePut_Default));
867:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEngineGet_C", VecMatlabEngineGet_Default));
868: #endif
869:   PetscFunctionReturn(PETSC_SUCCESS);
870: }

872: /*@
873:   VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
874:   where the user provides the array space to store the vector values.

876:   Collective

878:   Input Parameters:
879: + comm  - the communicator, should be `PETSC_COMM_SELF`
880: . bs    - the block size
881: . n     - the vector length
882: - array - memory where the vector elements are to be stored.

884:   Output Parameter:
885: . V - the vector

887:   Level: intermediate

889:   Notes:
890:   Use `VecDuplicate()` or `VecDuplicateVecs(`) to form additional vectors of the
891:   same type as an existing vector.

893:   If the user-provided array is` NULL`, then `VecPlaceArray()` can be used
894:   at a later stage to SET the array for storing the vector values.

896:   PETSc does NOT free the array when the vector is destroyed via `VecDestroy()`.
897:   The user should not free the array until the vector is destroyed.

899: .seealso: `VecCreateMPIWithArray()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`,
900:           `VecCreateGhost()`, `VecCreateSeq()`, `VecPlaceArray()`
901: @*/
902: PetscErrorCode VecCreateSeqWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar array[], Vec *V)
903: {
904:   PetscMPIInt size;

906:   PetscFunctionBegin;
907:   PetscCall(VecCreate(comm, V));
908:   PetscCall(VecSetSizes(*V, n, n));
909:   PetscCall(VecSetBlockSize(*V, bs));
910:   PetscCallMPI(MPI_Comm_size(comm, &size));
911:   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQ on more than one process");
912:   PetscCall(VecCreate_Seq_Private(*V, array));
913:   PetscFunctionReturn(PETSC_SUCCESS);
914: }