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