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