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