Actual source code: vector.c

  1: /*$Id: vector.c,v 1.238 2001/09/11 16:31:48 bsmith Exp $*/
  2: /*
  3:      Provides the interface functions for all vector operations.
  4:    These are the vector functions the user calls.
  5: */
 6:  #include src/vec/vecimpl.h

  8: /* Logging support */
  9: int VEC_COOKIE;
 10: int VEC_View, VEC_Max, VEC_Min, VEC_DotBarrier, VEC_Dot, VEC_MDotBarrier, VEC_MDot, VEC_TDot, VEC_MTDot, VEC_NormBarrier;
 11: int VEC_Norm, VEC_Scale, VEC_Copy, VEC_Set, VEC_AXPY, VEC_AYPX, VEC_WAXPY, VEC_MAXPY, VEC_Swap, VEC_AssemblyBegin;
 12: int VEC_AssemblyEnd, VEC_PointwiseMult, VEC_SetValues, VEC_Load, VEC_ScatterBarrier, VEC_ScatterBegin, VEC_ScatterEnd;
 13: int VEC_SetRandom, VEC_ReduceArithmetic, VEC_ReduceBarrier, VEC_ReduceCommunication;

 15: /*
 16:   VecSetTypeFromOptions_Private - Sets the type of vector from user options. Defaults to a PETSc sequential vector on one
 17:   processor and a PETSc MPI vector on more than one processor.

 19:   Collective on Vec

 21:   Input Parameter:
 22: . vec - The vector

 24:   Level: intermediate

 26: .keywords: Vec, set, options, database, type
 27: .seealso: VecSetFromOptions(), VecSetType()
 28: */
 29: static int VecSetTypeFromOptions_Private(Vec vec)
 30: {
 31:   PetscTruth opt;
 32:   char      *defaultType;
 33:   char       typeName[256];
 34:   int        numProcs;
 35:   int        ierr;

 38:   if (vec->type_name != PETSC_NULL) {
 39:     defaultType = vec->type_name;
 40:   } else {
 41:     MPI_Comm_size(vec->comm, &numProcs);
 42:     if (numProcs > 1) {
 43:       defaultType = VECMPI;
 44:     } else {
 45:       defaultType = VECSEQ;
 46:     }
 47:   }

 49:   if (!VecRegisterAllCalled) {
 50:     VecRegisterAll(PETSC_NULL);
 51:   }
 52:   PetscOptionsList("-vec_type", "Vector type"," VecSetType", VecList, defaultType, typeName, 256, &opt);
 53: 
 54:   if (opt == PETSC_TRUE) {
 55:     VecSetType(vec, typeName);
 56:   } else {
 57:     VecSetType(vec, defaultType);
 58:   }
 59:   return(0);
 60: }

 62: /*@C
 63:   VecSetFromOptions - Configures the vector from the options database.

 65:   Collective on Vec

 67:   Input Parameter:
 68: . vec - The vector

 70:   Notes:  To see all options, run your program with the -help option, or consult the users manual.
 71:           Must be called after VecCreate() but before the vector is used.

 73:   Level: intermediate

 75:   Concepts: vectors^setting options
 76:   Concepts: vectors^setting type

 78: .keywords: Vec, set, options, database
 79: .seealso: VecCreate(), VecPrintHelp(), VechSetOptionsPrefix()
 80: @*/
 81: int VecSetFromOptions(Vec vec)
 82: {
 83:   PetscTruth opt;
 84:   int        ierr;


 89:   PetscOptionsBegin(vec->comm, vec->prefix, "Vector options", "Vec");

 91:   /* Handle generic vector options */
 92:   PetscOptionsHasName(PETSC_NULL, "-help", &opt);
 93:   if (opt == PETSC_TRUE) {
 94:     VecPrintHelp(vec);
 95:   }

 97:   /* Handle vector type options */
 98:   VecSetTypeFromOptions_Private(vec);

100:   /* Handle specific vector options */
101:   if (vec->ops->setfromoptions != PETSC_NULL) {
102:     (*vec->ops->setfromoptions)(vec);
103:   }
104: #if defined(__cplusplus) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_MATSINGLE) && defined(PETSC_HAVE_CXX_NAMESPACE)
105:   VecESISetFromOptions(vec);
106: #endif

108:   PetscOptionsEnd();

110:   VecViewFromOptions(vec, vec->name);
111:   return(0);
112: }

114: /*@
115:   VecPrintHelp - Prints all options for the Vec.

117:   Input Parameter:
118: . vec - The vector

120:   Options Database Keys:
121: $  -help, -h

123:   Level: intermediate

125: .keywords: Vec, help
126: .seealso: VecSetFromOptions()
127: @*/
128: int VecPrintHelp(Vec vec)
129: {
130:   char p[64];
131:   int  ierr;


136:   PetscStrcpy(p, "-");
137:   if (vec->prefix != PETSC_NULL) {
138:     PetscStrcat(p, vec->prefix);
139:   }

141:   (*PetscHelpPrintf)(vec->comm, "Vec options ------------------------------------------------n");
142:   (*PetscHelpPrintf)(vec->comm,"   %svec_type <typename> : Sets the vector typen", p);
143:   return(0);
144: }

146: /*@
147:   VecSetSizes - Sets the local and global sizes, and checks to determine compatibility

149:   Collective on Vec

151:   Input Parameter:
152: + v - the vector
153: - bs - the blocksize

155:   Level: intermediate

157: .seealso: VecGetSize()
158: @*/
159: int VecSetSizes(Vec v, int n, int N)
160: {

165:   v->n = n;
166:   v->N = N;
167:   PetscSplitOwnership(v->comm, &v->n, &v->N);
168:   return(0);
169: }

171: /*@
172:    VecSetBlockSize - Sets the blocksize for future calls to VecSetValuesBlocked()
173:    and VecSetValuesBlockedLocal().

175:    Collective on Vec

177:    Input Parameter:
178: +  v - the vector
179: -  bs - the blocksize

181:    Notes:
182:    All vectors obtained by VecDuplicate() inherit the same blocksize.

184:    Level: advanced

186: .seealso: VecSetValuesBlocked(), VecSetLocalToGlobalMappingBlocked(), VecGetBlockSize()

188:   Concepts: block size^vectors
189: @*/
190: int VecSetBlockSize(Vec v,int bs)
191: {
194:   if (bs <= 0) bs = 1;
195:   if (bs == v->bs) return(0);
196:   if (v->bs != -1) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot reset blocksize. Current size %d new %d",v->bs,bs);
197:   if (v->N != -1 && v->N % bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Vector length not divisible by blocksize %d %d",v->N,bs);
198:   if (v->n != -1 && v->n % bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Local vector length not divisible by blocksize %d %dn
199:    Try setting blocksize before setting the vector type",v->n,bs);
200: 
201:   v->bs        = bs;
202:   v->bstash.bs = bs; /* use the same blocksize for the vec's block-stash */
203:   return(0);
204: }

206: /*@
207:    VecGetBlockSize - Gets the blocksize for the vector, i.e. what is used for VecSetValuesBlocked()
208:    and VecSetValuesBlockedLocal().

210:    Collective on Vec

212:    Input Parameter:
213: .  v - the vector

215:    Output Parameter:
216: .  bs - the blocksize

218:    Notes:
219:    All vectors obtained by VecDuplicate() inherit the same blocksize.

221:    Level: advanced

223: .seealso: VecSetValuesBlocked(), VecSetLocalToGlobalMappingBlocked(), VecSetBlockSize()

225:    Concepts: vector^block size
226:    Concepts: block^vector

228: @*/
229: int VecGetBlockSize(Vec v,int *bs)
230: {
233:   *bs = v->bs;
234:   return(0);
235: }

237: /*@
238:    VecValid - Checks whether a vector object is valid.

240:    Not Collective

242:    Input Parameter:
243: .  v - the object to check

245:    Output Parameter:
246:    flg - flag indicating vector status, either
247:    PETSC_TRUE if vector is valid, or PETSC_FALSE otherwise.

249:    Level: developer

251: @*/
252: int VecValid(Vec v,PetscTruth *flg)
253: {
256:   if (!v)                           *flg = PETSC_FALSE;
257:   else if (v->cookie != VEC_COOKIE) *flg = PETSC_FALSE;
258:   else                              *flg = PETSC_TRUE;
259:   return(0);
260: }

262: /*@
263:    VecDot - Computes the vector dot product.

265:    Collective on Vec

267:    Input Parameters:
268: .  x, y - the vectors

270:    Output Parameter:
271: .  alpha - the dot product

273:    Performance Issues:
274: +    per-processor memory bandwidth
275: .    interprocessor latency
276: -    work load inbalance that causes certain processes to arrive much earlier than
277:      others

279:    Notes for Users of Complex Numbers:
280:    For complex vectors, VecDot() computes 
281: $     val = (x,y) = y^H x,
282:    where y^H denotes the conjugate transpose of y.

284:    Use VecTDot() for the indefinite form
285: $     val = (x,y) = y^T x,
286:    where y^T denotes the transpose of y.

288:    Level: intermediate

290:    Concepts: inner product
291:    Concepts: vector^inner product

293: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd()
294: @*/
295: int VecDot(Vec x,Vec y,PetscScalar *val)
296: {

307:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
308:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

310:   PetscLogEventBarrierBegin(VEC_DotBarrier,x,y,0,0,x->comm);
311:   (*x->ops->dot)(x,y,val);
312:   PetscLogEventBarrierEnd(VEC_DotBarrier,x,y,0,0,x->comm);
313:   /*
314:      The next block is for incremental debugging
315:   */
316:   if (PetscCompare) {
317:     int flag;
318:     MPI_Comm_compare(PETSC_COMM_WORLD,x->comm,&flag);
319:     if (flag != MPI_UNEQUAL) {
320:       PetscCompareScalar(*val);
321:     }
322:   }
323:   return(0);
324: }

326: /*@
327:    VecNorm  - Computes the vector norm.

329:    Collective on Vec

331:    Input Parameters:
332: +  x - the vector
333: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
334:           NORM_1_AND_2, which computes both norms and stores them
335:           in a two element array.

337:    Output Parameter:
338: .  val - the norm 

340:    Notes:
341: $     NORM_1 denotes sum_i |x_i|
342: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
343: $     NORM_INFINITY denotes max_i |x_i|

345:    Level: intermediate

347:    Performance Issues:
348: +    per-processor memory bandwidth
349: .    interprocessor latency
350: -    work load inbalance that causes certain processes to arrive much earlier than
351:      others

353:    Compile Option:
354:    PETSC_HAVE_SLOW_NRM2 will cause a C (loop unrolled) version of the norm to be used, rather
355:  than the BLAS. This should probably only be used when one is using the FORTRAN BLAS routines 
356:  (as opposed to vendor provided) because the FORTRAN BLAS NRM2() routine is very slow. 

358:    Concepts: norm
359:    Concepts: vector^norm

361: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), 
362:           VecNormBegin(), VecNormEnd()

364: @*/
365: int VecNorm(Vec x,NormType type,PetscReal *val)
366: {

372:   PetscLogEventBarrierBegin(VEC_NormBarrier,x,0,0,0,x->comm);
373:   (*x->ops->norm)(x,type,val);
374:   PetscLogEventBarrierEnd(VEC_NormBarrier,x,0,0,0,x->comm);
375:   /*
376:      The next block is for incremental debugging
377:   */
378:   if (PetscCompare) {
379:     int flag;
380:     MPI_Comm_compare(PETSC_COMM_WORLD,x->comm,&flag);
381:     if (flag != MPI_UNEQUAL) {
382:       PetscCompareDouble(*val);
383:     }
384:   }
385:   return(0);
386: }

388: /*@C
389:    VecMax - Determines the maximum vector component and its location.

391:    Collective on Vec

393:    Input Parameter:
394: .  x - the vector

396:    Output Parameters:
397: +  val - the maximum component
398: -  p - the location of val

400:    Notes:
401:    Returns the value PETSC_MIN and p = -1 if the vector is of length 0.

403:    Level: intermediate

405:    Concepts: maximum^of vector
406:    Concepts: vector^maximum value

408: .seealso: VecNorm(), VecMin()
409: @*/
410: int VecMax(Vec x,int *p,PetscReal *val)
411: {

418:   PetscLogEventBegin(VEC_Max,x,0,0,0);
419:   (*x->ops->max)(x,p,val);
420:   PetscLogEventEnd(VEC_Max,x,0,0,0);
421:   return(0);
422: }

424: /*@
425:    VecMin - Determines the minimum vector component and its location.

427:    Collective on Vec

429:    Input Parameters:
430: .  x - the vector

432:    Output Parameter:
433: +  val - the minimum component
434: -  p - the location of val

436:    Level: intermediate

438:    Notes:
439:    Returns the value PETSC_MAX and p = -1 if the vector is of length 0.

441:    Concepts: minimum^of vector
442:    Concepts: vector^minimum entry

444: .seealso: VecMax()
445: @*/
446: int VecMin(Vec x,int *p,PetscReal *val)
447: {

454:   PetscLogEventBegin(VEC_Min,x,0,0,0);
455:   (*x->ops->min)(x,p,val);
456:   PetscLogEventEnd(VEC_Min,x,0,0,0);
457:   return(0);
458: }

460: /*@
461:    VecTDot - Computes an indefinite vector dot product. That is, this
462:    routine does NOT use the complex conjugate.

464:    Collective on Vec

466:    Input Parameters:
467: .  x, y - the vectors

469:    Output Parameter:
470: .  val - the dot product

472:    Notes for Users of Complex Numbers:
473:    For complex vectors, VecTDot() computes the indefinite form
474: $     val = (x,y) = y^T x,
475:    where y^T denotes the transpose of y.

477:    Use VecDot() for the inner product
478: $     val = (x,y) = y^H x,
479:    where y^H denotes the conjugate transpose of y.

481:    Level: intermediate

483:    Concepts: inner product^non-Hermitian
484:    Concepts: vector^inner product
485:    Concepts: non-Hermitian inner product

487: .seealso: VecDot(), VecMTDot()
488: @*/
489: int VecTDot(Vec x,Vec y,PetscScalar *val)
490: {

501:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
502:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

504:   PetscLogEventBegin(VEC_TDot,x,y,0,0);
505:   (*x->ops->tdot)(x,y,val);
506:   PetscLogEventEnd(VEC_TDot,x,y,0,0);
507:   return(0);
508: }

510: /*@
511:    VecScale - Scales a vector. 

513:    Collective on Vec

515:    Input Parameters:
516: +  x - the vector
517: -  alpha - the scalar

519:    Output Parameter:
520: .  x - the scaled vector

522:    Note:
523:    For a vector with n components, VecScale() computes 
524: $      x[i] = alpha * x[i], for i=1,...,n.

526:    Level: intermediate

528:    Concepts: vector^scaling
529:    Concepts: scaling^vector

531: @*/
532: int VecScale (const PetscScalar *alpha,Vec x)
533: {

540:   PetscLogEventBegin(VEC_Scale,x,0,0,0);
541:   (*x->ops->scale)(alpha,x);
542:   PetscLogEventEnd(VEC_Scale,x,0,0,0);
543:   return(0);
544: }

546: /*@
547:    VecCopy - Copies a vector. 

549:    Collective on Vec

551:    Input Parameter:
552: .  x - the vector

554:    Output Parameter:
555: .  y - the copy

557:    Notes:
558:    For default parallel PETSc vectors, both x and y must be distributed in
559:    the same manner; local copies are done.

561:    Level: beginner

563: .seealso: VecDuplicate()
564: @*/
565: int VecCopy(Vec x,Vec y)
566: {

575:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
576:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

578:   PetscLogEventBegin(VEC_Copy,x,y,0,0);
579:   (*x->ops->copy)(x,y);
580:   PetscLogEventEnd(VEC_Copy,x,y,0,0);
581:   return(0);
582: }

584: /*@
585:    VecSet - Sets all components of a vector to a single scalar value. 

587:    Collective on Vec

589:    Input Parameters:
590: +  alpha - the scalar
591: -  x  - the vector

593:    Output Parameter:
594: .  x  - the vector

596:    Note:
597:    For a vector of dimension n, VecSet() computes
598: $     x[i] = alpha, for i=1,...,n,
599:    so that all vector entries then equal the identical
600:    scalar value, alpha.  Use the more general routine
601:    VecSetValues() to set different vector entries.

603:    Level: beginner

605: .seealso VecSetValues(), VecSetValuesBlocked(), VecSetRandom()

607:    Concepts: vector^setting to constant

609: @*/
610: int VecSet(const PetscScalar *alpha,Vec x)
611: {


619:   PetscLogEventBegin(VEC_Set,x,0,0,0);
620:   (*x->ops->set)(alpha,x);
621:   PetscLogEventEnd(VEC_Set,x,0,0,0);
622:   return(0);
623: }

625: /*@C
626:    VecSetRandom - Sets all components of a vector to random numbers.

628:    Collective on Vec

630:    Input Parameters:
631: +  rctx - the random number context, formed by PetscRandomCreate(), or PETSC_NULL and
632:           it will create one internally.
633: -  x  - the vector

635:    Output Parameter:
636: .  x  - the vector

638:    Example of Usage:
639: .vb
640:      PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
641:      VecSetRandom(rctx,x);
642:      PetscRandomDestroy(rctx);
643: .ve

645:    Level: intermediate

647:    Concepts: vector^setting to random
648:    Concepts: random^vector

650: .seealso: VecSet(), VecSetValues(), PetscRandomCreate(), PetscRandomDestroy()
651: @*/
652: int VecSetRandom(PetscRandom rctx,Vec x)
653: {
654:   int         ierr;
655:   PetscRandom randObj = PETSC_NULL;


662:   if (!rctx) {
663:     MPI_Comm    comm;
664:     PetscObjectGetComm((PetscObject)x,&comm);
665:     PetscRandomCreate(comm,RANDOM_DEFAULT,&randObj);
666:     rctx = randObj;
667:   }

669:   PetscLogEventBegin(VEC_SetRandom,x,rctx,0,0);
670:   (*x->ops->setrandom)(rctx,x);
671:   PetscLogEventEnd(VEC_SetRandom,x,rctx,0,0);
672: 
673:   if (randObj) {
674:     PetscRandomDestroy(randObj);
675:   }
676:   return(0);
677: }

679: /*@
680:    VecAXPY - Computes y = alpha x + y. 

682:    Collective on Vec

684:    Input Parameters:
685: +  alpha - the scalar
686: -  x, y  - the vectors

688:    Output Parameter:
689: .  y - output vector

691:    Level: intermediate

693:    Concepts: vector^BLAS
694:    Concepts: BLAS

696: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
697: @*/
698: int VecAXPY(const PetscScalar *alpha,Vec x,Vec y)
699: {

710:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
711:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

713:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
714:   (*x->ops->axpy)(alpha,x,y);
715:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
716:   return(0);
717: }

719: /*@
720:    VecAXPBY - Computes y = alpha x + beta y. 

722:    Collective on Vec

724:    Input Parameters:
725: +  alpha,beta - the scalars
726: -  x, y  - the vectors

728:    Output Parameter:
729: .  y - output vector

731:    Level: intermediate

733:    Concepts: BLAS
734:    Concepts: vector^BLAS

736: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
737: @*/
738: int VecAXPBY(const PetscScalar *alpha,const PetscScalar *beta,Vec x,Vec y)
739: {

751:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
752:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

754:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
755:   (*x->ops->axpby)(alpha,beta,x,y);
756:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
757:   return(0);
758: }

760: /*@
761:    VecAYPX - Computes y = x + alpha y.

763:    Collective on Vec

765:    Input Parameters:
766: +  alpha - the scalar
767: -  x, y  - the vectors

769:    Output Parameter:
770: .  y - output vector

772:    Level: intermediate

774:    Concepts: vector^BLAS
775:    Concepts: BLAS

777: .seealso: VecAXPY(), VecWAXPY()
778: @*/
779: int VecAYPX(const PetscScalar *alpha,Vec x,Vec y)
780: {

791:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
792:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

794:   PetscLogEventBegin(VEC_AYPX,x,y,0,0);
795:    (*x->ops->aypx)(alpha,x,y);
796:   PetscLogEventEnd(VEC_AYPX,x,y,0,0);
797:   return(0);
798: }

800: /*@
801:    VecSwap - Swaps the vectors x and y.

803:    Collective on Vec

805:    Input Parameters:
806: .  x, y  - the vectors

808:    Level: advanced

810:    Concepts: vector^swapping values

812: @*/
813: int VecSwap(Vec x,Vec y)
814: {

824:   if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
825:   if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

827:   PetscLogEventBegin(VEC_Swap,x,y,0,0);
828:   (*x->ops->swap)(x,y);
829:   PetscLogEventEnd(VEC_Swap,x,y,0,0);
830:   return(0);
831: }

833: /*@
834:    VecWAXPY - Computes w = alpha x + y.

836:    Collective on Vec

838:    Input Parameters:
839: +  alpha - the scalar
840: -  x, y  - the vectors

842:    Output Parameter:
843: .  w - the result

845:    Level: intermediate

847:    Concepts: vector^BLAS
848:    Concepts: BLAS

850: .seealso: VecAXPY(), VecAYPX()
851: @*/
852: int VecWAXPY(const PetscScalar *alpha,Vec x,Vec y,Vec w)
853: {

866:   if (x->N != y->N || x->N != w->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
867:   if (x->n != y->n || x->n != w->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

869:   PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
870:    (*x->ops->waxpy)(alpha,x,y,w);
871:   PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
872:   return(0);
873: }

875: /*@
876:    VecPointwiseMult - Computes the componentwise multiplication w = x*y.

878:    Collective on Vec

880:    Input Parameters:
881: .  x, y  - the vectors

883:    Output Parameter:
884: .  w - the result

886:    Level: advanced

888:    Notes: any subset of the x, y, and w may be the same vector.

890:    Concepts: vector^pointwise multiply

892: .seealso: VecPointwiseDivide()
893: @*/
894: int VecPointwiseMult(Vec x,Vec y,Vec w)
895: {

907:   if (x->N != y->N || x->N != w->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
908:   if (x->n != y->n || x->n != w->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

910:   PetscLogEventBegin(VEC_PointwiseMult,x,y,w,0);
911:   (*x->ops->pointwisemult)(x,y,w);
912:   PetscLogEventEnd(VEC_PointwiseMult,x,y,w,0);
913:   return(0);
914: }

916: /*@
917:    VecPointwiseDivide - Computes the componentwise division w = x/y.

919:    Collective on Vec

921:    Input Parameters:
922: .  x, y  - the vectors

924:    Output Parameter:
925: .  w - the result

927:    Level: advanced

929:    Notes: any subset of the x, y, and w may be the same vector.

931:    Concepts: vector^pointwise divide

933: .seealso: VecPointwiseMult()
934: @*/
935: int VecPointwiseDivide(Vec x,Vec y,Vec w)
936: {

948:   if (x->N != y->N || x->N != w->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
949:   if (x->n != y->n || x->n != w->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

951:   (*x->ops->pointwisedivide)(x,y,w);
952:   return(0);
953: }

955: /*@C
956:    VecDuplicate - Creates a new vector of the same type as an existing vector.

958:    Collective on Vec

960:    Input Parameters:
961: .  v - a vector to mimic

963:    Output Parameter:
964: .  newv - location to put new vector

966:    Notes:
967:    VecDuplicate() does not copy the vector, but rather allocates storage
968:    for the new vector.  Use VecCopy() to copy a vector.

970:    Use VecDestroy() to free the space. Use VecDuplicateVecs() to get several
971:    vectors. 

973:    Level: beginner

975: .seealso: VecDestroy(), VecDuplicateVecs(), VecCreate(), VecCopy()
976: @*/
977: int VecDuplicate(Vec x,Vec *newv)
978: {

985:   (*x->ops->duplicate)(x,newv);
986:   return(0);
987: }

989: /*@C
990:    VecDestroy - Destroys a vector.

992:    Collective on Vec

994:    Input Parameters:
995: .  v  - the vector

997:    Level: beginner

999: .seealso: VecDuplicate()
1000: @*/
1001: int VecDestroy(Vec v)
1002: {

1007:   if (--v->refct > 0) return(0);
1008:   /* destroy the internal part */
1009:   if (v->ops->destroy) {
1010:     (*v->ops->destroy)(v);
1011:   }
1012:   /* destroy the external/common part */
1013:   if (v->mapping) {
1014:     ISLocalToGlobalMappingDestroy(v->mapping);
1015:   }
1016:   if (v->bmapping) {
1017:     ISLocalToGlobalMappingDestroy(v->bmapping);
1018:   }
1019:   if (v->map) {
1020:     PetscMapDestroy(v->map);
1021:   }
1022:   PetscLogObjectDestroy(v);
1023:   PetscHeaderDestroy(v);
1024:   return(0);
1025: }

1027: /*@C
1028:    VecDuplicateVecs - Creates several vectors of the same type as an existing vector.

1030:    Collective on Vec

1032:    Input Parameters:
1033: +  m - the number of vectors to obtain
1034: -  v - a vector to mimic

1036:    Output Parameter:
1037: .  V - location to put pointer to array of vectors

1039:    Notes:
1040:    Use VecDestroyVecs() to free the space. Use VecDuplicate() to form a single
1041:    vector.

1043:    Fortran Note:
1044:    The Fortran interface is slightly different from that given below, it 
1045:    requires one to pass in V a Vec (integer) array of size at least m.
1046:    See the Fortran chapter of the users manual and petsc/src/vec/examples for details.

1048:    Level: intermediate

1050: .seealso:  VecDestroyVecs(), VecDuplicate(), VecCreate(), VecDuplicateVecsF90()
1051: @*/
1052: int VecDuplicateVecs(Vec v,int m,Vec *V[])
1053: {

1060:   (*v->ops->duplicatevecs)(v, m,V);
1061:   return(0);
1062: }

1064: /*@C
1065:    VecDestroyVecs - Frees a block of vectors obtained with VecDuplicateVecs().

1067:    Collective on Vec

1069:    Input Parameters:
1070: +  vv - pointer to array of vector pointers
1071: -  m - the number of vectors previously obtained

1073:    Fortran Note:
1074:    The Fortran interface is slightly different from that given below.
1075:    See the Fortran chapter of the users manual and 
1076:    petsc/src/vec/examples for details.

1078:    Level: intermediate

1080: .seealso: VecDuplicateVecs(), VecDestroyVecsF90()
1081: @*/
1082: int VecDestroyVecs(const Vec vv[],int m)
1083: {

1087:   if (!vv) SETERRQ(PETSC_ERR_ARG_BADPTR,"Null vectors");
1090:   (*(*vv)->ops->destroyvecs)(vv,m);
1091:   return(0);
1092: }

1094: /*@
1095:    VecSetValues - Inserts or adds values into certain locations of a vector. 

1097:    Input Parameters:
1098:    Not Collective

1100: +  x - vector to insert in
1101: .  ni - number of elements to add
1102: .  ix - indices where to add
1103: .  y - array of values
1104: -  iora - either INSERT_VALUES or ADD_VALUES, where
1105:    ADD_VALUES adds values to any existing entries, and
1106:    INSERT_VALUES replaces existing entries with new values

1108:    Notes: 
1109:    VecSetValues() sets x[ix[i]] = y[i], for i=0,...,ni-1.

1111:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES 
1112:    options cannot be mixed without intervening calls to the assembly
1113:    routines.

1115:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
1116:    MUST be called after all calls to VecSetValues() have been completed.

1118:    VecSetValues() uses 0-based indices in Fortran as well as in C.

1120:    Negative indices may be passed in ix, these rows are 
1121:    simply ignored. This allows easily inserting element load matrices
1122:    with homogeneous Dirchlet boundary conditions that you don't want represented
1123:    in the vector.

1125:    Level: beginner

1127:    Concepts: vector^setting values

1129: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
1130:            VecSetValue(), VecSetValuesBlocked()
1131: @*/
1132: int VecSetValues(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1133: {

1141:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1142:   (*x->ops->setvalues)(x,ni,ix,y,iora);
1143:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1144:   return(0);
1145: }

1147: /*@
1148:    VecSetValuesBlocked - Inserts or adds blocks of values into certain locations of a vector. 

1150:    Not Collective

1152:    Input Parameters:
1153: +  x - vector to insert in
1154: .  ni - number of blocks to add
1155: .  ix - indices where to add in block count, rather than element count
1156: .  y - array of values
1157: -  iora - either INSERT_VALUES or ADD_VALUES, where
1158:    ADD_VALUES adds values to any existing entries, and
1159:    INSERT_VALUES replaces existing entries with new values

1161:    Notes: 
1162:    VecSetValuesBlocked() sets x[ix[bs*i]+j] = y[bs*i+j], 
1163:    for j=0,...,bs, for i=0,...,ni-1. where bs was set with VecSetBlockSize().

1165:    Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES 
1166:    options cannot be mixed without intervening calls to the assembly
1167:    routines.

1169:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
1170:    MUST be called after all calls to VecSetValuesBlocked() have been completed.

1172:    VecSetValuesBlocked() uses 0-based indices in Fortran as well as in C.

1174:    Negative indices may be passed in ix, these rows are 
1175:    simply ignored. This allows easily inserting element load matrices
1176:    with homogeneous Dirchlet boundary conditions that you don't want represented
1177:    in the vector.

1179:    Level: intermediate

1181:    Concepts: vector^setting values blocked

1183: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
1184:            VecSetValues()
1185: @*/
1186: int VecSetValuesBlocked(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1187: {

1195:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1196:   (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
1197:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1198:   return(0);
1199: }

1201: /*MC
1202:    VecSetValue - Set a single entry into a vector.

1204:    Synopsis:
1205:    int VecSetValue(Vec v,int row,PetscScalar value, InsertMode mode);

1207:    Not Collective

1209:    Input Parameters:
1210: +  v - the vector
1211: .  row - the row location of the entry
1212: .  value - the value to insert
1213: -  mode - either INSERT_VALUES or ADD_VALUES

1215:    Notes:
1216:    For efficiency one should use VecSetValues() and set several or 
1217:    many values simultaneously if possible.

1219:    Note that VecSetValue() does NOT return an error code (since this
1220:    is checked internally).

1222:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
1223:    MUST be called after all calls to VecSetValues() have been completed.

1225:    VecSetValues() uses 0-based indices in Fortran as well as in C.

1227:    Level: beginner

1229: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal()
1230: M*/

1232: /*@
1233:    VecSetLocalToGlobalMapping - Sets a local numbering to global numbering used
1234:    by the routine VecSetValuesLocal() to allow users to insert vector entries
1235:    using a local (per-processor) numbering.

1237:    Collective on Vec

1239:    Input Parameters:
1240: +  x - vector
1241: -  mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()

1243:    Notes: 
1244:    All vectors obtained with VecDuplicate() from this vector inherit the same mapping.

1246:    Level: intermediate

1248:    Concepts: vector^setting values with local numbering

1250: seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
1251:            VecSetLocalToGlobalMappingBlocked(), VecSetValuesBlockedLocal()
1252: @*/
1253: int VecSetLocalToGlobalMapping(Vec x,ISLocalToGlobalMapping mapping)
1254: {

1260:   if (x->mapping) {
1261:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for vector");
1262:   }

1264:   if (x->ops->setlocaltoglobalmapping) {
1265:     (*x->ops->setlocaltoglobalmapping)(x,mapping);
1266:   } else {
1267:     x->mapping = mapping;
1268:     PetscObjectReference((PetscObject)mapping);
1269:   }
1270:   return(0);
1271: }

1273: /*@
1274:    VecSetLocalToGlobalMappingBlock - Sets a local numbering to global numbering used
1275:    by the routine VecSetValuesBlockedLocal() to allow users to insert vector entries
1276:    using a local (per-processor) numbering.

1278:    Collective on Vec

1280:    Input Parameters:
1281: +  x - vector
1282: -  mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()

1284:    Notes: 
1285:    All vectors obtained with VecDuplicate() from this vector inherit the same mapping.

1287:    Level: intermediate

1289:    Concepts: vector^setting values blocked with local numbering

1291: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
1292:            VecSetLocalToGlobalMapping(), VecSetValuesBlockedLocal()
1293: @*/
1294: int VecSetLocalToGlobalMappingBlock(Vec x,ISLocalToGlobalMapping mapping)
1295: {

1301:   if (x->bmapping) {
1302:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for vector");
1303:   }
1304:   x->bmapping = mapping;
1305:   PetscObjectReference((PetscObject)mapping);
1306:   return(0);
1307: }

1309: /*@
1310:    VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
1311:    using a local ordering of the nodes. 

1313:    Not Collective

1315:    Input Parameters:
1316: +  x - vector to insert in
1317: .  ni - number of elements to add
1318: .  ix - indices where to add
1319: .  y - array of values
1320: -  iora - either INSERT_VALUES or ADD_VALUES, where
1321:    ADD_VALUES adds values to any existing entries, and
1322:    INSERT_VALUES replaces existing entries with new values

1324:    Level: intermediate

1326:    Notes: 
1327:    VecSetValuesLocal() sets x[ix[i]] = y[i], for i=0,...,ni-1.

1329:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES 
1330:    options cannot be mixed without intervening calls to the assembly
1331:    routines.

1333:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
1334:    MUST be called after all calls to VecSetValuesLocal() have been completed.

1336:    VecSetValuesLocal() uses 0-based indices in Fortran as well as in C.

1338:    Concepts: vector^setting values with local numbering

1340: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
1341:            VecSetValuesBlockedLocal()
1342: @*/
1343: int VecSetValuesLocal(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1344: {
1345:   int ierr,lixp[128],*lix = lixp;


1353:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1354:   if (!x->ops->setvalueslocal) {
1355:     if (!x->mapping) {
1356:       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
1357:     }
1358:     if (ni > 128) {
1359:       PetscMalloc(ni*sizeof(int),&lix);
1360:     }
1361:     ISLocalToGlobalMappingApply(x->mapping,ni,(int*)ix,lix);
1362:     (*x->ops->setvalues)(x,ni,lix,y,iora);
1363:     if (ni > 128) {
1364:       PetscFree(lix);
1365:     }
1366:   } else {
1367:     (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
1368:   }
1369:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1370:   return(0);
1371: }

1373: /*@
1374:    VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1375:    using a local ordering of the nodes. 

1377:    Not Collective

1379:    Input Parameters:
1380: +  x - vector to insert in
1381: .  ni - number of blocks to add
1382: .  ix - indices where to add in block count, not element count
1383: .  y - array of values
1384: -  iora - either INSERT_VALUES or ADD_VALUES, where
1385:    ADD_VALUES adds values to any existing entries, and
1386:    INSERT_VALUES replaces existing entries with new values

1388:    Level: intermediate

1390:    Notes: 
1391:    VecSetValuesBlockedLocal() sets x[bs*ix[i]+j] = y[bs*i+j], 
1392:    for j=0,..bs-1, for i=0,...,ni-1, where bs has been set with VecSetBlockSize().

1394:    Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES 
1395:    options cannot be mixed without intervening calls to the assembly
1396:    routines.

1398:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
1399:    MUST be called after all calls to VecSetValuesBlockedLocal() have been completed.

1401:    VecSetValuesBlockedLocal() uses 0-based indices in Fortran as well as in C.


1404:    Concepts: vector^setting values blocked with local numbering

1406: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(), 
1407:            VecSetLocalToGlobalMappingBlocked()
1408: @*/
1409: int VecSetValuesBlockedLocal(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1410: {
1411:   int ierr,lixp[128],*lix = lixp;

1418:   if (!x->bmapping) {
1419:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMappingBlocked()");
1420:   }
1421:   if (ni > 128) {
1422:     PetscMalloc(ni*sizeof(int),&lix);
1423:   }

1425:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1426:   ISLocalToGlobalMappingApply(x->bmapping,ni,(int*)ix,lix);
1427:   (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1428:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1429:   if (ni > 128) {
1430:     PetscFree(lix);
1431:   }
1432:   return(0);
1433: }

1435: /*@
1436:    VecAssemblyBegin - Begins assembling the vector.  This routine should
1437:    be called after completing all calls to VecSetValues().

1439:    Collective on Vec

1441:    Input Parameter:
1442: .  vec - the vector

1444:    Level: beginner

1446:    Concepts: assembly^vectors

1448: .seealso: VecAssemblyEnd(), VecSetValues()
1449: @*/
1450: int VecAssemblyBegin(Vec vec)
1451: {
1452:   int        ierr;
1453:   PetscTruth flg;


1459:   PetscOptionsHasName(vec->prefix,"-vec_view_stash",&flg);
1460:   if (flg) {
1461:     VecStashView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1462:   }

1464:   PetscLogEventBegin(VEC_AssemblyBegin,vec,0,0,0);
1465:   if (vec->ops->assemblybegin) {
1466:     (*vec->ops->assemblybegin)(vec);
1467:   }
1468:   PetscLogEventEnd(VEC_AssemblyBegin,vec,0,0,0);
1469:   return(0);
1470: }

1472: /*@
1473:    VecAssemblyEnd - Completes assembling the vector.  This routine should
1474:    be called after VecAssemblyBegin().

1476:    Collective on Vec

1478:    Input Parameter:
1479: .  vec - the vector

1481:    Options Database Keys:
1482: .  -vec_view - Prints vector in ASCII format
1483: .  -vec_view_matlab - Prints vector in Matlab format
1484: .  -vec_view_draw - Activates vector viewing using drawing tools
1485: .  -display <name> - Sets display name (default is host)
1486: .  -draw_pause <sec> - Sets number of seconds to pause after display
1487: .  -vec_view_socket - Activates vector viewing using a socket
1488: -  -vec_view_ams - Activates vector viewing using the ALICE Memory Snooper (AMS)
1489:  
1490:    Level: beginner

1492: .seealso: VecAssemblyBegin(), VecSetValues()
1493: @*/
1494: int VecAssemblyEnd(Vec vec)
1495: {
1496:   int        ierr;
1497:   PetscTruth flg;

1501:   PetscLogEventBegin(VEC_AssemblyEnd,vec,0,0,0);
1503:   if (vec->ops->assemblyend) {
1504:     (*vec->ops->assemblyend)(vec);
1505:   }
1506:   PetscLogEventEnd(VEC_AssemblyEnd,vec,0,0,0);
1507:   PetscOptionsHasName(vec->prefix,"-vec_view",&flg);
1508:   if (flg) {
1509:     VecView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1510:   }
1511:   PetscOptionsHasName(PETSC_NULL,"-vec_view_matlab",&flg);
1512:   if (flg) {
1513:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(vec->comm),PETSC_VIEWER_ASCII_MATLAB);
1514:     VecView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1515:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(vec->comm));
1516:   }
1517:   PetscOptionsHasName(PETSC_NULL,"-vec_view_draw",&flg);
1518:   if (flg) {
1519:     VecView(vec,PETSC_VIEWER_DRAW_(vec->comm));
1520:     PetscViewerFlush(PETSC_VIEWER_DRAW_(vec->comm));
1521:   }
1522:   PetscOptionsHasName(PETSC_NULL,"-vec_view_draw_lg",&flg);
1523:   if (flg) {
1524:     PetscViewerSetFormat(PETSC_VIEWER_DRAW_(vec->comm),PETSC_VIEWER_DRAW_LG);
1525:     VecView(vec,PETSC_VIEWER_DRAW_(vec->comm));
1526:     PetscViewerFlush(PETSC_VIEWER_DRAW_(vec->comm));
1527:   }
1528:   PetscOptionsHasName(PETSC_NULL,"-vec_view_socket",&flg);
1529:   if (flg) {
1530:     VecView(vec,PETSC_VIEWER_SOCKET_(vec->comm));
1531:     PetscViewerFlush(PETSC_VIEWER_SOCKET_(vec->comm));
1532:   }
1533:   PetscOptionsHasName(PETSC_NULL,"-vec_view_binary",&flg);
1534:   if (flg) {
1535:     VecView(vec,PETSC_VIEWER_BINARY_(vec->comm));
1536:     PetscViewerFlush(PETSC_VIEWER_BINARY_(vec->comm));
1537:   }
1538: #if defined(PETSC_HAVE_AMS)
1539:   PetscOptionsHasName(PETSC_NULL,"-vec_view_ams",&flg);
1540:   if (flg) {
1541:     VecView(vec,PETSC_VIEWER_AMS_(vec->comm));
1542:   }
1543: #endif
1544:   return(0);
1545: }


1548: /*@C
1549:    VecMTDot - Computes indefinite vector multiple dot products. 
1550:    That is, it does NOT use the complex conjugate.

1552:    Collective on Vec

1554:    Input Parameters:
1555: +  nv - number of vectors
1556: .  x - one vector
1557: -  y - array of vectors.  Note that vectors are pointers

1559:    Output Parameter:
1560: .  val - array of the dot products

1562:    Notes for Users of Complex Numbers:
1563:    For complex vectors, VecMTDot() computes the indefinite form
1564: $      val = (x,y) = y^T x,
1565:    where y^T denotes the transpose of y.

1567:    Use VecMDot() for the inner product
1568: $      val = (x,y) = y^H x,
1569:    where y^H denotes the conjugate transpose of y.

1571:    Level: intermediate

1573:    Concepts: inner product^multiple
1574:    Concepts: vector^multiple inner products

1576: .seealso: VecMDot(), VecTDot()
1577: @*/
1578: int VecMTDot(int nv,Vec x,const Vec y[],PetscScalar *val)
1579: {

1590:   if (x->N != (*y)->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1591:   if (x->n != (*y)->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

1593:   PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1594:   (*x->ops->mtdot)(nv,x,y,val);
1595:   PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1596:   return(0);
1597: }

1599: /*@C
1600:    VecMDot - Computes vector multiple dot products. 

1602:    Collective on Vec

1604:    Input Parameters:
1605: +  nv - number of vectors
1606: .  x - one vector
1607: -  y - array of vectors. 

1609:    Output Parameter:
1610: .  val - array of the dot products

1612:    Notes for Users of Complex Numbers:
1613:    For complex vectors, VecMDot() computes 
1614: $     val = (x,y) = y^H x,
1615:    where y^H denotes the conjugate transpose of y.

1617:    Use VecMTDot() for the indefinite form
1618: $     val = (x,y) = y^T x,
1619:    where y^T denotes the transpose of y.

1621:    Level: intermediate

1623:    Concepts: inner product^multiple
1624:    Concepts: vector^multiple inner products

1626: .seealso: VecMTDot(), VecDot()
1627: @*/
1628: int VecMDot(int nv,Vec x,const Vec y[],PetscScalar *val)
1629: {

1640:   if (x->N != (*y)->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1641:   if (x->n != (*y)->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

1643:   PetscLogEventBarrierBegin(VEC_MDotBarrier,x,*y,0,0,x->comm);
1644:   (*x->ops->mdot)(nv,x,y,val);
1645:   PetscLogEventBarrierEnd(VEC_MDotBarrier,x,*y,0,0,x->comm);
1646:   return(0);
1647: }

1649: /*@C
1650:    VecMAXPY - Computes y = y + sum alpha[j] x[j]

1652:    Collective on Vec

1654:    Input Parameters:
1655: +  nv - number of scalars and x-vectors
1656: .  alpha - array of scalars
1657: .  y - one vector
1658: -  x - array of vectors

1660:    Level: intermediate

1662:    Concepts: BLAS

1664: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1665: @*/
1666: int  VecMAXPY(int nv,const PetscScalar *alpha,Vec y,Vec *x)
1667: {

1678:   if (y->N != (*x)->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1679:   if (y->n != (*x)->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");

1681:   PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1682:   (*y->ops->maxpy)(nv,alpha,y,x);
1683:   PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1684:   return(0);
1685: }

1687: /*@C
1688:    VecGetArray - Returns a pointer to a contiguous array that contains this 
1689:    processor's portion of the vector data. For the standard PETSc
1690:    vectors, VecGetArray() returns a pointer to the local data array and
1691:    does not use any copies. If the underlying vector data is not stored
1692:    in a contiquous array this routine will copy the data to a contiquous
1693:    array and return a pointer to that. You MUST call VecRestoreArray() 
1694:    when you no longer need access to the array.

1696:    Not Collective

1698:    Input Parameter:
1699: .  x - the vector

1701:    Output Parameter:
1702: .  a - location to put pointer to the array

1704:    Fortran Note:
1705:    This routine is used differently from Fortran 77
1706: $    Vec         x
1707: $    PetscScalar x_array(1)
1708: $    PetscOffset i_x
1709: $    int         ierr
1710: $       call VecGetArray(x,x_array,i_x,ierr)
1711: $
1712: $   Access first local entry in vector with
1713: $      value = x_array(i_x + 1)
1714: $
1715: $      ...... other code
1716: $       call VecRestoreArray(x,x_array,i_x,ierr)
1717:    For Fortran 90 see VecGetArrayF90()

1719:    See the Fortran chapter of the users manual and 
1720:    petsc/src/snes/examples/tutorials/ex5f.F for details.

1722:    Level: beginner

1724:    Concepts: vector^accessing local values

1726: .seealso: VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(), VecGetArray2d()
1727: @*/
1728: int VecGetArray(Vec x,PetscScalar *a[])
1729: {

1736:   (*x->ops->getarray)(x,a);
1737:   return(0);
1738: }


1741: /*@C
1742:    VecGetArrays - Returns a pointer to the arrays in a set of vectors
1743:    that were created by a call to VecDuplicateVecs().  You MUST call
1744:    VecRestoreArrays() when you no longer need access to the array.

1746:    Not Collective

1748:    Input Parameter:
1749: +  x - the vectors
1750: -  n - the number of vectors

1752:    Output Parameter:
1753: .  a - location to put pointer to the array

1755:    Fortran Note:
1756:    This routine is not supported in Fortran.

1758:    Level: intermediate

1760: .seealso: VecGetArray(), VecRestoreArrays()
1761: @*/
1762: int VecGetArrays(const Vec x[],int n,PetscScalar **a[])
1763: {
1764:   int         i,ierr;
1765:   PetscScalar **q;

1770:   if (n <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %d",n);
1771:   PetscMalloc(n*sizeof(PetscScalar*),&q);
1772:   for (i=0; i<n; ++i) {
1773:     VecGetArray(x[i],&q[i]);
1774:   }
1775:   *a = q;
1776:   return(0);
1777: }

1779: /*@C
1780:    VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1781:    has been called.

1783:    Not Collective

1785:    Input Parameters:
1786: +  x - the vector
1787: .  n - the number of vectors
1788: -  a - location of pointer to arrays obtained from VecGetArrays()

1790:    Notes:
1791:    For regular PETSc vectors this routine does not involve any copies. For
1792:    any special vectors that do not store local vector data in a contiguous
1793:    array, this routine will copy the data back into the underlying 
1794:    vector data structure from the arrays obtained with VecGetArrays().

1796:    Fortran Note:
1797:    This routine is not supported in Fortran.

1799:    Level: intermediate

1801: .seealso: VecGetArrays(), VecRestoreArray()
1802: @*/
1803: int VecRestoreArrays(const Vec x[],int n,PetscScalar **a[])
1804: {
1805:   int         i,ierr;
1806:   PetscScalar **q = *a;


1812:   for(i=0;i<n;++i) {
1813:     VecRestoreArray(x[i],&q[i]);
1814:   }
1815:   PetscFree(q);
1816:   return(0);
1817: }

1819: /*@C
1820:    VecRestoreArray - Restores a vector after VecGetArray() has been called.

1822:    Not Collective

1824:    Input Parameters:
1825: +  x - the vector
1826: -  a - location of pointer to array obtained from VecGetArray()

1828:    Level: beginner

1830:    Notes:
1831:    For regular PETSc vectors this routine does not involve any copies. For
1832:    any special vectors that do not store local vector data in a contiguous
1833:    array, this routine will copy the data back into the underlying 
1834:    vector data structure from the array obtained with VecGetArray().

1836:    This routine actually zeros out the a pointer. This is to prevent accidental
1837:    us of the array after it has been restored. If you pass null for a it will 
1838:    not zero the array pointer a.

1840:    Fortran Note:
1841:    This routine is used differently from Fortran 77
1842: $    Vec         x
1843: $    PetscScalar x_array(1)
1844: $    PetscOffset i_x
1845: $    int         ierr
1846: $       call VecGetArray(x,x_array,i_x,ierr)
1847: $
1848: $   Access first local entry in vector with
1849: $      value = x_array(i_x + 1)
1850: $
1851: $      ...... other code
1852: $       call VecRestoreArray(x,x_array,i_x,ierr)

1854:    See the Fortran chapter of the users manual and 
1855:    petsc/src/snes/examples/tutorials/ex5f.F for details.
1856:    For Fortran 90 see VecRestoreArrayF90()

1858: .seealso: VecGetArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(), VecRestoreArray2d()
1859: @*/
1860: int VecRestoreArray(Vec x,PetscScalar *a[])
1861: {

1868:   if (x->ops->restorearray) {
1869:     (*x->ops->restorearray)(x,a);
1870:   }
1871:   return(0);
1872: }

1874: /*@
1875:   VecViewFromOptions - This function visualizes the vector based upon user options.

1877:   Collective on Vec

1879:   Input Parameters:
1880: . vec   - The vector
1881: . title - The title

1883:   Level: intermediate

1885: .keywords: Vec, view, options, database
1886: .seealso: VecSetFromOptions(), VecView()
1887: @*/
1888: int VecViewFromOptions(Vec vec, char *title)
1889: {
1890:   PetscViewer viewer;
1891:   PetscDraw   draw;
1892:   PetscTruth  opt;
1893:   char       *titleStr;
1894:   char        typeName[1024];
1895:   char        fileName[1024];
1896:   int         len;
1897:   int         ierr;

1900:   PetscOptionsHasName(vec->prefix, "-vec_view", &opt);
1901:   if (opt == PETSC_TRUE) {
1902:     PetscOptionsGetString(vec->prefix, "-vec_view", typeName, 1024, &opt);
1903:     PetscStrlen(typeName, &len);
1904:     if (len > 0) {
1905:       PetscViewerCreate(vec->comm, &viewer);
1906:       PetscViewerSetType(viewer, typeName);
1907:       PetscOptionsGetString(vec->prefix, "-vec_view_file", fileName, 1024, &opt);
1908:       if (opt == PETSC_TRUE) {
1909:         PetscViewerSetFilename(viewer, fileName);
1910:       } else {
1911:         PetscViewerSetFilename(viewer, vec->name);
1912:       }
1913:       VecView(vec, viewer);
1914:       PetscViewerFlush(viewer);
1915:       PetscViewerDestroy(viewer);
1916:     } else {
1917:       VecView(vec, PETSC_VIEWER_STDOUT_WORLD);
1918:     }
1919:   }
1920:   PetscOptionsHasName(vec->prefix, "-vec_view_draw", &opt);
1921:   if (opt == PETSC_TRUE) {
1922:     PetscViewerDrawOpen(vec->comm, 0, 0, 0, 0, 300, 300, &viewer);
1923:     PetscViewerDrawGetDraw(viewer, 0, &draw);
1924:     if (title != PETSC_NULL) {
1925:       titleStr = title;
1926:     } else {
1927:       PetscObjectName((PetscObject) vec);                                                          CHKERRQ(ierr) ;
1928:       titleStr = vec->name;
1929:     }
1930:     PetscDrawSetTitle(draw, titleStr);
1931:     VecView(vec, viewer);
1932:     PetscViewerFlush(viewer);
1933:     PetscDrawPause(draw);
1934:     PetscViewerDestroy(viewer);
1935:   }
1936:   return(0);
1937: }

1939: /*@C
1940:    VecView - Views a vector object. 

1942:    Collective on Vec

1944:    Input Parameters:
1945: +  v - the vector
1946: -  viewer - an optional visualization context

1948:    Notes:
1949:    The available visualization contexts include
1950: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1951: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1952:          output where only the first processor opens
1953:          the file.  All other processors send their 
1954:          data to the first processor to print. 

1956:    You can change the format the vector is printed using the 
1957:    option PetscViewerSetFormat().

1959:    The user can open alternative visualization contexts with
1960: +    PetscViewerASCIIOpen() - Outputs vector to a specified file
1961: .    PetscViewerBinaryOpen() - Outputs vector in binary to a
1962:          specified file; corresponding input uses VecLoad()
1963: .    PetscViewerDrawOpen() - Outputs vector to an X window display
1964: -    PetscViewerSocketOpen() - Outputs vector to Socket viewer

1966:    The user can call PetscViewerSetFormat() to specify the output
1967:    format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
1968:    PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen).  Available formats include
1969: +    PETSC_VIEWER_ASCII_DEFAULT - default, prints vector contents
1970: .    PETSC_VIEWER_ASCII_MATLAB - prints vector contents in Matlab format
1971: .    PETSC_VIEWER_ASCII_INDEX - prints vector contents, including indices of vector elements
1972: -    PETSC_VIEWER_ASCII_COMMON - prints vector contents, using a 
1973:          format common among all vector types

1975:    Level: beginner

1977:    Concepts: vector^printing
1978:    Concepts: vector^saving to disk

1980: .seealso: PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscDrawLGCreate(),
1981:           PetscViewerSocketOpen(), PetscViewerBinaryOpen(), VecLoad(), PetscViewerCreate(),
1982:           PetscRealView(), PetscScalarView(), PetscIntView()
1983: @*/
1984: int VecView(Vec vec,PetscViewer viewer)
1985: {
1986:   int               ierr;
1987:   PetscViewerFormat format;

1992:   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(vec->comm);
1995:   if (vec->stash.n || vec->bstash.n) SETERRQ(1,"Must call VecAssemblyBegin/End() before viewing this vector");

1997:   /*
1998:      Check if default viewer has been overridden, but user request it anyways
1999:   */
2000:   PetscViewerGetFormat(viewer,&format);
2001:   if (vec->ops->viewnative && format == PETSC_VIEWER_NATIVE) {
2002:     ierr   = PetscViewerPopFormat(viewer);
2003:     (*vec->ops->viewnative)(vec,viewer);
2004:     ierr   = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);
2005:   } else {
2006:     (*vec->ops->view)(vec,viewer);
2007:   }
2008:   return(0);
2009: }

2011: /*@
2012:    VecGetSize - Returns the global number of elements of the vector.

2014:    Not Collective

2016:    Input Parameter:
2017: .  x - the vector

2019:    Output Parameters:
2020: .  size - the global length of the vector

2022:    Level: beginner

2024:    Concepts: vector^local size

2026: .seealso: VecGetLocalSize()
2027: @*/
2028: int VecGetSize(Vec x,int *size)
2029: {

2036:   (*x->ops->getsize)(x,size);
2037:   return(0);
2038: }

2040: /*@
2041:    VecGetLocalSize - Returns the number of elements of the vector stored 
2042:    in local memory. This routine may be implementation dependent, so use 
2043:    with care.

2045:    Not Collective

2047:    Input Parameter:
2048: .  x - the vector

2050:    Output Parameter:
2051: .  size - the length of the local piece of the vector

2053:    Level: beginner

2055:    Concepts: vector^size

2057: .seealso: VecGetSize()
2058: @*/
2059: int VecGetLocalSize(Vec x,int *size)
2060: {

2067:   (*x->ops->getlocalsize)(x,size);
2068:   return(0);
2069: }

2071: /*@
2072:    VecGetOwnershipRange - Returns the range of indices owned by 
2073:    this processor, assuming that the vectors are laid out with the
2074:    first n1 elements on the first processor, next n2 elements on the
2075:    second, etc.  For certain parallel layouts this range may not be 
2076:    well defined. 

2078:    Not Collective

2080:    Input Parameter:
2081: .  x - the vector

2083:    Output Parameters:
2084: +  low - the first local element, pass in PETSC_NULL if not interested
2085: -  high - one more than the last local element, pass in PETSC_NULL if not interested

2087:    Note:
2088:    The high argument is one more than the last element stored locally.

2090:    Level: beginner

2092:    Concepts: ownership^of vectors
2093:    Concepts: vector^ownership of elements

2095: @*/
2096: int VecGetOwnershipRange(Vec x,int *low,int *high)
2097: {
2098:   int      ierr;

2105:   PetscMapGetLocalRange(x->map,low,high);
2106:   return(0);
2107: }

2109: /*@C
2110:    VecGetPetscMap - Returns the map associated with the vector

2112:    Not Collective

2114:    Input Parameter:
2115: .  x - the vector

2117:    Output Parameters:
2118: .  map - the map

2120:    Level: developer

2122: @*/
2123: int VecGetPetscMap(Vec x,PetscMap *map)
2124: {
2128:   *map = x->map;
2129:   return(0);
2130: }

2132: /*@
2133:    VecSetOption - Sets an option for controling a vector's behavior.

2135:    Collective on Vec

2137:    Input Parameter:
2138: +  x - the vector
2139: -  op - the option

2141:    Supported Options:
2142: +     VEC_IGNORE_OFF_PROC_ENTRIES, which causes VecSetValues() to ignore 
2143:       entries destined to be stored on a seperate processor. This can be used
2144:       to eliminate the global reduction in the VecAssemblyXXXX() if you know 
2145:       that you have only used VecSetValues() to set local elements
2146: -   VEC_TREAT_OFF_PROC_ENTRIES restores the treatment of off processor entries.

2148:    Level: intermediate

2150: @*/
2151: int VecSetOption(Vec x,VecOption op)
2152: {

2158:   if (x->ops->setoption) {
2159:     (*x->ops->setoption)(x,op);
2160:   }
2161:   return(0);
2162: }

2164: /* Default routines for obtaining and releasing; */
2165: /* may be used by any implementation */
2166: int VecDuplicateVecs_Default(Vec w,int m,Vec *V[])
2167: {
2168:   int  i,ierr;

2173:   if (m <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %d",m);
2174:   PetscMalloc(m*sizeof(Vec*),V);
2175:   for (i=0; i<m; i++) {VecDuplicate(w,*V+i);}
2176:   return(0);
2177: }

2179: int VecDestroyVecs_Default(const Vec v[], int m)
2180: {
2181:   int i,ierr;

2185:   if (m <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %d",m);
2186:   for (i=0; i<m; i++) {VecDestroy(v[i]);}
2187:   PetscFree((Vec*)v);
2188:   return(0);
2189: }

2191: /*@
2192:    VecPlaceArray - Allows one to replace the array in a vector with an
2193:    array provided by the user. This is useful to avoid copying an array
2194:    into a vector.

2196:    Not Collective

2198:    Input Parameters:
2199: +  vec - the vector
2200: -  array - the array

2202:    Notes:
2203:    You can return to the original array with a call to VecResetArray()

2205:    Level: developer

2207: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()

2209: @*/
2210: int VecPlaceArray(Vec vec,const PetscScalar array[])
2211: {

2217:   if (vec->ops->placearray) {
2218:     (*vec->ops->placearray)(vec,array);
2219:   } else {
2220:     SETERRQ(1,"Cannot place array in this type of vector");
2221:   }
2222:   return(0);
2223: }

2225: /*@
2226:    VecResetArray - Resets a vector to use its default memory. Call this 
2227:    after the use of VecPlaceArray().

2229:    Not Collective

2231:    Input Parameters:
2232: .  vec - the vector

2234:    Level: developer

2236: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecPlaceArray()

2238: @*/
2239: int VecResetArray(Vec vec)
2240: {

2246:   if (vec->ops->resetarray) {
2247:     (*vec->ops->resetarray)(vec);
2248:   } else {
2249:     SETERRQ(1,"Cannot reset array in this type of vector");
2250:   }
2251:   return(0);
2252: }

2254: /*@C
2255:    VecReplaceArray - Allows one to replace the array in a vector with an
2256:    array provided by the user. This is useful to avoid copying an array
2257:    into a vector.

2259:    Not Collective

2261:    Input Parameters:
2262: +  vec - the vector
2263: -  array - the array

2265:    Notes:
2266:    This permanently replaces the array and frees the memory associated
2267:    with the old array.

2269:    The memory passed in MUST be obtained with PetscMalloc() and CANNOT be
2270:    freed by the user. It will be freed when the vector is destroy. 

2272:    Not supported from Fortran

2274:    Level: developer

2276: .seealso: VecGetArray(), VecRestoreArray(), VecPlaceArray(), VecResetArray()

2278: @*/
2279: int VecReplaceArray(Vec vec,const PetscScalar array[])
2280: {

2286:   if (vec->ops->replacearray) {
2287:     (*vec->ops->replacearray)(vec,array);
2288:   } else {
2289:     SETERRQ(1,"Cannot replace array in this type of vector");
2290:   }
2291:   return(0);
2292: }

2294: /*MC
2295:     VecDuplicateVecsF90 - Creates several vectors of the same type as an existing vector
2296:     and makes them accessible via a Fortran90 pointer.

2298:     Synopsis:
2299:     VecDuplicateVecsF90(Vec x,int n,{Vec, pointer :: y(:)},integer ierr)

2301:     Collective on Vec

2303:     Input Parameters:
2304: +   x - a vector to mimic
2305: -   n - the number of vectors to obtain

2307:     Output Parameters:
2308: +   y - Fortran90 pointer to the array of vectors
2309: -   ierr - error code

2311:     Example of Usage: 
2312: .vb
2313:     Vec x
2314:     Vec, pointer :: y(:)
2315:     ....
2316:     call VecDuplicateVecsF90(x,2,y,ierr)
2317:     call VecSet(alpha,y(2),ierr)
2318:     call VecSet(alpha,y(2),ierr)
2319:     ....
2320:     call VecDestroyVecsF90(y,2,ierr)
2321: .ve

2323:     Notes:
2324:     Not yet supported for all F90 compilers

2326:     Use VecDestroyVecsF90() to free the space.

2328:     Level: beginner

2330: .seealso:  VecDestroyVecsF90(), VecDuplicateVecs()

2332: M*/

2334: /*MC
2335:     VecRestoreArrayF90 - Restores a vector to a usable state after a call to
2336:     VecGetArrayF90().

2338:     Synopsis:
2339:     VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2341:     Not collective

2343:     Input Parameters:
2344: +   x - vector
2345: -   xx_v - the Fortran90 pointer to the array

2347:     Output Parameter:
2348: .   ierr - error code

2350:     Example of Usage: 
2351: .vb
2352:     PetscScalar, pointer :: xx_v(:)
2353:     ....
2354:     call VecGetArrayF90(x,xx_v,ierr)
2355:     a = xx_v(3)
2356:     call VecRestoreArrayF90(x,xx_v,ierr)
2357: .ve
2358:    
2359:     Notes:
2360:     Not yet supported for all F90 compilers

2362:     Level: beginner

2364: .seealso:  VecGetArrayF90(), VecGetArray(), VecRestoreArray()

2366: M*/

2368: /*MC
2369:     VecDestroyVecsF90 - Frees a block of vectors obtained with VecDuplicateVecsF90().

2371:     Synopsis:
2372:     VecDestroyVecsF90({Vec, pointer :: x(:)},integer n,integer ierr)

2374:     Input Parameters:
2375: +   x - pointer to array of vector pointers
2376: -   n - the number of vectors previously obtained

2378:     Output Parameter:
2379: .   ierr - error code

2381:     Notes:
2382:     Not yet supported for all F90 compilers

2384:     Level: beginner

2386: .seealso:  VecDestroyVecs(), VecDuplicateVecsF90()

2388: M*/

2390: /*MC
2391:     VecGetArrayF90 - Accesses a vector array from Fortran90. For default PETSc
2392:     vectors, VecGetArrayF90() returns a pointer to the local data array. Otherwise,
2393:     this routine is implementation dependent. You MUST call VecRestoreArrayF90() 
2394:     when you no longer need access to the array.

2396:     Synopsis:
2397:     VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2399:     Not Collective 

2401:     Input Parameter:
2402: .   x - vector

2404:     Output Parameters:
2405: +   xx_v - the Fortran90 pointer to the array
2406: -   ierr - error code

2408:     Example of Usage: 
2409: .vb
2410:     PetscScalar, pointer :: xx_v(:)
2411:     ....
2412:     call VecGetArrayF90(x,xx_v,ierr)
2413:     a = xx_v(3)
2414:     call VecRestoreArrayF90(x,xx_v,ierr)
2415: .ve

2417:     Notes:
2418:     Not yet supported for all F90 compilers

2420:     Level: beginner

2422: .seealso:  VecRestoreArrayF90(), VecGetArray(), VecRestoreArray()

2424: M*/

2426: /*@C 
2427:   VecLoadIntoVector - Loads a vector that has been stored in binary format
2428:   with VecView().

2430:   Collective on PetscViewer 

2432:   Input Parameters:
2433: + viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
2434: - vec - vector to contain files values (must be of correct length)

2436:   Level: intermediate

2438:   Notes:
2439:   The input file must contain the full global vector, as
2440:   written by the routine VecView().

2442:   Use VecLoad() to create the vector as the values are read in

2444:   Notes for advanced users:
2445:   Most users should not need to know the details of the binary storage
2446:   format, since VecLoad() and VecView() completely hide these details.
2447:   But for anyone who's interested, the standard binary matrix storage
2448:   format is
2449: .vb
2450:      int    VEC_FILE_COOKIE
2451:      int    number of rows
2452:      PetscScalar *values of all nonzeros
2453: .ve

2455:    Note for Cray users, the int's stored in the binary file are 32 bit
2456: integers; not 64 as they are represented in the memory, so if you
2457: write your own routines to read/write these binary files from the Cray
2458: you need to adjust the integer sizes that you read in, see
2459: PetscReadBinary() and PetscWriteBinary() to see how this may be
2460: done.

2462:    In addition, PETSc automatically does the byte swapping for
2463: machines that store the bytes reversed, e.g.  DEC alpha, freebsd,
2464: linux, nt and the paragon; thus if you write your own binary
2465: read/write routines you have to swap the bytes; see PetscReadBinary()
2466: and PetscWriteBinary() to see how this may be done.

2468:    Concepts: vector^loading from file

2470: .seealso: PetscViewerBinaryOpen(), VecView(), MatLoad(), VecLoad() 
2471: @*/
2472: int VecLoadIntoVector(PetscViewer viewer,Vec vec)
2473: {

2480:   if (!vec->ops->loadintovector) {
2481:     SETERRQ(1,"Vector does not support load");
2482:   }
2483:   (*vec->ops->loadintovector)(viewer,vec);
2484:   return(0);
2485: }

2487: /*@
2488:    VecReciprocal - Replaces each component of a vector by its reciprocal.

2490:    Collective on Vec

2492:    Input Parameter:
2493: .  v - the vector 

2495:    Output Parameter:
2496: .  v - the vector reciprocal

2498:    Level: intermediate

2500:    Concepts: vector^reciprocal

2502: @*/
2503: int VecReciprocal(Vec vec)
2504: {
2505:   int    ierr;

2510:   if (!vec->ops->reciprocal) {
2511:     SETERRQ(1,"Vector does not support reciprocal operation");
2512:   }
2513:   (*vec->ops->reciprocal)(vec);
2514:   return(0);
2515: }

2517: int VecSetOperation(Vec vec,VecOperation op, void (*f)(void))
2518: {
2521:   /* save the native version of the viewer */
2522:   if (op == VECOP_VIEW && !vec->ops->viewnative) {
2523:     vec->ops->viewnative = vec->ops->view;
2524:   }
2525:   (((void(**)(void))vec->ops)[(int)op]) = f;
2526:   return(0);
2527: }

2529: /*@
2530:    VecSetStashInitialSize - sets the sizes of the vec-stash, that is
2531:    used during the assembly process to store values that belong to 
2532:    other processors.

2534:    Collective on Vec

2536:    Input Parameters:
2537: +  vec   - the vector
2538: .  size  - the initial size of the stash.
2539: -  bsize - the initial size of the block-stash(if used).

2541:    Options Database Keys:
2542: +   -vecstash_initial_size <size> or <size0,size1,...sizep-1>
2543: -   -vecstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1>

2545:    Level: intermediate

2547:    Notes: 
2548:      The block-stash is used for values set with VecSetValuesBlocked() while
2549:      the stash is used for values set with VecSetValues()

2551:      Run with the option -log_info and look for output of the form
2552:      VecAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
2553:      to determine the appropriate value, MM, to use for size and 
2554:      VecAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
2555:      to determine the value, BMM to use for bsize

2557:    Concepts: vector^stash
2558:    Concepts: stash^vector

2560: .seealso: VecSetBlockSize(), VecSetValues(), VecSetValuesBlocked(), VecStashView()

2562: @*/
2563: int VecSetStashInitialSize(Vec vec,int size,int bsize)
2564: {

2569:   VecStashSetInitialSize_Private(&vec->stash,size);
2570:   VecStashSetInitialSize_Private(&vec->bstash,bsize);
2571:   return(0);
2572: }

2574: /*@
2575:    VecStashView - Prints the entries in the vector stash and block stash.

2577:    Collective on Vec

2579:    Input Parameters:
2580: +  vec   - the vector
2581: -  viewer - the viewer

2583:    Level: advanced

2585:    Concepts: vector^stash
2586:    Concepts: stash^vector

2588: .seealso: VecSetBlockSize(), VecSetValues(), VecSetValuesBlocked()

2590: @*/
2591: int VecStashView(Vec v,PetscViewer viewer)
2592: {
2593:   int          ierr,rank,i,j;
2594:   PetscTruth   match;
2595:   VecStash     *s;
2596:   PetscScalar  val;


2603:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&match);
2604:   if (!match) SETERRQ1(1,"Stash viewer only works with ASCII viewer not %sn",((PetscObject)v)->type_name);
2605:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
2606:   MPI_Comm_rank(v->comm,&rank);
2607:   s = &v->bstash;

2609:   /* print block stash */
2610:   PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector Block stash size %d block size %dn",rank,s->n,s->bs);
2611:   for (i=0; i<s->n; i++) {
2612:     PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d ",rank,s->idx[i]);
2613:     for (j=0; j<s->bs; j++) {
2614:       val = s->array[i*s->bs+j];
2615: #if defined(PETSC_USE_COMPLEX)
2616:       PetscViewerASCIISynchronizedPrintf(viewer,"(%18.16e %18.16e) ",PetscRealPart(val),PetscImaginaryPart(val));
2617: #else
2618:       PetscViewerASCIISynchronizedPrintf(viewer,"%18.16e ",val);
2619: #endif
2620:     }
2621:     PetscViewerASCIISynchronizedPrintf(viewer,"n");
2622:   }
2623:   PetscViewerFlush(viewer);

2625:   s = &v->stash;

2627:   /* print basic stash */
2628:   PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector stash size %dn",rank,s->n);
2629:   for (i=0; i<s->n; i++) {
2630:     val = s->array[i];
2631: #if defined(PETSC_USE_COMPLEX)
2632:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d (%18.16e %18.16e) ",rank,s->idx[i],PetscRealPart(val),PetscImaginaryPart(val));
2633: #else
2634:     PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d %18.16en",rank,s->idx[i],val);
2635: #endif
2636:   }
2637:   PetscViewerFlush(viewer);

2639:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
2640:   return(0);
2641: }

2643: /*@C
2644:    VecGetArray2d - Returns a pointer to a 2d contiguous array that contains this 
2645:    processor's portion of the vector data.  You MUST call VecRestoreArray2d() 
2646:    when you no longer need access to the array.

2648:    Not Collective

2650:    Input Parameter:
2651: +  x - the vector
2652: .  m - first dimension of two dimensional array
2653: .  n - second dimension of two dimensional array
2654: .  mstart - first index you will use in first coordinate direction (often 0)
2655: -  nstart - first index in the second coordinate direction (often 0)

2657:    Output Parameter:
2658: .  a - location to put pointer to the array

2660:    Level: beginner

2662:   Notes:
2663:    For a vector obtained from DACreateLocalVector() mstart and nstart are likely
2664:    obtained from the corner indices obtained from DAGetGhostCorners() while for
2665:    DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
2666:    the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray2d().
2667:    
2668:    For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2670:    Concepts: vector^accessing local values as 2d array

2672: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2673:           VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2674:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2675: @*/
2676: int VecGetArray2d(Vec x,int m,int n,int mstart,int nstart,PetscScalar **a[])
2677: {
2678:   int         i,ierr,N;
2679:   PetscScalar *aa;

2685:   VecGetLocalSize(x,&N);
2686:   if (m*n != N) SETERRQ3(1,"Local array size %d does not match 2d array dimensions %d by %d",N,m,n);
2687:   VecGetArray(x,&aa);

2689:   PetscMalloc(m*sizeof(PetscScalar*),a);
2690:   for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
2691:   *a -= mstart;
2692:   return(0);
2693: }

2695: /*@C
2696:    VecRestoreArray2d - Restores a vector after VecGetArray2d() has been called.

2698:    Not Collective

2700:    Input Parameters:
2701: +  x - the vector
2702: .  m - first dimension of two dimensional array
2703: .  n - second dimension of the two dimensional array
2704: .  mstart - first index you will use in first coordinate direction (often 0)
2705: .  nstart - first index in the second coordinate direction (often 0)
2706: -  a - location of pointer to array obtained from VecGetArray2d()

2708:    Level: beginner

2710:    Notes:
2711:    For regular PETSc vectors this routine does not involve any copies. For
2712:    any special vectors that do not store local vector data in a contiguous
2713:    array, this routine will copy the data back into the underlying 
2714:    vector data structure from the array obtained with VecGetArray().

2716:    This routine actually zeros out the a pointer. 

2718: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2719:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2720:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2721: @*/
2722: int VecRestoreArray2d(Vec x,int m,int n,int mstart,int nstart,PetscScalar **a[])
2723: {

2729:   PetscFree(*a + mstart);
2730:   VecRestoreArray(x,PETSC_NULL);
2731:   return(0);
2732: }

2734: /*@C
2735:    VecGetArray1d - Returns a pointer to a 1d contiguous array that contains this 
2736:    processor's portion of the vector data.  You MUST call VecRestoreArray1d() 
2737:    when you no longer need access to the array.

2739:    Not Collective

2741:    Input Parameter:
2742: +  x - the vector
2743: .  m - first dimension of two dimensional array
2744: -  mstart - first index you will use in first coordinate direction (often 0)

2746:    Output Parameter:
2747: .  a - location to put pointer to the array

2749:    Level: beginner

2751:   Notes:
2752:    For a vector obtained from DACreateLocalVector() mstart are likely
2753:    obtained from the corner indices obtained from DAGetGhostCorners() while for
2754:    DACreateGlobalVector() they are the corner indices from DAGetCorners(). 
2755:    
2756:    For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2758: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2759:           VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2760:           VecGetarray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2761: @*/
2762: int VecGetArray1d(Vec x,int m,int mstart,PetscScalar *a[])
2763: {
2764:   int ierr,N;

2770:   VecGetLocalSize(x,&N);
2771:   if (m != N) SETERRQ2(1,"Local array size %d does not match 1d array dimensions %d",N,m);
2772:   VecGetArray(x,a);
2773:   *a  -= mstart;
2774:   return(0);
2775: }

2777: /*@C
2778:    VecRestoreArray1d - Restores a vector after VecGetArray1d() has been called.

2780:    Not Collective

2782:    Input Parameters:
2783: +  x - the vector
2784: .  m - first dimension of two dimensional array
2785: .  mstart - first index you will use in first coordinate direction (often 0)
2786: -  a - location of pointer to array obtained from VecGetArray21()

2788:    Level: beginner

2790:    Notes:
2791:    For regular PETSc vectors this routine does not involve any copies. For
2792:    any special vectors that do not store local vector data in a contiguous
2793:    array, this routine will copy the data back into the underlying 
2794:    vector data structure from the array obtained with VecGetArray1d().

2796:    This routine actually zeros out the a pointer. 

2798:    Concepts: vector^accessing local values as 1d array

2800: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2801:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2802:           VecGetarray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2803: @*/
2804: int VecRestoreArray1d(Vec x,int m,int mstart,PetscScalar *a[])
2805: {

2811:   VecRestoreArray(x,PETSC_NULL);
2812:   return(0);
2813: }

2815: /*@C
2816:    VecConjugate - Conjugates a vector.

2818:    Collective on Vec

2820:    Input Parameters:
2821: .  x - the vector

2823:    Level: intermediate

2825:    Concepts: vector^conjugate

2827: @*/
2828: int VecConjugate(Vec x)
2829: {
2830: #ifdef PETSC_USE_COMPLEX

2836:   (*x->ops->conjugate)(x);
2837:   return(0);
2838: #else
2839:   return(0);
2840: #endif
2841: }

2843: /*@C
2844:    VecGetArray3d - Returns a pointer to a 3d contiguous array that contains this 
2845:    processor's portion of the vector data.  You MUST call VecRestoreArray3d() 
2846:    when you no longer need access to the array.

2848:    Not Collective

2850:    Input Parameter:
2851: +  x - the vector
2852: .  m - first dimension of three dimensional array
2853: .  n - second dimension of three dimensional array
2854: .  p - third dimension of three dimensional array
2855: .  mstart - first index you will use in first coordinate direction (often 0)
2856: .  nstart - first index in the second coordinate direction (often 0)
2857: -  pstart - first index in the third coordinate direction (often 0)

2859:    Output Parameter:
2860: .  a - location to put pointer to the array

2862:    Level: beginner

2864:   Notes:
2865:    For a vector obtained from DACreateLocalVector() mstart, nstart, and pstart are likely
2866:    obtained from the corner indices obtained from DAGetGhostCorners() while for
2867:    DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
2868:    the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray3d().
2869:    
2870:    For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2872:    Concepts: vector^accessing local values as 3d array

2874: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2875:           VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2876:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2877: @*/
2878: int VecGetArray3d(Vec x,int m,int n,int p,int mstart,int nstart,int pstart,PetscScalar ***a[])
2879: {
2880:   int         i,ierr,N,j;
2881:   PetscScalar *aa,**b;

2887:   VecGetLocalSize(x,&N);
2888:   if (m*n*p != N) SETERRQ4(1,"Local array size %d does not match 3d array dimensions %d by %d by %d",N,m,n,p);
2889:   VecGetArray(x,&aa);

2891:   PetscMalloc(m*sizeof(PetscScalar**)+m*n*sizeof(PetscScalar*),a);
2892:   b    = (PetscScalar **)((*a) + m);
2893:   for (i=0; i<m; i++)   (*a)[i] = b + i*n - nstart;
2894:   for (i=0; i<m; i++) {
2895:     for (j=0; j<n; j++) {
2896:       b[i*n+j] = aa + i*n*p + j*p - pstart;
2897:     }
2898:   }
2899:   *a -= mstart;
2900:   return(0);
2901: }

2903: /*@C
2904:    VecRestoreArray3d - Restores a vector after VecGetArray3d() has been called.

2906:    Not Collective

2908:    Input Parameters:
2909: +  x - the vector
2910: .  m - first dimension of three dimensional array
2911: .  n - second dimension of the three dimensional array
2912: .  p - third dimension of the three dimensional array
2913: .  mstart - first index you will use in first coordinate direction (often 0)
2914: .  nstart - first index in the second coordinate direction (often 0)
2915: .  pstart - first index in the third coordinate direction (often 0)
2916: -  a - location of pointer to array obtained from VecGetArray3d()

2918:    Level: beginner

2920:    Notes:
2921:    For regular PETSc vectors this routine does not involve any copies. For
2922:    any special vectors that do not store local vector data in a contiguous
2923:    array, this routine will copy the data back into the underlying 
2924:    vector data structure from the array obtained with VecGetArray().

2926:    This routine actually zeros out the a pointer. 

2928: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2929:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2930:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2931: @*/
2932: int VecRestoreArray3d(Vec x,int m,int n,int p,int mstart,int nstart,int pstart,PetscScalar ***a[])
2933: {

2939:   PetscFree(*a + mstart);
2940:   VecRestoreArray(x,PETSC_NULL);
2941:   return(0);
2942: }