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: /*@
646: VecStrideGatherAll - Gathers all the single components from a multi-component vector into
647: separate vectors.
649: Collective
651: Input Parameters:
652: + v - the vector
653: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
655: Output Parameter:
656: . s - the location where the subvectors are stored
658: Level: advanced
660: Notes:
661: One must call `VecSetBlockSize()` before this routine to set the stride
662: information, or use a vector created from a multicomponent `DMDA`.
664: If x is the array representing the vector x then this gathers
665: the arrays (x[start],x[start+stride],x[start+2*stride], ....)
666: for start=0,1,2,...bs-1
668: The parallel layout of the vector and the subvector must be the same;
669: i.e., nlocal of v = stride*(nlocal of s)
671: Not optimized; could be easily
673: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGather()`,
674: `VecStrideScatterAll()`
675: @*/
676: PetscErrorCode VecStrideGatherAll(Vec v, Vec s[], InsertMode addv)
677: {
678: PetscInt i, n, n2, bs, j, k, *bss = NULL, nv, jj, nvc;
679: PetscScalar **y;
680: const PetscScalar *x;
682: PetscFunctionBegin;
684: PetscAssertPointer(s, 2);
686: PetscCall(VecGetLocalSize(v, &n));
687: PetscCall(VecGetLocalSize(s[0], &n2));
688: PetscCall(VecGetArrayRead(v, &x));
689: PetscCall(VecGetBlockSize(v, &bs));
690: PetscCheck(bs > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Input vector does not have a valid blocksize set");
692: PetscCall(PetscMalloc2(bs, &y, bs, &bss));
693: nv = 0;
694: nvc = 0;
695: for (i = 0; i < bs; i++) {
696: PetscCall(VecGetBlockSize(s[i], &bss[i]));
697: if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1 Re: [PETSC #8241] VecStrideGatherAll */
698: PetscCall(VecGetArray(s[i], &y[i]));
699: nvc += bss[i];
700: nv++;
701: PetscCheck(nvc <= bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of subvectors in subvectors > number of vectors in main vector");
702: if (nvc == bs) break;
703: }
705: n = n / bs;
707: jj = 0;
708: if (addv == INSERT_VALUES) {
709: for (j = 0; j < nv; j++) {
710: for (k = 0; k < bss[j]; k++) {
711: for (i = 0; i < n; i++) y[j][i * bss[j] + k] = x[bs * i + jj + k];
712: }
713: jj += bss[j];
714: }
715: } else if (addv == ADD_VALUES) {
716: for (j = 0; j < nv; j++) {
717: for (k = 0; k < bss[j]; k++) {
718: for (i = 0; i < n; i++) y[j][i * bss[j] + k] += x[bs * i + jj + k];
719: }
720: jj += bss[j];
721: }
722: #if !defined(PETSC_USE_COMPLEX)
723: } else if (addv == MAX_VALUES) {
724: for (j = 0; j < nv; j++) {
725: for (k = 0; k < bss[j]; k++) {
726: for (i = 0; i < n; i++) y[j][i * bss[j] + k] = PetscMax(y[j][i * bss[j] + k], x[bs * i + jj + k]);
727: }
728: jj += bss[j];
729: }
730: #endif
731: } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
733: PetscCall(VecRestoreArrayRead(v, &x));
734: for (i = 0; i < nv; i++) PetscCall(VecRestoreArray(s[i], &y[i]));
736: PetscCall(PetscFree2(y, bss));
737: PetscFunctionReturn(PETSC_SUCCESS);
738: }
740: /*@
741: VecStrideScatterAll - Scatters all the single components from separate vectors into
742: a multi-component vector.
744: Collective
746: Input Parameters:
747: + s - the location where the subvectors are stored
748: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
750: Output Parameter:
751: . v - the multicomponent vector
753: Level: advanced
755: Notes:
756: One must call `VecSetBlockSize()` before this routine to set the stride
757: information, or use a vector created from a multicomponent `DMDA`.
759: The parallel layout of the vector and the subvector must be the same;
760: i.e., nlocal of v = stride*(nlocal of s)
762: Not optimized; could be easily
764: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGather()`,
766: @*/
767: PetscErrorCode VecStrideScatterAll(Vec s[], Vec v, InsertMode addv)
768: {
769: PetscInt i, n, n2, bs, j, jj, k, *bss = NULL, nv, nvc;
770: PetscScalar *x;
771: PetscScalar const **y;
773: PetscFunctionBegin;
775: PetscAssertPointer(s, 1);
777: PetscCall(VecGetLocalSize(v, &n));
778: PetscCall(VecGetLocalSize(s[0], &n2));
779: PetscCall(VecGetArray(v, &x));
780: PetscCall(VecGetBlockSize(v, &bs));
781: PetscCheck(bs > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Input vector does not have a valid blocksize set");
783: PetscCall(PetscMalloc2(bs, (PetscScalar ***)&y, bs, &bss));
784: nv = 0;
785: nvc = 0;
786: for (i = 0; i < bs; i++) {
787: PetscCall(VecGetBlockSize(s[i], &bss[i]));
788: if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1 Re: [PETSC #8241] VecStrideGatherAll */
789: PetscCall(VecGetArrayRead(s[i], &y[i]));
790: nvc += bss[i];
791: nv++;
792: PetscCheck(nvc <= bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of subvectors in subvectors > number of vectors in main vector");
793: if (nvc == bs) break;
794: }
796: n = n / bs;
798: jj = 0;
799: if (addv == INSERT_VALUES) {
800: for (j = 0; j < nv; j++) {
801: for (k = 0; k < bss[j]; k++) {
802: for (i = 0; i < n; i++) x[bs * i + jj + k] = y[j][i * bss[j] + k];
803: }
804: jj += bss[j];
805: }
806: } else if (addv == ADD_VALUES) {
807: for (j = 0; j < nv; j++) {
808: for (k = 0; k < bss[j]; k++) {
809: for (i = 0; i < n; i++) x[bs * i + jj + k] += y[j][i * bss[j] + k];
810: }
811: jj += bss[j];
812: }
813: #if !defined(PETSC_USE_COMPLEX)
814: } else if (addv == MAX_VALUES) {
815: for (j = 0; j < nv; j++) {
816: for (k = 0; k < bss[j]; k++) {
817: for (i = 0; i < n; i++) x[bs * i + jj + k] = PetscMax(x[bs * i + jj + k], y[j][i * bss[j] + k]);
818: }
819: jj += bss[j];
820: }
821: #endif
822: } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
824: PetscCall(VecRestoreArray(v, &x));
825: for (i = 0; i < nv; i++) PetscCall(VecRestoreArrayRead(s[i], &y[i]));
826: PetscCall(PetscFree2(*(PetscScalar ***)&y, bss));
827: PetscFunctionReturn(PETSC_SUCCESS);
828: }
830: /*@
831: VecStrideGather - Gathers a single component from a multi-component vector into
832: another vector.
834: Collective
836: Input Parameters:
837: + v - the vector
838: . start - starting point of the subvector (defined by a stride)
839: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
841: Output Parameter:
842: . s - the location where the subvector is stored
844: Level: advanced
846: Notes:
847: One must call `VecSetBlockSize()` before this routine to set the stride
848: information, or use a vector created from a multicomponent `DMDA`.
850: If x is the array representing the vector x then this gathers
851: the array (x[start],x[start+stride],x[start+2*stride], ....)
853: The parallel layout of the vector and the subvector must be the same;
854: i.e., nlocal of v = stride*(nlocal of s)
856: Not optimized; could be easily
858: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
859: `VecStrideScatterAll()`
860: @*/
861: PetscErrorCode VecStrideGather(Vec v, PetscInt start, Vec s, InsertMode addv)
862: {
863: PetscFunctionBegin;
867: PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
868: 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,
869: v->map->bs);
870: PetscUseTypeMethod(v, stridegather, start, s, addv);
871: PetscFunctionReturn(PETSC_SUCCESS);
872: }
874: /*@
875: VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
877: Collective
879: Input Parameters:
880: + s - the single-component vector
881: . start - starting point of the subvector (defined by a stride)
882: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
884: Output Parameter:
885: . v - the location where the subvector is scattered (the multi-component vector)
887: Level: advanced
889: Notes:
890: One must call `VecSetBlockSize()` on the multi-component vector before this
891: routine to set the stride information, or use a vector created from a multicomponent `DMDA`.
893: The parallel layout of the vector and the subvector must be the same;
894: i.e., nlocal of v = stride*(nlocal of s)
896: Not optimized; could be easily
898: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
899: `VecStrideScatterAll()`, `VecStrideSubSetScatter()`, `VecStrideSubSetGather()`
900: @*/
901: PetscErrorCode VecStrideScatter(Vec s, PetscInt start, Vec v, InsertMode addv)
902: {
903: PetscFunctionBegin;
907: PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
908: 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,
909: v->map->bs);
910: PetscCall((*v->ops->stridescatter)(s, start, v, addv));
911: PetscFunctionReturn(PETSC_SUCCESS);
912: }
914: /*@
915: VecStrideSubSetGather - Gathers a subset of components from a multi-component vector into
916: another vector.
918: Collective
920: Input Parameters:
921: + v - the vector
922: . nidx - the number of indices
923: . idxv - the indices of the components 0 <= idxv[0] ...idxv[nidx-1] < bs(v), they need not be sorted
924: . 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`
925: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
927: Output Parameter:
928: . s - the location where the subvector is stored
930: Level: advanced
932: Notes:
933: One must call `VecSetBlockSize()` on both vectors before this routine to set the stride
934: information, or use a vector created from a multicomponent `DMDA`.
936: The parallel layout of the vector and the subvector must be the same;
938: Not optimized; could be easily
940: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideGather()`, `VecStrideSubSetScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
941: `VecStrideScatterAll()`
942: @*/
943: PetscErrorCode VecStrideSubSetGather(Vec v, PetscInt nidx, const PetscInt idxv[], const PetscInt idxs[], Vec s, InsertMode addv)
944: {
945: PetscFunctionBegin;
948: if (nidx == PETSC_DETERMINE) nidx = s->map->bs;
949: PetscUseTypeMethod(v, stridesubsetgather, nidx, idxv, idxs, s, addv);
950: PetscFunctionReturn(PETSC_SUCCESS);
951: }
953: /*@
954: VecStrideSubSetScatter - Scatters components from a vector into a subset of components of a multi-component vector.
956: Collective
958: Input Parameters:
959: + s - the smaller-component vector
960: . nidx - the number of indices in idx
961: . 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`
962: . idxv - the indices of the components in the larger-component vector, 0 <= idx[0] ...idx[nidx-1] < bs(v) they need not be sorted
963: - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
965: Output Parameter:
966: . v - the location where the subvector is into scattered (the multi-component vector)
968: Level: advanced
970: Notes:
971: One must call `VecSetBlockSize()` on the vectors before this
972: routine to set the stride information, or use a vector created from a multicomponent `DMDA`.
974: The parallel layout of the vector and the subvector must be the same;
976: Not optimized; could be easily
978: .seealso: `Vec`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideSubSetGather()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
979: `VecStrideScatterAll()`
980: @*/
981: PetscErrorCode VecStrideSubSetScatter(Vec s, PetscInt nidx, const PetscInt idxs[], const PetscInt idxv[], Vec v, InsertMode addv)
982: {
983: PetscFunctionBegin;
986: if (nidx == PETSC_DETERMINE) nidx = s->map->bs;
987: PetscCall((*v->ops->stridesubsetscatter)(s, nidx, idxs, idxv, v, addv));
988: PetscFunctionReturn(PETSC_SUCCESS);
989: }
991: PetscErrorCode VecStrideGather_Default(Vec v, PetscInt start, Vec s, InsertMode addv)
992: {
993: PetscInt i, n, bs, ns;
994: const PetscScalar *x;
995: PetscScalar *y;
997: PetscFunctionBegin;
998: PetscCall(VecGetLocalSize(v, &n));
999: PetscCall(VecGetLocalSize(s, &ns));
1000: PetscCall(VecGetArrayRead(v, &x));
1001: PetscCall(VecGetArray(s, &y));
1003: bs = v->map->bs;
1004: 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);
1005: x += start;
1006: n = n / bs;
1008: if (addv == INSERT_VALUES) {
1009: for (i = 0; i < n; i++) y[i] = x[bs * i];
1010: } else if (addv == ADD_VALUES) {
1011: for (i = 0; i < n; i++) y[i] += x[bs * i];
1012: #if !defined(PETSC_USE_COMPLEX)
1013: } else if (addv == MAX_VALUES) {
1014: for (i = 0; i < n; i++) y[i] = PetscMax(y[i], x[bs * i]);
1015: #endif
1016: } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1018: PetscCall(VecRestoreArrayRead(v, &x));
1019: PetscCall(VecRestoreArray(s, &y));
1020: PetscFunctionReturn(PETSC_SUCCESS);
1021: }
1023: PetscErrorCode VecStrideScatter_Default(Vec s, PetscInt start, Vec v, InsertMode addv)
1024: {
1025: PetscInt i, n, bs, ns;
1026: PetscScalar *x;
1027: const PetscScalar *y;
1029: PetscFunctionBegin;
1030: PetscCall(VecGetLocalSize(v, &n));
1031: PetscCall(VecGetLocalSize(s, &ns));
1032: PetscCall(VecGetArray(v, &x));
1033: PetscCall(VecGetArrayRead(s, &y));
1035: bs = v->map->bs;
1036: 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);
1037: x += start;
1038: n = n / bs;
1040: if (addv == INSERT_VALUES) {
1041: for (i = 0; i < n; i++) x[bs * i] = y[i];
1042: } else if (addv == ADD_VALUES) {
1043: for (i = 0; i < n; i++) x[bs * i] += y[i];
1044: #if !defined(PETSC_USE_COMPLEX)
1045: } else if (addv == MAX_VALUES) {
1046: for (i = 0; i < n; i++) x[bs * i] = PetscMax(y[i], x[bs * i]);
1047: #endif
1048: } else SETERRQ(PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1050: PetscCall(VecRestoreArray(v, &x));
1051: PetscCall(VecRestoreArrayRead(s, &y));
1052: PetscFunctionReturn(PETSC_SUCCESS);
1053: }
1055: PetscErrorCode VecStrideSubSetGather_Default(Vec v, PetscInt nidx, const PetscInt idxv[], const PetscInt idxs[], Vec s, InsertMode addv)
1056: {
1057: PetscInt i, j, n, bs, bss, ns;
1058: const PetscScalar *x;
1059: PetscScalar *y;
1061: PetscFunctionBegin;
1062: PetscCall(VecGetLocalSize(v, &n));
1063: PetscCall(VecGetLocalSize(s, &ns));
1064: PetscCall(VecGetArrayRead(v, &x));
1065: PetscCall(VecGetArray(s, &y));
1067: bs = v->map->bs;
1068: bss = s->map->bs;
1069: n = n / bs;
1071: if (PetscDefined(USE_DEBUG)) {
1072: PetscCheck(n == ns / bss, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout of vectors");
1073: for (j = 0; j < nidx; j++) {
1074: PetscCheck(idxv[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is negative", j, idxv[j]);
1075: 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);
1076: }
1077: PetscCheck(idxs || bss == nidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide idxs when not gathering into all locations");
1078: }
1080: if (addv == INSERT_VALUES) {
1081: if (!idxs) {
1082: for (i = 0; i < n; i++) {
1083: for (j = 0; j < bss; j++) y[bss * i + j] = x[bs * i + idxv[j]];
1084: }
1085: } else {
1086: for (i = 0; i < n; i++) {
1087: for (j = 0; j < bss; j++) y[bss * i + idxs[j]] = x[bs * i + idxv[j]];
1088: }
1089: }
1090: } else if (addv == ADD_VALUES) {
1091: if (!idxs) {
1092: for (i = 0; i < n; i++) {
1093: for (j = 0; j < bss; j++) y[bss * i + j] += x[bs * i + idxv[j]];
1094: }
1095: } else {
1096: for (i = 0; i < n; i++) {
1097: for (j = 0; j < bss; j++) y[bss * i + idxs[j]] += x[bs * i + idxv[j]];
1098: }
1099: }
1100: #if !defined(PETSC_USE_COMPLEX)
1101: } else if (addv == MAX_VALUES) {
1102: if (!idxs) {
1103: for (i = 0; i < n; i++) {
1104: for (j = 0; j < bss; j++) y[bss * i + j] = PetscMax(y[bss * i + j], x[bs * i + idxv[j]]);
1105: }
1106: } else {
1107: for (i = 0; i < n; i++) {
1108: for (j = 0; j < bss; j++) y[bss * i + idxs[j]] = PetscMax(y[bss * i + idxs[j]], x[bs * i + idxv[j]]);
1109: }
1110: }
1111: #endif
1112: } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1114: PetscCall(VecRestoreArrayRead(v, &x));
1115: PetscCall(VecRestoreArray(s, &y));
1116: PetscFunctionReturn(PETSC_SUCCESS);
1117: }
1119: PetscErrorCode VecStrideSubSetScatter_Default(Vec s, PetscInt nidx, const PetscInt idxs[], const PetscInt idxv[], Vec v, InsertMode addv)
1120: {
1121: PetscInt j, i, n, bs, ns, bss;
1122: PetscScalar *x;
1123: const PetscScalar *y;
1125: PetscFunctionBegin;
1126: PetscCall(VecGetLocalSize(v, &n));
1127: PetscCall(VecGetLocalSize(s, &ns));
1128: PetscCall(VecGetArray(v, &x));
1129: PetscCall(VecGetArrayRead(s, &y));
1131: bs = v->map->bs;
1132: bss = s->map->bs;
1133: n = n / bs;
1135: if (PetscDefined(USE_DEBUG)) {
1136: PetscCheck(n == ns / bss, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout of vectors");
1137: for (j = 0; j < bss; j++) {
1138: if (idxs) {
1139: PetscCheck(idxs[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is negative", j, idxs[j]);
1140: 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);
1141: }
1142: }
1143: PetscCheck(idxs || bss == nidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide idxs when not scattering from all locations");
1144: }
1146: if (addv == INSERT_VALUES) {
1147: if (!idxs) {
1148: for (i = 0; i < n; i++) {
1149: for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = y[bss * i + j];
1150: }
1151: } else {
1152: for (i = 0; i < n; i++) {
1153: for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = y[bss * i + idxs[j]];
1154: }
1155: }
1156: } else if (addv == ADD_VALUES) {
1157: if (!idxs) {
1158: for (i = 0; i < n; i++) {
1159: for (j = 0; j < bss; j++) x[bs * i + idxv[j]] += y[bss * i + j];
1160: }
1161: } else {
1162: for (i = 0; i < n; i++) {
1163: for (j = 0; j < bss; j++) x[bs * i + idxv[j]] += y[bss * i + idxs[j]];
1164: }
1165: }
1166: #if !defined(PETSC_USE_COMPLEX)
1167: } else if (addv == MAX_VALUES) {
1168: if (!idxs) {
1169: for (i = 0; i < n; i++) {
1170: for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = PetscMax(y[bss * i + j], x[bs * i + idxv[j]]);
1171: }
1172: } else {
1173: for (i = 0; i < n; i++) {
1174: for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = PetscMax(y[bss * i + idxs[j]], x[bs * i + idxv[j]]);
1175: }
1176: }
1177: #endif
1178: } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1180: PetscCall(VecRestoreArray(v, &x));
1181: PetscCall(VecRestoreArrayRead(s, &y));
1182: PetscFunctionReturn(PETSC_SUCCESS);
1183: }
1185: static PetscErrorCode VecApplyUnary_Private(Vec v, PetscDeviceContext dctx, const char async_name[], PetscErrorCode (*unary_op)(Vec), PetscScalar (*UnaryFunc)(PetscScalar))
1186: {
1187: PetscFunctionBegin;
1189: PetscCall(VecSetErrorIfLocked(v, 1));
1190: if (dctx) {
1191: PetscErrorCode (*unary_op_async)(Vec, PetscDeviceContext);
1193: PetscCall(PetscObjectQueryFunction((PetscObject)v, async_name, &unary_op_async));
1194: if (unary_op_async) {
1195: PetscCall((*unary_op_async)(v, dctx));
1196: PetscFunctionReturn(PETSC_SUCCESS);
1197: }
1198: }
1199: if (unary_op) {
1201: PetscCall((*unary_op)(v));
1202: } else {
1203: PetscInt n;
1204: PetscScalar *x;
1207: PetscCall(VecGetLocalSize(v, &n));
1208: PetscCall(VecGetArray(v, &x));
1209: for (PetscInt i = 0; i < n; ++i) x[i] = UnaryFunc(x[i]);
1210: PetscCall(VecRestoreArray(v, &x));
1211: }
1212: PetscFunctionReturn(PETSC_SUCCESS);
1213: }
1215: static PetscScalar ScalarReciprocal_Function(PetscScalar x)
1216: {
1217: const PetscScalar zero = 0.0;
1219: return x == zero ? zero : ((PetscScalar)1.0) / x;
1220: }
1222: PetscErrorCode VecReciprocalAsync_Private(Vec v, PetscDeviceContext dctx)
1223: {
1224: PetscFunctionBegin;
1225: PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Reciprocal), v->ops->reciprocal, ScalarReciprocal_Function));
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()`
1264: @*/
1265: PetscErrorCode VecExp(Vec v)
1266: {
1267: PetscFunctionBegin;
1268: PetscCall(VecExpAsync_Private(v, NULL));
1269: PetscFunctionReturn(PETSC_SUCCESS);
1270: }
1272: static PetscScalar ScalarLog_Function(PetscScalar x)
1273: {
1274: return PetscLogScalar(x);
1275: }
1277: PetscErrorCode VecLogAsync_Private(Vec v, PetscDeviceContext dctx)
1278: {
1279: PetscFunctionBegin;
1281: PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Log), v->ops->log, ScalarLog_Function));
1282: PetscFunctionReturn(PETSC_SUCCESS);
1283: }
1285: /*@
1286: VecLog - Replaces each component of a vector by log(x_i), the natural logarithm
1288: Not Collective
1290: Input Parameter:
1291: . v - The vector
1293: Output Parameter:
1294: . v - The vector of logs
1296: Level: beginner
1298: .seealso: `Vec`, `VecExp()`, `VecAbs()`, `VecSqrtAbs()`, `VecReciprocal()`
1300: @*/
1301: PetscErrorCode VecLog(Vec v)
1302: {
1303: PetscFunctionBegin;
1304: PetscCall(VecLogAsync_Private(v, NULL));
1305: PetscFunctionReturn(PETSC_SUCCESS);
1306: }
1308: static PetscScalar ScalarAbs_Function(PetscScalar x)
1309: {
1310: return PetscAbsScalar(x);
1311: }
1313: PetscErrorCode VecAbsAsync_Private(Vec v, PetscDeviceContext dctx)
1314: {
1315: PetscFunctionBegin;
1317: PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Abs), v->ops->abs, ScalarAbs_Function));
1318: PetscFunctionReturn(PETSC_SUCCESS);
1319: }
1321: /*@
1322: VecAbs - Replaces every element in a vector with its absolute value.
1324: Logically Collective
1326: Input Parameter:
1327: . v - the vector
1329: Level: intermediate
1331: .seealso: `Vec`, `VecExp()`, `VecSqrtAbs()`, `VecReciprocal()`, `VecLog()`
1332: @*/
1333: PetscErrorCode VecAbs(Vec v)
1334: {
1335: PetscFunctionBegin;
1336: PetscCall(VecAbsAsync_Private(v, NULL));
1337: PetscFunctionReturn(PETSC_SUCCESS);
1338: }
1340: static PetscScalar ScalarConjugate_Function(PetscScalar x)
1341: {
1342: return PetscConj(x);
1343: }
1345: PetscErrorCode VecConjugateAsync_Private(Vec v, PetscDeviceContext dctx)
1346: {
1347: PetscFunctionBegin;
1349: if (PetscDefined(USE_COMPLEX)) PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(Conjugate), v->ops->conjugate, ScalarConjugate_Function));
1350: PetscFunctionReturn(PETSC_SUCCESS);
1351: }
1353: /*@
1354: VecConjugate - Conjugates a vector. That is, replace every entry in a vector with its complex conjugate
1356: Logically Collective
1358: Input Parameter:
1359: . x - the vector
1361: Level: intermediate
1363: .seealso: [](ch_vectors), `Vec`, `VecSet()`
1364: @*/
1365: PetscErrorCode VecConjugate(Vec x)
1366: {
1367: PetscFunctionBegin;
1368: PetscCall(VecConjugateAsync_Private(x, NULL));
1369: PetscFunctionReturn(PETSC_SUCCESS);
1370: }
1372: static PetscScalar ScalarSqrtAbs_Function(PetscScalar x)
1373: {
1374: return PetscSqrtScalar(ScalarAbs_Function(x));
1375: }
1377: PetscErrorCode VecSqrtAbsAsync_Private(Vec v, PetscDeviceContext dctx)
1378: {
1379: PetscFunctionBegin;
1381: PetscCall(VecApplyUnary_Private(v, dctx, VecAsyncFnName(SqrtAbs), v->ops->sqrt, ScalarSqrtAbs_Function));
1382: PetscFunctionReturn(PETSC_SUCCESS);
1383: }
1385: /*@
1386: VecSqrtAbs - Replaces each component of a vector by the square root of its magnitude.
1388: Not Collective
1390: Input Parameter:
1391: . v - The vector
1393: Level: beginner
1395: Note:
1396: The actual function is sqrt(|x_i|)
1398: .seealso: `Vec`, `VecLog()`, `VecExp()`, `VecReciprocal()`, `VecAbs()`
1400: @*/
1401: PetscErrorCode VecSqrtAbs(Vec v)
1402: {
1403: PetscFunctionBegin;
1404: PetscCall(VecSqrtAbsAsync_Private(v, NULL));
1405: PetscFunctionReturn(PETSC_SUCCESS);
1406: }
1408: static PetscScalar ScalarImaginaryPart_Function(PetscScalar x)
1409: {
1410: const PetscReal imag = PetscImaginaryPart(x);
1412: #if PetscDefined(USE_COMPLEX)
1413: return PetscCMPLX(imag, 0.0);
1414: #else
1415: return imag;
1416: #endif
1417: }
1419: /*@
1420: VecImaginaryPart - Replaces a complex vector with its imginary part
1422: Collective
1424: Input Parameter:
1425: . v - the vector
1427: Level: beginner
1429: .seealso: `Vec`, `VecNorm()`, `VecRealPart()`
1430: @*/
1431: PetscErrorCode VecImaginaryPart(Vec v)
1432: {
1433: PetscFunctionBegin;
1435: PetscCall(VecApplyUnary_Private(v, NULL, NULL, NULL, ScalarImaginaryPart_Function));
1436: PetscFunctionReturn(PETSC_SUCCESS);
1437: }
1439: static PetscScalar ScalarRealPart_Function(PetscScalar x)
1440: {
1441: const PetscReal real = PetscRealPart(x);
1443: #if PetscDefined(USE_COMPLEX)
1444: return PetscCMPLX(real, 0.0);
1445: #else
1446: return real;
1447: #endif
1448: }
1450: /*@
1451: VecRealPart - Replaces a complex vector with its real part
1453: Collective
1455: Input Parameter:
1456: . v - the vector
1458: Level: beginner
1460: .seealso: `Vec`, `VecNorm()`, `VecImaginaryPart()`
1461: @*/
1462: PetscErrorCode VecRealPart(Vec v)
1463: {
1464: PetscFunctionBegin;
1466: PetscCall(VecApplyUnary_Private(v, NULL, NULL, NULL, ScalarRealPart_Function));
1467: PetscFunctionReturn(PETSC_SUCCESS);
1468: }
1470: /*@
1471: VecDotNorm2 - computes the inner product of two vectors and the 2-norm squared of the second vector
1473: Collective
1475: Input Parameters:
1476: + s - first vector
1477: - t - second vector
1479: Output Parameters:
1480: + dp - s'conj(t)
1481: - nm - t'conj(t)
1483: Level: advanced
1485: Note:
1486: conj(x) is the complex conjugate of x when x is complex
1488: .seealso: `Vec`, `VecDot()`, `VecNorm()`, `VecDotBegin()`, `VecNormBegin()`, `VecDotEnd()`, `VecNormEnd()`
1490: @*/
1491: PetscErrorCode VecDotNorm2(Vec s, Vec t, PetscScalar *dp, PetscReal *nm)
1492: {
1493: PetscScalar work[] = {0.0, 0.0};
1495: PetscFunctionBegin;
1498: PetscAssertPointer(dp, 3);
1499: PetscAssertPointer(nm, 4);
1502: PetscCheckSameTypeAndComm(s, 1, t, 2);
1503: PetscCheck(s->map->N == t->map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector global lengths");
1504: PetscCheck(s->map->n == t->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector local lengths");
1506: PetscCall(PetscLogEventBegin(VEC_DotNorm2, s, t, 0, 0));
1507: if (s->ops->dotnorm2) {
1508: PetscUseTypeMethod(s, dotnorm2, t, work, work + 1);
1509: } else {
1510: const PetscScalar *sx, *tx;
1511: PetscInt n;
1513: PetscCall(VecGetLocalSize(s, &n));
1514: PetscCall(VecGetArrayRead(s, &sx));
1515: PetscCall(VecGetArrayRead(t, &tx));
1516: for (PetscInt i = 0; i < n; ++i) {
1517: const PetscScalar txconj = PetscConj(tx[i]);
1519: work[0] += sx[i] * txconj;
1520: work[1] += tx[i] * txconj;
1521: }
1522: PetscCall(VecRestoreArrayRead(t, &tx));
1523: PetscCall(VecRestoreArrayRead(s, &sx));
1524: PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, work, 2, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)s)));
1525: PetscCall(PetscLogFlops(4.0 * n));
1526: }
1527: PetscCall(PetscLogEventEnd(VEC_DotNorm2, s, t, 0, 0));
1528: *dp = work[0];
1529: *nm = PetscRealPart(work[1]);
1530: PetscFunctionReturn(PETSC_SUCCESS);
1531: }
1533: /*@
1534: VecSum - Computes the sum of all the components of a vector.
1536: Collective
1538: Input Parameter:
1539: . v - the vector
1541: Output Parameter:
1542: . sum - the result
1544: Level: beginner
1546: .seealso: `Vec`, `VecMean()`, `VecNorm()`
1547: @*/
1548: PetscErrorCode VecSum(Vec v, PetscScalar *sum)
1549: {
1550: PetscScalar tmp = 0.0;
1552: PetscFunctionBegin;
1554: PetscAssertPointer(sum, 2);
1555: if (v->ops->sum) {
1556: PetscUseTypeMethod(v, sum, &tmp);
1557: } else {
1558: const PetscScalar *x;
1559: PetscInt n;
1561: PetscCall(VecGetLocalSize(v, &n));
1562: PetscCall(VecGetArrayRead(v, &x));
1563: for (PetscInt i = 0; i < n; ++i) tmp += x[i];
1564: PetscCall(VecRestoreArrayRead(v, &x));
1565: }
1566: PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &tmp, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)v)));
1567: *sum = tmp;
1568: PetscFunctionReturn(PETSC_SUCCESS);
1569: }
1571: /*@
1572: VecMean - Computes the arithmetic mean of all the components of a vector.
1574: Collective
1576: Input Parameter:
1577: . v - the vector
1579: Output Parameter:
1580: . mean - the result
1582: Level: beginner
1584: .seealso: `Vec`, `VecSum()`, `VecNorm()`
1585: @*/
1586: PetscErrorCode VecMean(Vec v, PetscScalar *mean)
1587: {
1588: PetscInt n;
1590: PetscFunctionBegin;
1592: PetscAssertPointer(mean, 2);
1593: PetscCall(VecGetSize(v, &n));
1594: PetscCall(VecSum(v, mean));
1595: *mean /= n;
1596: PetscFunctionReturn(PETSC_SUCCESS);
1597: }
1599: PetscErrorCode VecShiftAsync_Private(Vec v, PetscScalar shift, PetscDeviceContext dctx)
1600: {
1601: PetscErrorCode (*shift_async)(Vec, PetscScalar, PetscDeviceContext) = NULL;
1603: PetscFunctionBegin;
1604: if (dctx) {
1605: PetscErrorCode (*shift_async)(Vec, PetscScalar, PetscDeviceContext);
1607: PetscCall(PetscObjectQueryFunction((PetscObject)v, VecAsyncFnName(Shift), &shift_async));
1608: }
1609: if (shift_async) {
1610: PetscCall((*shift_async)(v, shift, dctx));
1611: } else if (v->ops->shift) {
1612: PetscUseTypeMethod(v, shift, shift);
1613: } else {
1614: PetscInt n;
1615: PetscScalar *x;
1617: PetscCall(VecGetLocalSize(v, &n));
1618: PetscCall(VecGetArray(v, &x));
1619: for (PetscInt i = 0; i < n; ++i) x[i] += shift;
1620: PetscCall(VecRestoreArray(v, &x));
1621: PetscCall(PetscLogFlops(n));
1622: }
1623: PetscFunctionReturn(PETSC_SUCCESS);
1624: }
1626: /*@
1627: VecShift - Shifts all of the components of a vector by computing
1628: `x[i] = x[i] + shift`.
1630: Logically Collective
1632: Input Parameters:
1633: + v - the vector
1634: - shift - the shift
1636: Level: intermediate
1638: .seealso: `Vec`, `VecISShift()`
1639: @*/
1640: PetscErrorCode VecShift(Vec v, PetscScalar shift)
1641: {
1642: PetscFunctionBegin;
1645: PetscCall(VecSetErrorIfLocked(v, 1));
1646: if (shift == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
1647: PetscCall(PetscLogEventBegin(VEC_Shift, v, 0, 0, 0));
1648: PetscCall(VecShiftAsync_Private(v, shift, NULL));
1649: PetscCall(PetscLogEventEnd(VEC_Shift, v, 0, 0, 0));
1650: PetscFunctionReturn(PETSC_SUCCESS);
1651: }
1653: /*@
1654: VecPermute - Permutes a vector in place using the given ordering.
1656: Input Parameters:
1657: + x - The vector
1658: . row - The ordering
1659: - inv - The flag for inverting the permutation
1661: Level: beginner
1663: Note:
1664: This function does not yet support parallel Index Sets with non-local permutations
1666: .seealso: `Vec`, `MatPermute()`
1667: @*/
1668: PetscErrorCode VecPermute(Vec x, IS row, PetscBool inv)
1669: {
1670: PetscScalar *array, *newArray;
1671: const PetscInt *idx;
1672: PetscInt i, rstart, rend;
1674: PetscFunctionBegin;
1677: PetscCall(VecSetErrorIfLocked(x, 1));
1678: PetscCall(VecGetOwnershipRange(x, &rstart, &rend));
1679: PetscCall(ISGetIndices(row, &idx));
1680: PetscCall(VecGetArray(x, &array));
1681: PetscCall(PetscMalloc1(x->map->n, &newArray));
1682: PetscCall(PetscArraycpy(newArray, array, x->map->n));
1683: if (PetscDefined(USE_DEBUG)) {
1684: 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]);
1685: }
1686: if (!inv) {
1687: for (i = 0; i < x->map->n; i++) array[i] = newArray[idx[i] - rstart];
1688: } else {
1689: for (i = 0; i < x->map->n; i++) array[idx[i] - rstart] = newArray[i];
1690: }
1691: PetscCall(VecRestoreArray(x, &array));
1692: PetscCall(ISRestoreIndices(row, &idx));
1693: PetscCall(PetscFree(newArray));
1694: PetscFunctionReturn(PETSC_SUCCESS);
1695: }
1697: /*@
1698: VecEqual - Compares two vectors. Returns true if the two vectors are either pointing to the same memory buffer,
1699: or if the two vectors have the same local and global layout as well as bitwise equality of all entries.
1700: Does NOT take round-off errors into account.
1702: Collective
1704: Input Parameters:
1705: + vec1 - the first vector
1706: - vec2 - the second vector
1708: Output Parameter:
1709: . flg - `PETSC_TRUE` if the vectors are equal; `PETSC_FALSE` otherwise.
1711: Level: intermediate
1713: .seealso: `Vec`
1714: @*/
1715: PetscErrorCode VecEqual(Vec vec1, Vec vec2, PetscBool *flg)
1716: {
1717: const PetscScalar *v1, *v2;
1718: PetscInt n1, n2, N1, N2;
1719: PetscBool flg1;
1721: PetscFunctionBegin;
1724: PetscAssertPointer(flg, 3);
1725: if (vec1 == vec2) *flg = PETSC_TRUE;
1726: else {
1727: PetscCall(VecGetSize(vec1, &N1));
1728: PetscCall(VecGetSize(vec2, &N2));
1729: if (N1 != N2) flg1 = PETSC_FALSE;
1730: else {
1731: PetscCall(VecGetLocalSize(vec1, &n1));
1732: PetscCall(VecGetLocalSize(vec2, &n2));
1733: if (n1 != n2) flg1 = PETSC_FALSE;
1734: else {
1735: PetscCall(VecGetArrayRead(vec1, &v1));
1736: PetscCall(VecGetArrayRead(vec2, &v2));
1737: PetscCall(PetscArraycmp(v1, v2, n1, &flg1));
1738: PetscCall(VecRestoreArrayRead(vec1, &v1));
1739: PetscCall(VecRestoreArrayRead(vec2, &v2));
1740: }
1741: }
1742: /* combine results from all processors */
1743: PetscCallMPI(MPIU_Allreduce(&flg1, flg, 1, MPIU_BOOL, MPI_MIN, PetscObjectComm((PetscObject)vec1)));
1744: }
1745: PetscFunctionReturn(PETSC_SUCCESS);
1746: }
1748: /*@
1749: VecUniqueEntries - Compute the number of unique entries, and those entries
1751: Collective
1753: Input Parameter:
1754: . vec - the vector
1756: Output Parameters:
1757: + n - The number of unique entries
1758: - e - The entries, each MPI process receives all the unique entries
1760: Level: intermediate
1762: .seealso: `Vec`
1763: @*/
1764: PetscErrorCode VecUniqueEntries(Vec vec, PetscInt *n, PetscScalar *e[])
1765: {
1766: const PetscScalar *v;
1767: PetscScalar *tmp, *vals;
1768: PetscMPIInt *N, *displs, l;
1769: PetscInt ng, m, i, j, p;
1770: PetscMPIInt size;
1772: PetscFunctionBegin;
1774: PetscAssertPointer(n, 2);
1775: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)vec), &size));
1776: PetscCall(VecGetLocalSize(vec, &m));
1777: PetscCall(VecGetArrayRead(vec, &v));
1778: PetscCall(PetscMalloc2(m, &tmp, size, &N));
1779: for (i = 0, l = 0; i < m; ++i) {
1780: /* Can speed this up with sorting */
1781: for (j = 0; j < l; ++j) {
1782: if (v[i] == tmp[j]) break;
1783: }
1784: if (j == l) {
1785: tmp[j] = v[i];
1786: ++l;
1787: }
1788: }
1789: PetscCall(VecRestoreArrayRead(vec, &v));
1790: /* Gather serial results */
1791: PetscCallMPI(MPI_Allgather(&l, 1, MPI_INT, N, 1, MPI_INT, PetscObjectComm((PetscObject)vec)));
1792: for (p = 0, ng = 0; p < size; ++p) ng += N[p];
1793: PetscCall(PetscMalloc2(ng, &vals, size + 1, &displs));
1794: for (p = 1, displs[0] = 0; p <= size; ++p) displs[p] = displs[p - 1] + N[p - 1];
1795: PetscCallMPI(MPI_Allgatherv(tmp, l, MPIU_SCALAR, vals, N, displs, MPIU_SCALAR, PetscObjectComm((PetscObject)vec)));
1796: /* Find unique entries */
1797: #ifdef PETSC_USE_COMPLEX
1798: SETERRQ(PetscObjectComm((PetscObject)vec), PETSC_ERR_SUP, "Does not work with complex numbers");
1799: #else
1800: *n = displs[size];
1801: PetscCall(PetscSortRemoveDupsReal(n, vals));
1802: if (e) {
1803: PetscAssertPointer(e, 3);
1804: PetscCall(PetscMalloc1(*n, e));
1805: for (i = 0; i < *n; ++i) (*e)[i] = vals[i];
1806: }
1807: PetscCall(PetscFree2(vals, displs));
1808: PetscCall(PetscFree2(tmp, N));
1809: PetscFunctionReturn(PETSC_SUCCESS);
1810: #endif
1811: }