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 (PetscScalar)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(VecGetArray(yin, &ya));
162:     PetscCall(PetscArraycpy(ya, xa, xin->map->n));
163:     PetscCall(VecRestoreArrayRead(xin, &xa));
164:     PetscCall(VecRestoreArray(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_VTK_DEPRECATED || format == PETSC_VIEWER_ASCII_VTK_CELL_DEPRECATED) {
275:     /*
276:        state 0: No header has been output
277:        state 1: Only POINT_DATA has been output
278:        state 2: Only CELL_DATA has been output
279:        state 3: Output both, POINT_DATA last
280:        state 4: Output both, CELL_DATA last
281:     */
282:     static PetscInt stateId     = -1;
283:     int             outputState = 0;
284:     PetscBool       hasState;
285:     int             doOutput = 0;
286:     PetscInt        bs, b;

288:     if (stateId < 0) PetscCall(PetscObjectComposedDataRegister(&stateId));
289:     PetscCall(PetscObjectComposedDataGetInt((PetscObject)viewer, stateId, outputState, hasState));
290:     if (!hasState) outputState = 0;
291:     PetscCall(PetscObjectGetName((PetscObject)xin, &name));
292:     PetscCall(VecGetBlockSize(xin, &bs));
293:     PetscCheck(bs >= 1 && bs <= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %" PetscInt_FMT, bs);
294:     if (format == PETSC_VIEWER_ASCII_VTK_DEPRECATED) {
295:       if (outputState == 0) {
296:         outputState = 1;
297:         doOutput    = 1;
298:       } else if (outputState == 1) doOutput = 0;
299:       else if (outputState == 2) {
300:         outputState = 3;
301:         doOutput    = 1;
302:       } else if (outputState == 3) doOutput = 0;
303:       else PetscCheck(outputState != 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Tried to output POINT_DATA again after intervening CELL_DATA");

305:       if (doOutput) PetscCall(PetscViewerASCIIPrintf(viewer, "POINT_DATA %" PetscInt_FMT "\n", n / bs));
306:     } else {
307:       if (outputState == 0) {
308:         outputState = 2;
309:         doOutput    = 1;
310:       } else if (outputState == 1) {
311:         outputState = 4;
312:         doOutput    = 1;
313:       } else if (outputState == 2) {
314:         doOutput = 0;
315:       } else {
316:         PetscCheck(outputState != 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Tried to output CELL_DATA again after intervening POINT_DATA");
317:         if (outputState == 4) doOutput = 0;
318:       }

320:       if (doOutput) PetscCall(PetscViewerASCIIPrintf(viewer, "CELL_DATA %" PetscInt_FMT "\n", n));
321:     }
322:     PetscCall(PetscObjectComposedDataSetInt((PetscObject)viewer, stateId, outputState));
323:     if (name) {
324:       if (bs == 3) {
325:         PetscCall(PetscViewerASCIIPrintf(viewer, "VECTORS %s double\n", name));
326:       } else {
327:         PetscCall(PetscViewerASCIIPrintf(viewer, "SCALARS %s double %" PetscInt_FMT "\n", name, bs));
328:       }
329:     } else {
330:       PetscCall(PetscViewerASCIIPrintf(viewer, "SCALARS scalars double %" PetscInt_FMT "\n", bs));
331:     }
332:     if (bs != 3) PetscCall(PetscViewerASCIIPrintf(viewer, "LOOKUP_TABLE default\n"));
333:     for (i = 0; i < n / bs; i++) {
334:       for (b = 0; b < bs; b++) {
335:         if (b > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " "));
336: #if !defined(PETSC_USE_COMPLEX)
337:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g", (double)xv[i * bs + b]));
338: #endif
339:       }
340:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
341:     }
342:   } else if (format == PETSC_VIEWER_ASCII_VTK_COORDS_DEPRECATED) {
343:     PetscInt bs, b;

345:     PetscCall(VecGetBlockSize(xin, &bs));
346:     PetscCheck(bs >= 1 && bs <= 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "VTK can only handle 3D objects, but vector dimension is %" PetscInt_FMT, bs);
347:     for (i = 0; i < n / bs; i++) {
348:       for (b = 0; b < bs; b++) {
349:         if (b > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " "));
350: #if !defined(PETSC_USE_COMPLEX)
351:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g", (double)xv[i * bs + b]));
352: #endif
353:       }
354:       for (b = bs; b < 3; b++) PetscCall(PetscViewerASCIIPrintf(viewer, " 0.0"));
355:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
356:     }
357:   } else if (format == PETSC_VIEWER_ASCII_PCICE) {
358:     PetscInt bs, b;

360:     PetscCall(VecGetBlockSize(xin, &bs));
361:     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);
362:     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", xin->map->N / bs));
363:     for (i = 0; i < n / bs; i++) {
364:       PetscCall(PetscViewerASCIIPrintf(viewer, "%7" PetscInt_FMT "   ", i + 1));
365:       for (b = 0; b < bs; b++) {
366:         if (b > 0) PetscCall(PetscViewerASCIIPrintf(viewer, " "));
367: #if !defined(PETSC_USE_COMPLEX)
368:         PetscCall(PetscViewerASCIIPrintf(viewer, "% 12.5E", (double)xv[i * bs + b]));
369: #endif
370:       }
371:       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
372:     }
373:   } else if (format == PETSC_VIEWER_ASCII_GLVIS) {
374:     /* GLVis ASCII visualization/dump: this function mimics mfem::GridFunction::Save() */
375:     const PetscScalar      *array;
376:     PetscInt                i, n, vdim, ordering = 1; /* mfem::FiniteElementSpace::Ordering::byVDIM */
377:     PetscContainer          glvis_container;
378:     PetscViewerGLVisVecInfo glvis_vec_info;
379:     PetscViewerGLVisInfo    glvis_info;

381:     /* mfem::FiniteElementSpace::Save() */
382:     PetscCall(VecGetBlockSize(xin, &vdim));
383:     PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
384:     PetscCall(PetscObjectQuery((PetscObject)xin, "_glvis_info_container", (PetscObject *)&glvis_container));
385:     PetscCheck(glvis_container, PetscObjectComm((PetscObject)xin), PETSC_ERR_PLIB, "Missing GLVis container");
386:     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_vec_info));
387:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", glvis_vec_info->fec_type));
388:     PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", vdim));
389:     PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: %" PetscInt_FMT "\n", ordering));
390:     PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
391:     /* mfem::Vector::Print() */
392:     PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
393:     PetscCheck(glvis_container, PetscObjectComm((PetscObject)viewer), PETSC_ERR_PLIB, "Missing GLVis container");
394:     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info));
395:     if (glvis_info->enabled) {
396:       PetscCall(VecGetLocalSize(xin, &n));
397:       PetscCall(VecGetArrayRead(xin, &array));
398:       for (i = 0; i < n; i++) {
399:         PetscCall(PetscViewerASCIIPrintf(viewer, glvis_info->fmt, (double)PetscRealPart(array[i])));
400:         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
401:       }
402:       PetscCall(VecRestoreArrayRead(xin, &array));
403:     }
404:   } else if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
405:     /* No info */
406:   } else {
407:     for (i = 0; i < n; i++) {
408:       if (format == PETSC_VIEWER_ASCII_INDEX) PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT ": ", i));
409: #if defined(PETSC_USE_COMPLEX)
410:       if (PetscImaginaryPart(xv[i]) > 0.0) {
411:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g + %g i\n", (double)PetscRealPart(xv[i]), (double)PetscImaginaryPart(xv[i])));
412:       } else if (PetscImaginaryPart(xv[i]) < 0.0) {
413:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g - %g i\n", (double)PetscRealPart(xv[i]), -(double)PetscImaginaryPart(xv[i])));
414:       } else {
415:         PetscCall(PetscViewerASCIIPrintf(viewer, "%g\n", (double)PetscRealPart(xv[i])));
416:       }
417: #else
418:       PetscCall(PetscViewerASCIIPrintf(viewer, "%g\n", (double)xv[i]));
419: #endif
420:     }
421:   }
422:   PetscCall(PetscViewerFlush(viewer));
423:   PetscCall(VecRestoreArrayRead(xin, &xv));
424:   PetscFunctionReturn(PETSC_SUCCESS);
425: }

427: #include <petscdraw.h>
428: static PetscErrorCode VecView_Seq_Draw_LG(Vec xin, PetscViewer v)
429: {
430:   PetscDraw          draw;
431:   PetscBool          isnull;
432:   PetscDrawLG        lg;
433:   PetscInt           i, c, bs = PetscAbs(xin->map->bs), n = xin->map->n / bs;
434:   const PetscScalar *xv;
435:   PetscReal         *xx, *yy, xmin, xmax, h;
436:   int                colors[] = {PETSC_DRAW_RED};
437:   PetscViewerFormat  format;
438:   PetscDrawAxis      axis;

440:   PetscFunctionBegin;
441:   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
442:   PetscCall(PetscDrawIsNull(draw, &isnull));
443:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

445:   PetscCall(PetscViewerGetFormat(v, &format));
446:   PetscCall(PetscMalloc2(n, &xx, n, &yy));
447:   PetscCall(VecGetArrayRead(xin, &xv));
448:   for (c = 0; c < bs; c++) {
449:     PetscCall(PetscViewerDrawGetDrawLG(v, c, &lg));
450:     PetscCall(PetscDrawLGReset(lg));
451:     PetscCall(PetscDrawLGSetDimension(lg, 1));
452:     PetscCall(PetscDrawLGSetColors(lg, colors));
453:     if (format == PETSC_VIEWER_DRAW_LG_XRANGE) {
454:       PetscCall(PetscDrawLGGetAxis(lg, &axis));
455:       PetscCall(PetscDrawAxisGetLimits(axis, &xmin, &xmax, NULL, NULL));
456:       h = (xmax - xmin) / n;
457:       for (i = 0; i < n; i++) xx[i] = i * h + 0.5 * h; /* cell center */
458:     } else {
459:       for (i = 0; i < n; i++) xx[i] = (PetscReal)i;
460:     }
461:     for (i = 0; i < n; i++) yy[i] = PetscRealPart(xv[c + i * bs]);

463:     PetscCall(PetscDrawLGAddPoints(lg, n, &xx, &yy));
464:     PetscCall(PetscDrawLGDraw(lg));
465:     PetscCall(PetscDrawLGSave(lg));
466:   }
467:   PetscCall(VecRestoreArrayRead(xin, &xv));
468:   PetscCall(PetscFree2(xx, yy));
469:   PetscFunctionReturn(PETSC_SUCCESS);
470: }

472: static PetscErrorCode VecView_Seq_Draw(Vec xin, PetscViewer v)
473: {
474:   PetscDraw draw;
475:   PetscBool isnull;

477:   PetscFunctionBegin;
478:   PetscCall(PetscViewerDrawGetDraw(v, 0, &draw));
479:   PetscCall(PetscDrawIsNull(draw, &isnull));
480:   if (isnull) PetscFunctionReturn(PETSC_SUCCESS);

482:   PetscCall(VecView_Seq_Draw_LG(xin, v));
483:   PetscFunctionReturn(PETSC_SUCCESS);
484: }

486: static PetscErrorCode VecView_Seq_Binary(Vec xin, PetscViewer viewer)
487: {
488:   return VecView_Binary(xin, viewer);
489: }

491: #if defined(PETSC_HAVE_MATLAB)
492: #include <petscmatlab.h>
493:   #include <mat.h> /* MATLAB include file */
494: PetscErrorCode VecView_Seq_Matlab(Vec vec, PetscViewer viewer)
495: {
496:   PetscInt           n;
497:   const PetscScalar *array;

499:   PetscFunctionBegin;
500:   PetscCall(VecGetLocalSize(vec, &n));
501:   PetscCall(PetscObjectName((PetscObject)vec));
502:   PetscCall(VecGetArrayRead(vec, &array));
503:   PetscCall(PetscViewerMatlabPutArray(viewer, n, 1, array, ((PetscObject)vec)->name));
504:   PetscCall(VecRestoreArrayRead(vec, &array));
505:   PetscFunctionReturn(PETSC_SUCCESS);
506: }
507: #endif

509: PetscErrorCode VecView_Seq(Vec xin, PetscViewer viewer)
510: {
511:   PetscBool isdraw, iascii, issocket, isbinary;
512: #if defined(PETSC_HAVE_MATHEMATICA)
513:   PetscBool ismathematica;
514: #endif
515: #if defined(PETSC_HAVE_MATLAB)
516:   PetscBool ismatlab;
517: #endif
518: #if defined(PETSC_HAVE_HDF5)
519:   PetscBool ishdf5;
520: #endif
521:   PetscBool isglvis;
522: #if defined(PETSC_HAVE_ADIOS)
523:   PetscBool isadios;
524: #endif

526:   PetscFunctionBegin;
527:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
528:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
529:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSOCKET, &issocket));
530:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
531: #if defined(PETSC_HAVE_MATHEMATICA)
532:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATHEMATICA, &ismathematica));
533: #endif
534: #if defined(PETSC_HAVE_HDF5)
535:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
536: #endif
537: #if defined(PETSC_HAVE_MATLAB)
538:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERMATLAB, &ismatlab));
539: #endif
540:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
541: #if defined(PETSC_HAVE_ADIOS)
542:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERADIOS, &isadios));
543: #endif

545:   if (isdraw) {
546:     PetscCall(VecView_Seq_Draw(xin, viewer));
547:   } else if (iascii) {
548:     PetscCall(VecView_Seq_ASCII(xin, viewer));
549:   } else if (isbinary) {
550:     PetscCall(VecView_Seq_Binary(xin, viewer));
551: #if defined(PETSC_HAVE_MATHEMATICA)
552:   } else if (ismathematica) {
553:     PetscCall(PetscViewerMathematicaPutVector(viewer, xin));
554: #endif
555: #if defined(PETSC_HAVE_HDF5)
556:   } else if (ishdf5) {
557:     PetscCall(VecView_MPI_HDF5(xin, viewer)); /* Reusing VecView_MPI_HDF5 ... don't want code duplication*/
558: #endif
559: #if defined(PETSC_HAVE_ADIOS)
560:   } else if (isadios) {
561:     PetscCall(VecView_MPI_ADIOS(xin, viewer)); /* Reusing VecView_MPI_ADIOS ... don't want code duplication*/
562: #endif
563: #if defined(PETSC_HAVE_MATLAB)
564:   } else if (ismatlab) {
565:     PetscCall(VecView_Seq_Matlab(xin, viewer));
566: #endif
567:   } else if (isglvis) PetscCall(VecView_GLVis(xin, viewer));
568:   PetscFunctionReturn(PETSC_SUCCESS);
569: }

571: PetscErrorCode VecGetValues_Seq(Vec xin, PetscInt ni, const PetscInt ix[], PetscScalar y[])
572: {
573:   const PetscBool    ignorenegidx = xin->stash.ignorenegidx;
574:   const PetscScalar *xx;

576:   PetscFunctionBegin;
577:   PetscCall(VecGetArrayRead(xin, &xx));
578:   for (PetscInt i = 0; i < ni; ++i) {
579:     if (ignorenegidx && (ix[i] < 0)) continue;
580:     if (PetscDefined(USE_DEBUG)) {
581:       PetscCheck(ix[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " cannot be negative", ix[i]);
582:       PetscCheck(ix[i] < xin->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " to large maximum allowed %" PetscInt_FMT, ix[i], xin->map->n);
583:     }
584:     y[i] = xx[ix[i]];
585:   }
586:   PetscCall(VecRestoreArrayRead(xin, &xx));
587:   PetscFunctionReturn(PETSC_SUCCESS);
588: }

590: PetscErrorCode VecSetValues_Seq(Vec xin, PetscInt ni, const PetscInt ix[], const PetscScalar y[], InsertMode m)
591: {
592:   const PetscBool ignorenegidx = xin->stash.ignorenegidx;
593:   PetscScalar    *xx;

595:   PetscFunctionBegin;
596:   // call to getarray (not e.g. getarraywrite() if m is INSERT_VALUES) is deliberate! If this
597:   // is secretly a VECSEQCUDA it may have values currently on the device, in which case --
598:   // unless we are replacing the entire array -- we need to copy them up
599:   PetscCall(VecGetArray(xin, &xx));
600:   for (PetscInt i = 0; i < ni; i++) {
601:     if (ignorenegidx && (ix[i] < 0)) continue;
602:     if (PetscDefined(USE_DEBUG)) {
603:       PetscCheck(ix[i] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " cannot be negative", ix[i]);
604:       PetscCheck(ix[i] < xin->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " maximum %" PetscInt_FMT, ix[i], xin->map->n);
605:     }
606:     if (m == INSERT_VALUES) {
607:       xx[ix[i]] = y[i];
608:     } else {
609:       xx[ix[i]] += y[i];
610:     }
611:   }
612:   PetscCall(VecRestoreArray(xin, &xx));
613:   PetscFunctionReturn(PETSC_SUCCESS);
614: }

616: PetscErrorCode VecSetValuesBlocked_Seq(Vec xin, PetscInt ni, const PetscInt ix[], const PetscScalar yin[], InsertMode m)
617: {
618:   PetscScalar *xx;
619:   PetscInt     bs;

621:   /* For optimization could treat bs = 2, 3, 4, 5 as special cases with loop unrolling */
622:   PetscFunctionBegin;
623:   PetscCall(VecGetBlockSize(xin, &bs));
624:   PetscCall(VecGetArray(xin, &xx));
625:   for (PetscInt i = 0; i < ni; ++i, yin += bs) {
626:     const PetscInt start = bs * ix[i];

628:     if (start < 0) continue;
629:     PetscCheck(start < xin->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Out of range index value %" PetscInt_FMT " maximum %" PetscInt_FMT, start, xin->map->n);
630:     for (PetscInt j = 0; j < bs; j++) {
631:       if (m == INSERT_VALUES) {
632:         xx[start + j] = yin[j];
633:       } else {
634:         xx[start + j] += yin[j];
635:       }
636:     }
637:   }
638:   PetscCall(VecRestoreArray(xin, &xx));
639:   PetscFunctionReturn(PETSC_SUCCESS);
640: }

642: static PetscErrorCode VecResetPreallocationCOO_Seq(Vec x)
643: {
644:   Vec_Seq *vs = (Vec_Seq *)x->data;

646:   PetscFunctionBegin;
647:   if (vs) {
648:     PetscCall(PetscFree(vs->jmap1)); /* Destroy old stuff */
649:     PetscCall(PetscFree(vs->perm1));
650:   }
651:   PetscFunctionReturn(PETSC_SUCCESS);
652: }

654: PetscErrorCode VecSetPreallocationCOO_Seq(Vec x, PetscCount coo_n, const PetscInt coo_i[])
655: {
656:   PetscInt    m, *i;
657:   PetscCount  k, nneg;
658:   PetscCount *perm1, *jmap1;
659:   Vec_Seq    *vs = (Vec_Seq *)x->data;

661:   PetscFunctionBegin;
662:   PetscCall(VecResetPreallocationCOO_Seq(x)); /* Destroy old stuff */
663:   PetscCall(PetscMalloc1(coo_n, &i));
664:   PetscCall(PetscArraycpy(i, coo_i, coo_n)); /* Make a copy since we'll modify it */
665:   PetscCall(PetscMalloc1(coo_n, &perm1));
666:   for (k = 0; k < coo_n; k++) perm1[k] = k;
667:   PetscCall(PetscSortIntWithCountArray(coo_n, i, perm1));
668:   for (k = 0; k < coo_n; k++) {
669:     if (i[k] >= 0) break;
670:   } /* Advance k to the first entry with a non-negative index */
671:   nneg = k;

673:   PetscCall(VecGetLocalSize(x, &m));
674:   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");
675:   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);

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

682:   if (nneg) { /* Discard leading negative indices */
683:     PetscCount *perm1_new;
684:     PetscCall(PetscMalloc1(coo_n - nneg, &perm1_new));
685:     PetscCall(PetscArraycpy(perm1_new, perm1 + nneg, coo_n - nneg));
686:     PetscCall(PetscFree(perm1));
687:     perm1 = perm1_new;
688:   }

690:   /* Record COO fields */
691:   vs->coo_n = coo_n;
692:   vs->tot1  = coo_n - nneg;
693:   vs->jmap1 = jmap1; /* [m+1] */
694:   vs->perm1 = perm1; /* [tot] */
695:   PetscFunctionReturn(PETSC_SUCCESS);
696: }

698: PetscErrorCode VecSetValuesCOO_Seq(Vec x, const PetscScalar coo_v[], InsertMode imode)
699: {
700:   Vec_Seq          *vs    = (Vec_Seq *)x->data;
701:   const PetscCount *perm1 = vs->perm1, *jmap1 = vs->jmap1;
702:   PetscScalar      *xv;
703:   PetscInt          m;

705:   PetscFunctionBegin;
706:   PetscCall(VecGetLocalSize(x, &m));
707:   PetscCall(VecGetArray(x, &xv));
708:   for (PetscInt i = 0; i < m; i++) {
709:     PetscScalar sum = 0.0;
710:     for (PetscCount j = jmap1[i]; j < jmap1[i + 1]; j++) sum += coo_v[perm1[j]];
711:     xv[i] = (imode == INSERT_VALUES ? 0.0 : xv[i]) + sum;
712:   }
713:   PetscCall(VecRestoreArray(x, &xv));
714:   PetscFunctionReturn(PETSC_SUCCESS);
715: }

717: PetscErrorCode VecDestroy_Seq(Vec v)
718: {
719:   Vec_Seq *vs = (Vec_Seq *)v->data;

721:   PetscFunctionBegin;
722:   PetscCall(PetscLogObjectState((PetscObject)v, "Length=%" PetscInt_FMT, v->map->n));
723:   if (vs) PetscCall(PetscFree(vs->array_allocated));
724:   PetscCall(VecResetPreallocationCOO_Seq(v));
725:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEnginePut_C", NULL));
726:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEngineGet_C", NULL));
727:   PetscCall(PetscFree(v->data));
728:   PetscFunctionReturn(PETSC_SUCCESS);
729: }

731: PetscErrorCode VecSetOption_Seq(Vec v, VecOption op, PetscBool flag)
732: {
733:   PetscFunctionBegin;
734:   if (op == VEC_IGNORE_NEGATIVE_INDICES) v->stash.ignorenegidx = flag;
735:   PetscFunctionReturn(PETSC_SUCCESS);
736: }

738: PetscErrorCode VecDuplicate_Seq(Vec win, Vec *V)
739: {
740:   PetscFunctionBegin;
741:   PetscCall(VecCreateWithLayout_Private(win->map, V));
742:   PetscCall(VecSetType(*V, ((PetscObject)win)->type_name));
743:   PetscCall(PetscObjectListDuplicate(((PetscObject)win)->olist, &((PetscObject)*V)->olist));
744:   PetscCall(PetscFunctionListDuplicate(((PetscObject)win)->qlist, &((PetscObject)*V)->qlist));

746:   (*V)->ops->view          = win->ops->view;
747:   (*V)->stash.ignorenegidx = win->stash.ignorenegidx;
748:   PetscFunctionReturn(PETSC_SUCCESS);
749: }

751: PetscErrorCode VecReplaceArray_Default_GEMV_Error(Vec v, const PetscScalar *a)
752: {
753:   PetscFunctionBegin;
754:   PetscCheck(PETSC_FALSE, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "VecReplaceArray() is not supported on the first Vec obtained from VecDuplicateVecs(). \
755: 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.");
756:   (void)a;
757:   PetscFunctionReturn(PETSC_SUCCESS);
758: }

760: static PetscErrorCode VecDuplicateVecs_Seq_GEMV(Vec w, PetscInt m, Vec *V[])
761: {
762:   PetscFunctionBegin;
763:   // This routine relies on the duplicate operation being VecDuplicate_Seq. If not, bail out to the default.
764:   if (w->ops->duplicate != VecDuplicate_Seq) {
765:     w->ops->duplicatevecs = VecDuplicateVecs_Default;
766:     PetscCall(VecDuplicateVecs(w, m, V));
767:   } else {
768:     PetscInt     nlocal;
769:     PetscScalar *array;
770:     PetscInt64   lda; // use 64-bit as we will do "m * lda"

772:     PetscCall(PetscMalloc1(m, V));
773:     PetscCall(VecGetLocalSize(w, &nlocal));
774:     lda = nlocal;
775:     lda = ((lda + 31) / 32) * 32; // make every vector 32-elements aligned

777:     PetscCall(PetscCalloc1(m * lda, &array));
778:     for (PetscInt i = 0; i < m; i++) {
779:       Vec v;
780:       PetscCall(VecCreateSeqWithLayoutAndArray_Private(w->map, PetscSafePointerPlusOffset(array, i * lda), &v));
781:       PetscCall(PetscObjectListDuplicate(((PetscObject)w)->olist, &((PetscObject)v)->olist));
782:       PetscCall(PetscFunctionListDuplicate(((PetscObject)w)->qlist, &((PetscObject)v)->qlist));
783:       v->ops->view          = w->ops->view;
784:       v->stash.ignorenegidx = w->stash.ignorenegidx;
785:       (*V)[i]               = v;
786:     }
787:     // so when the first vector is destroyed it will destroy the array
788:     if (m) ((Vec_Seq *)(*V)[0]->data)->array_allocated = array;
789:     // disable replacearray of the first vector, as freeing its memory also frees others in the group.
790:     // But replacearray of others is ok, as they don't own their array.
791:     if (m > 1) (*V)[0]->ops->replacearray = VecReplaceArray_Default_GEMV_Error;
792:   }
793:   PetscFunctionReturn(PETSC_SUCCESS);
794: }

796: static struct _VecOps DvOps = {
797:   PetscDesignatedInitializer(duplicate, VecDuplicate_Seq), /* 1 */
798:   PetscDesignatedInitializer(duplicatevecs, VecDuplicateVecs_Default),
799:   PetscDesignatedInitializer(destroyvecs, VecDestroyVecs_Default),
800:   PetscDesignatedInitializer(dot, VecDot_Seq),
801:   PetscDesignatedInitializer(mdot, VecMDot_Seq),
802:   PetscDesignatedInitializer(norm, VecNorm_Seq),
803:   PetscDesignatedInitializer(tdot, VecTDot_Seq),
804:   PetscDesignatedInitializer(mtdot, VecMTDot_Seq),
805:   PetscDesignatedInitializer(scale, VecScale_Seq),
806:   PetscDesignatedInitializer(copy, VecCopy_Seq), /* 10 */
807:   PetscDesignatedInitializer(set, VecSet_Seq),
808:   PetscDesignatedInitializer(swap, VecSwap_Seq),
809:   PetscDesignatedInitializer(axpy, VecAXPY_Seq),
810:   PetscDesignatedInitializer(axpby, VecAXPBY_Seq),
811:   PetscDesignatedInitializer(maxpy, VecMAXPY_Seq),
812:   PetscDesignatedInitializer(aypx, VecAYPX_Seq),
813:   PetscDesignatedInitializer(waxpy, VecWAXPY_Seq),
814:   PetscDesignatedInitializer(axpbypcz, VecAXPBYPCZ_Seq),
815:   PetscDesignatedInitializer(pointwisemult, VecPointwiseMult_Seq),
816:   PetscDesignatedInitializer(pointwisedivide, VecPointwiseDivide_Seq),
817:   PetscDesignatedInitializer(setvalues, VecSetValues_Seq), /* 20 */
818:   PetscDesignatedInitializer(assemblybegin, NULL),
819:   PetscDesignatedInitializer(assemblyend, NULL),
820:   PetscDesignatedInitializer(getarray, NULL),
821:   PetscDesignatedInitializer(getsize, VecGetSize_Seq),
822:   PetscDesignatedInitializer(getlocalsize, VecGetSize_Seq),
823:   PetscDesignatedInitializer(restorearray, NULL),
824:   PetscDesignatedInitializer(max, VecMax_Seq),
825:   PetscDesignatedInitializer(min, VecMin_Seq),
826:   PetscDesignatedInitializer(setrandom, VecSetRandom_Seq),
827:   PetscDesignatedInitializer(setoption, VecSetOption_Seq), /* 30 */
828:   PetscDesignatedInitializer(setvaluesblocked, VecSetValuesBlocked_Seq),
829:   PetscDesignatedInitializer(destroy, VecDestroy_Seq),
830:   PetscDesignatedInitializer(view, VecView_Seq),
831:   PetscDesignatedInitializer(placearray, VecPlaceArray_Seq),
832:   PetscDesignatedInitializer(replacearray, VecReplaceArray_Seq),
833:   PetscDesignatedInitializer(dot_local, VecDot_Seq),
834:   PetscDesignatedInitializer(tdot_local, VecTDot_Seq),
835:   PetscDesignatedInitializer(norm_local, VecNorm_Seq),
836:   PetscDesignatedInitializer(mdot_local, VecMDot_Seq),
837:   PetscDesignatedInitializer(mtdot_local, VecMTDot_Seq), /* 40 */
838:   PetscDesignatedInitializer(load, VecLoad_Default),
839:   PetscDesignatedInitializer(reciprocal, VecReciprocal_Default),
840:   PetscDesignatedInitializer(conjugate, VecConjugate_Seq),
841:   PetscDesignatedInitializer(setlocaltoglobalmapping, NULL),
842:   PetscDesignatedInitializer(getlocaltoglobalmapping, NULL),
843:   PetscDesignatedInitializer(setvalueslocal, NULL),
844:   PetscDesignatedInitializer(resetarray, VecResetArray_Seq),
845:   PetscDesignatedInitializer(setfromoptions, NULL),
846:   PetscDesignatedInitializer(maxpointwisedivide, VecMaxPointwiseDivide_Seq),
847:   PetscDesignatedInitializer(pointwisemax, VecPointwiseMax_Seq),
848:   PetscDesignatedInitializer(pointwisemaxabs, VecPointwiseMaxAbs_Seq),
849:   PetscDesignatedInitializer(pointwisemin, VecPointwiseMin_Seq),
850:   PetscDesignatedInitializer(getvalues, VecGetValues_Seq),
851:   PetscDesignatedInitializer(sqrt, NULL),
852:   PetscDesignatedInitializer(abs, NULL),
853:   PetscDesignatedInitializer(exp, NULL),
854:   PetscDesignatedInitializer(log, NULL),
855:   PetscDesignatedInitializer(shift, NULL),
856:   PetscDesignatedInitializer(create, NULL),
857:   PetscDesignatedInitializer(stridegather, VecStrideGather_Default),
858:   PetscDesignatedInitializer(stridescatter, VecStrideScatter_Default),
859:   PetscDesignatedInitializer(dotnorm2, NULL),
860:   PetscDesignatedInitializer(getsubvector, NULL),
861:   PetscDesignatedInitializer(restoresubvector, NULL),
862:   PetscDesignatedInitializer(getarrayread, NULL),
863:   PetscDesignatedInitializer(restorearrayread, NULL),
864:   PetscDesignatedInitializer(stridesubsetgather, VecStrideSubSetGather_Default),
865:   PetscDesignatedInitializer(stridesubsetscatter, VecStrideSubSetScatter_Default),
866:   PetscDesignatedInitializer(viewnative, VecView_Seq),
867:   PetscDesignatedInitializer(loadnative, NULL),
868:   PetscDesignatedInitializer(createlocalvector, NULL),
869:   PetscDesignatedInitializer(getlocalvector, NULL),
870:   PetscDesignatedInitializer(restorelocalvector, NULL),
871:   PetscDesignatedInitializer(getlocalvectorread, NULL),
872:   PetscDesignatedInitializer(restorelocalvectorread, NULL),
873:   PetscDesignatedInitializer(bindtocpu, NULL),
874:   PetscDesignatedInitializer(getarraywrite, NULL),
875:   PetscDesignatedInitializer(restorearraywrite, NULL),
876:   PetscDesignatedInitializer(getarrayandmemtype, NULL),
877:   PetscDesignatedInitializer(restorearrayandmemtype, NULL),
878:   PetscDesignatedInitializer(getarrayreadandmemtype, NULL),
879:   PetscDesignatedInitializer(restorearrayreadandmemtype, NULL),
880:   PetscDesignatedInitializer(getarraywriteandmemtype, NULL),
881:   PetscDesignatedInitializer(restorearraywriteandmemtype, NULL),
882:   PetscDesignatedInitializer(concatenate, NULL),
883:   PetscDesignatedInitializer(sum, NULL),
884:   PetscDesignatedInitializer(setpreallocationcoo, VecSetPreallocationCOO_Seq),
885:   PetscDesignatedInitializer(setvaluescoo, VecSetValuesCOO_Seq),
886:   PetscDesignatedInitializer(errorwnorm, NULL),
887: };

889: /*
890:   Create a VECSEQ with the given layout and array

892:   Input Parameter:
893: + map   - the layout
894: - array - the array on host

896:   Output Parameter:
897: . V  - The vector object
898: */
899: PetscErrorCode VecCreateSeqWithLayoutAndArray_Private(PetscLayout map, const PetscScalar array[], Vec *V)
900: {
901:   PetscMPIInt size;

903:   PetscFunctionBegin;
904:   PetscCall(VecCreateWithLayout_Private(map, V));
905:   PetscCallMPI(MPI_Comm_size(map->comm, &size));
906:   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQ on more than one process");
907:   PetscCall(VecCreate_Seq_Private(*V, array));
908:   PetscFunctionReturn(PETSC_SUCCESS);
909: }

911: /*
912:       This is called by VecCreate_Seq() (i.e. VecCreateSeq()) and VecCreateSeqWithArray()
913: */
914: PetscErrorCode VecCreate_Seq_Private(Vec v, const PetscScalar array[])
915: {
916:   Vec_Seq  *s;
917:   PetscBool mdot_use_gemv  = PETSC_TRUE;
918:   PetscBool maxpy_use_gemv = PETSC_FALSE; // default is false as we saw bad performance with vendors' GEMV with tall skinny matrices.

920:   PetscFunctionBegin;
921:   PetscCall(PetscNew(&s));
922:   v->ops[0] = DvOps;

924:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_mdot_use_gemv", &mdot_use_gemv, NULL));
925:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_maxpy_use_gemv", &maxpy_use_gemv, NULL));

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

930:   if (mdot_use_gemv) {
931:     v->ops[0].mdot        = VecMDot_Seq_GEMV;
932:     v->ops[0].mdot_local  = VecMDot_Seq_GEMV;
933:     v->ops[0].mtdot       = VecMTDot_Seq_GEMV;
934:     v->ops[0].mtdot_local = VecMTDot_Seq_GEMV;
935:   }
936:   if (maxpy_use_gemv) v->ops[0].maxpy = VecMAXPY_Seq_GEMV;

938:   v->data            = (void *)s;
939:   v->petscnative     = PETSC_TRUE;
940:   s->array           = (PetscScalar *)array;
941:   s->array_allocated = NULL;
942:   if (array) v->offloadmask = PETSC_OFFLOAD_CPU;

944:   PetscCall(PetscLayoutSetUp(v->map));
945:   PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQ));
946: #if defined(PETSC_HAVE_MATLAB)
947:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEnginePut_C", VecMatlabEnginePut_Default));
948:   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscMatlabEngineGet_C", VecMatlabEngineGet_Default));
949: #endif
950:   PetscFunctionReturn(PETSC_SUCCESS);
951: }

953: /*@C
954:   VecCreateSeqWithArray - Creates a standard,sequential array-style vector,
955:   where the user provides the array space to store the vector values.

957:   Collective

959:   Input Parameters:
960: + comm  - the communicator, should be `PETSC_COMM_SELF`
961: . bs    - the block size
962: . n     - the vector length
963: - array - memory where the vector elements are to be stored.

965:   Output Parameter:
966: . V - the vector

968:   Level: intermediate

970:   Notes:
971:   Use `VecDuplicate()` or `VecDuplicateVecs(`) to form additional vectors of the
972:   same type as an existing vector.

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

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

980: .seealso: `VecCreateMPIWithArray()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`,
981:           `VecCreateGhost()`, `VecCreateSeq()`, `VecPlaceArray()`
982: @*/
983: PetscErrorCode VecCreateSeqWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar array[], Vec *V)
984: {
985:   PetscMPIInt size;

987:   PetscFunctionBegin;
988:   PetscCall(VecCreate(comm, V));
989:   PetscCall(VecSetSizes(*V, n, n));
990:   PetscCall(VecSetBlockSize(*V, bs));
991:   PetscCallMPI(MPI_Comm_size(comm, &size));
992:   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQ on more than one process");
993:   PetscCall(VecCreate_Seq_Private(*V, array));
994:   PetscFunctionReturn(PETSC_SUCCESS);
995: }