xref: /petsc/src/vec/vec/utils/vinv.c (revision 3df6fe7bb68c28049dcc48161db5df4684fa05c1)
1 
2 /*
3      Some useful vector utility functions.
4 */
5 #include <../src/vec/vec/impls/mpi/pvecimpl.h> /*I "petscvec.h" I*/
6 
7 /*@
8   VecStrideSet - Sets a subvector of a vector defined
9   by a starting point and a stride with a given value
10 
11   Logically Collective
12 
13   Input Parameters:
14 + v     - the vector
15 . start - starting point of the subvector (defined by a stride)
16 - s     - value to set for each entry in that subvector
17 
18   Level: advanced
19 
20   Notes:
21   One must call `VecSetBlockSize()` before this routine to set the stride
22   information, or use a vector created from a multicomponent `DMDA`.
23 
24   This will only work if the desire subvector is a stride subvector
25 
26 .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideScale()`
27 @*/
28 PetscErrorCode VecStrideSet(Vec v, PetscInt start, PetscScalar s)
29 {
30   PetscInt     i, n, bs;
31   PetscScalar *x;
32 
33   PetscFunctionBegin;
34   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
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\n  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 }
44 
45 /*@
46   VecStrideScale - Scales a subvector of a vector defined
47   by a starting point and a stride.
48 
49   Logically Collective
50 
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
55 
56   Level: advanced
57 
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`.
61 
62   This will only work if the desire subvector is a stride subvector
63 
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;
70 
71   PetscFunctionBegin;
72   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
73   PetscCall(VecGetLocalSize(v, &n));
74   PetscCall(VecGetBlockSize(v, &bs));
75   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
76   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\n  Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
77   PetscCall(VecGetArray(v, &x));
78   for (i = start; i < n; i += bs) x[i] *= scale;
79   PetscCall(VecRestoreArray(v, &x));
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 
83 /*@
84   VecStrideNorm - Computes the norm of subvector of a vector defined
85   by a starting point and a stride.
86 
87   Collective
88 
89   Input Parameters:
90 + v     - the vector
91 . start - starting point of the subvector (defined by a stride)
92 - ntype - type of norm, one of `NORM_1`, `NORM_2`, `NORM_INFINITY`
93 
94   Output Parameter:
95 . nrm - the norm
96 
97   Level: advanced
98 
99   Notes:
100   One must call `VecSetBlockSize()` before this routine to set the stride
101   information, or use a vector created from a multicomponent `DMDA`.
102 
103   If x is the array representing the vector x then this computes the norm
104   of the array (x[start],x[start+stride],x[start+2*stride], ....)
105 
106   This is useful for computing, say the norm of the pressure variable when
107   the pressure is stored (interlaced) with other variables, say density etc.
108 
109   This will only work if the desire subvector is a stride subvector
110 
111 .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
112 @*/
113 PetscErrorCode VecStrideNorm(Vec v, PetscInt start, NormType ntype, PetscReal *nrm)
114 {
115   PetscInt           i, n, bs;
116   const PetscScalar *x;
117   PetscReal          tnorm;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
121   PetscValidLogicalCollectiveEnum(v, ntype, 3);
122   PetscAssertPointer(nrm, 4);
123   PetscCall(VecGetLocalSize(v, &n));
124   PetscCall(VecGetBlockSize(v, &bs));
125   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
126   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\n Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
127   PetscCall(VecGetArrayRead(v, &x));
128   if (ntype == NORM_2) {
129     PetscScalar sum = 0.0;
130     for (i = start; i < n; i += bs) sum += x[i] * (PetscConj(x[i]));
131     tnorm = PetscRealPart(sum);
132     PetscCall(MPIU_Allreduce(&tnorm, nrm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)v)));
133     *nrm = PetscSqrtReal(*nrm);
134   } else if (ntype == NORM_1) {
135     tnorm = 0.0;
136     for (i = start; i < n; i += bs) tnorm += PetscAbsScalar(x[i]);
137     PetscCall(MPIU_Allreduce(&tnorm, nrm, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)v)));
138   } else if (ntype == NORM_INFINITY) {
139     tnorm = 0.0;
140     for (i = start; i < n; i += bs) {
141       if (PetscAbsScalar(x[i]) > tnorm) tnorm = PetscAbsScalar(x[i]);
142     }
143     PetscCall(MPIU_Allreduce(&tnorm, nrm, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)v)));
144   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
145   PetscCall(VecRestoreArrayRead(v, &x));
146   PetscFunctionReturn(PETSC_SUCCESS);
147 }
148 
149 /*@
150   VecStrideMax - Computes the maximum of subvector of a vector defined
151   by a starting point and a stride and optionally its location.
152 
153   Collective
154 
155   Input Parameters:
156 + v     - the vector
157 - start - starting point of the subvector (defined by a stride)
158 
159   Output Parameters:
160 + idex - the location where the maximum occurred  (pass `NULL` if not required)
161 - nrm  - the maximum value in the subvector
162 
163   Level: advanced
164 
165   Notes:
166   One must call `VecSetBlockSize()` before this routine to set the stride
167   information, or use a vector created from a multicomponent `DMDA`.
168 
169   If xa is the array representing the vector x, then this computes the max
170   of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
171 
172   This is useful for computing, say the maximum of the pressure variable when
173   the pressure is stored (interlaced) with other variables, e.g., density, etc.
174   This will only work if the desire subvector is a stride subvector.
175 
176 .seealso: `Vec`, `VecMax()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`
177 @*/
178 PetscErrorCode VecStrideMax(Vec v, PetscInt start, PetscInt *idex, PetscReal *nrm)
179 {
180   PetscInt           i, n, bs, id = -1;
181   const PetscScalar *x;
182   PetscReal          max = PETSC_MIN_REAL;
183 
184   PetscFunctionBegin;
185   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
186   PetscAssertPointer(nrm, 4);
187   PetscCall(VecGetLocalSize(v, &n));
188   PetscCall(VecGetBlockSize(v, &bs));
189   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
190   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\n Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
191   PetscCall(VecGetArrayRead(v, &x));
192   for (i = start; i < n; i += bs) {
193     if (PetscRealPart(x[i]) > max) {
194       max = PetscRealPart(x[i]);
195       id  = i;
196     }
197   }
198   PetscCall(VecRestoreArrayRead(v, &x));
199 #if defined(PETSC_HAVE_MPIUNI)
200   *nrm = max;
201   if (idex) *idex = id;
202 #else
203   if (!idex) {
204     PetscCall(MPIU_Allreduce(&max, nrm, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)v)));
205   } else {
206     struct {
207       PetscReal v;
208       PetscInt  i;
209     } in, out;
210     PetscInt rstart;
211 
212     PetscCall(VecGetOwnershipRange(v, &rstart, NULL));
213     in.v = max;
214     in.i = rstart + id;
215     PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_REAL_INT, MPIU_MAXLOC, PetscObjectComm((PetscObject)v)));
216     *nrm  = out.v;
217     *idex = out.i;
218   }
219 #endif
220   PetscFunctionReturn(PETSC_SUCCESS);
221 }
222 
223 /*@
224   VecStrideMin - Computes the minimum of subvector of a vector defined
225   by a starting point and a stride and optionally its location.
226 
227   Collective
228 
229   Input Parameters:
230 + v     - the vector
231 - start - starting point of the subvector (defined by a stride)
232 
233   Output Parameters:
234 + idex - the location where the minimum occurred. (pass `NULL` if not required)
235 - nrm  - the minimum value in the subvector
236 
237   Level: advanced
238 
239   Notes:
240   One must call `VecSetBlockSize()` before this routine to set the stride
241   information, or use a vector created from a multicomponent `DMDA`.
242 
243   If xa is the array representing the vector x, then this computes the min
244   of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
245 
246   This is useful for computing, say the minimum of the pressure variable when
247   the pressure is stored (interlaced) with other variables, e.g., density, etc.
248   This will only work if the desire subvector is a stride subvector.
249 
250 .seealso: `Vec`, `VecMin()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMax()`
251 @*/
252 PetscErrorCode VecStrideMin(Vec v, PetscInt start, PetscInt *idex, PetscReal *nrm)
253 {
254   PetscInt           i, n, bs, id = -1;
255   const PetscScalar *x;
256   PetscReal          min = PETSC_MAX_REAL;
257 
258   PetscFunctionBegin;
259   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
260   PetscAssertPointer(nrm, 4);
261   PetscCall(VecGetLocalSize(v, &n));
262   PetscCall(VecGetBlockSize(v, &bs));
263   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
264   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\nHave you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
265   PetscCall(VecGetArrayRead(v, &x));
266   for (i = start; i < n; i += bs) {
267     if (PetscRealPart(x[i]) < min) {
268       min = PetscRealPart(x[i]);
269       id  = i;
270     }
271   }
272   PetscCall(VecRestoreArrayRead(v, &x));
273 #if defined(PETSC_HAVE_MPIUNI)
274   *nrm = min;
275   if (idex) *idex = id;
276 #else
277   if (!idex) {
278     PetscCall(MPIU_Allreduce(&min, nrm, 1, MPIU_REAL, MPIU_MIN, PetscObjectComm((PetscObject)v)));
279   } else {
280     struct {
281       PetscReal v;
282       PetscInt  i;
283     } in, out;
284     PetscInt rstart;
285 
286     PetscCall(VecGetOwnershipRange(v, &rstart, NULL));
287     in.v = min;
288     in.i = rstart + id;
289     PetscCall(MPIU_Allreduce(&in, &out, 1, MPIU_REAL_INT, MPIU_MINLOC, PetscObjectComm((PetscObject)v)));
290     *nrm  = out.v;
291     *idex = out.i;
292   }
293 #endif
294   PetscFunctionReturn(PETSC_SUCCESS);
295 }
296 
297 /*@
298   VecStrideSum - Computes the sum of subvector of a vector defined
299   by a starting point and a stride.
300 
301   Collective
302 
303   Input Parameters:
304 + v     - the vector
305 - start - starting point of the subvector (defined by a stride)
306 
307   Output Parameter:
308 . sum - the sum
309 
310   Level: advanced
311 
312   Notes:
313   One must call `VecSetBlockSize()` before this routine to set the stride
314   information, or use a vector created from a multicomponent `DMDA`.
315 
316   If x is the array representing the vector x then this computes the sum
317   of the array (x[start],x[start+stride],x[start+2*stride], ....)
318 
319 .seealso: `Vec`, `VecSum()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
320 @*/
321 PetscErrorCode VecStrideSum(Vec v, PetscInt start, PetscScalar *sum)
322 {
323   PetscInt           i, n, bs;
324   const PetscScalar *x;
325   PetscScalar        local_sum = 0.0;
326 
327   PetscFunctionBegin;
328   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
329   PetscAssertPointer(sum, 3);
330   PetscCall(VecGetLocalSize(v, &n));
331   PetscCall(VecGetBlockSize(v, &bs));
332   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
333   PetscCheck(start < bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\n Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start, bs);
334   PetscCall(VecGetArrayRead(v, &x));
335   for (i = start; i < n; i += bs) local_sum += x[i];
336   PetscCall(MPIU_Allreduce(&local_sum, sum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)v)));
337   PetscCall(VecRestoreArrayRead(v, &x));
338   PetscFunctionReturn(PETSC_SUCCESS);
339 }
340 
341 /*@
342   VecStrideScaleAll - Scales the subvectors of a vector defined
343   by a starting point and a stride.
344 
345   Logically Collective
346 
347   Input Parameters:
348 + v      - the vector
349 - scales - values to multiply each subvector entry by
350 
351   Level: advanced
352 
353   Notes:
354   One must call `VecSetBlockSize()` before this routine to set the stride
355   information, or use a vector created from a multicomponent `DMDA`.
356 
357   The dimension of scales must be the same as the vector block size
358 
359 .seealso: `Vec`, `VecNorm()`, `VecStrideScale()`, `VecScale()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
360 @*/
361 PetscErrorCode VecStrideScaleAll(Vec v, const PetscScalar *scales)
362 {
363   PetscInt     i, j, n, bs;
364   PetscScalar *x;
365 
366   PetscFunctionBegin;
367   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
368   PetscAssertPointer(scales, 2);
369   PetscCall(VecGetLocalSize(v, &n));
370   PetscCall(VecGetBlockSize(v, &bs));
371   PetscCall(VecGetArray(v, &x));
372   /* need to provide optimized code for each bs */
373   for (i = 0; i < n; i += bs) {
374     for (j = 0; j < bs; j++) x[i + j] *= scales[j];
375   }
376   PetscCall(VecRestoreArray(v, &x));
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 
380 /*@
381   VecStrideNormAll - Computes the norms of subvectors of a vector defined
382   by a starting point and a stride.
383 
384   Collective
385 
386   Input Parameters:
387 + v     - the vector
388 - ntype - type of norm, one of `NORM_1`, `NORM_2`, `NORM_INFINITY`
389 
390   Output Parameter:
391 . nrm - the norms
392 
393   Level: advanced
394 
395   Notes:
396   One must call `VecSetBlockSize()` before this routine to set the stride
397   information, or use a vector created from a multicomponent `DMDA`.
398 
399   If x is the array representing the vector x then this computes the norm
400   of the array (x[start],x[start+stride],x[start+2*stride], ....) for each start < stride
401 
402   The dimension of nrm must be the same as the vector block size
403 
404   This will only work if the desire subvector is a stride subvector
405 
406 .seealso: `Vec`, `VecNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
407 @*/
408 PetscErrorCode VecStrideNormAll(Vec v, NormType ntype, PetscReal nrm[])
409 {
410   PetscInt           i, j, n, bs;
411   const PetscScalar *x;
412   PetscReal          tnorm[128];
413   MPI_Comm           comm;
414 
415   PetscFunctionBegin;
416   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
417   PetscValidLogicalCollectiveEnum(v, ntype, 2);
418   PetscAssertPointer(nrm, 3);
419   PetscCall(VecGetLocalSize(v, &n));
420   PetscCall(VecGetArrayRead(v, &x));
421   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
422 
423   PetscCall(VecGetBlockSize(v, &bs));
424   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
425 
426   if (ntype == NORM_2) {
427     PetscScalar sum[128];
428     for (j = 0; j < bs; j++) sum[j] = 0.0;
429     for (i = 0; i < n; i += bs) {
430       for (j = 0; j < bs; j++) sum[j] += x[i + j] * (PetscConj(x[i + j]));
431     }
432     for (j = 0; j < bs; j++) tnorm[j] = PetscRealPart(sum[j]);
433 
434     PetscCall(MPIU_Allreduce(tnorm, nrm, bs, MPIU_REAL, MPIU_SUM, comm));
435     for (j = 0; j < bs; j++) nrm[j] = PetscSqrtReal(nrm[j]);
436   } else if (ntype == NORM_1) {
437     for (j = 0; j < bs; j++) tnorm[j] = 0.0;
438 
439     for (i = 0; i < n; i += bs) {
440       for (j = 0; j < bs; j++) tnorm[j] += PetscAbsScalar(x[i + j]);
441     }
442 
443     PetscCall(MPIU_Allreduce(tnorm, nrm, bs, MPIU_REAL, MPIU_SUM, comm));
444   } else if (ntype == NORM_INFINITY) {
445     PetscReal tmp;
446     for (j = 0; j < bs; j++) tnorm[j] = 0.0;
447 
448     for (i = 0; i < n; i += bs) {
449       for (j = 0; j < bs; j++) {
450         if ((tmp = PetscAbsScalar(x[i + j])) > tnorm[j]) tnorm[j] = tmp;
451         /* check special case of tmp == NaN */
452         if (tmp != tmp) {
453           tnorm[j] = tmp;
454           break;
455         }
456       }
457     }
458     PetscCall(MPIU_Allreduce(tnorm, nrm, bs, MPIU_REAL, MPIU_MAX, comm));
459   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
460   PetscCall(VecRestoreArrayRead(v, &x));
461   PetscFunctionReturn(PETSC_SUCCESS);
462 }
463 
464 /*@
465   VecStrideMaxAll - Computes the maximums of subvectors of a vector defined
466   by a starting point and a stride and optionally its location.
467 
468   Collective
469 
470   Input Parameter:
471 . v - the vector
472 
473   Output Parameters:
474 + idex - the location where the maximum occurred (not supported, pass `NULL`,
475            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
476 - nrm  - the maximum values of each subvector
477 
478   Level: advanced
479 
480   Notes:
481   One must call `VecSetBlockSize()` before this routine to set the stride
482   information, or use a vector created from a multicomponent `DMDA`.
483 
484   The dimension of nrm must be the same as the vector block size
485 
486 .seealso: `Vec`, `VecMax()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`
487 @*/
488 PetscErrorCode VecStrideMaxAll(Vec v, PetscInt idex[], PetscReal nrm[])
489 {
490   PetscInt           i, j, n, bs;
491   const PetscScalar *x;
492   PetscReal          max[128], tmp;
493   MPI_Comm           comm;
494 
495   PetscFunctionBegin;
496   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
497   PetscAssertPointer(nrm, 3);
498   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");
499   PetscCall(VecGetLocalSize(v, &n));
500   PetscCall(VecGetArrayRead(v, &x));
501   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
502 
503   PetscCall(VecGetBlockSize(v, &bs));
504   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
505 
506   if (!n) {
507     for (j = 0; j < bs; j++) max[j] = PETSC_MIN_REAL;
508   } else {
509     for (j = 0; j < bs; j++) max[j] = PetscRealPart(x[j]);
510 
511     for (i = bs; i < n; i += bs) {
512       for (j = 0; j < bs; j++) {
513         if ((tmp = PetscRealPart(x[i + j])) > max[j]) max[j] = tmp;
514       }
515     }
516   }
517   PetscCall(MPIU_Allreduce(max, nrm, bs, MPIU_REAL, MPIU_MAX, comm));
518 
519   PetscCall(VecRestoreArrayRead(v, &x));
520   PetscFunctionReturn(PETSC_SUCCESS);
521 }
522 
523 /*@
524   VecStrideMinAll - Computes the minimum of subvector of a vector defined
525   by a starting point and a stride and optionally its location.
526 
527   Collective
528 
529   Input Parameter:
530 . v - the vector
531 
532   Output Parameters:
533 + idex - the location where the minimum occurred (not supported, pass `NULL`,
534            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
535 - nrm  - the minimums of each subvector
536 
537   Level: advanced
538 
539   Notes:
540   One must call `VecSetBlockSize()` before this routine to set the stride
541   information, or use a vector created from a multicomponent `DMDA`.
542 
543   The dimension of `nrm` must be the same as the vector block size
544 
545 .seealso: `Vec`, `VecMin()`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMax()`
546 @*/
547 PetscErrorCode VecStrideMinAll(Vec v, PetscInt idex[], PetscReal nrm[])
548 {
549   PetscInt           i, n, bs, j;
550   const PetscScalar *x;
551   PetscReal          min[128], tmp;
552   MPI_Comm           comm;
553 
554   PetscFunctionBegin;
555   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
556   PetscAssertPointer(nrm, 3);
557   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");
558   PetscCall(VecGetLocalSize(v, &n));
559   PetscCall(VecGetArrayRead(v, &x));
560   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
561 
562   PetscCall(VecGetBlockSize(v, &bs));
563   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
564 
565   if (!n) {
566     for (j = 0; j < bs; j++) min[j] = PETSC_MAX_REAL;
567   } else {
568     for (j = 0; j < bs; j++) min[j] = PetscRealPart(x[j]);
569 
570     for (i = bs; i < n; i += bs) {
571       for (j = 0; j < bs; j++) {
572         if ((tmp = PetscRealPart(x[i + j])) < min[j]) min[j] = tmp;
573       }
574     }
575   }
576   PetscCall(MPIU_Allreduce(min, nrm, bs, MPIU_REAL, MPIU_MIN, comm));
577 
578   PetscCall(VecRestoreArrayRead(v, &x));
579   PetscFunctionReturn(PETSC_SUCCESS);
580 }
581 
582 /*@
583   VecStrideSumAll - Computes the sums of subvectors of a vector defined by a stride.
584 
585   Collective
586 
587   Input Parameter:
588 . v - the vector
589 
590   Output Parameter:
591 . sums - the sums
592 
593   Level: advanced
594 
595   Notes:
596   One must call `VecSetBlockSize()` before this routine to set the stride
597   information, or use a vector created from a multicomponent `DMDA`.
598 
599   If x is the array representing the vector x then this computes the sum
600   of the array (x[start],x[start+stride],x[start+2*stride], ....) for each start < stride
601 
602 .seealso: `Vec`, `VecSum()`, `VecStrideGather()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`
603 @*/
604 PetscErrorCode VecStrideSumAll(Vec v, PetscScalar sums[])
605 {
606   PetscInt           i, j, n, bs;
607   const PetscScalar *x;
608   PetscScalar        local_sums[128];
609   MPI_Comm           comm;
610 
611   PetscFunctionBegin;
612   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
613   PetscAssertPointer(sums, 2);
614   PetscCall(VecGetLocalSize(v, &n));
615   PetscCall(VecGetArrayRead(v, &x));
616   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
617 
618   PetscCall(VecGetBlockSize(v, &bs));
619   PetscCheck(bs <= 128, comm, PETSC_ERR_SUP, "Currently supports only blocksize up to 128");
620 
621   for (j = 0; j < bs; j++) local_sums[j] = 0.0;
622   for (i = 0; i < n; i += bs) {
623     for (j = 0; j < bs; j++) local_sums[j] += x[i + j];
624   }
625   PetscCall(MPIU_Allreduce(local_sums, sums, bs, MPIU_SCALAR, MPIU_SUM, comm));
626 
627   PetscCall(VecRestoreArrayRead(v, &x));
628   PetscFunctionReturn(PETSC_SUCCESS);
629 }
630 
631 /*----------------------------------------------------------------------------------------------*/
632 /*@
633   VecStrideGatherAll - Gathers all the single components from a multi-component vector into
634   separate vectors.
635 
636   Collective
637 
638   Input Parameters:
639 + v    - the vector
640 - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
641 
642   Output Parameter:
643 . s - the location where the subvectors are stored
644 
645   Level: advanced
646 
647   Notes:
648   One must call `VecSetBlockSize()` before this routine to set the stride
649   information, or use a vector created from a multicomponent `DMDA`.
650 
651   If x is the array representing the vector x then this gathers
652   the arrays (x[start],x[start+stride],x[start+2*stride], ....)
653   for start=0,1,2,...bs-1
654 
655   The parallel layout of the vector and the subvector must be the same;
656   i.e., nlocal of v = stride*(nlocal of s)
657 
658   Not optimized; could be easily
659 
660 .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGather()`,
661           `VecStrideScatterAll()`
662 @*/
663 PetscErrorCode VecStrideGatherAll(Vec v, Vec s[], InsertMode addv)
664 {
665   PetscInt           i, n, n2, bs, j, k, *bss = NULL, nv, jj, nvc;
666   PetscScalar      **y;
667   const PetscScalar *x;
668 
669   PetscFunctionBegin;
670   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
671   PetscAssertPointer(s, 2);
672   PetscValidHeaderSpecific(*s, VEC_CLASSID, 2);
673   PetscCall(VecGetLocalSize(v, &n));
674   PetscCall(VecGetLocalSize(s[0], &n2));
675   PetscCall(VecGetArrayRead(v, &x));
676   PetscCall(VecGetBlockSize(v, &bs));
677   PetscCheck(bs > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Input vector does not have a valid blocksize set");
678 
679   PetscCall(PetscMalloc2(bs, &y, bs, &bss));
680   nv  = 0;
681   nvc = 0;
682   for (i = 0; i < bs; i++) {
683     PetscCall(VecGetBlockSize(s[i], &bss[i]));
684     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
685     PetscCall(VecGetArray(s[i], &y[i]));
686     nvc += bss[i];
687     nv++;
688     PetscCheck(nvc <= bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of subvectors in subvectors > number of vectors in main vector");
689     if (nvc == bs) break;
690   }
691 
692   n = n / bs;
693 
694   jj = 0;
695   if (addv == INSERT_VALUES) {
696     for (j = 0; j < nv; j++) {
697       for (k = 0; k < bss[j]; k++) {
698         for (i = 0; i < n; i++) y[j][i * bss[j] + k] = x[bs * i + jj + k];
699       }
700       jj += bss[j];
701     }
702   } else if (addv == ADD_VALUES) {
703     for (j = 0; j < nv; j++) {
704       for (k = 0; k < bss[j]; k++) {
705         for (i = 0; i < n; i++) y[j][i * bss[j] + k] += x[bs * i + jj + k];
706       }
707       jj += bss[j];
708     }
709 #if !defined(PETSC_USE_COMPLEX)
710   } else if (addv == MAX_VALUES) {
711     for (j = 0; j < nv; j++) {
712       for (k = 0; k < bss[j]; k++) {
713         for (i = 0; i < n; i++) y[j][i * bss[j] + k] = PetscMax(y[j][i * bss[j] + k], x[bs * i + jj + k]);
714       }
715       jj += bss[j];
716     }
717 #endif
718   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
719 
720   PetscCall(VecRestoreArrayRead(v, &x));
721   for (i = 0; i < nv; i++) PetscCall(VecRestoreArray(s[i], &y[i]));
722 
723   PetscCall(PetscFree2(y, bss));
724   PetscFunctionReturn(PETSC_SUCCESS);
725 }
726 
727 /*@
728   VecStrideScatterAll - Scatters all the single components from separate vectors into
729   a multi-component vector.
730 
731   Collective
732 
733   Input Parameters:
734 + s    - the location where the subvectors are stored
735 - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
736 
737   Output Parameter:
738 . v - the multicomponent vector
739 
740   Level: advanced
741 
742   Notes:
743   One must call `VecSetBlockSize()` before this routine to set the stride
744   information, or use a vector created from a multicomponent `DMDA`.
745 
746   The parallel layout of the vector and the subvector must be the same;
747   i.e., nlocal of v = stride*(nlocal of s)
748 
749   Not optimized; could be easily
750 
751 .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGather()`,
752 
753 @*/
754 PetscErrorCode VecStrideScatterAll(Vec s[], Vec v, InsertMode addv)
755 {
756   PetscInt            i, n, n2, bs, j, jj, k, *bss = NULL, nv, nvc;
757   PetscScalar        *x;
758   PetscScalar const **y;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(v, VEC_CLASSID, 2);
762   PetscAssertPointer(s, 1);
763   PetscValidHeaderSpecific(*s, VEC_CLASSID, 1);
764   PetscCall(VecGetLocalSize(v, &n));
765   PetscCall(VecGetLocalSize(s[0], &n2));
766   PetscCall(VecGetArray(v, &x));
767   PetscCall(VecGetBlockSize(v, &bs));
768   PetscCheck(bs > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Input vector does not have a valid blocksize set");
769 
770   PetscCall(PetscMalloc2(bs, (PetscScalar ***)&y, bs, &bss));
771   nv  = 0;
772   nvc = 0;
773   for (i = 0; i < bs; i++) {
774     PetscCall(VecGetBlockSize(s[i], &bss[i]));
775     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
776     PetscCall(VecGetArrayRead(s[i], &y[i]));
777     nvc += bss[i];
778     nv++;
779     PetscCheck(nvc <= bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of subvectors in subvectors > number of vectors in main vector");
780     if (nvc == bs) break;
781   }
782 
783   n = n / bs;
784 
785   jj = 0;
786   if (addv == INSERT_VALUES) {
787     for (j = 0; j < nv; j++) {
788       for (k = 0; k < bss[j]; k++) {
789         for (i = 0; i < n; i++) x[bs * i + jj + k] = y[j][i * bss[j] + k];
790       }
791       jj += bss[j];
792     }
793   } else if (addv == ADD_VALUES) {
794     for (j = 0; j < nv; j++) {
795       for (k = 0; k < bss[j]; k++) {
796         for (i = 0; i < n; i++) x[bs * i + jj + k] += y[j][i * bss[j] + k];
797       }
798       jj += bss[j];
799     }
800 #if !defined(PETSC_USE_COMPLEX)
801   } else if (addv == MAX_VALUES) {
802     for (j = 0; j < nv; j++) {
803       for (k = 0; k < bss[j]; k++) {
804         for (i = 0; i < n; i++) x[bs * i + jj + k] = PetscMax(x[bs * i + jj + k], y[j][i * bss[j] + k]);
805       }
806       jj += bss[j];
807     }
808 #endif
809   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
810 
811   PetscCall(VecRestoreArray(v, &x));
812   for (i = 0; i < nv; i++) PetscCall(VecRestoreArrayRead(s[i], &y[i]));
813   PetscCall(PetscFree2(*(PetscScalar ***)&y, bss));
814   PetscFunctionReturn(PETSC_SUCCESS);
815 }
816 
817 /*@
818   VecStrideGather - Gathers a single component from a multi-component vector into
819   another vector.
820 
821   Collective
822 
823   Input Parameters:
824 + v     - the vector
825 . start - starting point of the subvector (defined by a stride)
826 - addv  - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
827 
828   Output Parameter:
829 . s - the location where the subvector is stored
830 
831   Level: advanced
832 
833   Notes:
834   One must call `VecSetBlockSize()` before this routine to set the stride
835   information, or use a vector created from a multicomponent `DMDA`.
836 
837   If x is the array representing the vector x then this gathers
838   the array (x[start],x[start+stride],x[start+2*stride], ....)
839 
840   The parallel layout of the vector and the subvector must be the same;
841   i.e., nlocal of v = stride*(nlocal of s)
842 
843   Not optimized; could be easily
844 
845 .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
846           `VecStrideScatterAll()`
847 @*/
848 PetscErrorCode VecStrideGather(Vec v, PetscInt start, Vec s, InsertMode addv)
849 {
850   PetscFunctionBegin;
851   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
852   PetscValidHeaderSpecific(s, VEC_CLASSID, 3);
853   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
854   PetscCheck(start < v->map->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\n Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start,
855              v->map->bs);
856   PetscUseTypeMethod(v, stridegather, start, s, addv);
857   PetscFunctionReturn(PETSC_SUCCESS);
858 }
859 
860 /*@
861   VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
862 
863   Collective
864 
865   Input Parameters:
866 + s     - the single-component vector
867 . start - starting point of the subvector (defined by a stride)
868 - addv  - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
869 
870   Output Parameter:
871 . v - the location where the subvector is scattered (the multi-component vector)
872 
873   Level: advanced
874 
875   Notes:
876   One must call `VecSetBlockSize()` on the multi-component vector before this
877   routine to set the stride  information, or use a vector created from a multicomponent `DMDA`.
878 
879   The parallel layout of the vector and the subvector must be the same;
880   i.e., nlocal of v = stride*(nlocal of s)
881 
882   Not optimized; could be easily
883 
884 .seealso: `Vec`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
885           `VecStrideScatterAll()`, `VecStrideSubSetScatter()`, `VecStrideSubSetGather()`
886 @*/
887 PetscErrorCode VecStrideScatter(Vec s, PetscInt start, Vec v, InsertMode addv)
888 {
889   PetscFunctionBegin;
890   PetscValidHeaderSpecific(s, VEC_CLASSID, 1);
891   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
892   PetscCheck(start >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative start %" PetscInt_FMT, start);
893   PetscCheck(start < v->map->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Start of stride subvector (%" PetscInt_FMT ") is too large for stride\n Have you set the vector blocksize (%" PetscInt_FMT ") correctly with VecSetBlockSize()?", start,
894              v->map->bs);
895   PetscCall((*v->ops->stridescatter)(s, start, v, addv));
896   PetscFunctionReturn(PETSC_SUCCESS);
897 }
898 
899 /*@
900   VecStrideSubSetGather - Gathers a subset of components from a multi-component vector into
901   another vector.
902 
903   Collective
904 
905   Input Parameters:
906 + v    - the vector
907 . nidx - the number of indices
908 . idxv - the indices of the components 0 <= idxv[0] ...idxv[nidx-1] < bs(v), they need not be sorted
909 . 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`
910 - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
911 
912   Output Parameter:
913 . s - the location where the subvector is stored
914 
915   Level: advanced
916 
917   Notes:
918   One must call `VecSetBlockSize()` on both vectors before this routine to set the stride
919   information, or use a vector created from a multicomponent `DMDA`.
920 
921   The parallel layout of the vector and the subvector must be the same;
922 
923   Not optimized; could be easily
924 
925 .seealso: `Vec`, `VecStrideNorm()`, `VecStrideScatter()`, `VecStrideGather()`, `VecStrideSubSetScatter()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
926           `VecStrideScatterAll()`
927 @*/
928 PetscErrorCode VecStrideSubSetGather(Vec v, PetscInt nidx, const PetscInt idxv[], const PetscInt idxs[], Vec s, InsertMode addv)
929 {
930   PetscFunctionBegin;
931   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
932   PetscValidHeaderSpecific(s, VEC_CLASSID, 5);
933   if (nidx == PETSC_DETERMINE) nidx = s->map->bs;
934   PetscUseTypeMethod(v, stridesubsetgather, nidx, idxv, idxs, s, addv);
935   PetscFunctionReturn(PETSC_SUCCESS);
936 }
937 
938 /*@
939   VecStrideSubSetScatter - Scatters components from a vector into a subset of components of a multi-component vector.
940 
941   Collective
942 
943   Input Parameters:
944 + s    - the smaller-component vector
945 . nidx - the number of indices in idx
946 . 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`
947 . idxv - the indices of the components in the larger-component vector, 0 <= idx[0] ...idx[nidx-1] < bs(v) they need not be sorted
948 - addv - one of `ADD_VALUES`, `INSERT_VALUES`, `MAX_VALUES`
949 
950   Output Parameter:
951 . v - the location where the subvector is into scattered (the multi-component vector)
952 
953   Level: advanced
954 
955   Notes:
956   One must call `VecSetBlockSize()` on the vectors before this
957   routine to set the stride  information, or use a vector created from a multicomponent `DMDA`.
958 
959   The parallel layout of the vector and the subvector must be the same;
960 
961   Not optimized; could be easily
962 
963 .seealso: `Vec`, `VecStrideNorm()`, `VecStrideGather()`, `VecStrideSubSetGather()`, `VecStrideMin()`, `VecStrideMax()`, `VecStrideGatherAll()`,
964           `VecStrideScatterAll()`
965 @*/
966 PetscErrorCode VecStrideSubSetScatter(Vec s, PetscInt nidx, const PetscInt idxs[], const PetscInt idxv[], Vec v, InsertMode addv)
967 {
968   PetscFunctionBegin;
969   PetscValidHeaderSpecific(s, VEC_CLASSID, 1);
970   PetscValidHeaderSpecific(v, VEC_CLASSID, 5);
971   if (nidx == PETSC_DETERMINE) nidx = s->map->bs;
972   PetscCall((*v->ops->stridesubsetscatter)(s, nidx, idxs, idxv, v, addv));
973   PetscFunctionReturn(PETSC_SUCCESS);
974 }
975 
976 PetscErrorCode VecStrideGather_Default(Vec v, PetscInt start, Vec s, InsertMode addv)
977 {
978   PetscInt           i, n, bs, ns;
979   const PetscScalar *x;
980   PetscScalar       *y;
981 
982   PetscFunctionBegin;
983   PetscCall(VecGetLocalSize(v, &n));
984   PetscCall(VecGetLocalSize(s, &ns));
985   PetscCall(VecGetArrayRead(v, &x));
986   PetscCall(VecGetArray(s, &y));
987 
988   bs = v->map->bs;
989   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);
990   x += start;
991   n = n / bs;
992 
993   if (addv == INSERT_VALUES) {
994     for (i = 0; i < n; i++) y[i] = x[bs * i];
995   } else if (addv == ADD_VALUES) {
996     for (i = 0; i < n; i++) y[i] += x[bs * i];
997 #if !defined(PETSC_USE_COMPLEX)
998   } else if (addv == MAX_VALUES) {
999     for (i = 0; i < n; i++) y[i] = PetscMax(y[i], x[bs * i]);
1000 #endif
1001   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1002 
1003   PetscCall(VecRestoreArrayRead(v, &x));
1004   PetscCall(VecRestoreArray(s, &y));
1005   PetscFunctionReturn(PETSC_SUCCESS);
1006 }
1007 
1008 PetscErrorCode VecStrideScatter_Default(Vec s, PetscInt start, Vec v, InsertMode addv)
1009 {
1010   PetscInt           i, n, bs, ns;
1011   PetscScalar       *x;
1012   const PetscScalar *y;
1013 
1014   PetscFunctionBegin;
1015   PetscCall(VecGetLocalSize(v, &n));
1016   PetscCall(VecGetLocalSize(s, &ns));
1017   PetscCall(VecGetArray(v, &x));
1018   PetscCall(VecGetArrayRead(s, &y));
1019 
1020   bs = v->map->bs;
1021   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);
1022   x += start;
1023   n = n / bs;
1024 
1025   if (addv == INSERT_VALUES) {
1026     for (i = 0; i < n; i++) x[bs * i] = y[i];
1027   } else if (addv == ADD_VALUES) {
1028     for (i = 0; i < n; i++) x[bs * i] += y[i];
1029 #if !defined(PETSC_USE_COMPLEX)
1030   } else if (addv == MAX_VALUES) {
1031     for (i = 0; i < n; i++) x[bs * i] = PetscMax(y[i], x[bs * i]);
1032 #endif
1033   } else SETERRQ(PetscObjectComm((PetscObject)s), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1034 
1035   PetscCall(VecRestoreArray(v, &x));
1036   PetscCall(VecRestoreArrayRead(s, &y));
1037   PetscFunctionReturn(PETSC_SUCCESS);
1038 }
1039 
1040 PetscErrorCode VecStrideSubSetGather_Default(Vec v, PetscInt nidx, const PetscInt idxv[], const PetscInt idxs[], Vec s, InsertMode addv)
1041 {
1042   PetscInt           i, j, n, bs, bss, ns;
1043   const PetscScalar *x;
1044   PetscScalar       *y;
1045 
1046   PetscFunctionBegin;
1047   PetscCall(VecGetLocalSize(v, &n));
1048   PetscCall(VecGetLocalSize(s, &ns));
1049   PetscCall(VecGetArrayRead(v, &x));
1050   PetscCall(VecGetArray(s, &y));
1051 
1052   bs  = v->map->bs;
1053   bss = s->map->bs;
1054   n   = n / bs;
1055 
1056   if (PetscDefined(USE_DEBUG)) {
1057     PetscCheck(n == ns / bss, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout of vectors");
1058     for (j = 0; j < nidx; j++) {
1059       PetscCheck(idxv[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is negative", j, idxv[j]);
1060       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);
1061     }
1062     PetscCheck(idxs || bss == nidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide idxs when not gathering into all locations");
1063   }
1064 
1065   if (addv == INSERT_VALUES) {
1066     if (!idxs) {
1067       for (i = 0; i < n; i++) {
1068         for (j = 0; j < bss; j++) y[bss * i + j] = x[bs * i + idxv[j]];
1069       }
1070     } else {
1071       for (i = 0; i < n; i++) {
1072         for (j = 0; j < bss; j++) y[bss * i + idxs[j]] = x[bs * i + idxv[j]];
1073       }
1074     }
1075   } else if (addv == ADD_VALUES) {
1076     if (!idxs) {
1077       for (i = 0; i < n; i++) {
1078         for (j = 0; j < bss; j++) y[bss * i + j] += x[bs * i + idxv[j]];
1079       }
1080     } else {
1081       for (i = 0; i < n; i++) {
1082         for (j = 0; j < bss; j++) y[bss * i + idxs[j]] += x[bs * i + idxv[j]];
1083       }
1084     }
1085 #if !defined(PETSC_USE_COMPLEX)
1086   } else if (addv == MAX_VALUES) {
1087     if (!idxs) {
1088       for (i = 0; i < n; i++) {
1089         for (j = 0; j < bss; j++) y[bss * i + j] = PetscMax(y[bss * i + j], x[bs * i + idxv[j]]);
1090       }
1091     } else {
1092       for (i = 0; i < n; i++) {
1093         for (j = 0; j < bss; j++) y[bss * i + idxs[j]] = PetscMax(y[bss * i + idxs[j]], x[bs * i + idxv[j]]);
1094       }
1095     }
1096 #endif
1097   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1098 
1099   PetscCall(VecRestoreArrayRead(v, &x));
1100   PetscCall(VecRestoreArray(s, &y));
1101   PetscFunctionReturn(PETSC_SUCCESS);
1102 }
1103 
1104 PetscErrorCode VecStrideSubSetScatter_Default(Vec s, PetscInt nidx, const PetscInt idxs[], const PetscInt idxv[], Vec v, InsertMode addv)
1105 {
1106   PetscInt           j, i, n, bs, ns, bss;
1107   PetscScalar       *x;
1108   const PetscScalar *y;
1109 
1110   PetscFunctionBegin;
1111   PetscCall(VecGetLocalSize(v, &n));
1112   PetscCall(VecGetLocalSize(s, &ns));
1113   PetscCall(VecGetArray(v, &x));
1114   PetscCall(VecGetArrayRead(s, &y));
1115 
1116   bs  = v->map->bs;
1117   bss = s->map->bs;
1118   n   = n / bs;
1119 
1120   if (PetscDefined(USE_DEBUG)) {
1121     PetscCheck(n == ns / bss, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout of vectors");
1122     for (j = 0; j < bss; j++) {
1123       if (idxs) {
1124         PetscCheck(idxs[j] >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "idx[%" PetscInt_FMT "] %" PetscInt_FMT " is negative", j, idxs[j]);
1125         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);
1126       }
1127     }
1128     PetscCheck(idxs || bss == nidx, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must provide idxs when not scattering from all locations");
1129   }
1130 
1131   if (addv == INSERT_VALUES) {
1132     if (!idxs) {
1133       for (i = 0; i < n; i++) {
1134         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = y[bss * i + j];
1135       }
1136     } else {
1137       for (i = 0; i < n; i++) {
1138         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = y[bss * i + idxs[j]];
1139       }
1140     }
1141   } else if (addv == ADD_VALUES) {
1142     if (!idxs) {
1143       for (i = 0; i < n; i++) {
1144         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] += y[bss * i + j];
1145       }
1146     } else {
1147       for (i = 0; i < n; i++) {
1148         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] += y[bss * i + idxs[j]];
1149       }
1150     }
1151 #if !defined(PETSC_USE_COMPLEX)
1152   } else if (addv == MAX_VALUES) {
1153     if (!idxs) {
1154       for (i = 0; i < n; i++) {
1155         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = PetscMax(y[bss * i + j], x[bs * i + idxv[j]]);
1156       }
1157     } else {
1158       for (i = 0; i < n; i++) {
1159         for (j = 0; j < bss; j++) x[bs * i + idxv[j]] = PetscMax(y[bss * i + idxs[j]], x[bs * i + idxv[j]]);
1160       }
1161     }
1162 #endif
1163   } else SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown norm type");
1164 
1165   PetscCall(VecRestoreArray(v, &x));
1166   PetscCall(VecRestoreArrayRead(s, &y));
1167   PetscFunctionReturn(PETSC_SUCCESS);
1168 }
1169 
1170 static PetscErrorCode VecApplyUnary_Private(Vec v, PetscErrorCode (*unary_op)(Vec), PetscScalar (*UnaryFunc)(PetscScalar))
1171 {
1172   PetscFunctionBegin;
1173   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1174   PetscCall(VecSetErrorIfLocked(v, 1));
1175   if (unary_op) {
1176     PetscValidFunction(unary_op, 2);
1177     PetscCall((*unary_op)(v));
1178   } else {
1179     PetscInt     n;
1180     PetscScalar *x;
1181 
1182     PetscValidFunction(UnaryFunc, 3);
1183     PetscCall(VecGetLocalSize(v, &n));
1184     PetscCall(VecGetArray(v, &x));
1185     for (PetscInt i = 0; i < n; ++i) x[i] = UnaryFunc(x[i]);
1186     PetscCall(VecRestoreArray(v, &x));
1187   }
1188   PetscFunctionReturn(PETSC_SUCCESS);
1189 }
1190 
1191 static PetscScalar ScalarReciprocal_Fn(PetscScalar x)
1192 {
1193   const PetscScalar zero = 0.0;
1194 
1195   return x == zero ? zero : ((PetscScalar)1.0) / x;
1196 }
1197 
1198 PetscErrorCode VecReciprocal_Default(Vec v)
1199 {
1200   PetscFunctionBegin;
1201   PetscCall(VecApplyUnary_Private(v, NULL, ScalarReciprocal_Fn));
1202   PetscFunctionReturn(PETSC_SUCCESS);
1203 }
1204 
1205 static PetscScalar ScalarExp_Fn(PetscScalar x)
1206 {
1207   return PetscExpScalar(x);
1208 }
1209 
1210 /*@
1211   VecExp - Replaces each component of a vector by e^x_i
1212 
1213   Not Collective
1214 
1215   Input Parameter:
1216 . v - The vector
1217 
1218   Output Parameter:
1219 . v - The vector of exponents
1220 
1221   Level: beginner
1222 
1223 .seealso: `Vec`, `VecLog()`, `VecAbs()`, `VecSqrtAbs()`, `VecReciprocal()`
1224 
1225 @*/
1226 PetscErrorCode VecExp(Vec v)
1227 {
1228   PetscFunctionBegin;
1229   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1230   PetscCall(VecApplyUnary_Private(v, v->ops->exp, ScalarExp_Fn));
1231   PetscFunctionReturn(PETSC_SUCCESS);
1232 }
1233 
1234 static PetscScalar ScalarLog_Fn(PetscScalar x)
1235 {
1236   return PetscLogScalar(x);
1237 }
1238 
1239 /*@
1240   VecLog - Replaces each component of a vector by log(x_i), the natural logarithm
1241 
1242   Not Collective
1243 
1244   Input Parameter:
1245 . v - The vector
1246 
1247   Output Parameter:
1248 . v - The vector of logs
1249 
1250   Level: beginner
1251 
1252 .seealso: `Vec`, `VecExp()`, `VecAbs()`, `VecSqrtAbs()`, `VecReciprocal()`
1253 
1254 @*/
1255 PetscErrorCode VecLog(Vec v)
1256 {
1257   PetscFunctionBegin;
1258   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1259   PetscCall(VecApplyUnary_Private(v, v->ops->log, ScalarLog_Fn));
1260   PetscFunctionReturn(PETSC_SUCCESS);
1261 }
1262 
1263 static PetscScalar ScalarAbs_Fn(PetscScalar x)
1264 {
1265   return PetscAbsScalar(x);
1266 }
1267 
1268 /*@
1269   VecAbs - Replaces every element in a vector with its absolute value.
1270 
1271   Logically Collective
1272 
1273   Input Parameter:
1274 . v - the vector
1275 
1276   Level: intermediate
1277 
1278 .seealso: `Vec`, `VecExp()`, `VecSqrtAbs()`, `VecReciprocal()`, `VecLog()`
1279 @*/
1280 PetscErrorCode VecAbs(Vec v)
1281 {
1282   PetscFunctionBegin;
1283   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1284   PetscCall(VecApplyUnary_Private(v, v->ops->abs, ScalarAbs_Fn));
1285   PetscFunctionReturn(PETSC_SUCCESS);
1286 }
1287 
1288 static PetscScalar ScalarSqrtAbs_Fn(PetscScalar x)
1289 {
1290   return PetscSqrtScalar(ScalarAbs_Fn(x));
1291 }
1292 
1293 /*@
1294   VecSqrtAbs - Replaces each component of a vector by the square root of its magnitude.
1295 
1296   Not Collective
1297 
1298   Input Parameter:
1299 . v - The vector
1300 
1301   Level: beginner
1302 
1303   Note:
1304   The actual function is sqrt(|x_i|)
1305 
1306 .seealso: `Vec`, `VecLog()`, `VecExp()`, `VecReciprocal()`, `VecAbs()`
1307 
1308 @*/
1309 PetscErrorCode VecSqrtAbs(Vec v)
1310 {
1311   PetscFunctionBegin;
1312   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1313   PetscCall(VecApplyUnary_Private(v, v->ops->sqrt, ScalarSqrtAbs_Fn));
1314   PetscFunctionReturn(PETSC_SUCCESS);
1315 }
1316 
1317 static PetscScalar ScalarImaginaryPart_Fn(PetscScalar x)
1318 {
1319   const PetscReal imag = PetscImaginaryPart(x);
1320 
1321 #if PetscDefined(USE_COMPLEX)
1322   return PetscCMPLX(imag, 0.0);
1323 #else
1324   return imag;
1325 #endif
1326 }
1327 
1328 /*@
1329   VecImaginaryPart - Replaces a complex vector with its imginary part
1330 
1331   Collective
1332 
1333   Input Parameter:
1334 . v - the vector
1335 
1336   Level: beginner
1337 
1338 .seealso: `Vec`, `VecNorm()`, `VecRealPart()`
1339 @*/
1340 PetscErrorCode VecImaginaryPart(Vec v)
1341 {
1342   PetscFunctionBegin;
1343   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1344   PetscCall(VecApplyUnary_Private(v, NULL, ScalarImaginaryPart_Fn));
1345   PetscFunctionReturn(PETSC_SUCCESS);
1346 }
1347 
1348 static PetscScalar ScalarRealPart_Fn(PetscScalar x)
1349 {
1350   const PetscReal real = PetscRealPart(x);
1351 
1352 #if PetscDefined(USE_COMPLEX)
1353   return PetscCMPLX(real, 0.0);
1354 #else
1355   return real;
1356 #endif
1357 }
1358 
1359 /*@
1360   VecRealPart - Replaces a complex vector with its real part
1361 
1362   Collective
1363 
1364   Input Parameter:
1365 . v - the vector
1366 
1367   Level: beginner
1368 
1369 .seealso: `Vec`, `VecNorm()`, `VecImaginaryPart()`
1370 @*/
1371 PetscErrorCode VecRealPart(Vec v)
1372 {
1373   PetscFunctionBegin;
1374   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1375   PetscCall(VecApplyUnary_Private(v, NULL, ScalarRealPart_Fn));
1376   PetscFunctionReturn(PETSC_SUCCESS);
1377 }
1378 
1379 /*@
1380   VecDotNorm2 - computes the inner product of two vectors and the 2-norm squared of the second vector
1381 
1382   Collective
1383 
1384   Input Parameters:
1385 + s - first vector
1386 - t - second vector
1387 
1388   Output Parameters:
1389 + dp - s'conj(t)
1390 - nm - t'conj(t)
1391 
1392   Level: advanced
1393 
1394   Note:
1395   conj(x) is the complex conjugate of x when x is complex
1396 
1397 .seealso: `Vec`, `VecDot()`, `VecNorm()`, `VecDotBegin()`, `VecNormBegin()`, `VecDotEnd()`, `VecNormEnd()`
1398 
1399 @*/
1400 PetscErrorCode VecDotNorm2(Vec s, Vec t, PetscScalar *dp, PetscReal *nm)
1401 {
1402   PetscScalar work[] = {0.0, 0.0};
1403 
1404   PetscFunctionBegin;
1405   PetscValidHeaderSpecific(s, VEC_CLASSID, 1);
1406   PetscValidHeaderSpecific(t, VEC_CLASSID, 2);
1407   PetscAssertPointer(dp, 3);
1408   PetscAssertPointer(nm, 4);
1409   PetscValidType(s, 1);
1410   PetscValidType(t, 2);
1411   PetscCheckSameTypeAndComm(s, 1, t, 2);
1412   PetscCheck(s->map->N == t->map->N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector global lengths");
1413   PetscCheck(s->map->n == t->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector local lengths");
1414 
1415   PetscCall(PetscLogEventBegin(VEC_DotNorm2, s, t, 0, 0));
1416   if (s->ops->dotnorm2) {
1417     PetscUseTypeMethod(s, dotnorm2, t, work, work + 1);
1418   } else {
1419     const PetscScalar *sx, *tx;
1420     PetscInt           n;
1421 
1422     PetscCall(VecGetLocalSize(s, &n));
1423     PetscCall(VecGetArrayRead(s, &sx));
1424     PetscCall(VecGetArrayRead(t, &tx));
1425     for (PetscInt i = 0; i < n; ++i) {
1426       const PetscScalar txconj = PetscConj(tx[i]);
1427 
1428       work[0] += sx[i] * txconj;
1429       work[1] += tx[i] * txconj;
1430     }
1431     PetscCall(VecRestoreArrayRead(t, &tx));
1432     PetscCall(VecRestoreArrayRead(s, &sx));
1433     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, work, 2, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)s)));
1434     PetscCall(PetscLogFlops(4.0 * n));
1435   }
1436   PetscCall(PetscLogEventEnd(VEC_DotNorm2, s, t, 0, 0));
1437   *dp = work[0];
1438   *nm = PetscRealPart(work[1]);
1439   PetscFunctionReturn(PETSC_SUCCESS);
1440 }
1441 
1442 /*@
1443   VecSum - Computes the sum of all the components of a vector.
1444 
1445   Collective
1446 
1447   Input Parameter:
1448 . v - the vector
1449 
1450   Output Parameter:
1451 . sum - the result
1452 
1453   Level: beginner
1454 
1455 .seealso: `Vec`, `VecMean()`, `VecNorm()`
1456 @*/
1457 PetscErrorCode VecSum(Vec v, PetscScalar *sum)
1458 {
1459   PetscScalar tmp = 0.0;
1460 
1461   PetscFunctionBegin;
1462   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1463   PetscAssertPointer(sum, 2);
1464   if (v->ops->sum) {
1465     PetscUseTypeMethod(v, sum, &tmp);
1466   } else {
1467     const PetscScalar *x;
1468     PetscInt           n;
1469 
1470     PetscCall(VecGetLocalSize(v, &n));
1471     PetscCall(VecGetArrayRead(v, &x));
1472     for (PetscInt i = 0; i < n; ++i) tmp += x[i];
1473     PetscCall(VecRestoreArrayRead(v, &x));
1474   }
1475   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &tmp, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)v)));
1476   *sum = tmp;
1477   PetscFunctionReturn(PETSC_SUCCESS);
1478 }
1479 
1480 /*@
1481   VecMean - Computes the arithmetic mean of all the components of a vector.
1482 
1483   Collective
1484 
1485   Input Parameter:
1486 . v - the vector
1487 
1488   Output Parameter:
1489 . mean - the result
1490 
1491   Level: beginner
1492 
1493 .seealso: `Vec`, `VecSum()`, `VecNorm()`
1494 @*/
1495 PetscErrorCode VecMean(Vec v, PetscScalar *mean)
1496 {
1497   PetscInt n;
1498 
1499   PetscFunctionBegin;
1500   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1501   PetscAssertPointer(mean, 2);
1502   PetscCall(VecGetSize(v, &n));
1503   PetscCall(VecSum(v, mean));
1504   *mean /= n;
1505   PetscFunctionReturn(PETSC_SUCCESS);
1506 }
1507 
1508 /*@
1509   VecShift - Shifts all of the components of a vector by computing
1510   `x[i] = x[i] + shift`.
1511 
1512   Logically Collective
1513 
1514   Input Parameters:
1515 + v     - the vector
1516 - shift - the shift
1517 
1518   Level: intermediate
1519 
1520 .seealso: `Vec`
1521 @*/
1522 PetscErrorCode VecShift(Vec v, PetscScalar shift)
1523 {
1524   PetscFunctionBegin;
1525   PetscValidHeaderSpecific(v, VEC_CLASSID, 1);
1526   PetscValidLogicalCollectiveScalar(v, shift, 2);
1527   PetscCall(VecSetErrorIfLocked(v, 1));
1528   if (shift == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
1529 
1530   if (v->ops->shift) {
1531     PetscUseTypeMethod(v, shift, shift);
1532   } else {
1533     PetscInt     n;
1534     PetscScalar *x;
1535 
1536     PetscCall(VecGetLocalSize(v, &n));
1537     PetscCall(VecGetArray(v, &x));
1538     for (PetscInt i = 0; i < n; ++i) x[i] += shift;
1539     PetscCall(VecRestoreArray(v, &x));
1540   }
1541   PetscFunctionReturn(PETSC_SUCCESS);
1542 }
1543 
1544 /*@
1545   VecPermute - Permutes a vector in place using the given ordering.
1546 
1547   Input Parameters:
1548 + x   - The vector
1549 . row - The ordering
1550 - inv - The flag for inverting the permutation
1551 
1552   Level: beginner
1553 
1554   Note:
1555   This function does not yet support parallel Index Sets with non-local permutations
1556 
1557 .seealso: `Vec`, `MatPermute()`
1558 @*/
1559 PetscErrorCode VecPermute(Vec x, IS row, PetscBool inv)
1560 {
1561   const PetscScalar *array;
1562   PetscScalar       *newArray;
1563   const PetscInt    *idx;
1564   PetscInt           i, rstart, rend;
1565 
1566   PetscFunctionBegin;
1567   PetscValidHeaderSpecific(x, VEC_CLASSID, 1);
1568   PetscValidHeaderSpecific(row, IS_CLASSID, 2);
1569   PetscCall(VecSetErrorIfLocked(x, 1));
1570   PetscCall(VecGetOwnershipRange(x, &rstart, &rend));
1571   PetscCall(ISGetIndices(row, &idx));
1572   PetscCall(VecGetArrayRead(x, &array));
1573   PetscCall(PetscMalloc1(x->map->n, &newArray));
1574   if (PetscDefined(USE_DEBUG)) {
1575     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]);
1576   }
1577   if (!inv) {
1578     for (i = 0; i < x->map->n; i++) newArray[i] = array[idx[i] - rstart];
1579   } else {
1580     for (i = 0; i < x->map->n; i++) newArray[idx[i] - rstart] = array[i];
1581   }
1582   PetscCall(VecRestoreArrayRead(x, &array));
1583   PetscCall(ISRestoreIndices(row, &idx));
1584   PetscCall(VecReplaceArray(x, newArray));
1585   PetscFunctionReturn(PETSC_SUCCESS);
1586 }
1587 
1588 /*@
1589   VecEqual - Compares two vectors. Returns true if the two vectors are either pointing to the same memory buffer,
1590   or if the two vectors have the same local and global layout as well as bitwise equality of all entries.
1591   Does NOT take round-off errors into account.
1592 
1593   Collective
1594 
1595   Input Parameters:
1596 + vec1 - the first vector
1597 - vec2 - the second vector
1598 
1599   Output Parameter:
1600 . flg - `PETSC_TRUE` if the vectors are equal; `PETSC_FALSE` otherwise.
1601 
1602   Level: intermediate
1603 
1604 .seealso: `Vec`
1605 @*/
1606 PetscErrorCode VecEqual(Vec vec1, Vec vec2, PetscBool *flg)
1607 {
1608   const PetscScalar *v1, *v2;
1609   PetscInt           n1, n2, N1, N2;
1610   PetscBool          flg1;
1611 
1612   PetscFunctionBegin;
1613   PetscValidHeaderSpecific(vec1, VEC_CLASSID, 1);
1614   PetscValidHeaderSpecific(vec2, VEC_CLASSID, 2);
1615   PetscAssertPointer(flg, 3);
1616   if (vec1 == vec2) *flg = PETSC_TRUE;
1617   else {
1618     PetscCall(VecGetSize(vec1, &N1));
1619     PetscCall(VecGetSize(vec2, &N2));
1620     if (N1 != N2) flg1 = PETSC_FALSE;
1621     else {
1622       PetscCall(VecGetLocalSize(vec1, &n1));
1623       PetscCall(VecGetLocalSize(vec2, &n2));
1624       if (n1 != n2) flg1 = PETSC_FALSE;
1625       else {
1626         PetscCall(VecGetArrayRead(vec1, &v1));
1627         PetscCall(VecGetArrayRead(vec2, &v2));
1628         PetscCall(PetscArraycmp(v1, v2, n1, &flg1));
1629         PetscCall(VecRestoreArrayRead(vec1, &v1));
1630         PetscCall(VecRestoreArrayRead(vec2, &v2));
1631       }
1632     }
1633     /* combine results from all processors */
1634     PetscCall(MPIU_Allreduce(&flg1, flg, 1, MPIU_BOOL, MPI_MIN, PetscObjectComm((PetscObject)vec1)));
1635   }
1636   PetscFunctionReturn(PETSC_SUCCESS);
1637 }
1638 
1639 /*@
1640   VecUniqueEntries - Compute the number of unique entries, and those entries
1641 
1642   Collective
1643 
1644   Input Parameter:
1645 . vec - the vector
1646 
1647   Output Parameters:
1648 + n - The number of unique entries
1649 - e - The entries
1650 
1651   Level: intermediate
1652 
1653 .seealso: `Vec`
1654 @*/
1655 PetscErrorCode VecUniqueEntries(Vec vec, PetscInt *n, PetscScalar **e)
1656 {
1657   const PetscScalar *v;
1658   PetscScalar       *tmp, *vals;
1659   PetscMPIInt       *N, *displs, l;
1660   PetscInt           ng, m, i, j, p;
1661   PetscMPIInt        size;
1662 
1663   PetscFunctionBegin;
1664   PetscValidHeaderSpecific(vec, VEC_CLASSID, 1);
1665   PetscAssertPointer(n, 2);
1666   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)vec), &size));
1667   PetscCall(VecGetLocalSize(vec, &m));
1668   PetscCall(VecGetArrayRead(vec, &v));
1669   PetscCall(PetscMalloc2(m, &tmp, size, &N));
1670   for (i = 0, j = 0, l = 0; i < m; ++i) {
1671     /* Can speed this up with sorting */
1672     for (j = 0; j < l; ++j) {
1673       if (v[i] == tmp[j]) break;
1674     }
1675     if (j == l) {
1676       tmp[j] = v[i];
1677       ++l;
1678     }
1679   }
1680   PetscCall(VecRestoreArrayRead(vec, &v));
1681   /* Gather serial results */
1682   PetscCallMPI(MPI_Allgather(&l, 1, MPI_INT, N, 1, MPI_INT, PetscObjectComm((PetscObject)vec)));
1683   for (p = 0, ng = 0; p < size; ++p) ng += N[p];
1684   PetscCall(PetscMalloc2(ng, &vals, size + 1, &displs));
1685   for (p = 1, displs[0] = 0; p <= size; ++p) displs[p] = displs[p - 1] + N[p - 1];
1686   PetscCallMPI(MPI_Allgatherv(tmp, l, MPIU_SCALAR, vals, N, displs, MPIU_SCALAR, PetscObjectComm((PetscObject)vec)));
1687   /* Find unique entries */
1688 #ifdef PETSC_USE_COMPLEX
1689   SETERRQ(PetscObjectComm((PetscObject)vec), PETSC_ERR_SUP, "Does not work with complex numbers");
1690 #else
1691   *n = displs[size];
1692   PetscCall(PetscSortRemoveDupsReal(n, (PetscReal *)vals));
1693   if (e) {
1694     PetscAssertPointer(e, 3);
1695     PetscCall(PetscMalloc1(*n, e));
1696     for (i = 0; i < *n; ++i) (*e)[i] = vals[i];
1697   }
1698   PetscCall(PetscFree2(vals, displs));
1699   PetscCall(PetscFree2(tmp, N));
1700   PetscFunctionReturn(PETSC_SUCCESS);
1701 #endif
1702 }
1703