xref: /petsc/src/vec/vec/utils/vinv.c (revision ca9f406c45a43819c7a962370166b42afc3e7ff4)
1 #define PETSCVEC_DLL
2 
3 /*
4      Some useful vector utility functions.
5 */
6 #include "private/vecimpl.h"          /*I "petscvec.h" I*/
7 
8 extern MPI_Op VecMax_Local_Op;
9 extern MPI_Op VecMin_Local_Op;
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "VecStrideScale"
13 /*@
14    VecStrideScale - Scales a subvector of a vector defined
15    by a starting point and a stride.
16 
17    Collective on Vec
18 
19    Input Parameter:
20 +  v - the vector
21 .  start - starting point of the subvector (defined by a stride)
22 -  scale - value to multiply each subvector entry by
23 
24    Notes:
25    One must call VecSetBlockSize() before this routine to set the stride
26    information, or use a vector created from a multicomponent DA.
27 
28    This will only work if the desire subvector is a stride subvector
29 
30    Level: advanced
31 
32    Concepts: scale^on stride of vector
33    Concepts: stride^scale
34 
35 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideScale()
36 @*/
37 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScale(Vec v,PetscInt start,PetscScalar scale)
38 {
39   PetscErrorCode ierr;
40   PetscInt       i,n,bs;
41   PetscScalar    *x;
42 
43   PetscFunctionBegin;
44   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
45   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
46   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
47 
48   bs   = v->map.bs;
49   if (start < 0) {
50     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
51   } else if (start >= bs) {
52     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
53             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
54   }
55   x += start;
56 
57   for (i=0; i<n; i+=bs) {
58     x[i] *= scale;
59   }
60   x -= start;
61 
62   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "VecStrideNorm"
68 /*@
69    VecStrideNorm - Computes the norm of subvector of a vector defined
70    by a starting point and a stride.
71 
72    Collective on Vec
73 
74    Input Parameter:
75 +  v - the vector
76 .  start - starting point of the subvector (defined by a stride)
77 -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
78 
79    Output Parameter:
80 .  norm - the norm
81 
82    Notes:
83    One must call VecSetBlockSize() before this routine to set the stride
84    information, or use a vector created from a multicomponent DA.
85 
86    If x is the array representing the vector x then this computes the norm
87    of the array (x[start],x[start+stride],x[start+2*stride], ....)
88 
89    This is useful for computing, say the norm of the pressure variable when
90    the pressure is stored (interlaced) with other variables, say density etc.
91 
92    This will only work if the desire subvector is a stride subvector
93 
94    Level: advanced
95 
96    Concepts: norm^on stride of vector
97    Concepts: stride^norm
98 
99 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
100 @*/
101 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideNorm(Vec v,PetscInt start,NormType ntype,PetscReal *nrm)
102 {
103   PetscErrorCode ierr;
104   PetscInt       i,n,bs;
105   PetscScalar    *x;
106   PetscReal      tnorm;
107   MPI_Comm       comm;
108 
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
111   PetscValidDoublePointer(nrm,3);
112   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
113   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
114   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
115 
116   bs   = v->map.bs;
117   if (start < 0) {
118     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
119   } else if (start >= bs) {
120     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
121             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
122   }
123   x += start;
124 
125   if (ntype == NORM_2) {
126     PetscScalar sum = 0.0;
127     for (i=0; i<n; i+=bs) {
128       sum += x[i]*(PetscConj(x[i]));
129     }
130     tnorm  = PetscRealPart(sum);
131     ierr   = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);CHKERRQ(ierr);
132     *nrm = sqrt(*nrm);
133   } else if (ntype == NORM_1) {
134     tnorm = 0.0;
135     for (i=0; i<n; i+=bs) {
136       tnorm += PetscAbsScalar(x[i]);
137     }
138     ierr   = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);CHKERRQ(ierr);
139   } else if (ntype == NORM_INFINITY) {
140     PetscReal tmp;
141     tnorm = 0.0;
142 
143     for (i=0; i<n; i+=bs) {
144       if ((tmp = PetscAbsScalar(x[i])) > tnorm) tnorm = tmp;
145       /* check special case of tmp == NaN */
146       if (tmp != tmp) {tnorm = tmp; break;}
147     }
148     ierr   = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_MAX,comm);CHKERRQ(ierr);
149   } else {
150     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown norm type");
151   }
152 
153   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
154   PetscFunctionReturn(0);
155 }
156 
157 #undef __FUNCT__
158 #define __FUNCT__ "VecStrideMax"
159 /*@
160    VecStrideMax - Computes the maximum of subvector of a vector defined
161    by a starting point and a stride and optionally its location.
162 
163    Collective on Vec
164 
165    Input Parameter:
166 +  v - the vector
167 -  start - starting point of the subvector (defined by a stride)
168 
169    Output Parameter:
170 +  index - the location where the maximum occurred  (pass PETSC_NULL if not required)
171 -  nrm - the max
172 
173    Notes:
174    One must call VecSetBlockSize() before this routine to set the stride
175    information, or use a vector created from a multicomponent DA.
176 
177    If xa is the array representing the vector x, then this computes the max
178    of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
179 
180    This is useful for computing, say the maximum of the pressure variable when
181    the pressure is stored (interlaced) with other variables, e.g., density, etc.
182    This will only work if the desire subvector is a stride subvector.
183 
184    Level: advanced
185 
186    Concepts: maximum^on stride of vector
187    Concepts: stride^maximum
188 
189 .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
190 @*/
191 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideMax(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
192 {
193   PetscErrorCode ierr;
194   PetscInt       i,n,bs,id;
195   PetscScalar    *x;
196   PetscReal      max,tmp;
197   MPI_Comm       comm;
198 
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
201   PetscValidDoublePointer(nrm,3);
202 
203   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
204   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
205   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
206 
207   bs   = v->map.bs;
208   if (start < 0) {
209     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
210   } else if (start >= bs) {
211     SETERRQ2(PETSC_ERR_ARG_WRONG,"Start of stride subvector (%D) is too large for stride\n\
212             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
213   }
214   x += start;
215 
216   id = -1;
217   if (!n) {
218     max = PETSC_MIN;
219   } else {
220     id  = 0;
221 #if defined(PETSC_USE_COMPLEX)
222     max = PetscRealPart(x[0]);
223 #else
224     max = x[0];
225 #endif
226     for (i=bs; i<n; i+=bs) {
227 #if defined(PETSC_USE_COMPLEX)
228       if ((tmp = PetscRealPart(x[i])) > max) { max = tmp; id = i;}
229 #else
230       if ((tmp = x[i]) > max) { max = tmp; id = i;}
231 #endif
232     }
233   }
234   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
235 
236   if (!idex) {
237     ierr   = MPI_Allreduce(&max,nrm,1,MPIU_REAL,MPI_MAX,comm);CHKERRQ(ierr);
238   } else {
239     PetscReal in[2],out[2];
240     PetscInt  rstart;
241 
242     ierr  = VecGetOwnershipRange(v,&rstart,PETSC_NULL);CHKERRQ(ierr);
243     in[0] = max;
244     in[1] = rstart+id;
245     ierr  = MPI_Allreduce(in,out,2,MPIU_REAL,VecMax_Local_Op,v->comm);CHKERRQ(ierr);
246     *nrm  = out[0];
247     *idex = (PetscInt)out[1];
248   }
249 
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "VecStrideMin"
255 /*@
256    VecStrideMin - Computes the minimum of subvector of a vector defined
257    by a starting point and a stride and optionally its location.
258 
259    Collective on Vec
260 
261    Input Parameter:
262 +  v - the vector
263 -  start - starting point of the subvector (defined by a stride)
264 
265    Output Parameter:
266 +  idex - the location where the minimum occurred. (pass PETSC_NULL if not required)
267 -  nrm - the min
268 
269    Level: advanced
270 
271    Notes:
272    One must call VecSetBlockSize() before this routine to set the stride
273    information, or use a vector created from a multicomponent DA.
274 
275    If xa is the array representing the vector x, then this computes the min
276    of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)
277 
278    This is useful for computing, say the minimum of the pressure variable when
279    the pressure is stored (interlaced) with other variables, e.g., density, etc.
280    This will only work if the desire subvector is a stride subvector.
281 
282    Concepts: minimum^on stride of vector
283    Concepts: stride^minimum
284 
285 .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
286 @*/
287 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideMin(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
288 {
289   PetscErrorCode ierr;
290   PetscInt       i,n,bs,id;
291   PetscScalar    *x;
292   PetscReal      min,tmp;
293   MPI_Comm       comm;
294 
295   PetscFunctionBegin;
296   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
297   PetscValidDoublePointer(nrm,4);
298 
299   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
300   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
301   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
302 
303   bs   = v->map.bs;
304   if (start < 0) {
305     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
306   } else if (start >= bs) {
307     SETERRQ2(PETSC_ERR_ARG_WRONG,"Start of stride subvector (%D) is too large for stride\n\
308             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
309   }
310   x += start;
311 
312   id = -1;
313   if (!n) {
314     min = PETSC_MAX;
315   } else {
316     id = 0;
317 #if defined(PETSC_USE_COMPLEX)
318     min = PetscRealPart(x[0]);
319 #else
320     min = x[0];
321 #endif
322     for (i=bs; i<n; i+=bs) {
323 #if defined(PETSC_USE_COMPLEX)
324       if ((tmp = PetscRealPart(x[i])) < min) { min = tmp; id = i;}
325 #else
326       if ((tmp = x[i]) < min) { min = tmp; id = i;}
327 #endif
328     }
329   }
330   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
331 
332   if (!idex) {
333     ierr   = MPI_Allreduce(&min,nrm,1,MPIU_REAL,MPI_MIN,comm);CHKERRQ(ierr);
334   } else {
335     PetscReal in[2],out[2];
336     PetscInt  rstart;
337 
338     ierr  = VecGetOwnershipRange(v,&rstart,PETSC_NULL);CHKERRQ(ierr);
339     in[0] = min;
340     in[1] = rstart+id;
341     ierr  = MPI_Allreduce(in,out,2,MPIU_REAL,VecMin_Local_Op,v->comm);CHKERRQ(ierr);
342     *nrm  = out[0];
343     *idex = (PetscInt)out[1];
344   }
345 
346   PetscFunctionReturn(0);
347 }
348 
349 #undef __FUNCT__
350 #define __FUNCT__ "VecStrideScaleAll"
351 /*@
352    VecStrideScaleAll - Scales the subvectors of a vector defined
353    by a starting point and a stride.
354 
355    Collective on Vec
356 
357    Input Parameter:
358 +  v - the vector
359 -  scales - values to multiply each subvector entry by
360 
361    Notes:
362    One must call VecSetBlockSize() before this routine to set the stride
363    information, or use a vector created from a multicomponent DA.
364 
365 
366    Level: advanced
367 
368    Concepts: scale^on stride of vector
369    Concepts: stride^scale
370 
371 .seealso: VecNorm(), VecStrideScale(), VecScale(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
372 @*/
373 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScaleAll(Vec v,PetscScalar *scales)
374 {
375   PetscErrorCode ierr;
376   PetscInt       i,j,n,bs;
377   PetscScalar    *x;
378 
379   PetscFunctionBegin;
380   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
381   PetscValidScalarPointer(scales,2);
382   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
383   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
384 
385   bs   = v->map.bs;
386 
387   /* need to provide optimized code for each bs */
388   for (i=0; i<n; i+=bs) {
389     for (j=0; j<bs; j++) {
390       x[i+j] *= scales[j];
391     }
392   }
393   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
394   PetscFunctionReturn(0);
395 }
396 
397 #undef __FUNCT__
398 #define __FUNCT__ "VecStrideNormAll"
399 /*@
400    VecStrideNormAll - Computes the norms  subvectors of a vector defined
401    by a starting point and a stride.
402 
403    Collective on Vec
404 
405    Input Parameter:
406 +  v - the vector
407 -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY
408 
409    Output Parameter:
410 .  nrm - the norms
411 
412    Notes:
413    One must call VecSetBlockSize() before this routine to set the stride
414    information, or use a vector created from a multicomponent DA.
415 
416    If x is the array representing the vector x then this computes the norm
417    of the array (x[start],x[start+stride],x[start+2*stride], ....)
418 
419    This is useful for computing, say the norm of the pressure variable when
420    the pressure is stored (interlaced) with other variables, say density etc.
421 
422    This will only work if the desire subvector is a stride subvector
423 
424    Level: advanced
425 
426    Concepts: norm^on stride of vector
427    Concepts: stride^norm
428 
429 .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
430 @*/
431 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideNormAll(Vec v,NormType ntype,PetscReal *nrm)
432 {
433   PetscErrorCode ierr;
434   PetscInt       i,j,n,bs;
435   PetscScalar    *x;
436   PetscReal      tnorm[128];
437   MPI_Comm       comm;
438 
439   PetscFunctionBegin;
440   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
441   PetscValidDoublePointer(nrm,3);
442   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
443   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
444   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
445 
446   bs   = v->map.bs;
447   if (bs > 128) SETERRQ(PETSC_ERR_SUP,"Currently supports only blocksize up to 128");
448 
449   if (ntype == NORM_2) {
450     PetscScalar sum[128];
451     for (j=0; j<bs; j++) sum[j] = 0.0;
452     for (i=0; i<n; i+=bs) {
453       for (j=0; j<bs; j++) {
454 	sum[j] += x[i+j]*(PetscConj(x[i+j]));
455       }
456     }
457     for (j=0; j<bs; j++) {
458       tnorm[j]  = PetscRealPart(sum[j]);
459     }
460     ierr   = MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_SUM,comm);CHKERRQ(ierr);
461     for (j=0; j<bs; j++) {
462       nrm[j] = sqrt(nrm[j]);
463     }
464   } else if (ntype == NORM_1) {
465     for (j=0; j<bs; j++) {
466       tnorm[j] = 0.0;
467     }
468     for (i=0; i<n; i+=bs) {
469       for (j=0; j<bs; j++) {
470 	tnorm[j] += PetscAbsScalar(x[i+j]);
471       }
472     }
473     ierr   = MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_SUM,comm);CHKERRQ(ierr);
474   } else if (ntype == NORM_INFINITY) {
475     PetscReal tmp;
476     for (j=0; j<bs; j++) {
477       tnorm[j] = 0.0;
478     }
479 
480     for (i=0; i<n; i+=bs) {
481       for (j=0; j<bs; j++) {
482 	if ((tmp = PetscAbsScalar(x[i+j])) > tnorm[j]) tnorm[j] = tmp;
483 	/* check special case of tmp == NaN */
484 	if (tmp != tmp) {tnorm[j] = tmp; break;}
485       }
486     }
487     ierr   = MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_MAX,comm);CHKERRQ(ierr);
488   } else {
489     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown norm type");
490   }
491 
492   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
493   PetscFunctionReturn(0);
494 }
495 
496 #undef __FUNCT__
497 #define __FUNCT__ "VecStrideMaxAll"
498 /*@
499    VecStrideMaxAll - Computes the maximums of subvectors of a vector defined
500    by a starting point and a stride and optionally its location.
501 
502    Collective on Vec
503 
504    Input Parameter:
505 .  v - the vector
506 
507    Output Parameter:
508 +  index - the location where the maximum occurred (not supported, pass PETSC_NULL,
509            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
510 -  nrm - the maximums
511 
512    Notes:
513    One must call VecSetBlockSize() before this routine to set the stride
514    information, or use a vector created from a multicomponent DA.
515 
516    This is useful for computing, say the maximum of the pressure variable when
517    the pressure is stored (interlaced) with other variables, e.g., density, etc.
518    This will only work if the desire subvector is a stride subvector.
519 
520    Level: advanced
521 
522    Concepts: maximum^on stride of vector
523    Concepts: stride^maximum
524 
525 .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
526 @*/
527 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideMaxAll(Vec v,PetscInt *idex,PetscReal *nrm)
528 {
529   PetscErrorCode ierr;
530   PetscInt       i,j,n,bs;
531   PetscScalar    *x;
532   PetscReal      max[128],tmp;
533   MPI_Comm       comm;
534 
535   PetscFunctionBegin;
536   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
537   PetscValidDoublePointer(nrm,3);
538   if (idex) {
539     SETERRQ(PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
540   }
541   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
542   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
543   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
544 
545   bs   = v->map.bs;
546   if (bs > 128) SETERRQ(PETSC_ERR_SUP,"Currently supports only blocksize up to 128");
547 
548   if (!n) {
549     for (j=0; j<bs; j++) {
550       max[j] = PETSC_MIN;
551     }
552   } else {
553     for (j=0; j<bs; j++) {
554 #if defined(PETSC_USE_COMPLEX)
555       max[j] = PetscRealPart(x[j]);
556 #else
557       max[j] = x[j];
558 #endif
559     }
560     for (i=bs; i<n; i+=bs) {
561       for (j=0; j<bs; j++) {
562 #if defined(PETSC_USE_COMPLEX)
563 	if ((tmp = PetscRealPart(x[i+j])) > max[j]) { max[j] = tmp;}
564 #else
565 	if ((tmp = x[i+j]) > max[j]) { max[j] = tmp; }
566 #endif
567       }
568     }
569   }
570   ierr   = MPI_Allreduce(max,nrm,bs,MPIU_REAL,MPI_MAX,comm);CHKERRQ(ierr);
571 
572   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
573   PetscFunctionReturn(0);
574 }
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "VecStrideMinAll"
578 /*@
579    VecStrideMinAll - Computes the minimum of subvector of a vector defined
580    by a starting point and a stride and optionally its location.
581 
582    Collective on Vec
583 
584    Input Parameter:
585 .  v - the vector
586 
587    Output Parameter:
588 +  idex - the location where the minimum occurred (not supported, pass PETSC_NULL,
589            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
590 -  nrm - the minimums
591 
592    Level: advanced
593 
594    Notes:
595    One must call VecSetBlockSize() before this routine to set the stride
596    information, or use a vector created from a multicomponent DA.
597 
598    This is useful for computing, say the minimum of the pressure variable when
599    the pressure is stored (interlaced) with other variables, e.g., density, etc.
600    This will only work if the desire subvector is a stride subvector.
601 
602    Concepts: minimum^on stride of vector
603    Concepts: stride^minimum
604 
605 .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
606 @*/
607 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideMinAll(Vec v,PetscInt *idex,PetscReal *nrm)
608 {
609   PetscErrorCode ierr;
610   PetscInt       i,n,bs,j;
611   PetscScalar    *x;
612   PetscReal      min[128],tmp;
613   MPI_Comm       comm;
614 
615   PetscFunctionBegin;
616   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
617   PetscValidDoublePointer(nrm,3);
618   if (idex) {
619     SETERRQ(PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
620   }
621   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
622   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
623   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
624 
625   bs   = v->map.bs;
626   if (bs > 128) SETERRQ(PETSC_ERR_SUP,"Currently supports only blocksize up to 128");
627 
628   if (!n) {
629     for (j=0; j<bs; j++) {
630       min[j] = PETSC_MAX;
631     }
632   } else {
633     for (j=0; j<bs; j++) {
634 #if defined(PETSC_USE_COMPLEX)
635       min[j] = PetscRealPart(x[j]);
636 #else
637       min[j] = x[j];
638 #endif
639     }
640     for (i=bs; i<n; i+=bs) {
641       for (j=0; j<bs; j++) {
642 #if defined(PETSC_USE_COMPLEX)
643 	if ((tmp = PetscRealPart(x[i+j])) < min[j]) { min[j] = tmp;}
644 #else
645 	if ((tmp = x[i+j]) < min[j]) { min[j] = tmp; }
646 #endif
647       }
648     }
649   }
650   ierr   = MPI_Allreduce(min,nrm,bs,MPIU_REAL,MPI_MIN,comm);CHKERRQ(ierr);
651 
652   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
653   PetscFunctionReturn(0);
654 }
655 
656 /*----------------------------------------------------------------------------------------------*/
657 #undef __FUNCT__
658 #define __FUNCT__ "VecStrideGatherAll"
659 /*@
660    VecStrideGatherAll - Gathers all the single components from a multi-component vector into
661    separate vectors.
662 
663    Collective on Vec
664 
665    Input Parameter:
666 +  v - the vector
667 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
668 
669    Output Parameter:
670 .  s - the location where the subvectors are stored
671 
672    Notes:
673    One must call VecSetBlockSize() before this routine to set the stride
674    information, or use a vector created from a multicomponent DA.
675 
676    If x is the array representing the vector x then this gathers
677    the arrays (x[start],x[start+stride],x[start+2*stride], ....)
678    for start=0,1,2,...bs-1
679 
680    The parallel layout of the vector and the subvector must be the same;
681    i.e., nlocal of v = stride*(nlocal of s)
682 
683    Not optimized; could be easily
684 
685    Level: advanced
686 
687    Concepts: gather^into strided vector
688 
689 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
690           VecStrideScatterAll()
691 @*/
692 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideGatherAll(Vec v,Vec *s,InsertMode addv)
693 {
694   PetscErrorCode ierr;
695   PetscInt       i,n,bs,j,k,*bss = PETSC_NULL,nv,jj,nvc;
696   PetscScalar    *x,**y;
697 
698   PetscFunctionBegin;
699   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
700   PetscValidPointer(s,2);
701   PetscValidHeaderSpecific(*s,VEC_COOKIE,2);
702   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
703   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
704   bs   = v->map.bs;
705   if (bs < 0) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set");
706   ierr = PetscMalloc2(bs,PetscReal*,&y,bs,PetscInt,&bss);CHKERRQ(ierr);
707   nv   = 0;
708   nvc  = 0;
709   for (i=0; i<bs; i++) {
710     ierr = VecGetBlockSize(s[i],&bss[i]);CHKERRQ(ierr);
711     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
712     ierr = VecGetArray(s[i],&y[i]);CHKERRQ(ierr);
713     nvc  += bss[i];
714     nv++;
715     if (nvc > bs)  SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector");
716     if (nvc == bs) break;
717   }
718 
719   n =  n/bs;
720 
721   jj = 0;
722   if (addv == INSERT_VALUES) {
723     for (j=0; j<nv; j++) {
724       for (k=0; k<bss[j]; k++) {
725 	for (i=0; i<n; i++) {
726 	  y[j][i*bss[j] + k] = x[bs*i+jj+k];
727         }
728       }
729       jj += bss[j];
730     }
731   } else if (addv == ADD_VALUES) {
732     for (j=0; j<nv; j++) {
733       for (k=0; k<bss[j]; k++) {
734 	for (i=0; i<n; i++) {
735 	  y[j][i*bss[j] + k] += x[bs*i+jj+k];
736         }
737       }
738       jj += bss[j];
739     }
740 #if !defined(PETSC_USE_COMPLEX)
741   } else if (addv == MAX_VALUES) {
742     for (j=0; j<nv; j++) {
743       for (k=0; k<bss[j]; k++) {
744 	for (i=0; i<n; i++) {
745 	  y[j][i*bss[j] + k] = PetscMax(y[j][i*bss[j] + k],x[bs*i+jj+k]);
746         }
747       }
748       jj += bss[j];
749     }
750 #endif
751   } else {
752     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
753   }
754 
755   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
756   for (i=0; i<nv; i++) {
757     ierr = VecRestoreArray(s[i],&y[i]);CHKERRQ(ierr);
758   }
759 
760   ierr = PetscFree2(y,bss);CHKERRQ(ierr);
761   PetscFunctionReturn(0);
762 }
763 
764 #undef __FUNCT__
765 #define __FUNCT__ "VecStrideScatterAll"
766 /*@
767    VecStrideScatterAll - Scatters all the single components from separate vectors into
768      a multi-component vector.
769 
770    Collective on Vec
771 
772    Input Parameter:
773 +  s - the location where the subvectors are stored
774 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
775 
776    Output Parameter:
777 .  v - the multicomponent vector
778 
779    Notes:
780    One must call VecSetBlockSize() before this routine to set the stride
781    information, or use a vector created from a multicomponent DA.
782 
783    The parallel layout of the vector and the subvector must be the same;
784    i.e., nlocal of v = stride*(nlocal of s)
785 
786    Not optimized; could be easily
787 
788    Level: advanced
789 
790    Concepts:  scatter^into strided vector
791 
792 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
793           VecStrideScatterAll()
794 @*/
795 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScatterAll(Vec *s,Vec v,InsertMode addv)
796 {
797   PetscErrorCode ierr;
798   PetscInt        i,n,bs,j,jj,k,*bss = PETSC_NULL,nv,nvc;
799   PetscScalar     *x,**y;
800 
801   PetscFunctionBegin;
802   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
803   PetscValidPointer(s,2);
804   PetscValidHeaderSpecific(*s,VEC_COOKIE,2);
805   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
806   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
807   bs   = v->map.bs;
808   if (bs < 0) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set");
809 
810   ierr = PetscMalloc2(bs,PetscScalar**,&y,bs,PetscInt,&bss);CHKERRQ(ierr);
811   nv  = 0;
812   nvc = 0;
813   for (i=0; i<bs; i++) {
814     ierr = VecGetBlockSize(s[i],&bss[i]);CHKERRQ(ierr);
815     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
816     ierr = VecGetArray(s[i],&y[i]);CHKERRQ(ierr);
817     nvc  += bss[i];
818     nv++;
819     if (nvc > bs)  SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector");
820     if (nvc == bs) break;
821   }
822 
823   n =  n/bs;
824 
825   jj = 0;
826   if (addv == INSERT_VALUES) {
827     for (j=0; j<nv; j++) {
828       for (k=0; k<bss[j]; k++) {
829 	for (i=0; i<n; i++) {
830 	  x[bs*i+jj+k] = y[j][i*bss[j] + k];
831         }
832       }
833       jj += bss[j];
834     }
835   } else if (addv == ADD_VALUES) {
836     for (j=0; j<nv; j++) {
837       for (k=0; k<bss[j]; k++) {
838 	for (i=0; i<n; i++) {
839 	  x[bs*i+jj+k] += y[j][i*bss[j] + k];
840         }
841       }
842       jj += bss[j];
843     }
844 #if !defined(PETSC_USE_COMPLEX)
845   } else if (addv == MAX_VALUES) {
846     for (j=0; j<nv; j++) {
847       for (k=0; k<bss[j]; k++) {
848 	for (i=0; i<n; i++) {
849 	  x[bs*i+jj+k] = PetscMax(x[bs*i+jj+k],y[j][i*bss[j] + k]);
850         }
851       }
852       jj += bss[j];
853     }
854 #endif
855   } else {
856     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
857   }
858 
859   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
860   for (i=0; i<nv; i++) {
861     ierr = VecRestoreArray(s[i],&y[i]);CHKERRQ(ierr);
862   }
863   ierr = PetscFree2(y,bss);CHKERRQ(ierr);
864   PetscFunctionReturn(0);
865 }
866 
867 #undef __FUNCT__
868 #define __FUNCT__ "VecStrideGather"
869 /*@
870    VecStrideGather - Gathers a single component from a multi-component vector into
871    another vector.
872 
873    Collective on Vec
874 
875    Input Parameter:
876 +  v - the vector
877 .  start - starting point of the subvector (defined by a stride)
878 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
879 
880    Output Parameter:
881 .  s - the location where the subvector is stored
882 
883    Notes:
884    One must call VecSetBlockSize() before this routine to set the stride
885    information, or use a vector created from a multicomponent DA.
886 
887    If x is the array representing the vector x then this gathers
888    the array (x[start],x[start+stride],x[start+2*stride], ....)
889 
890    The parallel layout of the vector and the subvector must be the same;
891    i.e., nlocal of v = stride*(nlocal of s)
892 
893    Not optimized; could be easily
894 
895    Level: advanced
896 
897    Concepts: gather^into strided vector
898 
899 .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
900           VecStrideScatterAll()
901 @*/
902 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideGather(Vec v,PetscInt start,Vec s,InsertMode addv)
903 {
904   PetscErrorCode ierr;
905   PetscInt       i,n,bs,ns;
906   PetscScalar    *x,*y;
907 
908   PetscFunctionBegin;
909   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
910   PetscValidHeaderSpecific(s,VEC_COOKIE,3);
911   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
912   ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr);
913   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
914   ierr = VecGetArray(s,&y);CHKERRQ(ierr);
915 
916   bs   = v->map.bs;
917   if (start < 0) {
918     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
919   } else if (start >= bs) {
920     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
921             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
922   }
923   if (n != ns*bs) {
924     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Subvector length * blocksize %D not correct for gather from original vector %D",ns*bs,n);
925   }
926   x += start;
927   n =  n/bs;
928 
929   if (addv == INSERT_VALUES) {
930     for (i=0; i<n; i++) {
931       y[i] = x[bs*i];
932     }
933   } else if (addv == ADD_VALUES) {
934     for (i=0; i<n; i++) {
935       y[i] += x[bs*i];
936     }
937 #if !defined(PETSC_USE_COMPLEX)
938   } else if (addv == MAX_VALUES) {
939     for (i=0; i<n; i++) {
940       y[i] = PetscMax(y[i],x[bs*i]);
941     }
942 #endif
943   } else {
944     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
945   }
946 
947   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
948   ierr = VecRestoreArray(s,&y);CHKERRQ(ierr);
949   PetscFunctionReturn(0);
950 }
951 
952 #undef __FUNCT__
953 #define __FUNCT__ "VecStrideScatter"
954 /*@
955    VecStrideScatter - Scatters a single component from a vector into a multi-component vector.
956 
957    Collective on Vec
958 
959    Input Parameter:
960 +  s - the single-component vector
961 .  start - starting point of the subvector (defined by a stride)
962 -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES
963 
964    Output Parameter:
965 .  v - the location where the subvector is scattered (the multi-component vector)
966 
967    Notes:
968    One must call VecSetBlockSize() on the multi-component vector before this
969    routine to set the stride  information, or use a vector created from a multicomponent DA.
970 
971    The parallel layout of the vector and the subvector must be the same;
972    i.e., nlocal of v = stride*(nlocal of s)
973 
974    Not optimized; could be easily
975 
976    Level: advanced
977 
978    Concepts: scatter^into strided vector
979 
980 .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
981           VecStrideScatterAll()
982 @*/
983 PetscErrorCode PETSCVEC_DLLEXPORT VecStrideScatter(Vec s,PetscInt start,Vec v,InsertMode addv)
984 {
985   PetscErrorCode ierr;
986   PetscInt       i,n,bs,ns;
987   PetscScalar    *x,*y;
988 
989   PetscFunctionBegin;
990   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
991   PetscValidHeaderSpecific(s,VEC_COOKIE,3);
992   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
993   ierr = VecGetLocalSize(s,&ns);CHKERRQ(ierr);
994   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
995   ierr = VecGetArray(s,&y);CHKERRQ(ierr);
996 
997   bs   = v->map.bs;
998   if (start < 0) {
999     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative start %D",start);
1000   } else if (start >= bs) {
1001     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
1002             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
1003   }
1004   if (n != ns*bs) {
1005     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Subvector length * blocksize %D not correct for scatter to multicomponent vector %D",ns*bs,n);
1006   }
1007   x += start;
1008   n =  n/bs;
1009 
1010 
1011   if (addv == INSERT_VALUES) {
1012     for (i=0; i<n; i++) {
1013       x[bs*i] = y[i];
1014     }
1015   } else if (addv == ADD_VALUES) {
1016     for (i=0; i<n; i++) {
1017       x[bs*i] += y[i];
1018     }
1019 #if !defined(PETSC_USE_COMPLEX)
1020   } else if (addv == MAX_VALUES) {
1021     for (i=0; i<n; i++) {
1022       x[bs*i] = PetscMax(y[i],x[bs*i]);
1023     }
1024 #endif
1025   } else {
1026     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
1027   }
1028 
1029 
1030   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1031   ierr = VecRestoreArray(s,&y);CHKERRQ(ierr);
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 #undef __FUNCT__
1036 #define __FUNCT__ "VecReciprocal_Default"
1037 PetscErrorCode VecReciprocal_Default(Vec v)
1038 {
1039   PetscErrorCode ierr;
1040   PetscInt         i,n;
1041   PetscScalar *x;
1042 
1043   PetscFunctionBegin;
1044   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1045   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1046   for (i=0; i<n; i++) {
1047     if (x[i] != 0.0) x[i] = 1.0/x[i];
1048   }
1049   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 #undef __FUNCT__
1054 #define __FUNCT__ "VecSqrt"
1055 /*@
1056   VecSqrt - Replaces each component of a vector by the square root of its magnitude.
1057 
1058   Not collective
1059 
1060   Input Parameter:
1061 . v - The vector
1062 
1063   Output Parameter:
1064 . v - The vector square root
1065 
1066   Level: beginner
1067 
1068   Note: The actual function is sqrt(|x_i|)
1069 
1070 .keywords: vector, sqrt, square root
1071 @*/
1072 PetscErrorCode PETSCVEC_DLLEXPORT VecSqrt(Vec v)
1073 {
1074   PetscScalar *x;
1075   PetscInt         i, n;
1076   PetscErrorCode ierr;
1077 
1078   PetscFunctionBegin;
1079   PetscValidHeaderSpecific(v, VEC_COOKIE,1);
1080   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
1081   ierr = VecGetArray(v, &x);CHKERRQ(ierr);
1082   for(i = 0; i < n; i++) {
1083     x[i] = sqrt(PetscAbsScalar(x[i]));
1084   }
1085   ierr = VecRestoreArray(v, &x);CHKERRQ(ierr);
1086   PetscFunctionReturn(0);
1087 }
1088 
1089 #undef __FUNCT__
1090 #define __FUNCT__ "VecSum"
1091 /*@
1092    VecSum - Computes the sum of all the components of a vector.
1093 
1094    Collective on Vec
1095 
1096    Input Parameter:
1097 .  v - the vector
1098 
1099    Output Parameter:
1100 .  sum - the result
1101 
1102    Level: beginner
1103 
1104    Concepts: sum^of vector entries
1105 
1106 .seealso: VecNorm()
1107 @*/
1108 PetscErrorCode PETSCVEC_DLLEXPORT VecSum(Vec v,PetscScalar *sum)
1109 {
1110   PetscErrorCode ierr;
1111   PetscInt       i,n;
1112   PetscScalar    *x,lsum = 0.0;
1113 
1114   PetscFunctionBegin;
1115   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
1116   PetscValidScalarPointer(sum,2);
1117   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1118   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1119   for (i=0; i<n; i++) {
1120     lsum += x[i];
1121   }
1122   ierr = MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,PetscSum_Op,v->comm);CHKERRQ(ierr);
1123   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNCT__
1128 #define __FUNCT__ "VecShift"
1129 /*@
1130    VecShift - Shifts all of the components of a vector by computing
1131    x[i] = x[i] + shift.
1132 
1133    Collective on Vec
1134 
1135    Input Parameters:
1136 +  v - the vector
1137 -  shift - the shift
1138 
1139    Output Parameter:
1140 .  v - the shifted vector
1141 
1142    Level: intermediate
1143 
1144    Concepts: vector^adding constant
1145 
1146 @*/
1147 PetscErrorCode PETSCVEC_DLLEXPORT VecShift(Vec v,PetscScalar shift)
1148 {
1149   PetscErrorCode ierr;
1150   PetscInt       i,n;
1151   PetscScalar    *x;
1152 
1153   PetscFunctionBegin;
1154   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
1155   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1156   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1157   for (i=0; i<n; i++) {
1158     x[i] += shift;
1159   }
1160   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1161   PetscFunctionReturn(0);
1162 }
1163 
1164 #undef __FUNCT__
1165 #define __FUNCT__ "VecAbs"
1166 /*@
1167    VecAbs - Replaces every element in a vector with its absolute value.
1168 
1169    Collective on Vec
1170 
1171    Input Parameters:
1172 .  v - the vector
1173 
1174    Level: intermediate
1175 
1176    Concepts: vector^absolute value
1177 
1178 @*/
1179 PetscErrorCode PETSCVEC_DLLEXPORT VecAbs(Vec v)
1180 {
1181   PetscErrorCode ierr;
1182   PetscInt       i,n;
1183   PetscScalar    *x;
1184 
1185   PetscFunctionBegin;
1186   PetscValidHeaderSpecific(v,VEC_COOKIE,1);
1187   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1188   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1189   for (i=0; i<n; i++) {
1190     x[i] = PetscAbsScalar(x[i]);
1191   }
1192   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 #undef __FUNCT__
1197 #define __FUNCT__ "VecPermute"
1198 /*@
1199   VecPermute - Permutes a vector in place using the given ordering.
1200 
1201   Input Parameters:
1202 + vec   - The vector
1203 . order - The ordering
1204 - inv   - The flag for inverting the permutation
1205 
1206   Level: beginner
1207 
1208   Note: This function does not yet support parallel Index Sets
1209 
1210 .seealso: MatPermute()
1211 .keywords: vec, permute
1212 @*/
1213 PetscErrorCode PETSCVEC_DLLEXPORT VecPermute(Vec x, IS row, PetscTruth inv)
1214 {
1215   PetscScalar    *array, *newArray;
1216   PetscInt       *idx;
1217   PetscInt       i;
1218   PetscErrorCode ierr;
1219 
1220   PetscFunctionBegin;
1221   ierr = ISGetIndices(row, &idx);CHKERRQ(ierr);
1222   ierr = VecGetArray(x, &array);CHKERRQ(ierr);
1223   ierr = PetscMalloc(x->map.n*sizeof(PetscScalar), &newArray);CHKERRQ(ierr);
1224 #ifdef PETSC_USE_DEBUG
1225   for(i = 0; i < x->map.n; i++) {
1226     if ((idx[i] < 0) || (idx[i] >= x->map.n)) {
1227       SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Permutation index %D is out of bounds: %D", i, idx[i]);
1228     }
1229   }
1230 #endif
1231   if (!inv) {
1232     for(i = 0; i < x->map.n; i++) newArray[i]      = array[idx[i]];
1233   } else {
1234     for(i = 0; i < x->map.n; i++) newArray[idx[i]] = array[i];
1235   }
1236   ierr = VecRestoreArray(x, &array);CHKERRQ(ierr);
1237   ierr = ISRestoreIndices(row, &idx);CHKERRQ(ierr);
1238   ierr = VecReplaceArray(x, newArray);CHKERRQ(ierr);
1239   PetscFunctionReturn(0);
1240 }
1241 
1242 #undef __FUNCT__
1243 #define __FUNCT__ "VecEqual"
1244 /*@
1245    VecEqual - Compares two vectors.
1246 
1247    Collective on Vec
1248 
1249    Input Parameters:
1250 +  vec1 - the first vector
1251 -  vec2 - the second vector
1252 
1253    Output Parameter:
1254 .  flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.
1255 
1256    Level: intermediate
1257 
1258    Concepts: equal^two vectors
1259    Concepts: vector^equality
1260 
1261 @*/
1262 PetscErrorCode PETSCVEC_DLLEXPORT VecEqual(Vec vec1,Vec vec2,PetscTruth *flg)
1263 {
1264   PetscScalar    *v1,*v2;
1265   PetscErrorCode ierr;
1266   PetscInt       n1,n2,N1,N2;
1267   PetscInt       state1,state2;
1268   PetscTruth     flg1;
1269 
1270   PetscFunctionBegin;
1271   PetscValidHeaderSpecific(vec1,VEC_COOKIE,1);
1272   PetscValidHeaderSpecific(vec2,VEC_COOKIE,2);
1273   PetscValidPointer(flg,3);
1274   if (vec1 == vec2) {
1275     *flg = PETSC_TRUE;
1276   } else {
1277     ierr = VecGetSize(vec1,&N1);CHKERRQ(ierr);
1278     ierr = VecGetSize(vec2,&N2);CHKERRQ(ierr);
1279     if (N1 != N2) {
1280       flg1 = PETSC_FALSE;
1281     } else {
1282       ierr = VecGetLocalSize(vec1,&n1);CHKERRQ(ierr);
1283       ierr = VecGetLocalSize(vec2,&n2);CHKERRQ(ierr);
1284       if (n1 != n2) {
1285         flg1 = PETSC_FALSE;
1286       } else {
1287         ierr = PetscObjectStateQuery((PetscObject) vec1,&state1);CHKERRQ(ierr);
1288         ierr = PetscObjectStateQuery((PetscObject) vec2,&state2);CHKERRQ(ierr);
1289         ierr = VecGetArray(vec1,&v1);CHKERRQ(ierr);
1290         ierr = VecGetArray(vec2,&v2);CHKERRQ(ierr);
1291         ierr = PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);CHKERRQ(ierr);
1292         ierr = VecRestoreArray(vec1,&v1);CHKERRQ(ierr);
1293         ierr = VecRestoreArray(vec2,&v2);CHKERRQ(ierr);
1294         ierr = PetscObjectSetState((PetscObject) vec1,state1);CHKERRQ(ierr);
1295         ierr = PetscObjectSetState((PetscObject) vec2,state2);CHKERRQ(ierr);
1296       }
1297     }
1298     /* combine results from all processors */
1299     ierr = MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,vec1->comm);CHKERRQ(ierr);
1300   }
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 
1305 
1306