Actual source code: vinv.c

  1: /*$Id: vinv.c,v 1.71 2001/09/11 16:31:37 bsmith Exp $*/
  2: /*
  3:      Some useful vector utility functions.
  4: */
 5:  #include src/vec/vecimpl.h

  7: /*@C
  8:    VecStrideNorm - Computes the norm of subvector of a vector defined 
  9:    by a starting point and a stride.

 11:    Collective on Vec

 13:    Input Parameter:
 14: +  v - the vector 
 15: .  start - starting point of the subvector (defined by a stride)
 16: -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY

 18:    Output Parameter:
 19: .  norm - the norm

 21:    Notes:
 22:    One must call VecSetBlockSize() before this routine to set the stride 
 23:    information, or use a vector created from a multicomponent DA.

 25:    If x is the array representing the vector x then this computes the norm 
 26:    of the array (x[start],x[start+stride],x[start+2*stride], ....)

 28:    This is useful for computing, say the norm of the pressure variable when
 29:    the pressure is stored (interlaced) with other variables, say density etc.

 31:    This will only work if the desire subvector is a stride subvector

 33:    Level: advanced

 35:    Concepts: norm^on stride of vector
 36:    Concepts: stride^norm

 38: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
 39: @*/
 40: int VecStrideNorm(Vec v,int start,NormType ntype,PetscReal *nrm)
 41: {
 42:   int         i,n,ierr,bs;
 43:   PetscScalar *x;
 44:   PetscReal   tnorm;
 45:   MPI_Comm    comm;

 49:   VecGetLocalSize(v,&n);
 50:   VecGetArray(v,&x);
 51:   PetscObjectGetComm((PetscObject)v,&comm);

 53:   bs   = v->bs;
 54:   if (start >= bs) {
 55:     SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
 56:             Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
 57:   }
 58:   x += start;

 60:   if (ntype == NORM_2) {
 61:     PetscScalar sum = 0.0;
 62:     for (i=0; i<n; i+=bs) {
 63:       sum += x[i]*(PetscConj(x[i]));
 64:     }
 65:     tnorm  = PetscRealPart(sum);
 66:     ierr   = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);
 67:     *nrm = sqrt(*nrm);
 68:   } else if (ntype == NORM_1) {
 69:     tnorm = 0.0;
 70:     for (i=0; i<n; i+=bs) {
 71:       tnorm += PetscAbsScalar(x[i]);
 72:     }
 73:     ierr   = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);
 74:   } else if (ntype == NORM_INFINITY) {
 75:     PetscReal tmp;
 76:     tnorm = 0.0;

 78:     for (i=0; i<n; i+=bs) {
 79:       if ((tmp = PetscAbsScalar(x[i])) > tnorm) tnorm = tmp;
 80:       /* check special case of tmp == NaN */
 81:       if (tmp != tmp) {tnorm = tmp; break;}
 82:     }
 83:     ierr   = MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_MAX,comm);
 84:   } else {
 85:     SETERRQ(1,"Unknown norm type");
 86:   }

 88:   VecRestoreArray(v,&x);
 89:   return(0);
 90: }

 92: /*@C
 93:    VecStrideMax - Computes the maximum of subvector of a vector defined 
 94:    by a starting point and a stride and optionally its location.

 96:    Collective on Vec

 98:    Input Parameter:
 99: +  v - the vector 
100: -  start - starting point of the subvector (defined by a stride)

102:    Output Parameter:
103: +  index - the location where the maximum occurred (not supported, pass PETSC_NULL,
104:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
105: -  nrm - the max

107:    Notes:
108:    One must call VecSetBlockSize() before this routine to set the stride 
109:    information, or use a vector created from a multicomponent DA.

111:    If xa is the array representing the vector x, then this computes the max
112:    of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)

114:    This is useful for computing, say the maximum of the pressure variable when
115:    the pressure is stored (interlaced) with other variables, e.g., density, etc.
116:    This will only work if the desire subvector is a stride subvector.

118:    Level: advanced

120:    Concepts: maximum^on stride of vector
121:    Concepts: stride^maximum

123: .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
124: @*/
125: int VecStrideMax(Vec v,int start,int *index,PetscReal *nrm)
126: {
127:   int         i,n,ierr,bs;
128:   PetscScalar *x;
129:   PetscReal   max,tmp;
130:   MPI_Comm    comm;

134:   if (index) {
135:     SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
136:   }
137:   VecGetLocalSize(v,&n);
138:   VecGetArray(v,&x);
139:   PetscObjectGetComm((PetscObject)v,&comm);

141:   bs   = v->bs;
142:   if (start >= bs) {
143:     SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
144:             Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
145:   }
146:   x += start;

148:   if (!n) {
149:     max = PETSC_MIN;
150:   } else {
151: #if defined(PETSC_USE_COMPLEX)
152:     max = PetscRealPart(x[0]);
153: #else
154:     max = x[0];
155: #endif
156:     for (i=bs; i<n; i+=bs) {
157: #if defined(PETSC_USE_COMPLEX)
158:       if ((tmp = PetscRealPart(x[i])) > max) { max = tmp;}
159: #else
160:       if ((tmp = x[i]) > max) { max = tmp; }
161: #endif
162:     }
163:   }
164:   ierr   = MPI_Allreduce(&max,nrm,1,MPIU_REAL,MPI_MAX,comm);

166:   VecRestoreArray(v,&x);
167:   return(0);
168: }

170: /*@C
171:    VecStrideMin - Computes the minimum of subvector of a vector defined 
172:    by a starting point and a stride and optionally its location.

174:    Collective on Vec

176:    Input Parameter:
177: +  v - the vector 
178: -  start - starting point of the subvector (defined by a stride)

180:    Output Parameter:
181: +  index - the location where the minimum occurred (not supported, pass PETSC_NULL,
182:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
183: -  nrm - the min

185:    Level: advanced

187:    Notes:
188:    One must call VecSetBlockSize() before this routine to set the stride 
189:    information, or use a vector created from a multicomponent DA.

191:    If xa is the array representing the vector x, then this computes the min
192:    of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)

194:    This is useful for computing, say the minimum of the pressure variable when
195:    the pressure is stored (interlaced) with other variables, e.g., density, etc.
196:    This will only work if the desire subvector is a stride subvector.

198:    Concepts: minimum^on stride of vector
199:    Concepts: stride^minimum

201: .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
202: @*/
203: int VecStrideMin(Vec v,int start,int *index,PetscReal *nrm)
204: {
205:   int         i,n,ierr,bs;
206:   PetscScalar *x;
207:   PetscReal   min,tmp;
208:   MPI_Comm    comm;

212:   if (index) {
213:     SETERRQ(1,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
214:   }
215:   VecGetLocalSize(v,&n);
216:   VecGetArray(v,&x);
217:   PetscObjectGetComm((PetscObject)v,&comm);

219:   bs   = v->bs;
220:   if (start >= bs) {
221:     SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
222:             Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
223:   }
224:   x += start;

226:   if (!n) {
227:     min = PETSC_MAX;
228:   } else {
229: #if defined(PETSC_USE_COMPLEX)
230:     min = PetscRealPart(x[0]);
231: #else
232:     min = x[0];
233: #endif
234:     for (i=bs; i<n; i+=bs) {
235: #if defined(PETSC_USE_COMPLEX)
236:       if ((tmp = PetscRealPart(x[i])) < min) { min = tmp;}
237: #else
238:       if ((tmp = x[i]) < min) { min = tmp; }
239: #endif
240:     }
241:   }
242:   ierr   = MPI_Allreduce(&min,nrm,1,MPIU_REAL,MPI_MIN,comm);

244:   VecRestoreArray(v,&x);
245:   return(0);
246: }

248: /*@
249:    VecStrideGatherAll - Gathers all the single components from a multi-component vector into
250:    seperate vectors.

252:    Collective on Vec

254:    Input Parameter:
255: +  v - the vector 
256: -  addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES

258:    Output Parameter:
259: .  s - the location where the subvectors are stored

261:    Notes:
262:    One must call VecSetBlockSize() before this routine to set the stride 
263:    information, or use a vector created from a multicomponent DA.

265:    If x is the array representing the vector x then this gathers
266:    the arrays (x[start],x[start+stride],x[start+2*stride], ....)
267:    for start=0,1,2,...bs-1

269:    The parallel layout of the vector and the subvector must be the same;
270:    i.e., nlocal of v = stride*(nlocal of s) 

272:    Level: advanced

274:    Concepts: gather^into strided vector

276: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
277:           VecStrideScatterAll()
278: @*/
279: int VecStrideGatherAll(Vec v,Vec *s,InsertMode addv)
280: {
281:   int          i,n,ierr,bs,ns,j;
282:   PetscScalar  *x,**y;

287:   VecGetLocalSize(v,&n);
288:   VecGetLocalSize(*s,&ns);
289:   VecGetArray(v,&x);
290:   bs   = v->bs;

292:   PetscMalloc(bs*sizeof(PetscReal*),&y);
293:   for (i=0; i<bs; i++) {
294:     VecGetArray(s[i],&y[i]);
295:   }

297:   if (n != ns*bs) {
298:     SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
299:   }
300:   n =  n/bs;

302:   if (addv == INSERT_VALUES) {
303:     for (j=0; j<bs; j++) {
304:       for (i=0; i<n; i++) {
305:         y[j][i] = x[bs*i+j];
306:       }
307:     }
308:   } else if (addv == ADD_VALUES) {
309:     for (j=0; j<bs; j++) {
310:       for (i=0; i<n; i++) {
311:         y[j][i] += x[bs*i+j];
312:       }
313:     }
314: #if !defined(PETSC_USE_COMPLEX)
315:   } else if (addv == MAX_VALUES) {
316:     for (j=0; j<bs; j++) {
317:       for (i=0; i<n; i++) {
318:         y[j][i] = PetscMax(y[j][i],x[bs*i+j]);
319:       }
320:     }
321: #endif
322:   } else {
323:     SETERRQ(1,"Unknown insert type");
324:   }

326:   VecRestoreArray(v,&x);
327:   for (i=0; i<bs; i++) {
328:     VecRestoreArray(s[i],&y[i]);
329:   }
330:   PetscFree(y);
331:   return(0);
332: }
333: /*@
334:    VecStrideScatterAll - Scatters all the single components from seperate vectors into 
335:      a multi-component vector.

337:    Collective on Vec

339:    Input Parameter:
340: +  s - the location where the subvectors are stored
341: -  addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES

343:    Output Parameter:
344: .  v - the multicomponent vector 

346:    Notes:
347:    One must call VecSetBlockSize() before this routine to set the stride 
348:    information, or use a vector created from a multicomponent DA.

350:    The parallel layout of the vector and the subvector must be the same;
351:    i.e., nlocal of v = stride*(nlocal of s) 

353:    Level: advanced

355:    Concepts:  scatter^into strided vector

357: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
358:           VecStrideScatterAll()
359: @*/
360: int VecStrideScatterAll(Vec *s,Vec v,InsertMode addv)
361: {
362:   int          i,n,ierr,bs,ns,j;
363:   PetscScalar  *x,**y;

368:   VecGetLocalSize(v,&n);
369:   VecGetLocalSize(*s,&ns);
370:   VecGetArray(v,&x);
371:   bs   = v->bs;

373:   PetscMalloc(bs*sizeof(PetscReal*),&y);
374:   for (i=0; i<bs; i++) {
375:     VecGetArray(s[i],&y[i]);
376:   }

378:   if (n != ns*bs) {
379:     SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
380:   }
381:   n =  n/bs;

383:   if (addv == INSERT_VALUES) {
384:     for (j=0; j<bs; j++) {
385:       for (i=0; i<n; i++) {
386:         x[bs*i+j] = y[j][i];
387:       }
388:     }
389:   } else if (addv == ADD_VALUES) {
390:     for (j=0; j<bs; j++) {
391:       for (i=0; i<n; i++) {
392:         x[bs*i+j] += y[j][i];
393:       }
394:     }
395: #if !defined(PETSC_USE_COMPLEX)
396:   } else if (addv == MAX_VALUES) {
397:     for (j=0; j<bs; j++) {
398:       for (i=0; i<n; i++) {
399:         x[bs*i+j] = PetscMax(y[j][i],x[bs*i+j]);
400:       }
401:     }
402: #endif
403:   } else {
404:     SETERRQ(1,"Unknown insert type");
405:   }

407:   VecRestoreArray(v,&x);
408:   for (i=0; i<bs; i++) {
409:     VecRestoreArray(s[i],&y[i]);
410:   }
411:   PetscFree(y);
412:   return(0);
413: }

415: /*@
416:    VecStrideGather - Gathers a single component from a multi-component vector into
417:    another vector.

419:    Collective on Vec

421:    Input Parameter:
422: +  v - the vector 
423: .  start - starting point of the subvector (defined by a stride)
424: -  addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES

426:    Output Parameter:
427: .  s - the location where the subvector is stored

429:    Notes:
430:    One must call VecSetBlockSize() before this routine to set the stride 
431:    information, or use a vector created from a multicomponent DA.

433:    If x is the array representing the vector x then this gathers
434:    the array (x[start],x[start+stride],x[start+2*stride], ....)

436:    The parallel layout of the vector and the subvector must be the same;
437:    i.e., nlocal of v = stride*(nlocal of s) 

439:    Level: advanced

441:    Concepts: gather^into strided vector

443: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
444:           VecStrideScatterAll()
445: @*/
446: int VecStrideGather(Vec v,int start,Vec s,InsertMode addv)
447: {
448:   int          i,n,ierr,bs,ns;
449:   PetscScalar  *x,*y;

454:   VecGetLocalSize(v,&n);
455:   VecGetLocalSize(s,&ns);
456:   VecGetArray(v,&x);
457:   VecGetArray(s,&y);

459:   bs   = v->bs;
460:   if (start >= bs) {
461:     SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
462:             Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
463:   }
464:   if (n != ns*bs) {
465:     SETERRQ2(1,"Subvector length * blocksize %d not correct for gather from original vector %d",ns*bs,n);
466:   }
467:   x += start;
468:   n =  n/bs;

470:   if (addv == INSERT_VALUES) {
471:     for (i=0; i<n; i++) {
472:       y[i] = x[bs*i];
473:     }
474:   } else if (addv == ADD_VALUES) {
475:     for (i=0; i<n; i++) {
476:       y[i] += x[bs*i];
477:     }
478: #if !defined(PETSC_USE_COMPLEX)
479:   } else if (addv == MAX_VALUES) {
480:     for (i=0; i<n; i++) {
481:       y[i] = PetscMax(y[i],x[bs*i]);
482:     }
483: #endif
484:   } else {
485:     SETERRQ(1,"Unknown insert type");
486:   }

488:   VecRestoreArray(v,&x);
489:   VecRestoreArray(s,&y);
490:   return(0);
491: }

493: /*@
494:    VecStrideScatter - Scatters a single component from a vector into a multi-component vector.

496:    Collective on Vec

498:    Input Parameter:
499: +  s - the single-component vector 
500: .  start - starting point of the subvector (defined by a stride)
501: -  addv - one of ADD_VALUES,SET_VALUES,MAX_VALUES

503:    Output Parameter:
504: .  v - the location where the subvector is scattered (the multi-component vector)

506:    Notes:
507:    One must call VecSetBlockSize() on the multi-component vector before this
508:    routine to set the stride  information, or use a vector created from a multicomponent DA.

510:    The parallel layout of the vector and the subvector must be the same;
511:    i.e., nlocal of v = stride*(nlocal of s) 

513:    Level: advanced

515:    Concepts: scatter^into strided vector

517: .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
518:           VecStrideScatterAll()
519: @*/
520: int VecStrideScatter(Vec s,int start,Vec v,InsertMode addv)
521: {
522:   int          i,n,ierr,bs,ns;
523:   PetscScalar  *x,*y;

528:   VecGetLocalSize(v,&n);
529:   VecGetLocalSize(s,&ns);
530:   VecGetArray(v,&x);
531:   VecGetArray(s,&y);

533:   bs   = v->bs;
534:   if (start >= bs) {
535:     SETERRQ2(1,"Start of stride subvector (%d) is too large for striden
536:             Have you set the vector blocksize (%d) correctly with VecSetBlockSize()?",start,bs);
537:   }
538:   if (n != ns*bs) {
539:     SETERRQ2(1,"Subvector length * blocksize %d not correct for scatter to multicomponent vector %d",ns*bs,n);
540:   }
541:   x += start;
542:   n =  n/bs;


545:   if (addv == INSERT_VALUES) {
546:     for (i=0; i<n; i++) {
547:       x[bs*i] = y[i];
548:     }
549:   } else if (addv == ADD_VALUES) {
550:     for (i=0; i<n; i++) {
551:       x[bs*i] += y[i];
552:     }
553: #if !defined(PETSC_USE_COMPLEX)
554:   } else if (addv == MAX_VALUES) {
555:     for (i=0; i<n; i++) {
556:       x[bs*i] = PetscMax(y[i],x[bs*i]);
557:     }
558: #endif
559:   } else {
560:     SETERRQ(1,"Unknown insert type");
561:   }


564:   VecRestoreArray(v,&x);
565:   VecRestoreArray(s,&y);
566:   return(0);
567: }

569: int VecReciprocal_Default(Vec v)
570: {
571:   int         i,n,ierr;
572:   PetscScalar *x;

575:   VecGetLocalSize(v,&n);
576:   VecGetArrayFast(v,&x);
577:   for (i=0; i<n; i++) {
578:     if (x[i] != 0.0) x[i] = 1.0/x[i];
579:   }
580:   VecRestoreArrayFast(v,&x);
581:   return(0);
582: }

584: /*@
585:   VecSqrt - Replaces each component of a vector by the square root of its magnitude.

587:   Not collective

589:   Input Parameter:
590: . v - The vector

592:   Output Parameter:
593: . v - The vector square root

595:   Level: beginner

597:   Note: The actual function is sqrt(|x_i|)

599: .keywords: vector, sqrt, square root
600: @*/
601: int VecSqrt(Vec v)
602: {
603:   PetscScalar *x;
604:   int         i, n;
605:   int         ierr;

609:   VecGetLocalSize(v, &n);
610:   VecGetArray(v, &x);
611:   for(i = 0; i < n; i++) {
612:     x[i] = sqrt(PetscAbsScalar(x[i]));
613:   }
614:   VecRestoreArray(v, &x);
615:   return(0);
616: }

618: /*@
619:    VecSum - Computes the sum of all the components of a vector.

621:    Collective on Vec

623:    Input Parameter:
624: .  v - the vector 

626:    Output Parameter:
627: .  sum - the result

629:    Level: beginner

631:    Concepts: sum^of vector entries

633: .seealso: VecNorm()
634: @*/
635: int VecSum(Vec v,PetscScalar *sum)
636: {
637:   int         i,n,ierr;
638:   PetscScalar *x,lsum = 0.0;

642:   VecGetLocalSize(v,&n);
643:   VecGetArray(v,&x);
644:   for (i=0; i<n; i++) {
645:     lsum += x[i];
646:   }
647:   MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,PetscSum_Op,v->comm);
648:   VecRestoreArray(v,&x);
649:   return(0);
650: }

652: /*@
653:    VecShift - Shifts all of the components of a vector by computing
654:    x[i] = x[i] + shift.

656:    Collective on Vec

658:    Input Parameters:
659: +  v - the vector 
660: -  shift - the shift

662:    Output Parameter:
663: .  v - the shifted vector 

665:    Level: intermediate

667:    Concepts: vector^adding constant

669: @*/
670: int VecShift(const PetscScalar *shift,Vec v)
671: {
672:   int         i,n,ierr;
673:   PetscScalar *x,lsum = *shift;

677:   VecGetLocalSize(v,&n);
678:   VecGetArray(v,&x);
679:   for (i=0; i<n; i++) {
680:     x[i] += lsum;
681:   }
682:   VecRestoreArray(v,&x);
683:   return(0);
684: }

686: /*@
687:    VecAbs - Replaces every element in a vector with its absolute value.

689:    Collective on Vec

691:    Input Parameters:
692: .  v - the vector 

694:    Level: intermediate

696:    Concepts: vector^absolute value

698: @*/
699: int VecAbs(Vec v)
700: {
701:   int         i,n,ierr;
702:   PetscScalar *x;

706:   VecGetLocalSize(v,&n);
707:   VecGetArray(v,&x);
708:   for (i=0; i<n; i++) {
709:     x[i] = PetscAbsScalar(x[i]);
710:   }
711:   VecRestoreArray(v,&x);
712:   return(0);
713: }

715: /*@
716:   VecPermute - Permutes a vector in place using the given ordering.

718:   Input Parameters:
719: + vec   - The vector
720: . order - The ordering
721: - inv   - The flag for inverting the permutation

723:   Level: beginner

725:   Note: This function does not yet support parallel Index Sets

727: .seealso: MatPermute()
728: .keywords: vec, permute
729: @*/
730: int VecPermute(Vec x, IS row, PetscTruth inv)
731: {
732:   PetscScalar *array, *newArray;
733:   int         *idx;
734:   int          i;
735:   int          ierr;

738:   ISGetIndices(row, &idx);
739:   VecGetArray(x, &array);
740:   PetscMalloc((x->n+1) * sizeof(PetscScalar), &newArray);
741: #ifdef PETSC_USE_BOPT_g
742:   for(i = 0; i < x->n; i++) {
743:     if ((idx[i] < 0) || (idx[i] >= x->n)) {
744:       SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Permutation index %d is out of bounds: %d", i, idx[i]);
745:     }
746:   }
747: #endif
748:   if (inv == PETSC_FALSE) {
749:     for(i = 0; i < x->n; i++) newArray[i]      = array[idx[i]];
750:   } else {
751:     for(i = 0; i < x->n; i++) newArray[idx[i]] = array[i];
752:   }
753:   VecRestoreArray(x, &array);
754:   ISRestoreIndices(row, &idx);
755:   VecReplaceArray(x, newArray);
756:   return(0);
757: }

759: /*@
760:    VecEqual - Compares two vectors.

762:    Collective on Vec

764:    Input Parameters:
765: +  vec1 - the first matrix
766: -  vec2 - the second matrix

768:    Output Parameter:
769: .  flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.

771:    Level: intermediate

773:    Concepts: equal^two vectors
774:    Concepts: vector^equality

776: @*/
777: int VecEqual(Vec vec1,Vec vec2,PetscTruth *flg)
778: {
779:   PetscScalar  *v1,*v2;
780:   int          n1,n2,ierr;
781:   PetscTruth   flg1;

784:   VecGetSize(vec1,&n1);
785:   VecGetSize(vec2,&n2);
786:   if (n1 != n2) {
787:     flg1 = PETSC_FALSE;
788:   } else {
789:     VecGetArray(vec1,&v1);
790:     VecGetArray(vec2,&v2);
791:     PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);
792:     VecRestoreArray(vec1,&v1);
793:     VecRestoreArray(vec2,&v2);
794:   }

796:   /* combine results from all processors */
797:   MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,vec1->comm);
798: 

800:   return(0);
801: }