Actual source code: vinv.c

  1: /*
  2:      Some useful vector utility functions.
  3: */
 4:  #include vecimpl.h

  8: /*@C
  9:    VecStrideScale - Scales a subvector of a vector defined 
 10:    by a starting point and a stride.

 12:    Collective on Vec

 14:    Input Parameter:
 15: +  v - the vector 
 16: .  start - starting point of the subvector (defined by a stride)
 17: -  scale - value to multiply each subvector entry by

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

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

 25:    Level: advanced

 27:    Concepts: scale^on stride of vector
 28:    Concepts: stride^scale

 30: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideScale()
 31: @*/
 32: PetscErrorCode VecStrideScale(Vec v,PetscInt start,PetscScalar *scale)
 33: {
 35:   PetscInt       i,n,bs;
 36:   PetscScalar    *x,xscale = *scale;

 41:   VecGetLocalSize(v,&n);
 42:   VecGetArray(v,&x);

 44:   bs   = v->bs;
 45:   if (start >= bs) {
 46:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
 47:             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
 48:   }
 49:   x += start;

 51:   for (i=0; i<n; i+=bs) {
 52:     x[i] *= xscale;
 53:   }
 54:   x -= start;
 55:   VecRestoreArray(v,&x);
 56:   return(0);
 57: }

 61: /*@C
 62:    VecStrideNorm - Computes the norm of subvector of a vector defined 
 63:    by a starting point and a stride.

 65:    Collective on Vec

 67:    Input Parameter:
 68: +  v - the vector 
 69: .  start - starting point of the subvector (defined by a stride)
 70: -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY

 72:    Output Parameter:
 73: .  norm - the norm

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

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

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

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

 87:    Level: advanced

 89:    Concepts: norm^on stride of vector
 90:    Concepts: stride^norm

 92: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
 93: @*/
 94: PetscErrorCode VecStrideNorm(Vec v,PetscInt start,NormType ntype,PetscReal *nrm)
 95: {
 97:   PetscInt       i,n,bs;
 98:   PetscScalar    *x;
 99:   PetscReal      tnorm;
100:   MPI_Comm       comm;

105:   VecGetLocalSize(v,&n);
106:   VecGetArray(v,&x);
107:   PetscObjectGetComm((PetscObject)v,&comm);

109:   bs   = v->bs;
110:   if (start >= bs) {
111:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
112:             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
113:   }
114:   x += start;

116:   if (ntype == NORM_2) {
117:     PetscScalar sum = 0.0;
118:     for (i=0; i<n; i+=bs) {
119:       sum += x[i]*(PetscConj(x[i]));
120:     }
121:     tnorm  = PetscRealPart(sum);
122:     MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);
123:     *nrm = sqrt(*nrm);
124:   } else if (ntype == NORM_1) {
125:     tnorm = 0.0;
126:     for (i=0; i<n; i+=bs) {
127:       tnorm += PetscAbsScalar(x[i]);
128:     }
129:     MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_SUM,comm);
130:   } else if (ntype == NORM_INFINITY) {
131:     PetscReal tmp;
132:     tnorm = 0.0;

134:     for (i=0; i<n; i+=bs) {
135:       if ((tmp = PetscAbsScalar(x[i])) > tnorm) tnorm = tmp;
136:       /* check special case of tmp == NaN */
137:       if (tmp != tmp) {tnorm = tmp; break;}
138:     }
139:     MPI_Allreduce(&tnorm,nrm,1,MPIU_REAL,MPI_MAX,comm);
140:   } else {
141:     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown norm type");
142:   }

144:   VecRestoreArray(v,&x);
145:   return(0);
146: }

150: /*@
151:    VecStrideMax - Computes the maximum of subvector of a vector defined 
152:    by a starting point and a stride and optionally its location.

154:    Collective on Vec

156:    Input Parameter:
157: +  v - the vector 
158: -  start - starting point of the subvector (defined by a stride)

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

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

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

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

176:    Level: advanced

178:    Concepts: maximum^on stride of vector
179:    Concepts: stride^maximum

181: .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
182: @*/
183: PetscErrorCode VecStrideMax(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
184: {
186:   PetscInt       i,n,bs;
187:   PetscScalar    *x;
188:   PetscReal      max,tmp;
189:   MPI_Comm       comm;

194:   if (idex) {
195:     SETERRQ(PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
196:   }
197:   VecGetLocalSize(v,&n);
198:   VecGetArray(v,&x);
199:   PetscObjectGetComm((PetscObject)v,&comm);

201:   bs   = v->bs;
202:   if (start >= bs) {
203:     SETERRQ2(PETSC_ERR_ARG_WRONG,"Start of stride subvector (%D) is too large for stride\n\
204:             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
205:   }
206:   x += start;

208:   if (!n) {
209:     max = PETSC_MIN;
210:   } else {
211: #if defined(PETSC_USE_COMPLEX)
212:     max = PetscRealPart(x[0]);
213: #else
214:     max = x[0];
215: #endif
216:     for (i=bs; i<n; i+=bs) {
217: #if defined(PETSC_USE_COMPLEX)
218:       if ((tmp = PetscRealPart(x[i])) > max) { max = tmp;}
219: #else
220:       if ((tmp = x[i]) > max) { max = tmp; }
221: #endif
222:     }
223:   }
224:   MPI_Allreduce(&max,nrm,1,MPIU_REAL,MPI_MAX,comm);

226:   VecRestoreArray(v,&x);
227:   return(0);
228: }

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

236:    Collective on Vec

238:    Input Parameter:
239: +  v - the vector 
240: -  start - starting point of the subvector (defined by a stride)

242:    Output Parameter:
243: +  idex - the location where the minimum occurred (not supported, pass PETSC_NULL,
244:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
245: -  nrm - the min

247:    Level: advanced

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

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

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

260:    Concepts: minimum^on stride of vector
261:    Concepts: stride^minimum

263: .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
264: @*/
265: PetscErrorCode VecStrideMin(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
266: {
268:   PetscInt       i,n,bs;
269:   PetscScalar    *x;
270:   PetscReal      min,tmp;
271:   MPI_Comm       comm;

276:   if (idex) {
277:     SETERRQ(PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
278:   }
279:   VecGetLocalSize(v,&n);
280:   VecGetArray(v,&x);
281:   PetscObjectGetComm((PetscObject)v,&comm);

283:   bs   = v->bs;
284:   if (start >= bs) {
285:     SETERRQ2(PETSC_ERR_ARG_WRONG,"Start of stride subvector (%D) is too large for stride\n\
286:             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
287:   }
288:   x += start;

290:   if (!n) {
291:     min = PETSC_MAX;
292:   } else {
293: #if defined(PETSC_USE_COMPLEX)
294:     min = PetscRealPart(x[0]);
295: #else
296:     min = x[0];
297: #endif
298:     for (i=bs; i<n; i+=bs) {
299: #if defined(PETSC_USE_COMPLEX)
300:       if ((tmp = PetscRealPart(x[i])) < min) { min = tmp;}
301: #else
302:       if ((tmp = x[i]) < min) { min = tmp; }
303: #endif
304:     }
305:   }
306:   MPI_Allreduce(&min,nrm,1,MPIU_REAL,MPI_MIN,comm);

308:   VecRestoreArray(v,&x);
309:   return(0);
310: }

314: /*@C
315:    VecStrideScaleAll - Scales the subvectors of a vector defined 
316:    by a starting point and a stride.

318:    Collective on Vec

320:    Input Parameter:
321: +  v - the vector 
322: -  scales - values to multiply each subvector entry by

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


329:    Level: advanced

331:    Concepts: scale^on stride of vector
332:    Concepts: stride^scale

334: .seealso: VecNorm(), VecStrideScale(), VecScale(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
335: @*/
336: PetscErrorCode VecStrideScaleAll(Vec v,PetscScalar *scales)
337: {
339:   PetscInt         i,j,n,bs;
340:   PetscScalar *x;

345:   VecGetLocalSize(v,&n);
346:   VecGetArray(v,&x);

348:   bs   = v->bs;

350:   /* need to provide optimized code for each bs */
351:   for (i=0; i<n; i+=bs) {
352:     for (j=0; j<bs; j++) {
353:       x[i+j] *= scales[j];
354:     }
355:   }
356:   VecRestoreArray(v,&x);
357:   return(0);
358: }

362: /*@C
363:    VecStrideNormAll - Computes the norms  subvectors of a vector defined 
364:    by a starting point and a stride.

366:    Collective on Vec

368:    Input Parameter:
369: +  v - the vector 
370: -  ntype - type of norm, one of NORM_1, NORM_2, NORM_INFINITY

372:    Output Parameter:
373: .  nrm - the norms

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

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

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

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

387:    Level: advanced

389:    Concepts: norm^on stride of vector
390:    Concepts: stride^norm

392: .seealso: VecNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin(), VecStrideMax()
393: @*/
394: PetscErrorCode VecStrideNormAll(Vec v,NormType ntype,PetscReal *nrm)
395: {
397:   PetscInt         i,j,n,bs;
398:   PetscScalar *x;
399:   PetscReal   tnorm[128];
400:   MPI_Comm    comm;

405:   VecGetLocalSize(v,&n);
406:   VecGetArray(v,&x);
407:   PetscObjectGetComm((PetscObject)v,&comm);

409:   bs   = v->bs;
410:   if (bs > 128) SETERRQ(PETSC_ERR_SUP,"Currently supports only blocksize up to 128");

412:   if (ntype == NORM_2) {
413:     PetscScalar sum[128];
414:     for (j=0; j<bs; j++) sum[j] = 0.0;
415:     for (i=0; i<n; i+=bs) {
416:       for (j=0; j<bs; j++) {
417:         sum[j] += x[i+j]*(PetscConj(x[i+j]));
418:       }
419:     }
420:     for (j=0; j<bs; j++) {
421:       tnorm[j]  = PetscRealPart(sum[j]);
422:     }
423:     MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_SUM,comm);
424:     for (j=0; j<bs; j++) {
425:       nrm[j] = sqrt(nrm[j]);
426:     }
427:   } else if (ntype == NORM_1) {
428:     for (j=0; j<bs; j++) {
429:       tnorm[j] = 0.0;
430:     }
431:     for (i=0; i<n; i+=bs) {
432:       for (j=0; j<bs; j++) {
433:         tnorm[j] += PetscAbsScalar(x[i+j]);
434:       }
435:     }
436:     MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_SUM,comm);
437:   } else if (ntype == NORM_INFINITY) {
438:     PetscReal tmp;
439:     for (j=0; j<bs; j++) {
440:       tnorm[j] = 0.0;
441:     }

443:     for (i=0; i<n; i+=bs) {
444:       for (j=0; j<bs; j++) {
445:         if ((tmp = PetscAbsScalar(x[i+j])) > tnorm[j]) tnorm[j] = tmp;
446:         /* check special case of tmp == NaN */
447:         if (tmp != tmp) {tnorm[j] = tmp; break;}
448:       }
449:     }
450:     MPI_Allreduce(tnorm,nrm,bs,MPIU_REAL,MPI_MAX,comm);
451:   } else {
452:     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown norm type");
453:   }

455:   VecRestoreArray(v,&x);
456:   return(0);
457: }

461: /*@C
462:    VecStrideMaxAll - Computes the maximums of subvectors of a vector defined 
463:    by a starting point and a stride and optionally its location.

465:    Collective on Vec

467:    Input Parameter:
468: .  v - the vector 

470:    Output Parameter:
471: +  index - the location where the maximum occurred (not supported, pass PETSC_NULL,
472:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
473: -  nrm - the maximums

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

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

483:    Level: advanced

485:    Concepts: maximum^on stride of vector
486:    Concepts: stride^maximum

488: .seealso: VecMax(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMin()
489: @*/
490: PetscErrorCode VecStrideMaxAll(Vec v,PetscInt *idex,PetscReal *nrm)
491: {
493:   PetscInt         i,j,n,bs;
494:   PetscScalar *x;
495:   PetscReal   max[128],tmp;
496:   MPI_Comm    comm;

501:   if (idex) {
502:     SETERRQ(PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
503:   }
504:   VecGetLocalSize(v,&n);
505:   VecGetArray(v,&x);
506:   PetscObjectGetComm((PetscObject)v,&comm);

508:   bs   = v->bs;
509:   if (bs > 128) SETERRQ(PETSC_ERR_SUP,"Currently supports only blocksize up to 128");

511:   if (!n) {
512:     for (j=0; j<bs; j++) {
513:       max[j] = PETSC_MIN;
514:     }
515:   } else {
516:     for (j=0; j<bs; j++) {
517: #if defined(PETSC_USE_COMPLEX)
518:       max[j] = PetscRealPart(x[j]);
519: #else
520:       max[j] = x[j];
521: #endif
522:     }
523:     for (i=bs; i<n; i+=bs) {
524:       for (j=0; j<bs; j++) {
525: #if defined(PETSC_USE_COMPLEX)
526:         if ((tmp = PetscRealPart(x[i+j])) > max[j]) { max[j] = tmp;}
527: #else
528:         if ((tmp = x[i+j]) > max[j]) { max[j] = tmp; }
529: #endif
530:       }
531:     }
532:   }
533:   MPI_Allreduce(max,nrm,bs,MPIU_REAL,MPI_MAX,comm);

535:   VecRestoreArray(v,&x);
536:   return(0);
537: }

541: /*@C
542:    VecStrideMinAll - Computes the minimum of subvector of a vector defined 
543:    by a starting point and a stride and optionally its location.

545:    Collective on Vec

547:    Input Parameter:
548: .  v - the vector 

550:    Output Parameter:
551: +  idex - the location where the minimum occurred (not supported, pass PETSC_NULL,
552:            if you need this, send mail to petsc-maint@mcs.anl.gov to request it)
553: -  nrm - the minimums

555:    Level: advanced

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

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

565:    Concepts: minimum^on stride of vector
566:    Concepts: stride^minimum

568: .seealso: VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()
569: @*/
570: PetscErrorCode VecStrideMinAll(Vec v,PetscInt *idex,PetscReal *nrm)
571: {
573:   PetscInt         i,n,bs,j;
574:   PetscScalar *x;
575:   PetscReal   min[128],tmp;
576:   MPI_Comm    comm;

581:   if (idex) {
582:     SETERRQ(PETSC_ERR_SUP,"No support yet for returning index; send mail to petsc-maint@mcs.anl.gov asking for it");
583:   }
584:   VecGetLocalSize(v,&n);
585:   VecGetArray(v,&x);
586:   PetscObjectGetComm((PetscObject)v,&comm);

588:   bs   = v->bs;
589:   if (bs > 128) SETERRQ(PETSC_ERR_SUP,"Currently supports only blocksize up to 128");

591:   if (!n) {
592:     for (j=0; j<bs; j++) {
593:       min[j] = PETSC_MAX;
594:     }
595:   } else {
596:     for (j=0; j<bs; j++) {
597: #if defined(PETSC_USE_COMPLEX)
598:       min[j] = PetscRealPart(x[j]);
599: #else
600:       min[j] = x[j];
601: #endif
602:     }
603:     for (i=bs; i<n; i+=bs) {
604:       for (j=0; j<bs; j++) {
605: #if defined(PETSC_USE_COMPLEX)
606:         if ((tmp = PetscRealPart(x[i+j])) < min[j]) { min[j] = tmp;}
607: #else
608:         if ((tmp = x[i+j]) < min[j]) { min[j] = tmp; }
609: #endif
610:       }
611:     }
612:   }
613:   MPI_Allreduce(min,nrm,bs,MPIU_REAL,MPI_MIN,comm);

615:   VecRestoreArray(v,&x);
616:   return(0);
617: }

619: /*----------------------------------------------------------------------------------------------*/
622: /*@
623:    VecStrideGatherAll - Gathers all the single components from a multi-component vector into
624:    seperate vectors.

626:    Collective on Vec

628:    Input Parameter:
629: +  v - the vector 
630: -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES

632:    Output Parameter:
633: .  s - the location where the subvectors are stored

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

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

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

646:    Not optimized; could be easily

648:    Level: advanced

650:    Concepts: gather^into strided vector

652: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
653:           VecStrideScatterAll()
654: @*/
655: PetscErrorCode VecStrideGatherAll(Vec v,Vec *s,InsertMode addv)
656: {
658:   PetscInt       i,n,bs,j,k,*bss,nv,jj,nvc;
659:   PetscScalar    *x,**y;

665:   VecGetLocalSize(v,&n);
666:   VecGetArray(v,&x);
667:   bs   = v->bs;
668:   if (bs < 0) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set");

670:   PetscMalloc2(bs,PetscReal*,&y,bs,PetscInt,&bss);
671:   nv   = 0;
672:   nvc  = 0;
673:   for (i=0; i<bs; i++) {
674:     VecGetBlockSize(s[i],&bss[i]);
675:     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
676:     VecGetArray(s[i],&y[i]);
677:     nvc  += bss[i];
678:     nv++;
679:     if (nvc > bs)  SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector");
680:     if (nvc == bs) break;
681:   }

683:   n =  n/bs;

685:   jj = 0;
686:   if (addv == INSERT_VALUES) {
687:     for (j=0; j<nv; j++) {
688:       for (k=0; k<bss[j]; k++) {
689:         for (i=0; i<n; i++) {
690:           y[j][i*bss[j] + k] = x[bs*i+jj+k];
691:         }
692:       }
693:       jj += bss[j];
694:     }
695:   } else if (addv == ADD_VALUES) {
696:     for (j=0; j<nv; j++) {
697:       for (k=0; k<bss[j]; k++) {
698:         for (i=0; i<n; i++) {
699:           y[j][i*bss[j] + k] += x[bs*i+jj+k];
700:         }
701:       }
702:       jj += bss[j];
703:     }
704: #if !defined(PETSC_USE_COMPLEX)
705:   } else if (addv == MAX_VALUES) {
706:     for (j=0; j<nv; j++) {
707:       for (k=0; k<bss[j]; k++) {
708:         for (i=0; i<n; i++) {
709:           y[j][i*bss[j] + k] = PetscMax(y[j][i*bss[j] + k],x[bs*i+jj+k]);
710:         }
711:       }
712:       jj += bss[j];
713:     }
714: #endif
715:   } else {
716:     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
717:   }

719:   VecRestoreArray(v,&x);
720:   for (i=0; i<nv; i++) {
721:     VecRestoreArray(s[i],&y[i]);
722:   }
723:   PetscFree2(y,bss);
724:   return(0);
725: }

729: /*@
730:    VecStrideScatterAll - Scatters all the single components from seperate vectors into 
731:      a multi-component vector.

733:    Collective on Vec

735:    Input Parameter:
736: +  s - the location where the subvectors are stored
737: -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES

739:    Output Parameter:
740: .  v - the multicomponent vector 

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

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

749:    Not optimized; could be easily

751:    Level: advanced

753:    Concepts:  scatter^into strided vector

755: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGather(),
756:           VecStrideScatterAll()
757: @*/
758: PetscErrorCode VecStrideScatterAll(Vec *s,Vec v,InsertMode addv)
759: {
761:   PetscInt        i,n,bs,j,jj,k,*bss,nv,nvc;
762:   PetscScalar     *x,**y;

768:   VecGetLocalSize(v,&n);
769:   VecGetArray(v,&x);
770:   bs   = v->bs;
771:   if (bs < 0) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Input vector does not have a valid blocksize set");

773:   PetscMalloc2(bs,PetscScalar**,&y,bs,PetscInt,&bss);
774:   nv  = 0;
775:   nvc = 0;
776:   for (i=0; i<bs; i++) {
777:     VecGetBlockSize(s[i],&bss[i]);
778:     if (bss[i] < 1) bss[i] = 1; /* if user never set it then assume 1  Re: [PETSC #8241] VecStrideGatherAll */
779:     VecGetArray(s[i],&y[i]);
780:     nvc  += bss[i];
781:     nv++;
782:     if (nvc > bs)  SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of subvectors in subvectors > number of vectors in main vector");
783:     if (nvc == bs) break;
784:   }

786:   n =  n/bs;

788:   jj = 0;
789:   if (addv == INSERT_VALUES) {
790:     for (j=0; j<nv; j++) {
791:       for (k=0; k<bss[j]; k++) {
792:         for (i=0; i<n; i++) {
793:           x[bs*i+jj+k] = y[j][i*bss[j] + k];
794:         }
795:       }
796:       jj += bss[j];
797:     }
798:   } else if (addv == ADD_VALUES) {
799:     for (j=0; j<nv; j++) {
800:       for (k=0; k<bss[j]; k++) {
801:         for (i=0; i<n; i++) {
802:           x[bs*i+jj+k] += y[j][i*bss[j] + k];
803:         }
804:       }
805:       jj += bss[j];
806:     }
807: #if !defined(PETSC_USE_COMPLEX)
808:   } else if (addv == MAX_VALUES) {
809:     for (j=0; j<nv; j++) {
810:       for (k=0; k<bss[j]; k++) {
811:         for (i=0; i<n; i++) {
812:           x[bs*i+jj+k] = PetscMax(x[bs*i+jj+k],y[j][i*bss[j] + k]);
813:         }
814:       }
815:       jj += bss[j];
816:     }
817: #endif
818:   } else {
819:     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
820:   }

822:   VecRestoreArray(v,&x);
823:   for (i=0; i<nv; i++) {
824:     VecRestoreArray(s[i],&y[i]);
825:   }
826:   PetscFree2(y,bss);
827:   return(0);
828: }

832: /*@
833:    VecStrideGather - Gathers a single component from a multi-component vector into
834:    another vector.

836:    Collective on Vec

838:    Input Parameter:
839: +  v - the vector 
840: .  start - starting point of the subvector (defined by a stride)
841: -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES

843:    Output Parameter:
844: .  s - the location where the subvector is stored

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

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

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

856:    Not optimized; could be easily

858:    Level: advanced

860:    Concepts: gather^into strided vector

862: .seealso: VecStrideNorm(), VecStrideScatter(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
863:           VecStrideScatterAll()
864: @*/
865: PetscErrorCode VecStrideGather(Vec v,PetscInt start,Vec s,InsertMode addv)
866: {
868:   PetscInt       i,n,bs,ns;
869:   PetscScalar    *x,*y;

874:   VecGetLocalSize(v,&n);
875:   VecGetLocalSize(s,&ns);
876:   VecGetArray(v,&x);
877:   VecGetArray(s,&y);

879:   bs   = v->bs;
880:   if (start >= bs) {
881:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
882:             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
883:   }
884:   if (n != ns*bs) {
885:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Subvector length * blocksize %D not correct for gather from original vector %D",ns*bs,n);
886:   }
887:   x += start;
888:   n =  n/bs;

890:   if (addv == INSERT_VALUES) {
891:     for (i=0; i<n; i++) {
892:       y[i] = x[bs*i];
893:     }
894:   } else if (addv == ADD_VALUES) {
895:     for (i=0; i<n; i++) {
896:       y[i] += x[bs*i];
897:     }
898: #if !defined(PETSC_USE_COMPLEX)
899:   } else if (addv == MAX_VALUES) {
900:     for (i=0; i<n; i++) {
901:       y[i] = PetscMax(y[i],x[bs*i]);
902:     }
903: #endif
904:   } else {
905:     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
906:   }

908:   VecRestoreArray(v,&x);
909:   VecRestoreArray(s,&y);
910:   return(0);
911: }

915: /*@
916:    VecStrideScatter - Scatters a single component from a vector into a multi-component vector.

918:    Collective on Vec

920:    Input Parameter:
921: +  s - the single-component vector 
922: .  start - starting point of the subvector (defined by a stride)
923: -  addv - one of ADD_VALUES,INSERT_VALUES,MAX_VALUES

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

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

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

935:    Not optimized; could be easily

937:    Level: advanced

939:    Concepts: scatter^into strided vector

941: .seealso: VecStrideNorm(), VecStrideGather(), VecStrideMin(), VecStrideMax(), VecStrideGatherAll(),
942:           VecStrideScatterAll()
943: @*/
944: PetscErrorCode VecStrideScatter(Vec s,PetscInt start,Vec v,InsertMode addv)
945: {
947:   PetscInt       i,n,bs,ns;
948:   PetscScalar    *x,*y;

953:   VecGetLocalSize(v,&n);
954:   VecGetLocalSize(s,&ns);
955:   VecGetArray(v,&x);
956:   VecGetArray(s,&y);

958:   bs   = v->bs;
959:   if (start >= bs) {
960:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Start of stride subvector (%D) is too large for stride\n\
961:             Have you set the vector blocksize (%D) correctly with VecSetBlockSize()?",start,bs);
962:   }
963:   if (n != ns*bs) {
964:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Subvector length * blocksize %D not correct for scatter to multicomponent vector %D",ns*bs,n);
965:   }
966:   x += start;
967:   n =  n/bs;


970:   if (addv == INSERT_VALUES) {
971:     for (i=0; i<n; i++) {
972:       x[bs*i] = y[i];
973:     }
974:   } else if (addv == ADD_VALUES) {
975:     for (i=0; i<n; i++) {
976:       x[bs*i] += y[i];
977:     }
978: #if !defined(PETSC_USE_COMPLEX)
979:   } else if (addv == MAX_VALUES) {
980:     for (i=0; i<n; i++) {
981:       x[bs*i] = PetscMax(y[i],x[bs*i]);
982:     }
983: #endif
984:   } else {
985:     SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown insert type");
986:   }


989:   VecRestoreArray(v,&x);
990:   VecRestoreArray(s,&y);
991:   return(0);
992: }

996: PetscErrorCode VecReciprocal_Default(Vec v)
997: {
999:   PetscInt         i,n;
1000:   PetscScalar *x;

1003:   VecGetLocalSize(v,&n);
1004:   VecGetArray(v,&x);
1005:   for (i=0; i<n; i++) {
1006:     if (x[i] != 0.0) x[i] = 1.0/x[i];
1007:   }
1008:   VecRestoreArray(v,&x);
1009:   return(0);
1010: }

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

1017:   Not collective

1019:   Input Parameter:
1020: . v - The vector

1022:   Output Parameter:
1023: . v - The vector square root

1025:   Level: beginner

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

1029: .keywords: vector, sqrt, square root
1030: @*/
1031: PetscErrorCode VecSqrt(Vec v)
1032: {
1033:   PetscScalar *x;
1034:   PetscInt         i, n;

1039:   VecGetLocalSize(v, &n);
1040:   VecGetArray(v, &x);
1041:   for(i = 0; i < n; i++) {
1042:     x[i] = sqrt(PetscAbsScalar(x[i]));
1043:   }
1044:   VecRestoreArray(v, &x);
1045:   return(0);
1046: }

1050: /*@
1051:    VecSum - Computes the sum of all the components of a vector.

1053:    Collective on Vec

1055:    Input Parameter:
1056: .  v - the vector 

1058:    Output Parameter:
1059: .  sum - the result

1061:    Level: beginner

1063:    Concepts: sum^of vector entries

1065: .seealso: VecNorm()
1066: @*/
1067: PetscErrorCode VecSum(Vec v,PetscScalar *sum)
1068: {
1070:   PetscInt         i,n;
1071:   PetscScalar *x,lsum = 0.0;

1076:   VecGetLocalSize(v,&n);
1077:   VecGetArray(v,&x);
1078:   for (i=0; i<n; i++) {
1079:     lsum += x[i];
1080:   }
1081:   MPI_Allreduce(&lsum,sum,1,MPIU_SCALAR,PetscSum_Op,v->comm);
1082:   VecRestoreArray(v,&x);
1083:   return(0);
1084: }

1088: /*@
1089:    VecShift - Shifts all of the components of a vector by computing
1090:    x[i] = x[i] + shift.

1092:    Collective on Vec

1094:    Input Parameters:
1095: +  v - the vector 
1096: -  shift - the shift

1098:    Output Parameter:
1099: .  v - the shifted vector 

1101:    Level: intermediate

1103:    Concepts: vector^adding constant

1105: @*/
1106: PetscErrorCode VecShift(const PetscScalar *shift,Vec v)
1107: {
1109:   PetscInt         i,n;
1110:   PetscScalar *x,lsum = *shift;

1115:   VecGetLocalSize(v,&n);
1116:   VecGetArray(v,&x);
1117:   for (i=0; i<n; i++) {
1118:     x[i] += lsum;
1119:   }
1120:   VecRestoreArray(v,&x);
1121:   return(0);
1122: }

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

1129:    Collective on Vec

1131:    Input Parameters:
1132: .  v - the vector 

1134:    Level: intermediate

1136:    Concepts: vector^absolute value

1138: @*/
1139: PetscErrorCode VecAbs(Vec v)
1140: {
1142:   PetscInt       i,n;
1143:   PetscScalar    *x;

1147:   VecGetLocalSize(v,&n);
1148:   VecGetArray(v,&x);
1149:   for (i=0; i<n; i++) {
1150:     x[i] = PetscAbsScalar(x[i]);
1151:   }
1152:   VecRestoreArray(v,&x);
1153:   return(0);
1154: }

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

1161:   Input Parameters:
1162: + vec   - The vector
1163: . order - The ordering
1164: - inv   - The flag for inverting the permutation

1166:   Level: beginner

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

1170: .seealso: MatPermute()
1171: .keywords: vec, permute
1172: @*/
1173: PetscErrorCode VecPermute(Vec x, IS row, PetscTruth inv)
1174: {
1175:   PetscScalar    *array, *newArray;
1176:   PetscInt       *idx;
1177:   PetscInt       i;

1181:   ISGetIndices(row, &idx);
1182:   VecGetArray(x, &array);
1183:   PetscMalloc(x->n*sizeof(PetscScalar), &newArray);
1184: #ifdef PETSC_USE_BOPT_g
1185:   for(i = 0; i < x->n; i++) {
1186:     if ((idx[i] < 0) || (idx[i] >= x->n)) {
1187:       SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Permutation index %D is out of bounds: %D", i, idx[i]);
1188:     }
1189:   }
1190: #endif
1191:   if (inv == PETSC_FALSE) {
1192:     for(i = 0; i < x->n; i++) newArray[i]      = array[idx[i]];
1193:   } else {
1194:     for(i = 0; i < x->n; i++) newArray[idx[i]] = array[i];
1195:   }
1196:   VecRestoreArray(x, &array);
1197:   ISRestoreIndices(row, &idx);
1198:   VecReplaceArray(x, newArray);
1199:   return(0);
1200: }

1204: /*@
1205:    VecEqual - Compares two vectors.

1207:    Collective on Vec

1209:    Input Parameters:
1210: +  vec1 - the first matrix
1211: -  vec2 - the second matrix

1213:    Output Parameter:
1214: .  flg - PETSC_TRUE if the vectors are equal; PETSC_FALSE otherwise.

1216:    Level: intermediate

1218:    Concepts: equal^two vectors
1219:    Concepts: vector^equality

1221: @*/
1222: PetscErrorCode VecEqual(Vec vec1,Vec vec2,PetscTruth *flg)
1223: {
1224:   PetscScalar    *v1,*v2;
1226:   PetscInt       n1,n2;
1227:   PetscTruth     flg1;

1230:   VecGetSize(vec1,&n1);
1231:   VecGetSize(vec2,&n2);
1232:   if (vec1 == vec2) {
1233:     flg1 = PETSC_TRUE;
1234:   } else if (n1 != n2) {
1235:     flg1 = PETSC_FALSE;
1236:   } else {
1237:     VecGetArray(vec1,&v1);
1238:     VecGetArray(vec2,&v2);
1239:     PetscMemcmp(v1,v2,n1*sizeof(PetscScalar),&flg1);
1240:     VecRestoreArray(vec1,&v1);
1241:     VecRestoreArray(vec2,&v2);
1242:   }

1244:   /* combine results from all processors */
1245:   MPI_Allreduce(&flg1,flg,1,MPI_INT,MPI_MIN,vec1->comm);
1246:   return(0);
1247: }