Actual source code: vinv.c

  1: /*
  2:      Some useful vector utility functions.
  3: */
  4: #include <../src/vec/vec/impls/mpi/pvecimpl.h>

  6: /*@
  7:   VecStrideSet - Sets a subvector of a vector defined
  8:   by a starting point and a stride with a given value

 10:   Logically Collective

 12:   Input Parameters:
 13: + v     - the vector
 14: . start - starting point of the subvector (defined by a stride)
 15: - s     - value to set for each entry in that subvector

 17:   Level: advanced

 19:   Notes:
 20:   One must call `VecSetBlockSize()` before this routine to set the stride
 21:   information, or use a vector created from a multicomponent `DMDA`.

 23:   This will only work if the desire subvector is a stride subvector

 25: .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideScale()`
 26: @*/
 27: PetscErrorCode VecStrideSet(Vec v, PetscInt start, PetscScalar s)
 28: {
 29:   PetscInt     i, n, bs;
 30:   PetscScalar *x;

 32:   PetscFunctionBegin;
 35:   PetscCall(VecGetLocalSize(v, &n));
 36:   PetscCall(VecGetBlockSize(v, &bs));
 37:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
 38:   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
 39:   PetscCall(VecGetArray(v, &x));
 40:   for (i = start; i < n; i += bs) x[i] = s;
 41:   PetscCall(VecRestoreArray(v, &x));
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: /*@
 46:   VecStrideScale - Scales a subvector of a vector defined
 47:   by a starting point and a stride.

 49:   Logically Collective

 51:   Input Parameters:
 52: + v     - the vector
 53: . start - starting point of the subvector (defined by a stride)
 54: - scale - value to multiply each subvector entry by

 56:   Level: advanced

 58:   Notes:
 59:   One must call `VecSetBlockSize()` before this routine to set the stride
 60:   information, or use a vector created from a multicomponent `DMDA`.

 62:   This will only work if the desire subvector is a stride subvector

 64: .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
 65: @*/
 66: PetscErrorCode VecStrideScale(Vec v, PetscInt start, PetscScalar scale)
 67: {
 68:   PetscInt     i, n, bs;
 69:   PetscScalar *x;

 71:   PetscFunctionBegin;
 75:   PetscCall(VecGetLocalSize(v, &n));
 76:   PetscCall(VecGetBlockSize(v, &bs));
 77:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
 78:   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
 79:   PetscCall(VecGetArray(v, &x));
 80:   for (i = start; i < n; i += bs) x[i] *= scale;
 81:   PetscCall(VecRestoreArray(v, &x));
 82:   PetscFunctionReturn(PETSC_SUCCESS);
 83: }

 85: /*@
 86:   VecStrideNorm - Computes the norm of subvector of a vector defined
 87:   by a starting point and a stride.

 89:   Collective

 91:   Input Parameters:
 92: + v     - the vector
 93: . start - starting point of the subvector (defined by a stride)
 94: - ntype - type of norm, one of `NORM_1`, `NORM_2`, `NORM_INFINITY`

 96:   Output Parameter:
 97: . nrm - the norm

 99:   Level: advanced

101:   Notes:
102:   One must call `VecSetBlockSize()` before this routine to set the stride
103:   information, or use a vector created from a multicomponent `DMDA`.

105:   If x is the array representing the vector x then this computes the norm
106:   of the array (x[start],x[start+stride],x[start+2*stride], ....)

108:   This is useful for computing, say the norm of the pressure variable when
109:   the pressure is stored (interlaced) with other variables, say density etc.

111:   This will only work if the desire subvector is a stride subvector

113: .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
114: @*/
115: PetscErrorCode VecStrideNorm(Vec v, PetscInt start, NormType ntype, PetscReal *nrm)
116: {
117:   PetscInt           i, n, bs;
118:   const PetscScalar *x;
119:   PetscReal          tnorm;

121:   PetscFunctionBegin;
125:   PetscAssertPointer(nrm, 4);
126:   PetscCall(VecGetLocalSize(v, &n));
127:   PetscCall(VecGetBlockSize(v, &bs));
128:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
129:   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
130:   PetscCall(VecGetArrayRead(v, &x));
131:   if (ntype == NORM_2) {
132:     PetscScalar sum = 0.0;
133:     for (i = start; i < n; i += bs) sum += x[i] * (PetscConj(x[i]));
134:     tnorm = PetscRealPart(sum);
135:     PetscCallMPI(MPIU_Allreduce(&tnorm, nrm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)v)));
136:     *nrm = PetscSqrtReal(*nrm);
137:   } else if (ntype == NORM_1) {
138:     tnorm = 0.0;
139:     for (i = start; i < n; i += bs) tnorm += PetscAbsScalar(x[i]);
140:     PetscCallMPI(MPIU_Allreduce(&tnorm, nrm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)v)));
141:   } else if (ntype == NORM_INFINITY) {
142:     tnorm = 0.0;
143:     for (i = start; i < n; i += bs) {
144:       if (PetscAbsScalar(x[i]) > tnorm) tnorm = PetscAbsScalar(x[i]);
145:     }
146:     PetscCallMPI(MPIU_Allreduce(&tnorm, nrm, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)v)));
147:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
148:   PetscCall(VecRestoreArrayRead(v, &x));
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: /*@
153:   VecStrideMax - Computes the maximum of subvector of a vector defined
154:   by a starting point and a stride and optionally its location.

156:   Collective

158:   Input Parameters:
159: + v     - the vector
160: - start - starting point of the subvector (defined by a stride)

162:   Output Parameters:
163: + idex - the location where the maximum occurred  (pass `NULL` if not required)
164: - nrm  - the maximum value in the subvector

166:   Level: advanced

168:   Notes:
169:   One must call `VecSetBlockSize()` before this routine to set the stride
170:   information, or use a vector created from a multicomponent `DMDA`.

172:   If xa is the array representing the vector x, then this computes the max
173:   of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)

175:   This is useful for computing, say the maximum of the pressure variable when
176:   the pressure is stored (interlaced) with other variables, e.g., density, etc.
177:   This will only work if the desire subvector is a stride subvector.

179: .seealso: `Vec`, `VecMax()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`
180: @*/
181: PetscErrorCode VecStrideMax(Vec v, PetscInt start, PetscInt *idex, PetscReal *nrm)
182: {
183:   PetscInt           i, n, bs, id = -1;
184:   const PetscScalar *x;
185:   PetscReal          max = PETSC_MIN_REAL;

187:   PetscFunctionBegin;
190:   PetscAssertPointer(nrm, 4);
191:   PetscCall(VecGetLocalSize(v, &n));
192:   PetscCall(VecGetBlockSize(v, &bs));
193:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
194:   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
195:   PetscCall(VecGetArrayRead(v, &x));
196:   for (i = start; i < n; i += bs) {
197:     if (PetscRealPart(x[i]) > max) {
198:       max = PetscRealPart(x[i]);
199:       id  = i;
200:     }
201:   }
202:   PetscCall(VecRestoreArrayRead(v, &x));
203: #if defined(PETSC_HAVE_MPIUNI)
204:   *nrm = max;
205:   if (idex) *idex = id;
206: #else
207:   if (!idex) {
208:     PetscCallMPI(MPIU_Allreduce(&max, nrm, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)v)));
209:   } else {
210:     struct {
211:       PetscReal v;
212:       PetscInt  i;
213:     } in, out;
214:     PetscInt rstart;

216:     PetscCall(VecGetOwnershipRange(v, &rstart, NULL));
217:     in.v = max;
218:     in.i = rstart + id;
219:     PetscCallMPI(MPIU_Allreduce(&in, &out, 1, MPIU_REAL_INT, MPIU_MAXLOC, PetscObjectComm((PetscObject)v)));
220:     *nrm  = out.v;
221:     *idex = out.i;
222:   }
223: #endif
224:   PetscFunctionReturn(PETSC_SUCCESS);
225: }

227: /*@
228:   VecStrideMin - Computes the minimum of subvector of a vector defined
229:   by a starting point and a stride and optionally its location.

231:   Collective

233:   Input Parameters:
234: + v     - the vector
235: - start - starting point of the subvector (defined by a stride)

237:   Output Parameters:
238: + idex - the location where the minimum occurred. (pass `NULL` if not required)
239: - nrm  - the minimum value in the subvector

241:   Level: advanced

243:   Notes:
244:   One must call `VecSetBlockSize()` before this routine to set the stride
245:   information, or use a vector created from a multicomponent `DMDA`.

247:   If xa is the array representing the vector x, then this computes the min
248:   of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)

250:   This is useful for computing, say the minimum of the pressure variable when
251:   the pressure is stored (interlaced) with other variables, e.g., density, etc.
252:   This will only work if the desire subvector is a stride subvector.

254: .seealso: `Vec`, `VecMin()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMax()`
255: @*/
256: PetscErrorCode VecStrideMin(Vec v, PetscInt start, PetscInt *idex, PetscReal *nrm)
257: {
258:   PetscInt           i, n, bs, id = -1;
259:   const PetscScalar *x;
260:   PetscReal          min = PETSC_MAX_REAL;

262:   PetscFunctionBegin;
265:   PetscAssertPointer(nrm, 4);
266:   PetscCall(VecGetLocalSize(v, &n));
267:   PetscCall(VecGetBlockSize(v, &bs));
268:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
269:   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
270:   PetscCall(VecGetArrayRead(v, &x));
271:   for (i = start; i < n; i += bs) {
272:     if (PetscRealPart(x[i]) < min) {
273:       min = PetscRealPart(x[i]);
274:       id  = i;
275:     }
276:   }
277:   PetscCall(VecRestoreArrayRead(v, &x));
278: #if defined(PETSC_HAVE_MPIUNI)
279:   *nrm = min;
280:   if (idex) *idex = id;
281: #else
282:   if (!idex) {
283:     PetscCallMPI(MPIU_Allreduce(&min, nrm, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)v)));
284:   } else {
285:     struct {
286:       PetscReal v;
287:       PetscInt  i;
288:     } in, out;
289:     PetscInt rstart;

291:     PetscCall(VecGetOwnershipRange(v, &rstart, NULL));
292:     in.v = min;
293:     in.i = rstart + id;
294:     PetscCallMPI(MPIU_Allreduce(&in, &out, 1, MPIU_REAL_INT, MPIU_MINLOC, PetscObjectComm((PetscObject)v)));
295:     *nrm  = out.v;
296:     *idex = out.i;
297:   }
298: #endif
299:   PetscFunctionReturn(PETSC_SUCCESS);
300: }

302: /*@
303:   VecStrideSum - Computes the sum of subvector of a vector defined
304:   by a starting point and a stride.

306:   Collective

308:   Input Parameters:
309: + v     - the vector
310: - start - starting point of the subvector (defined by a stride)

312:   Output Parameter:
313: . sum - the sum

315:   Level: advanced

317:   Notes:
318:   One must call `VecSetBlockSize()` before this routine to set the stride
319:   information, or use a vector created from a multicomponent `DMDA`.

321:   If x is the array representing the vector x then this computes the sum
322:   of the array (x[start],x[start+stride],x[start+2*stride], ....)

324: .seealso: `Vec`, `VecSum()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
325: @*/
326: PetscErrorCode VecStrideSum(Vec v, PetscInt start, PetscScalar *sum)
327: {
328:   PetscInt           i, n, bs;
329:   const PetscScalar *x;
330:   PetscScalar        local_sum = 0.0;

332:   PetscFunctionBegin;
335:   PetscAssertPointer(sum, 3);
336:   PetscCall(VecGetLocalSize(v, &n));
337:   PetscCall(VecGetBlockSize(v, &bs));
338:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
339:   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
340:   PetscCall(VecGetArrayRead(v, &x));
341:   for (i = start; i < n; i += bs) local_sum += x[i];
342:   PetscCallMPI(MPIU_Allreduce(&local_sum, sum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)v)));
343:   PetscCall(VecRestoreArrayRead(v, &x));
344:   PetscFunctionReturn(PETSC_SUCCESS);
345: }

347: /*@
348:   VecStrideScaleAll - Scales the subvectors of a vector defined
349:   by a starting point and a stride.

351:   Logically Collective

353:   Input Parameters:
354: + v      - the vector
355: - scales - values to multiply each subvector entry by

357:   Level: advanced

359:   Notes:
360:   One must call `VecSetBlockSize()` before this routine to set the stride
361:   information, or use a vector created from a multicomponent `DMDA`.

363:   The dimension of scales must be the same as the vector block size

365: .seealso: `Vec`, `VecNorm()`, `VecStrideScale()`, `VecScale()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
366: @*/
367: PetscErrorCode VecStrideScaleAll(Vec v, const PetscScalar *scales)
368: {
369:   PetscInt     i, j, n, bs;
370:   PetscScalar *x;

372:   PetscFunctionBegin;
374:   PetscAssertPointer(scales, 2);
375:   PetscCall(VecGetLocalSize(v, &n));
376:   PetscCall(VecGetBlockSize(v, &bs));
377:   PetscCall(VecGetArray(v, &x));
378:   /* need to provide optimized code for each bs */
379:   for (i = 0; i < n; i += bs) {
380:     for (j = 0; j < bs; j++) x[i + j] *= scales[j];
381:   }
382:   PetscCall(VecRestoreArray(v, &x));
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }

386: /*@
387:   VecStrideNormAll - Computes the norms of subvectors of a vector defined
388:   by a starting point and a stride.

390:   Collective

392:   Input Parameters:
393: + v     - the vector
394: - ntype - type of norm, one of `NORM_1`, `NORM_2`, `NORM_INFINITY`

396:   Output Parameter:
397: . nrm - the norms

399:   Level: advanced

401:   Notes:
402:   One must call `VecSetBlockSize()` before this routine to set the stride
403:   information, or use a vector created from a multicomponent `DMDA`.

405:   If x is the array representing the vector x then this computes the norm
406:   of the array (x[start],x[start+stride],x[start+2*stride], ....) for each start < stride

408:   The dimension of nrm must be the same as the vector block size

410:   This will only work if the desire subvector is a stride subvector

412: .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
413: @*/
414: PetscErrorCode VecStrideNormAll(Vec v, NormType ntype, PetscReal nrm[])
415: {
416:   PetscInt           i, j, n, bs;
417:   const PetscScalar *x;
418:   PetscReal          tnorm[128];
419:   MPI_Comm           comm;
420:   PetscMPIInt        ibs;

422:   PetscFunctionBegin;
425:   PetscAssertPointer(nrm, 3);
426:   PetscCall(VecGetLocalSize(v, &n));
427:   PetscCall(VecGetArrayRead(v, &x));
428:   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));

430:   PetscCall(VecGetBlockSize(v, &bs));
431:   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
432:   PetscCall(PetscMPIIntCast(bs, &ibs));
433:   if (ntype == NORM_2) {
434:     PetscScalar sum[128];
435:     for (j = 0; j < bs; j++) sum[j] = 0.0;
436:     for (i = 0; i < n; i += bs) {
437:       for (j = 0; j < bs; j++) sum[j] += x[i + j] * (PetscConj(x[i + j]));
438:     }
439:     for (j = 0; j < bs; j++) tnorm[j] = PetscRealPart(sum[j]);

441:     PetscCallMPI(MPIU_Allreduce(tnorm, nrm, ibs, MPIU_REAL, MPIU_SUM, comm));
442:     for (j = 0; j < bs; j++) nrm[j] = PetscSqrtReal(nrm[j]);
443:   } else if (ntype == NORM_1) {
444:     for (j = 0; j < bs; j++) tnorm[j] = 0.0;

446:     for (i = 0; i < n; i += bs) {
447:       for (j = 0; j < bs; j++) tnorm[j] += PetscAbsScalar(x[i + j]);
448:     }

450:     PetscCallMPI(MPIU_Allreduce(tnorm, nrm, ibs, MPIU_REAL, MPIU_SUM, comm));
451:   } else if (ntype == NORM_INFINITY) {
452:     PetscReal tmp;
453:     for (j = 0; j < bs; j++) tnorm[j] = 0.0;

455:     for (i = 0; i < n; i += bs) {
456:       for (j = 0; j < bs; j++) {
457:         if ((tmp = PetscAbsScalar(x[i + j])) > tnorm[j]) tnorm[j] = tmp;
458:         /* check special case of tmp == NaN */
459:         if (tmp != tmp) {
460:           tnorm[j] = tmp;
461:           break;
462:         }
463:       }
464:     }
465:     PetscCallMPI(MPIU_Allreduce(tnorm, nrm, ibs, MPIU_REAL, MPIU_MAX, comm));
466:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
467:   PetscCall(VecRestoreArrayRead(v, &x));
468:   PetscFunctionReturn(PETSC_SUCCESS);
469: }

471: /*@
472:   VecStrideMaxAll - Computes the maximums of subvectors of a vector defined
473:   by a starting point and a stride and optionally its location.

475:   Collective

477:   Input Parameter:
478: . v - the vector

480:   Output Parameters:
481: + idex - the location where the maximum occurred (not supported, pass `NULL`,
482:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
483: - nrm  - the maximum values of each subvector

485:   Level: advanced

487:   Notes:
488:   One must call `VecSetBlockSize()` before this routine to set the stride
489:   information, or use a vector created from a multicomponent `DMDA`.

491:   The dimension of nrm must be the same as the vector block size

493: .seealso: `Vec`, `VecMax()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`
494: @*/
495: PetscErrorCode VecStrideMaxAll(Vec v, PetscInt idex[], PetscReal nrm[])
496: {
497:   PetscInt           i, j, n, bs;
498:   const PetscScalar *x;
499:   PetscReal          max[128], tmp;
500:   MPI_Comm           comm;
501:   PetscMPIInt        ibs;

503:   PetscFunctionBegin;
505:   PetscAssertPointer(nrm, 3);
506:   PetscCheck(!idex, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
507:   PetscCall(VecGetLocalSize(v, &n));
508:   PetscCall(VecGetArrayRead(v, &x));
509:   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));

511:   PetscCall(VecGetBlockSize(v, &bs));
512:   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
513:   PetscCall(PetscMPIIntCast(bs, &ibs));

515:   if (!n) {
516:     for (j = 0; j < bs; j++) max[j] = PETSC_MIN_REAL;
517:   } else {
518:     for (j = 0; j < bs; j++) max[j] = PetscRealPart(x[j]);

520:     for (i = bs; i < n; i += bs) {
521:       for (j = 0; j < bs; j++) {
522:         if ((tmp = PetscRealPart(x[i + j])) > max[j]) max[j] = tmp;
523:       }
524:     }
525:   }
526:   PetscCallMPI(MPIU_Allreduce(max, nrm, ibs, MPIU_REAL, MPIU_MAX, comm));

528:   PetscCall(VecRestoreArrayRead(v, &x));
529:   PetscFunctionReturn(PETSC_SUCCESS);
530: }

532: /*@
533:   VecStrideMinAll - Computes the minimum of subvector of a vector defined
534:   by a starting point and a stride and optionally its location.

536:   Collective

538:   Input Parameter:
539: . v - the vector

541:   Output Parameters:
542: + idex - the location where the minimum occurred (not supported, pass `NULL`,
543:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
544: - nrm  - the minimums of each subvector

546:   Level: advanced

548:   Notes:
549:   One must call `VecSetBlockSize()` before this routine to set the stride
550:   information, or use a vector created from a multicomponent `DMDA`.

552:   The dimension of `nrm` must be the same as the vector block size

554: .seealso: `Vec`, `VecMin()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMax()`
555: @*/
556: PetscErrorCode VecStrideMinAll(Vec v, PetscInt idex[], PetscReal nrm[])
557: {
558:   PetscInt           i, n, bs, j;
559:   const PetscScalar *x;
560:   PetscReal          min[128], tmp;
561:   MPI_Comm           comm;
562:   PetscMPIInt        ibs;

564:   PetscFunctionBegin;
566:   PetscAssertPointer(nrm, 3);
567:   PetscCheck(!idex, PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
568:   PetscCall(VecGetLocalSize(v, &n));
569:   PetscCall(VecGetArrayRead(v, &x));
570:   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));

572:   PetscCall(VecGetBlockSize(v, &bs));
573:   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
574:   PetscCall(PetscMPIIntCast(bs, &ibs));

576:   if (!n) {
577:     for (j = 0; j < bs; j++) min[j] = PETSC_MAX_REAL;
578:   } else {
579:     for (j = 0; j < bs; j++) min[j] = PetscRealPart(x[j]);

581:     for (i = bs; i < n; i += bs) {
582:       for (j = 0; j < bs; j++) {
583:         if ((tmp = PetscRealPart(x[i + j])) < min[j]) min[j] = tmp;
584:       }
585:     }
586:   }
587:   PetscCallMPI(MPIU_Allreduce(min, nrm, ibs, MPIU_REAL, MPIU_MIN, comm));

589:   PetscCall(VecRestoreArrayRead(v, &x));
590:   PetscFunctionReturn(PETSC_SUCCESS);
591: }

593: /*@
594:   VecStrideSumAll - Computes the sums of subvectors of a vector defined by a stride.

596:   Collective

598:   Input Parameter:
599: . v - the vector

601:   Output Parameter:
602: . sums - the sums

604:   Level: advanced

606:   Notes:
607:   One must call `VecSetBlockSize()` before this routine to set the stride
608:   information, or use a vector created from a multicomponent `DMDA`.

610:   If x is the array representing the vector x then this computes the sum
611:   of the array (x[start],x[start+stride],x[start+2*stride], ....) for each start < stride

613: .seealso: `Vec`, `VecSum()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
614: @*/
615: PetscErrorCode VecStrideSumAll(Vec v, PetscScalar sums[])
616: {
617:   PetscInt           i, j, n, bs;
618:   const PetscScalar *x;
619:   PetscScalar        local_sums[128];
620:   MPI_Comm           comm;
621:   PetscMPIInt        ibs;

623:   PetscFunctionBegin;
625:   PetscAssertPointer(sums, 2);
626:   PetscCall(VecGetLocalSize(v, &n));
627:   PetscCall(VecGetArrayRead(v, &x));
628:   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));

630:   PetscCall(VecGetBlockSize(v, &bs));
631:   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
632:   PetscCall(PetscMPIIntCast(bs, &ibs));

634:   for (j = 0; j < bs; j++) local_sums[j] = 0.0;
635:   for (i = 0; i < n; i += bs) {
636:     for (j = 0; j < bs; j++) local_sums[j] += x[i + j];
637:   }
638:   PetscCallMPI(MPIU_Allreduce(local_sums, sums, ibs, MPIU_SCALAR, MPIU_SUM, comm));

640:   PetscCall(VecRestoreArrayRead(v, &x));
641:   PetscFunctionReturn(PETSC_SUCCESS);
642: }

644: /*@
645:   VecStrideGatherAll - Gathers all the single components from a multi-component vector into
646:   separate vectors.

648:   Collective

650:   Input Parameters:
651: + v    - the vector
652: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`

654:   Output Parameter:
655: . s - the location where the subvectors are stored

657:   Level: advanced

659:   Notes:
660:   One must call `VecSetBlockSize()` before this routine to set the stride
661:   information, or use a vector created from a multicomponent `DMDA`.

663:   If x is the array representing the vector x then this gathers
664:   the arrays (x[start],x[start+stride],x[start+2*stride], ....)
665:   for start=0,1,2,...bs-1

667:   The parallel layout of the vector and the subvector must be the same;
668:   i.e., nlocal of v = stride*(nlocal of s)

670:   Not optimized; could be easily

672: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGather()`,
673:           `VecStrideScatterAll()`
674: @*/
675: PetscErrorCode VecStrideGatherAll(Vec v, Vec s[], InsertMode addv)
676: {
677:   PetscInt           i, n, n2, bs, j, k, *bss = NULL, nv, jj, nvc;
678:   PetscScalar      **y;
679:   const PetscScalar *x;

681:   PetscFunctionBegin;
683:   PetscAssertPointer(s, 2);
685:   PetscCall(VecGetLocalSize(v, &n));
686:   PetscCall(VecGetLocalSize(s[0], &n2));
687:   PetscCall(VecGetArrayRead(v, &x));
688:   PetscCall(VecGetBlockSize(v, &bs));
689:   PetscCheck(bs > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Input vector does not have a valid blocksize set");

691:   PetscCall(PetscMalloc2(bs, &y, bs, &bss));
692:   nv  = 0;
693:   nvc = 0;
694:   for (i = 0; i < bs; i++) {
695:     PetscCall(VecGetBlockSize(s[i], &bss[i]));
696:     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
697:     PetscCall(VecGetArray(s[i], &y[i]));
698:     nvc += bss[i];
699:     nv++;
700:     PetscCheck(nvc <= bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of subvectors in subvectors > number of vectors in main vector");
701:     if (nvc == bs) break;
702:   }

704:   n = n / bs;

706:   jj = 0;
707:   if (addv == INSERT_VALUES) {
708:     for (j = 0; j < nv; j++) {
709:       for (k = 0; k < bss[j]; k++) {
710:         for (i = 0; i < n; i++) y[j][i * bss[j] + k] = x[bs * i + jj + k];
711:       }
712:       jj += bss[j];
713:     }
714:   } else if (addv == ADD_VALUES) {
715:     for (j = 0; j < nv; j++) {
716:       for (k = 0; k < bss[j]; k++) {
717:         for (i = 0; i < n; i++) y[j][i * bss[j] + k] += x[bs * i + jj + k];
718:       }
719:       jj += bss[j];
720:     }
721: #if !defined(PETSC_USE_COMPLEX)
722:   } else if (addv == MAX_VALUES) {
723:     for (j = 0; j < nv; j++) {
724:       for (k = 0; k < bss[j]; k++) {
725:         for (i = 0; i < n; i++) y[j][i * bss[j] + k] = PetscMax(y[j][i * bss[j] + k], x[bs * i + jj + k]);
726:       }
727:       jj += bss[j];
728:     }
729: #endif
730:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");

732:   PetscCall(VecRestoreArrayRead(v, &x));
733:   for (i = 0; i < nv; i++) PetscCall(VecRestoreArray(s[i], &y[i]));

735:   PetscCall(PetscFree2(y, bss));
736:   PetscFunctionReturn(PETSC_SUCCESS);
737: }

739: /*@
740:   VecStrideScatterAll - Scatters all the single components from separate vectors into
741:   a multi-component vector.

743:   Collective

745:   Input Parameters:
746: + s    - the location where the subvectors are stored
747: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`

749:   Output Parameter:
750: . v - the multicomponent vector

752:   Level: advanced

754:   Notes:
755:   One must call `VecSetBlockSize()` before this routine to set the stride
756:   information, or use a vector created from a multicomponent `DMDA`.

758:   The parallel layout of the vector and the subvector must be the same;
759:   i.e., nlocal of v = stride*(nlocal of s)

761:   Not optimized; could be easily

763: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGather()`
764: @*/
765: PetscErrorCode VecStrideScatterAll(Vec s[], Vec v, InsertMode addv)
766: {
767:   PetscInt            i, n, n2, bs, j, jj, k, *bss = NULL, nv, nvc;
768:   PetscScalar        *x;
769:   PetscScalar const **y;

771:   PetscFunctionBegin;
773:   PetscAssertPointer(s, 1);
775:   PetscCall(VecGetLocalSize(v, &n));
776:   PetscCall(VecGetLocalSize(s[0], &n2));
777:   PetscCall(VecGetArray(v, &x));
778:   PetscCall(VecGetBlockSize(v, &bs));
779:   PetscCheck(bs > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Input vector does not have a valid blocksize set");

781:   PetscCall(PetscMalloc2(bs, (PetscScalar ***)&y, bs, &bss));
782:   nv  = 0;
783:   nvc = 0;
784:   for (i = 0; i < bs; i++) {
785:     PetscCall(VecGetBlockSize(s[i], &bss[i]));
786:     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
787:     PetscCall(VecGetArrayRead(s[i], &y[i]));
788:     nvc += bss[i];
789:     nv++;
790:     PetscCheck(nvc <= bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of subvectors in subvectors > number of vectors in main vector");
791:     if (nvc == bs) break;
792:   }

794:   n = n / bs;

796:   jj = 0;
797:   if (addv == INSERT_VALUES) {
798:     for (j = 0; j < nv; j++) {
799:       for (k = 0; k < bss[j]; k++) {
800:         for (i = 0; i < n; i++) x[bs * i + jj + k] = y[j][i * bss[j] + k];
801:       }
802:       jj += bss[j];
803:     }
804:   } else if (addv == ADD_VALUES) {
805:     for (j = 0; j < nv; j++) {
806:       for (k = 0; k < bss[j]; k++) {
807:         for (i = 0; i < n; i++) x[bs * i + jj + k] += y[j][i * bss[j] + k];
808:       }
809:       jj += bss[j];
810:     }
811: #if !defined(PETSC_USE_COMPLEX)
812:   } else if (addv == MAX_VALUES) {
813:     for (j = 0; j < nv; j++) {
814:       for (k = 0; k < bss[j]; k++) {
815:         for (i = 0; i < n; i++) x[bs * i + jj + k] = PetscMax(x[bs * i + jj + k], y[j][i * bss[j] + k]);
816:       }
817:       jj += bss[j];
818:     }
819: #endif
820:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");

822:   PetscCall(VecRestoreArray(v, &x));
823:   for (i = 0; i < nv; i++) PetscCall(VecRestoreArrayRead(s[i], &y[i]));
824:   PetscCall(PetscFree2(*(PetscScalar ***)&y, bss));
825:   PetscFunctionReturn(PETSC_SUCCESS);
826: }

828: /*@
829:   VecStrideGather - Gathers a single component from a multi-component vector into
830:   another vector.

832:   Collective

834:   Input Parameters:
835: + v     - the vector
836: . start - starting point of the subvector (defined by a stride)
837: - addv  - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`

839:   Output Parameter:
840: . s - the location where the subvector is stored

842:   Level: advanced

844:   Notes:
845:   One must call `VecSetBlockSize()` before this routine to set the stride
846:   information, or use a vector created from a multicomponent `DMDA`.

848:   If x is the array representing the vector x then this gathers
849:   the array (x[start],x[start+stride],x[start+2*stride], ....)

851:   The parallel layout of the vector and the subvector must be the same;
852:   i.e., nlocal of v = stride*(nlocal of s)

854:   Not optimized; could be easily

856: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
857:           `VecStrideScatterAll()`
858: @*/
859: PetscErrorCode VecStrideGather(Vec v, PetscInt start, Vec s, InsertMode addv)
860: {
861:   PetscFunctionBegin;
865:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
866:   PetscCheck(start < v->map->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start,
867:              v->map->bs);
868:   PetscUseTypeMethod(v, stridegather, start, s, addv);
869:   PetscFunctionReturn(PETSC_SUCCESS);
870: }

872: /*@
873:   VecStrideScatter - Scatters a single component from a vector into a multi-component vector.

875:   Collective

877:   Input Parameters:
878: + s     - the single-component vector
879: . start - starting point of the subvector (defined by a stride)
880: - addv  - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`

882:   Output Parameter:
883: . v - the location where the subvector is scattered (the multi-component vector)

885:   Level: advanced

887:   Notes:
888:   One must call `VecSetBlockSize()` on the multi-component vector before this
889:   routine to set the stride  information, or use a vector created from a multicomponent `DMDA`.

891:   The parallel layout of the vector and the subvector must be the same;
892:   i.e., nlocal of v = stride*(nlocal of s)

894:   Not optimized; could be easily

896: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
897:           `VecStrideScatterAll()`, `VecStrideSubSetScatter()`, `VecStrideSubSetGather()`
898: @*/
899: PetscErrorCode VecStrideScatter(Vec s, PetscInt start, Vec v, InsertMode addv)
900: {
901:   PetscFunctionBegin;
905:   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
906:   PetscCheck(start < v->map->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride. Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start,
907:              v->map->bs);
908:   PetscCall((*v->ops->stridescatter)(s, start, v, addv));
909:   PetscFunctionReturn(PETSC_SUCCESS);
910: }

912: /*@
913:   VecStrideSubSetGather - Gathers a subset of components from a multi-component vector into
914:   another vector.

916:   Collective

918:   Input Parameters:
919: + v    - the vector
920: . nidx - the number of indices
921: . idxv - the indices of the components 0 <= idxv[0] ...idxv[nidx-1] < bs(v), they need not be sorted
922: . idxs - the indices of the components 0 <= idxs[0] ...idxs[nidx-1] < bs(s), they need not be sorted, may be null if nidx == bs(s) or is `PETSC_DETERMINE`
923: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`

925:   Output Parameter:
926: . s - the location where the subvector is stored

928:   Level: advanced

930:   Notes:
931:   One must call `VecSetBlockSize()` on both vectors before this routine to set the stride
932:   information, or use a vector created from a multicomponent `DMDA`.

934:   The parallel layout of the vector and the subvector must be the same;

936:   Not optimized; could be easily

938: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideGather()`, `VecStrideSubSetScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
939:           `VecStrideScatterAll()`
940: @*/
941: PetscErrorCode VecStrideSubSetGather(Vec v, PetscInt nidx, const PetscInt idxv[], const PetscInt idxs[], Vec s, InsertMode addv)
942: {
943:   PetscFunctionBegin;
946:   if (nidx == PETSC_DETERMINE) nidx = s->map->bs;
947:   PetscUseTypeMethod(v, stridesubsetgather, nidx, idxv, idxs, s, addv);
948:   PetscFunctionReturn(PETSC_SUCCESS);
949: }

951: /*@
952:   VecStrideSubSetScatter - Scatters components from a vector into a subset of components of a multi-component vector.

954:   Collective

956:   Input Parameters:
957: + s    - the smaller-component vector
958: . nidx - the number of indices in idx
959: . idxs - the indices of the components in the smaller-component vector, 0 <= idxs[0] ...idxs[nidx-1] < bs(s) they need not be sorted, may be null if nidx == bs(s) or is `PETSC_DETERMINE`
960: . idxv - the indices of the components in the larger-component vector, 0 <= idx[0] ...idx[nidx-1] < bs(v) they need not be sorted
961: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`

963:   Output Parameter:
964: . v - the location where the subvector is into scattered (the multi-component vector)

966:   Level: advanced

968:   Notes:
969:   One must call `VecSetBlockSize()` on the vectors before this
970:   routine to set the stride  information, or use a vector created from a multicomponent `DMDA`.

972:   The parallel layout of the vector and the subvector must be the same;

974:   Not optimized; could be easily

976: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideSubSetGather()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
977:           `VecStrideScatterAll()`
978: @*/
979: PetscErrorCode VecStrideSubSetScatter(Vec s, PetscInt nidx, const PetscInt idxs[], const PetscInt idxv[], Vec v, InsertMode addv)
980: {
981:   PetscFunctionBegin;
984:   if (nidx == PETSC_DETERMINE) nidx = s->map->bs;
985:   PetscCall((*v->ops->stridesubsetscatter)(s, nidx, idxs, idxv, v, addv));
986:   PetscFunctionReturn(PETSC_SUCCESS);
987: }

989: PetscErrorCode VecStrideGather_Default(Vec v, PetscInt start, Vec s, InsertMode addv)
990: {
991:   PetscInt           i, n, bs, ns;
992:   const PetscScalar *x;
993:   PetscScalar       *y;

995:   PetscFunctionBegin;
996:   PetscCall(VecGetLocalSize(v, &n));
997:   PetscCall(VecGetLocalSize(s, &ns));
998:   PetscCall(VecGetArrayRead(v, &x));
999:   PetscCall(VecGetArray(s, &y));

1001:   bs = v->map->bs;
1002:   PetscCheck(n == ns * bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Subvector length * blocksize %" PetscInt_FMT " not correct for gather from original vector %" PetscInt_FMT, ns * bs, n);
1003:   x += start;
1004:   n = n / bs;

1006:   if (addv == INSERT_VALUES) {
1007:     for (i = 0; i < n; i++) y[i] = x[bs * i];
1008:   } else if (addv == ADD_VALUES) {
1009:     for (i = 0; i < n; i++) y[i] += x[bs * i];
1010: #if !defined(PETSC_USE_COMPLEX)
1011:   } else if (addv == MAX_VALUES) {
1012:     for (i = 0; i < n; i++) y[i] = PetscMax(y[i], x[bs * i]);
1013: #endif
1014:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");

1016:   PetscCall(VecRestoreArrayRead(v, &x));
1017:   PetscCall(VecRestoreArray(s, &y));
1018:   PetscFunctionReturn(PETSC_SUCCESS);
1019: }

1021: PetscErrorCode VecStrideScatter_Default(Vec s, PetscInt start, Vec v, InsertMode addv)
1022: {
1023:   PetscInt           i, n, bs, ns;
1024:   PetscScalar       *x;
1025:   const PetscScalar *y;

1027:   PetscFunctionBegin;
1028:   PetscCall(VecGetLocalSize(v, &n));
1029:   PetscCall(VecGetLocalSize(s, &ns));
1030:   PetscCall(VecGetArray(v, &x));
1031:   PetscCall(VecGetArrayRead(s, &y));

1033:   bs = v->map->bs;
1034:   PetscCheck(n == ns * bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Subvector length * blocksize %" PetscInt_FMT " not correct for scatter to multicomponent vector %" PetscInt_FMT, ns * bs, n);
1035:   x += start;
1036:   n = n / bs;

1038:   if (addv == INSERT_VALUES) {
1039:     for (i = 0; i < n; i++) x[bs * i] = y[i];
1040:   } else if (addv == ADD_VALUES) {
1041:     for (i = 0; i < n; i++) x[bs * i] += y[i];
1042: #if !defined(PETSC_USE_COMPLEX)
1043:   } else if (addv == MAX_VALUES) {
1044:     for (i = 0; i < n; i++) x[bs * i] = PetscMax(y[i], x[bs * i]);
1045: #endif
1046:   } else SETERRQ(PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");

1048:   PetscCall(VecRestoreArray(v, &x));
1049:   PetscCall(VecRestoreArrayRead(s, &y));
1050:   PetscFunctionReturn(PETSC_SUCCESS);
1051: }

1053: PetscErrorCode VecStrideSubSetGather_Default(Vec v, PetscInt nidx, const PetscInt idxv[], const PetscInt idxs[], Vec s, InsertMode addv)
1054: {
1055:   PetscInt           i, j, n, bs, bss, ns;
1056:   const PetscScalar *x;
1057:   PetscScalar       *y;

1059:   PetscFunctionBegin;
1060:   PetscCall(VecGetLocalSize(v, &n));
1061:   PetscCall(VecGetLocalSize(s, &ns));
1062:   PetscCall(VecGetArrayRead(v, &x));
1063:   PetscCall(VecGetArray(s, &y));

1065:   bs  = v->map->bs;
1066:   bss = s->map->bs;
1067:   n   = n / bs;

1069:   if (PetscDefined(USE_DEBUG)) {
1070:     PetscCheck(n == ns / bss, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout of vectors");
1071:     for (j = 0; j < nidx; j++) {
1072:       PetscCheck(idxv[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is negative", j, idxv[j]);
1073:       PetscCheck(idxv[j] < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is greater than or equal to vector blocksize %" PetscInt_FMT, j, idxv[j], bs);
1074:     }
1075:     PetscCheck(idxs || bss == nidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide idxs when not gathering into all locations");
1076:   }

1078:   if (addv == INSERT_VALUES) {
1079:     if (!idxs) {
1080:       for (i = 0; i < n; i++) {
1081:         for (j = 0; j < bss; j++) y[bss * i + j] = x[bs * i + idxv[j]];
1082:       }
1083:     } else {
1084:       for (i = 0; i < n; i++) {
1085:         for (j = 0; j < bss; j++) y[bss * i + idxs[j]] = x[bs * i + idxv[j]];
1086:       }
1087:     }
1088:   } else if (addv == ADD_VALUES) {
1089:     if (!idxs) {
1090:       for (i = 0; i < n; i++) {
1091:         for (j = 0; j < bss; j++) y[bss * i + j] += x[bs * i + idxv[j]];
1092:       }
1093:     } else {
1094:       for (i = 0; i < n; i++) {
1095:         for (j = 0; j < bss; j++) y[bss * i + idxs[j]] += x[bs * i + idxv[j]];
1096:       }
1097:     }
1098: #if !defined(PETSC_USE_COMPLEX)
1099:   } else if (addv == MAX_VALUES) {
1100:     if (!idxs) {
1101:       for (i = 0; i < n; i++) {
1102:         for (j = 0; j < bss; j++) y[bss * i + j] = PetscMax(y[bss * i + j], x[bs * i + idxv[j]]);
1103:       }
1104:     } else {
1105:       for (i = 0; i < n; i++) {
1106:         for (j = 0; j < bss; j++) y[bss * i + idxs[j]] = PetscMax(y[bss * i + idxs[j]], x[bs * i + idxv[j]]);
1107:       }
1108:     }
1109: #endif
1110:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");

1112:   PetscCall(VecRestoreArrayRead(v, &x));
1113:   PetscCall(VecRestoreArray(s, &y));
1114:   PetscFunctionReturn(PETSC_SUCCESS);
1115: }

1117: PetscErrorCode VecStrideSubSetScatter_Default(Vec s, PetscInt nidx, const PetscInt idxs[], const PetscInt idxv[], Vec v, InsertMode addv)
1118: {
1119:   PetscInt           j, i, n, bs, ns, bss;
1120:   PetscScalar       *x;
1121:   const PetscScalar *y;

1123:   PetscFunctionBegin;
1124:   PetscCall(VecGetLocalSize(v, &n));
1125:   PetscCall(VecGetLocalSize(s, &ns));
1126:   PetscCall(VecGetArray(v, &x));
1127:   PetscCall(VecGetArrayRead(s, &y));

1129:   bs  = v->map->bs;
1130:   bss = s->map->bs;
1131:   n   = n / bs;

1133:   if (PetscDefined(USE_DEBUG)) {
1134:     PetscCheck(n == ns / bss, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout of vectors");
1135:     for (j = 0; j < bss; j++) {
1136:       if (idxs) {
1137:         PetscCheck(idxs[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is negative", j, idxs[j]);
1138:         PetscCheck(idxs[j] < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is greater than or equal to vector blocksize %" PetscInt_FMT, j, idxs[j], bs);
1139:       }
1140:     }
1141:     PetscCheck(idxs || bss == nidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide idxs when not scattering from all locations");
1142:   }

1144:   if (addv == INSERT_VALUES) {
1145:     if (!idxs) {
1146:       for (i = 0; i < n; i++) {
1147:         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = y[bss * i + j];
1148:       }
1149:     } else {
1150:       for (i = 0; i < n; i++) {
1151:         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = y[bss * i + idxs[j]];
1152:       }
1153:     }
1154:   } else if (addv == ADD_VALUES) {
1155:     if (!idxs) {
1156:       for (i = 0; i < n; i++) {
1157:         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] += y[bss * i + j];
1158:       }
1159:     } else {
1160:       for (i = 0; i < n; i++) {
1161:         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] += y[bss * i + idxs[j]];
1162:       }
1163:     }
1164: #if !defined(PETSC_USE_COMPLEX)
1165:   } else if (addv == MAX_VALUES) {
1166:     if (!idxs) {
1167:       for (i = 0; i < n; i++) {
1168:         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = PetscMax(y[bss * i + j], x[bs * i + idxv[j]]);
1169:       }
1170:     } else {
1171:       for (i = 0; i < n; i++) {
1172:         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = PetscMax(y[bss * i + idxs[j]], x[bs * i + idxv[j]]);
1173:       }
1174:     }
1175: #endif
1176:   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");

1178:   PetscCall(VecRestoreArray(v, &x));
1179:   PetscCall(VecRestoreArrayRead(s, &y));
1180:   PetscFunctionReturn(PETSC_SUCCESS);
1181: }

1183: static PetscErrorCode VecApplyUnary_Private(Vec v, PetscDeviceContext dctx, const char async_name[], PetscErrorCode (*unary_op)(Vec), PetscScalar (*UnaryFunc)(PetscScalar))
1184: {
1185:   PetscFunctionBegin;
1187:   PetscCall(VecSetErrorIfLocked(v, 1));
1188:   if (dctx) {
1189:     PetscErrorCode (*unary_op_async)(Vec, PetscDeviceContext);

1191:     PetscCall(PetscObjectQueryFunction((PetscObject)v, async_name, &unary_op_async));
1192:     if (unary_op_async) {
1193:       PetscCall((*unary_op_async)(v, dctx));
1194:       PetscFunctionReturn(PETSC_SUCCESS);
1195:     }
1196:   }
1197:   if (unary_op) {
1199:     PetscCall((*unary_op)(v));
1200:   } else {
1201:     PetscInt     n;
1202:     PetscScalar *x;

1205:     PetscCall(VecGetLocalSize(v, &n));
1206:     PetscCall(VecGetArray(v, &x));
1207:     for (PetscInt i = 0; i < n; ++i) x[i] = UnaryFunc(x[i]);
1208:     PetscCall(VecRestoreArray(v, &x));
1209:   }
1210:   PetscFunctionReturn(PETSC_SUCCESS);
1211: }

1213: static PetscScalar ScalarReciprocal_Function(PetscScalar x)
1214: {
1215:   const PetscScalar zero = 0.0;

1217:   return x == zero ? zero : ((PetscScalar)1.0) / x;
1218: }

1220: PetscErrorCode VecReciprocalAsync_Private(Vec v, PetscDeviceContext dctx)
1221: {
1222:   PetscFunctionBegin;
1223:   PetscCall(PetscLogEventBegin(VEC_Reciprocal, v, NULL, NULL, NULL));
1224:   PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Reciprocal), v->ops->reciprocal, ScalarReciprocal_Function));
1225:   PetscCall(PetscLogEventEnd(VEC_Reciprocal, v, NULL, NULL, NULL));
1226:   PetscFunctionReturn(PETSC_SUCCESS);
1227: }

1229: PetscErrorCode VecReciprocal_Default(Vec v)
1230: {
1231:   PetscFunctionBegin;
1232:   PetscCall(VecApplyUnary_Private(v, NULL, NULL, NULL, ScalarReciprocal_Function));
1233:   PetscFunctionReturn(PETSC_SUCCESS);
1234: }

1236: static PetscScalar ScalarExp_Function(PetscScalar x)
1237: {
1238:   return PetscExpScalar(x);
1239: }

1241: PetscErrorCode VecExpAsync_Private(Vec v, PetscDeviceContext dctx)
1242: {
1243:   PetscFunctionBegin;
1245:   PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Exp), v->ops->exp, ScalarExp_Function));
1246:   PetscFunctionReturn(PETSC_SUCCESS);
1247: }

1249: /*@
1250:   VecExp - Replaces each component of a vector by e^x_i

1252:   Not Collective

1254:   Input Parameter:
1255: . v - The vector

1257:   Output Parameter:
1258: . v - The vector of exponents

1260:   Level: beginner

1262: .seealso: `Vec`, `VecLog()`, `VecAbs()`, `VecSqrtAbs()`, `VecReciprocal()`
1263: @*/
1264: PetscErrorCode VecExp(Vec v)
1265: {
1266:   PetscFunctionBegin;
1267:   PetscCall(VecExpAsync_Private(v, NULL));
1268:   PetscFunctionReturn(PETSC_SUCCESS);
1269: }

1271: static PetscScalar ScalarLog_Function(PetscScalar x)
1272: {
1273:   return PetscLogScalar(x);
1274: }

1276: PetscErrorCode VecLogAsync_Private(Vec v, PetscDeviceContext dctx)
1277: {
1278:   PetscFunctionBegin;
1280:   PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Log), v->ops->log, ScalarLog_Function));
1281:   PetscFunctionReturn(PETSC_SUCCESS);
1282: }

1284: /*@
1285:   VecLog - Replaces each component of a vector by log(x_i), the natural logarithm

1287:   Not Collective

1289:   Input Parameter:
1290: . v - The vector

1292:   Output Parameter:
1293: . v - The vector of logs

1295:   Level: beginner

1297: .seealso: `Vec`, `VecExp()`, `VecAbs()`, `VecSqrtAbs()`, `VecReciprocal()`
1298: @*/
1299: PetscErrorCode VecLog(Vec v)
1300: {
1301:   PetscFunctionBegin;
1302:   PetscCall(VecLogAsync_Private(v, NULL));
1303:   PetscFunctionReturn(PETSC_SUCCESS);
1304: }

1306: static PetscScalar ScalarAbs_Function(PetscScalar x)
1307: {
1308:   return PetscAbsScalar(x);
1309: }

1311: PetscErrorCode VecAbsAsync_Private(Vec v, PetscDeviceContext dctx)
1312: {
1313:   PetscFunctionBegin;
1315:   PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Abs), v->ops->abs, ScalarAbs_Function));
1316:   PetscFunctionReturn(PETSC_SUCCESS);
1317: }

1319: /*@
1320:   VecAbs - Replaces every element in a vector with its absolute value.

1322:   Logically Collective

1324:   Input Parameter:
1325: . v - the vector

1327:   Level: intermediate

1329: .seealso: `Vec`, `VecExp()`, `VecSqrtAbs()`, `VecReciprocal()`, `VecLog()`, `VecPointwiseSign()`
1330: @*/
1331: PetscErrorCode VecAbs(Vec v)
1332: {
1333:   PetscFunctionBegin;
1334:   PetscCall(VecAbsAsync_Private(v, NULL));
1335:   PetscFunctionReturn(PETSC_SUCCESS);
1336: }

1338: static PetscScalar ScalarConjugate_Function(PetscScalar x)
1339: {
1340:   return PetscConj(x);
1341: }

1343: PetscErrorCode VecConjugateAsync_Private(Vec v, PetscDeviceContext dctx)
1344: {
1345:   PetscFunctionBegin;
1347:   if (PetscDefined(USE_COMPLEX)) PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Conjugate), v->ops->conjugate, ScalarConjugate_Function));
1348:   PetscFunctionReturn(PETSC_SUCCESS);
1349: }

1351: /*@
1352:   VecConjugate - Conjugates a vector. That is, replace every entry in a vector with its complex conjugate

1354:   Logically Collective

1356:   Input Parameter:
1357: . x - the vector

1359:   Level: intermediate

1361: .seealso: [](ch_vectors), `Vec`, `VecSet()`
1362: @*/
1363: PetscErrorCode VecConjugate(Vec x)
1364: {
1365:   PetscFunctionBegin;
1366:   PetscCall(VecConjugateAsync_Private(x, NULL));
1367:   PetscFunctionReturn(PETSC_SUCCESS);
1368: }

1370: static PetscScalar ScalarSqrtAbs_Function(PetscScalar x)
1371: {
1372:   return PetscSqrtScalar(ScalarAbs_Function(x));
1373: }

1375: PetscErrorCode VecSqrtAbsAsync_Private(Vec v, PetscDeviceContext dctx)
1376: {
1377:   PetscFunctionBegin;
1379:   PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(SqrtAbs), v->ops->sqrt, ScalarSqrtAbs_Function));
1380:   PetscFunctionReturn(PETSC_SUCCESS);
1381: }

1383: /*@
1384:   VecSqrtAbs - Replaces each component of a vector by the square root of its magnitude.

1386:   Not Collective

1388:   Input Parameter:
1389: . v - The vector

1391:   Level: beginner

1393:   Note:
1394:   The actual function is sqrt(|x_i|)

1396: .seealso: `Vec`, `VecLog()`, `VecExp()`, `VecReciprocal()`, `VecAbs()`
1397: @*/
1398: PetscErrorCode VecSqrtAbs(Vec v)
1399: {
1400:   PetscFunctionBegin;
1401:   PetscCall(VecSqrtAbsAsync_Private(v, NULL));
1402:   PetscFunctionReturn(PETSC_SUCCESS);
1403: }

1405: static PetscScalar ScalarImaginaryPart_Function(PetscScalar x)
1406: {
1407:   const PetscReal imag = PetscImaginaryPart(x);

1409: #if PetscDefined(USE_COMPLEX)
1410:   return PetscCMPLX(imag, 0.0);
1411: #else
1412:   return imag;
1413: #endif
1414: }

1416: /*@
1417:   VecImaginaryPart - Replaces a complex vector with its imaginary part

1419:   Collective

1421:   Input Parameter:
1422: . v - the vector

1424:   Level: beginner

1426: .seealso: `Vec`, `VecNorm()`, `VecRealPart()`
1427: @*/
1428: PetscErrorCode VecImaginaryPart(Vec v)
1429: {
1430:   PetscFunctionBegin;
1432:   PetscCall(VecApplyUnary_Private(v, NULL, NULL, NULL, ScalarImaginaryPart_Function));
1433:   PetscFunctionReturn(PETSC_SUCCESS);
1434: }

1436: static PetscScalar ScalarRealPart_Function(PetscScalar x)
1437: {
1438:   const PetscReal real = PetscRealPart(x);

1440: #if PetscDefined(USE_COMPLEX)
1441:   return PetscCMPLX(real, 0.0);
1442: #else
1443:   return real;
1444: #endif
1445: }

1447: /*@
1448:   VecRealPart - Replaces a complex vector with its real part

1450:   Collective

1452:   Input Parameter:
1453: . v - the vector

1455:   Level: beginner

1457: .seealso: `Vec`, `VecNorm()`, `VecImaginaryPart()`
1458: @*/
1459: PetscErrorCode VecRealPart(Vec v)
1460: {
1461:   PetscFunctionBegin;
1463:   PetscCall(VecApplyUnary_Private(v, NULL, NULL, NULL, ScalarRealPart_Function));
1464:   PetscFunctionReturn(PETSC_SUCCESS);
1465: }

1467: /*@
1468:   VecDotNorm2 - computes the inner product of two vectors and the 2-norm squared of the second vector

1470:   Collective

1472:   Input Parameters:
1473: + s - first vector
1474: - t - second vector

1476:   Output Parameters:
1477: + dp - s'conj(t)
1478: - nm - t'conj(t)

1480:   Level: advanced

1482:   Note:
1483:   conj(x) is the complex conjugate of x when x is complex

1485: .seealso: `Vec`, `VecDot()`, `VecNorm()`, `VecDotBegin()`, `VecNormBegin()`, `VecDotEnd()`, `VecNormEnd()`
1486: @*/
1487: PetscErrorCode VecDotNorm2(Vec s, Vec t, PetscScalar *dp, PetscReal *nm)
1488: {
1489:   PetscScalar work[] = {0.0, 0.0};

1491:   PetscFunctionBegin;
1494:   PetscAssertPointer(dp, 3);
1495:   PetscAssertPointer(nm, 4);
1498:   PetscCheckSameTypeAndComm(s, 1, t, 2);
1499:   PetscCheck(s->map->N == t->map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector global lengths");
1500:   PetscCheck(s->map->n == t->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector local lengths");

1502:   PetscCall(PetscLogEventBegin(VEC_DotNorm2, s, t, 0, 0));
1503:   if (s->ops->dotnorm2) {
1504:     PetscUseTypeMethod(s, dotnorm2, t, work, work + 1);
1505:   } else {
1506:     const PetscScalar *sx, *tx;
1507:     PetscInt           n;

1509:     PetscCall(VecGetLocalSize(s, &n));
1510:     PetscCall(VecGetArrayRead(s, &sx));
1511:     PetscCall(VecGetArrayRead(t, &tx));
1512:     for (PetscInt i = 0; i < n; ++i) {
1513:       const PetscScalar txconj = PetscConj(tx[i]);

1515:       work[0] += sx[i] * txconj;
1516:       work[1] += tx[i] * txconj;
1517:     }
1518:     PetscCall(VecRestoreArrayRead(t, &tx));
1519:     PetscCall(VecRestoreArrayRead(s, &sx));
1520:     PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, work, 2, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)s)));
1521:     PetscCall(PetscLogFlops(4.0 * n));
1522:   }
1523:   PetscCall(PetscLogEventEnd(VEC_DotNorm2, s, t, 0, 0));
1524:   *dp = work[0];
1525:   *nm = PetscRealPart(work[1]);
1526:   PetscFunctionReturn(PETSC_SUCCESS);
1527: }

1529: /*@
1530:   VecSum - Computes the sum of all the components of a vector.

1532:   Collective

1534:   Input Parameter:
1535: . v - the vector

1537:   Output Parameter:
1538: . sum - the result

1540:   Level: beginner

1542: .seealso: `Vec`, `VecMean()`, `VecNorm()`
1543: @*/
1544: PetscErrorCode VecSum(Vec v, PetscScalar *sum)
1545: {
1546:   PetscScalar tmp = 0.0;

1548:   PetscFunctionBegin;
1550:   PetscAssertPointer(sum, 2);
1551:   if (v->ops->sum) PetscUseTypeMethod(v, sum, &tmp);
1552:   else {
1553:     const PetscScalar *x;
1554:     PetscInt           n;

1556:     PetscCall(VecGetLocalSize(v, &n));
1557:     PetscCall(VecGetArrayRead(v, &x));
1558:     for (PetscInt i = 0; i < n; ++i) tmp += x[i];
1559:     PetscCall(VecRestoreArrayRead(v, &x));
1560:   }
1561:   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &tmp, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)v)));
1562:   *sum = tmp;
1563:   PetscFunctionReturn(PETSC_SUCCESS);
1564: }

1566: /*@
1567:   VecMean - Computes the arithmetic mean of all the components of a vector.

1569:   Collective

1571:   Input Parameter:
1572: . v - the vector

1574:   Output Parameter:
1575: . mean - the result

1577:   Level: beginner

1579: .seealso: `Vec`, `VecSum()`, `VecNorm()`
1580: @*/
1581: PetscErrorCode VecMean(Vec v, PetscScalar *mean)
1582: {
1583:   PetscInt n;

1585:   PetscFunctionBegin;
1587:   PetscAssertPointer(mean, 2);
1588:   PetscCall(VecGetSize(v, &n));
1589:   PetscCall(VecSum(v, mean));
1590:   *mean /= n;
1591:   PetscFunctionReturn(PETSC_SUCCESS);
1592: }

1594: PetscErrorCode VecShiftAsync_Private(Vec v, PetscScalar shift, PetscDeviceContext dctx)
1595: {
1596:   PetscErrorCode (*shift_async)(Vec, PetscScalar, PetscDeviceContext) = NULL;

1598:   PetscFunctionBegin;
1599:   if (dctx) {
1600:     PetscErrorCode (*shift_async)(Vec, PetscScalar, PetscDeviceContext);

1602:     PetscCall(PetscObjectQueryFunction((PetscObject)v, VecAsyncFnName(Shift), &shift_async));
1603:   }
1604:   if (shift_async) PetscCall((*shift_async)(v, shift, dctx));
1605:   else if (v->ops->shift) PetscUseTypeMethod(v, shift, shift);
1606:   else {
1607:     PetscInt     n;
1608:     PetscScalar *x;

1610:     PetscCall(VecGetLocalSize(v, &n));
1611:     PetscCall(VecGetArray(v, &x));
1612:     for (PetscInt i = 0; i < n; ++i) x[i] += shift;
1613:     PetscCall(VecRestoreArray(v, &x));
1614:     PetscCall(PetscLogFlops(n));
1615:   }
1616:   PetscFunctionReturn(PETSC_SUCCESS);
1617: }

1619: /*@
1620:   VecShift - Shifts all of the components of a vector by computing
1621:   `x[i] = x[i] + shift`.

1623:   Logically Collective

1625:   Input Parameters:
1626: + v     - the vector
1627: - shift - the shift

1629:   Level: intermediate

1631: .seealso: `Vec`, `VecISShift()`
1632: @*/
1633: PetscErrorCode VecShift(Vec v, PetscScalar shift)
1634: {
1635:   PetscFunctionBegin;
1638:   PetscCall(VecSetErrorIfLocked(v, 1));
1639:   if (shift == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
1640:   PetscCall(PetscLogEventBegin(VEC_Shift, v, 0, 0, 0));
1641:   PetscCall(VecShiftAsync_Private(v, shift, NULL));
1642:   PetscCall(PetscLogEventEnd(VEC_Shift, v, 0, 0, 0));
1643:   PetscFunctionReturn(PETSC_SUCCESS);
1644: }

1646: /*@
1647:   VecPermute - Permutes a vector in place using the given ordering.

1649:   Input Parameters:
1650: + x   - The vector
1651: . row - The ordering
1652: - inv - The flag for inverting the permutation

1654:   Level: beginner

1656:   Note:
1657:   This function does not yet support parallel Index Sets with non-local permutations

1659: .seealso: `Vec`, `MatPermute()`
1660: @*/
1661: PetscErrorCode VecPermute(Vec x, IS row, PetscBool inv)
1662: {
1663:   PetscScalar    *array, *newArray;
1664:   const PetscInt *idx;
1665:   PetscInt        i, rstart, rend;

1667:   PetscFunctionBegin;
1670:   PetscCall(VecSetErrorIfLocked(x, 1));
1671:   PetscCall(VecGetOwnershipRange(x, &rstart, &rend));
1672:   PetscCall(ISGetIndices(row, &idx));
1673:   PetscCall(VecGetArray(x, &array));
1674:   PetscCall(PetscMalloc1(x->map->n, &newArray));
1675:   PetscCall(PetscArraycpy(newArray, array, x->map->n));
1676:   if (PetscDefined(USE_DEBUG)) {
1677:     for (i = 0; i < x->map->n; i++) PetscCheck(!(idx[i] < rstart) && !(idx[i] >= rend), PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "Permutation index %" PetscInt_FMT " is out of bounds: %" PetscInt_FMT, i, idx[i]);
1678:   }
1679:   if (!inv) {
1680:     for (i = 0; i < x->map->n; i++) array[i] = newArray[idx[i] - rstart];
1681:   } else {
1682:     for (i = 0; i < x->map->n; i++) array[idx[i] - rstart] = newArray[i];
1683:   }
1684:   PetscCall(VecRestoreArray(x, &array));
1685:   PetscCall(ISRestoreIndices(row, &idx));
1686:   PetscCall(PetscFree(newArray));
1687:   PetscFunctionReturn(PETSC_SUCCESS);
1688: }

1690: /*@
1691:   VecEqual - Compares two vectors. Returns true if the two vectors are either pointing to the same memory buffer,
1692:   or if the two vectors have the same local and global layout as well as bitwise equality of all entries.
1693:   Does NOT take round-off errors into account.

1695:   Collective

1697:   Input Parameters:
1698: + vec1 - the first vector
1699: - vec2 - the second vector

1701:   Output Parameter:
1702: . flg - `PETSC_TRUE` if the vectors are equal; `PETSC_FALSE` otherwise.

1704:   Level: intermediate

1706: .seealso: `Vec`
1707: @*/
1708: PetscErrorCode VecEqual(Vec vec1, Vec vec2, PetscBool *flg)
1709: {
1710:   const PetscScalar *v1, *v2;
1711:   PetscInt           n1, n2, N1, N2;
1712:   PetscBool          flg1;

1714:   PetscFunctionBegin;
1717:   PetscAssertPointer(flg, 3);
1718:   if (vec1 == vec2) *flg = PETSC_TRUE;
1719:   else {
1720:     PetscCall(VecGetSize(vec1, &N1));
1721:     PetscCall(VecGetSize(vec2, &N2));
1722:     if (N1 != N2) flg1 = PETSC_FALSE;
1723:     else {
1724:       PetscCall(VecGetLocalSize(vec1, &n1));
1725:       PetscCall(VecGetLocalSize(vec2, &n2));
1726:       if (n1 != n2) flg1 = PETSC_FALSE;
1727:       else {
1728:         PetscCall(VecGetArrayRead(vec1, &v1));
1729:         PetscCall(VecGetArrayRead(vec2, &v2));
1730:         PetscCall(PetscArraycmp(v1, v2, n1, &flg1));
1731:         PetscCall(VecRestoreArrayRead(vec1, &v1));
1732:         PetscCall(VecRestoreArrayRead(vec2, &v2));
1733:       }
1734:     }
1735:     /* combine results from all processors */
1736:     PetscCallMPI(MPIU_Allreduce(&flg1, flg, 1, MPI_C_BOOL, MPI_LAND, PetscObjectComm((PetscObject)vec1)));
1737:   }
1738:   PetscFunctionReturn(PETSC_SUCCESS);
1739: }

1741: /*@
1742:   VecUniqueEntries - Compute the number of unique entries, and those entries

1744:   Collective

1746:   Input Parameter:
1747: . vec - the vector

1749:   Output Parameters:
1750: + n - The number of unique entries
1751: - e - The entries, each MPI process receives all the unique entries

1753:   Level: intermediate

1755: .seealso: `Vec`
1756: @*/
1757: PetscErrorCode VecUniqueEntries(Vec vec, PetscInt *n, PetscScalar *e[])
1758: {
1759:   const PetscScalar *v;
1760:   PetscScalar       *tmp, *vals;
1761:   PetscMPIInt       *N, *displs, l;
1762:   PetscInt           ng, m, i, j, p;
1763:   PetscMPIInt        size;

1765:   PetscFunctionBegin;
1767:   PetscAssertPointer(n, 2);
1768:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)vec), &size));
1769:   PetscCall(VecGetLocalSize(vec, &m));
1770:   PetscCall(VecGetArrayRead(vec, &v));
1771:   PetscCall(PetscMalloc2(m, &tmp, size, &N));
1772:   for (i = 0, l = 0; i < m; ++i) {
1773:     /* Can speed this up with sorting */
1774:     for (j = 0; j < l; ++j) {
1775:       if (v[i] == tmp[j]) break;
1776:     }
1777:     if (j == l) {
1778:       tmp[j] = v[i];
1779:       ++l;
1780:     }
1781:   }
1782:   PetscCall(VecRestoreArrayRead(vec, &v));
1783:   /* Gather serial results */
1784:   PetscCallMPI(MPI_Allgather(&l, 1, MPI_INT, N, 1, MPI_INT, PetscObjectComm((PetscObject)vec)));
1785:   for (p = 0, ng = 0; p < size; ++p) ng += N[p];
1786:   PetscCall(PetscMalloc2(ng, &vals, size + 1, &displs));
1787:   for (p = 1, displs[0] = 0; p <= size; ++p) displs[p] = displs[p - 1] + N[p - 1];
1788:   PetscCallMPI(MPI_Allgatherv(tmp, l, MPIU_SCALAR, vals, N, displs, MPIU_SCALAR, PetscObjectComm((PetscObject)vec)));
1789:   /* Find unique entries */
1790: #ifdef PETSC_USE_COMPLEX
1791:   SETERRQ(PetscObjectComm((PetscObject)vec), PETSC_ERR_SUP, "Does not work with complex numbers");
1792: #else
1793:   *n = displs[size];
1794:   PetscCall(PetscSortRemoveDupsReal(n, vals));
1795:   if (e) {
1796:     PetscAssertPointer(e, 3);
1797:     PetscCall(PetscMalloc1(*n, e));
1798:     for (i = 0; i < *n; ++i) (*e)[i] = vals[i];
1799:   }
1800:   PetscCall(PetscFree2(vals, displs));
1801:   PetscCall(PetscFree2(tmp, N));
1802:   PetscFunctionReturn(PETSC_SUCCESS);
1803: #endif
1804: }