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: beginner

 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:   PetscOptionsEnd();

106: #if defined(__cplusplus) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_MATSINGLE) && defined(PETSC_HAVE_CXX_NAMESPACE)
107:   VecESISetFromOptions(vec);
108: #endif

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 Parameters:
152: + v - the vector
153: . n - the local size (or PETSC_DECIDE to have it set)
154: - N - the global size (or PETSC_DECIDE)

156:   Notes:
157:   n and N cannot be both PETSC_DECIDE
158:   If one processor calls this with N of PETSC_DECIDE then all processors must, otherwise the program will hang.

160:   Level: intermediate

162: .seealso: VecGetSize(), PetscSplitOwnership()
163: @*/
164: int VecSetSizes(Vec v, int n, int N)
165: {

170:   v->n = n;
171:   v->N = N;
172:   PetscSplitOwnership(v->comm, &v->n, &v->N);
173:   return(0);
174: }

176: /*@
177:    VecSetBlockSize - Sets the blocksize for future calls to VecSetValuesBlocked()
178:    and VecSetValuesBlockedLocal().

180:    Collective on Vec

182:    Input Parameter:
183: +  v - the vector
184: -  bs - the blocksize

186:    Notes:
187:    All vectors obtained by VecDuplicate() inherit the same blocksize.

189:    Level: advanced

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

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

211: /*@
212:    VecGetBlockSize - Gets the blocksize for the vector, i.e. what is used for VecSetValuesBlocked()
213:    and VecSetValuesBlockedLocal().

215:    Collective on Vec

217:    Input Parameter:
218: .  v - the vector

220:    Output Parameter:
221: .  bs - the blocksize

223:    Notes:
224:    All vectors obtained by VecDuplicate() inherit the same blocksize.

226:    Level: advanced

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

230:    Concepts: vector^block size
231:    Concepts: block^vector

233: @*/
234: int VecGetBlockSize(Vec v,int *bs)
235: {
238:   *bs = v->bs;
239:   return(0);
240: }

242: /*@
243:    VecValid - Checks whether a vector object is valid.

245:    Not Collective

247:    Input Parameter:
248: .  v - the object to check

250:    Output Parameter:
251:    flg - flag indicating vector status, either
252:    PETSC_TRUE if vector is valid, or PETSC_FALSE otherwise.

254:    Level: developer

256: @*/
257: int VecValid(Vec v,PetscTruth *flg)
258: {
261:   if (!v)                           *flg = PETSC_FALSE;
262:   else if (v->cookie != VEC_COOKIE) *flg = PETSC_FALSE;
263:   else                              *flg = PETSC_TRUE;
264:   return(0);
265: }

267: /*@
268:    VecDot - Computes the vector dot product.

270:    Collective on Vec

272:    Input Parameters:
273: .  x, y - the vectors

275:    Output Parameter:
276: .  alpha - the dot product

278:    Performance Issues:
279: +    per-processor memory bandwidth
280: .    interprocessor latency
281: -    work load inbalance that causes certain processes to arrive much earlier than
282:      others

284:    Notes for Users of Complex Numbers:
285:    For complex vectors, VecDot() computes 
286: $     val = (x,y) = y^H x,
287:    where y^H denotes the conjugate transpose of y.

289:    Use VecTDot() for the indefinite form
290: $     val = (x,y) = y^T x,
291:    where y^T denotes the transpose of y.

293:    Level: intermediate

295:    Concepts: inner product
296:    Concepts: vector^inner product

298: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd()
299: @*/
300: int VecDot(Vec x,Vec y,PetscScalar *val)
301: {

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

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

331: /*@
332:    VecNorm  - Computes the vector norm.

334:    Collective on Vec

336:    Input Parameters:
337: +  x - the vector
338: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
339:           NORM_1_AND_2, which computes both norms and stores them
340:           in a two element array.

342:    Output Parameter:
343: .  val - the norm 

345:    Notes:
346: $     NORM_1 denotes sum_i |x_i|
347: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
348: $     NORM_INFINITY denotes max_i |x_i|

350:    Level: intermediate

352:    Performance Issues:
353: +    per-processor memory bandwidth
354: .    interprocessor latency
355: -    work load inbalance that causes certain processes to arrive much earlier than
356:      others

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

363:    Concepts: norm
364:    Concepts: vector^norm

366: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), 
367:           VecNormBegin(), VecNormEnd()

369: @*/
370: int VecNorm(Vec x,NormType type,PetscReal *val)
371: {

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

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

396:    Collective on Vec

398:    Input Parameter:
399: .  x - the vector

401:    Output Parameters:
402: +  val - the maximum component
403: -  p - the location of val

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

408:    Level: intermediate

410:    Concepts: maximum^of vector
411:    Concepts: vector^maximum value

413: .seealso: VecNorm(), VecMin()
414: @*/
415: int VecMax(Vec x,int *p,PetscReal *val)
416: {

423:   PetscLogEventBegin(VEC_Max,x,0,0,0);
424:   (*x->ops->max)(x,p,val);
425:   PetscLogEventEnd(VEC_Max,x,0,0,0);
426:   return(0);
427: }

429: /*@
430:    VecMin - Determines the minimum vector component and its location.

432:    Collective on Vec

434:    Input Parameters:
435: .  x - the vector

437:    Output Parameter:
438: +  val - the minimum component
439: -  p - the location of val

441:    Level: intermediate

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

446:    Concepts: minimum^of vector
447:    Concepts: vector^minimum entry

449: .seealso: VecMax()
450: @*/
451: int VecMin(Vec x,int *p,PetscReal *val)
452: {

459:   PetscLogEventBegin(VEC_Min,x,0,0,0);
460:   (*x->ops->min)(x,p,val);
461:   PetscLogEventEnd(VEC_Min,x,0,0,0);
462:   return(0);
463: }

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

469:    Collective on Vec

471:    Input Parameters:
472: .  x, y - the vectors

474:    Output Parameter:
475: .  val - the dot product

477:    Notes for Users of Complex Numbers:
478:    For complex vectors, VecTDot() computes the indefinite form
479: $     val = (x,y) = y^T x,
480:    where y^T denotes the transpose of y.

482:    Use VecDot() for the inner product
483: $     val = (x,y) = y^H x,
484:    where y^H denotes the conjugate transpose of y.

486:    Level: intermediate

488:    Concepts: inner product^non-Hermitian
489:    Concepts: vector^inner product
490:    Concepts: non-Hermitian inner product

492: .seealso: VecDot(), VecMTDot()
493: @*/
494: int VecTDot(Vec x,Vec y,PetscScalar *val)
495: {

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

509:   PetscLogEventBegin(VEC_TDot,x,y,0,0);
510:   (*x->ops->tdot)(x,y,val);
511:   PetscLogEventEnd(VEC_TDot,x,y,0,0);
512:   return(0);
513: }

515: /*@
516:    VecScale - Scales a vector. 

518:    Collective on Vec

520:    Input Parameters:
521: +  x - the vector
522: -  alpha - the scalar

524:    Output Parameter:
525: .  x - the scaled vector

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

531:    Level: intermediate

533:    Concepts: vector^scaling
534:    Concepts: scaling^vector

536: @*/
537: int VecScale (const PetscScalar *alpha,Vec x)
538: {

545:   PetscLogEventBegin(VEC_Scale,x,0,0,0);
546:   (*x->ops->scale)(alpha,x);
547:   PetscLogEventEnd(VEC_Scale,x,0,0,0);
548:   return(0);
549: }

551: /*@
552:    VecCopy - Copies a vector. 

554:    Collective on Vec

556:    Input Parameter:
557: .  x - the vector

559:    Output Parameter:
560: .  y - the copy

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

566:    Level: beginner

568: .seealso: VecDuplicate()
569: @*/
570: int VecCopy(Vec x,Vec y)
571: {

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

583:   PetscLogEventBegin(VEC_Copy,x,y,0,0);
584:   (*x->ops->copy)(x,y);
585:   PetscLogEventEnd(VEC_Copy,x,y,0,0);
586:   return(0);
587: }

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

592:    Collective on Vec

594:    Input Parameters:
595: +  alpha - the scalar
596: -  x  - the vector

598:    Output Parameter:
599: .  x  - the vector

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

608:    Level: beginner

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

612:    Concepts: vector^setting to constant

614: @*/
615: int VecSet(const PetscScalar *alpha,Vec x)
616: {


624:   PetscLogEventBegin(VEC_Set,x,0,0,0);
625:   (*x->ops->set)(alpha,x);
626:   PetscLogEventEnd(VEC_Set,x,0,0,0);
627:   return(0);
628: }

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

633:    Collective on Vec

635:    Input Parameters:
636: +  rctx - the random number context, formed by PetscRandomCreate(), or PETSC_NULL and
637:           it will create one internally.
638: -  x  - the vector

640:    Output Parameter:
641: .  x  - the vector

643:    Example of Usage:
644: .vb
645:      PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
646:      VecSetRandom(rctx,x);
647:      PetscRandomDestroy(rctx);
648: .ve

650:    Level: intermediate

652:    Concepts: vector^setting to random
653:    Concepts: random^vector

655: .seealso: VecSet(), VecSetValues(), PetscRandomCreate(), PetscRandomDestroy()
656: @*/
657: int VecSetRandom(PetscRandom rctx,Vec x)
658: {
659:   int         ierr;
660:   PetscRandom randObj = PETSC_NULL;


667:   if (!rctx) {
668:     MPI_Comm    comm;
669:     PetscObjectGetComm((PetscObject)x,&comm);
670:     PetscRandomCreate(comm,RANDOM_DEFAULT,&randObj);
671:     rctx = randObj;
672:   }

674:   PetscLogEventBegin(VEC_SetRandom,x,rctx,0,0);
675:   (*x->ops->setrandom)(rctx,x);
676:   PetscLogEventEnd(VEC_SetRandom,x,rctx,0,0);
677: 
678:   if (randObj) {
679:     PetscRandomDestroy(randObj);
680:   }
681:   return(0);
682: }

684: /*@
685:    VecAXPY - Computes y = alpha x + y. 

687:    Collective on Vec

689:    Input Parameters:
690: +  alpha - the scalar
691: -  x, y  - the vectors

693:    Output Parameter:
694: .  y - output vector

696:    Level: intermediate

698:    Concepts: vector^BLAS
699:    Concepts: BLAS

701: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
702: @*/
703: int VecAXPY(const PetscScalar *alpha,Vec x,Vec y)
704: {

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

718:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
719:   (*x->ops->axpy)(alpha,x,y);
720:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
721:   return(0);
722: }

724: /*@
725:    VecAXPBY - Computes y = alpha x + beta y. 

727:    Collective on Vec

729:    Input Parameters:
730: +  alpha,beta - the scalars
731: -  x, y  - the vectors

733:    Output Parameter:
734: .  y - output vector

736:    Level: intermediate

738:    Concepts: BLAS
739:    Concepts: vector^BLAS

741: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
742: @*/
743: int VecAXPBY(const PetscScalar *alpha,const PetscScalar *beta,Vec x,Vec y)
744: {

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

759:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
760:   (*x->ops->axpby)(alpha,beta,x,y);
761:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
762:   return(0);
763: }

765: /*@
766:    VecAYPX - Computes y = x + alpha y.

768:    Collective on Vec

770:    Input Parameters:
771: +  alpha - the scalar
772: -  x, y  - the vectors

774:    Output Parameter:
775: .  y - output vector

777:    Level: intermediate

779:    Concepts: vector^BLAS
780:    Concepts: BLAS

782: .seealso: VecAXPY(), VecWAXPY()
783: @*/
784: int VecAYPX(const PetscScalar *alpha,Vec x,Vec y)
785: {

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

799:   PetscLogEventBegin(VEC_AYPX,x,y,0,0);
800:    (*x->ops->aypx)(alpha,x,y);
801:   PetscLogEventEnd(VEC_AYPX,x,y,0,0);
802:   return(0);
803: }

805: /*@
806:    VecSwap - Swaps the vectors x and y.

808:    Collective on Vec

810:    Input Parameters:
811: .  x, y  - the vectors

813:    Level: advanced

815:    Concepts: vector^swapping values

817: @*/
818: int VecSwap(Vec x,Vec y)
819: {

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

832:   PetscLogEventBegin(VEC_Swap,x,y,0,0);
833:   (*x->ops->swap)(x,y);
834:   PetscLogEventEnd(VEC_Swap,x,y,0,0);
835:   return(0);
836: }

838: /*@
839:    VecWAXPY - Computes w = alpha x + y.

841:    Collective on Vec

843:    Input Parameters:
844: +  alpha - the scalar
845: -  x, y  - the vectors

847:    Output Parameter:
848: .  w - the result

850:    Level: intermediate

852:    Concepts: vector^BLAS
853:    Concepts: BLAS

855: .seealso: VecAXPY(), VecAYPX()
856: @*/
857: int VecWAXPY(const PetscScalar *alpha,Vec x,Vec y,Vec w)
858: {

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

874:   PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
875:    (*x->ops->waxpy)(alpha,x,y,w);
876:   PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
877:   return(0);
878: }

880: /*@
881:    VecPointwiseMult - Computes the componentwise multiplication w = x*y.

883:    Collective on Vec

885:    Input Parameters:
886: .  x, y  - the vectors

888:    Output Parameter:
889: .  w - the result

891:    Level: advanced

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

895:    Concepts: vector^pointwise multiply

897: .seealso: VecPointwiseDivide()
898: @*/
899: int VecPointwiseMult(Vec x,Vec y,Vec w)
900: {

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

915:   PetscLogEventBegin(VEC_PointwiseMult,x,y,w,0);
916:   (*x->ops->pointwisemult)(x,y,w);
917:   PetscLogEventEnd(VEC_PointwiseMult,x,y,w,0);
918:   return(0);
919: }

921: /*@
922:    VecPointwiseDivide - Computes the componentwise division w = x/y.

924:    Collective on Vec

926:    Input Parameters:
927: .  x, y  - the vectors

929:    Output Parameter:
930: .  w - the result

932:    Level: advanced

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

936:    Concepts: vector^pointwise divide

938: .seealso: VecPointwiseMult()
939: @*/
940: int VecPointwiseDivide(Vec x,Vec y,Vec w)
941: {

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

956:   (*x->ops->pointwisedivide)(x,y,w);
957:   return(0);
958: }

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

963:    Collective on Vec

965:    Input Parameters:
966: .  v - a vector to mimic

968:    Output Parameter:
969: .  newv - location to put new vector

971:    Notes:
972:    VecDuplicate() does not copy the vector, but rather allocates storage
973:    for the new vector.  Use VecCopy() to copy a vector.

975:    Use VecDestroy() to free the space. Use VecDuplicateVecs() to get several
976:    vectors. 

978:    Level: beginner

980: .seealso: VecDestroy(), VecDuplicateVecs(), VecCreate(), VecCopy()
981: @*/
982: int VecDuplicate(Vec x,Vec *newv)
983: {

990:   (*x->ops->duplicate)(x,newv);
991:   return(0);
992: }

994: /*@C
995:    VecDestroy - Destroys a vector.

997:    Collective on Vec

999:    Input Parameters:
1000: .  v  - the vector

1002:    Level: beginner

1004: .seealso: VecDuplicate()
1005: @*/
1006: int VecDestroy(Vec v)
1007: {

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

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

1035:    Collective on Vec

1037:    Input Parameters:
1038: +  m - the number of vectors to obtain
1039: -  v - a vector to mimic

1041:    Output Parameter:
1042: .  V - location to put pointer to array of vectors

1044:    Notes:
1045:    Use VecDestroyVecs() to free the space. Use VecDuplicate() to form a single
1046:    vector.

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

1053:    Level: intermediate

1055: .seealso:  VecDestroyVecs(), VecDuplicate(), VecCreate(), VecDuplicateVecsF90()
1056: @*/
1057: int VecDuplicateVecs(Vec v,int m,Vec *V[])
1058: {

1065:   (*v->ops->duplicatevecs)(v, m,V);
1066:   return(0);
1067: }

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

1072:    Collective on Vec

1074:    Input Parameters:
1075: +  vv - pointer to array of vector pointers
1076: -  m - the number of vectors previously obtained

1078:    Fortran Note:
1079:    The Fortran interface is slightly different from that given below.
1080:    See the Fortran chapter of the users manual and 
1081:    petsc/src/vec/examples for details.

1083:    Level: intermediate

1085: .seealso: VecDuplicateVecs(), VecDestroyVecsF90()
1086: @*/
1087: int VecDestroyVecs(const Vec vv[],int m)
1088: {

1092:   if (!vv) SETERRQ(PETSC_ERR_ARG_BADPTR,"Null vectors");
1095:   (*(*vv)->ops->destroyvecs)(vv,m);
1096:   return(0);
1097: }

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

1102:    Input Parameters:
1103:    Not Collective

1105: +  x - vector to insert in
1106: .  ni - number of elements to add
1107: .  ix - indices where to add
1108: .  y - array of values
1109: -  iora - either INSERT_VALUES or ADD_VALUES, where
1110:    ADD_VALUES adds values to any existing entries, and
1111:    INSERT_VALUES replaces existing entries with new values

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

1116:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES 
1117:    options cannot be mixed without intervening calls to the assembly
1118:    routines.

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

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

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

1130:    Level: beginner

1132:    Concepts: vector^setting values

1134: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
1135:            VecSetValue(), VecSetValuesBlocked()
1136: @*/
1137: int VecSetValues(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1138: {

1146:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1147:   (*x->ops->setvalues)(x,ni,ix,y,iora);
1148:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1149:   return(0);
1150: }

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

1155:    Not Collective

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

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

1170:    Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES 
1171:    options cannot be mixed without intervening calls to the assembly
1172:    routines.

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

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

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

1184:    Level: intermediate

1186:    Concepts: vector^setting values blocked

1188: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
1189:            VecSetValues()
1190: @*/
1191: int VecSetValuesBlocked(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1192: {

1200:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1201:   (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
1202:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1203:   return(0);
1204: }

1206: /*MC
1207:    VecSetValue - Set a single entry into a vector.

1209:    Synopsis:
1210:    int VecSetValue(Vec v,int row,PetscScalar value, InsertMode mode);

1212:    Not Collective

1214:    Input Parameters:
1215: +  v - the vector
1216: .  row - the row location of the entry
1217: .  value - the value to insert
1218: -  mode - either INSERT_VALUES or ADD_VALUES

1220:    Notes:
1221:    For efficiency one should use VecSetValues() and set several or 
1222:    many values simultaneously if possible.

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

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

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

1232:    Level: beginner

1234: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal()
1235: M*/

1237: /*@
1238:    VecSetLocalToGlobalMapping - Sets a local numbering to global numbering used
1239:    by the routine VecSetValuesLocal() to allow users to insert vector entries
1240:    using a local (per-processor) numbering.

1242:    Collective on Vec

1244:    Input Parameters:
1245: +  x - vector
1246: -  mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()

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

1251:    Level: intermediate

1253:    Concepts: vector^setting values with local numbering

1255: seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
1256:            VecSetLocalToGlobalMappingBlocked(), VecSetValuesBlockedLocal()
1257: @*/
1258: int VecSetLocalToGlobalMapping(Vec x,ISLocalToGlobalMapping mapping)
1259: {

1265:   if (x->mapping) {
1266:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for vector");
1267:   }

1269:   if (x->ops->setlocaltoglobalmapping) {
1270:     (*x->ops->setlocaltoglobalmapping)(x,mapping);
1271:   } else {
1272:     x->mapping = mapping;
1273:     PetscObjectReference((PetscObject)mapping);
1274:   }
1275:   return(0);
1276: }

1278: /*@
1279:    VecSetLocalToGlobalMappingBlock - Sets a local numbering to global numbering used
1280:    by the routine VecSetValuesBlockedLocal() to allow users to insert vector entries
1281:    using a local (per-processor) numbering.

1283:    Collective on Vec

1285:    Input Parameters:
1286: +  x - vector
1287: -  mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()

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

1292:    Level: intermediate

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

1296: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
1297:            VecSetLocalToGlobalMapping(), VecSetValuesBlockedLocal()
1298: @*/
1299: int VecSetLocalToGlobalMappingBlock(Vec x,ISLocalToGlobalMapping mapping)
1300: {

1306:   if (x->bmapping) {
1307:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for vector");
1308:   }
1309:   x->bmapping = mapping;
1310:   PetscObjectReference((PetscObject)mapping);
1311:   return(0);
1312: }

1314: /*@
1315:    VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
1316:    using a local ordering of the nodes. 

1318:    Not Collective

1320:    Input Parameters:
1321: +  x - vector to insert in
1322: .  ni - number of elements to add
1323: .  ix - indices where to add
1324: .  y - array of values
1325: -  iora - either INSERT_VALUES or ADD_VALUES, where
1326:    ADD_VALUES adds values to any existing entries, and
1327:    INSERT_VALUES replaces existing entries with new values

1329:    Level: intermediate

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

1334:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES 
1335:    options cannot be mixed without intervening calls to the assembly
1336:    routines.

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

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

1343:    Concepts: vector^setting values with local numbering

1345: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
1346:            VecSetValuesBlockedLocal()
1347: @*/
1348: int VecSetValuesLocal(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1349: {
1350:   int ierr,lixp[128],*lix = lixp;


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

1378: /*@
1379:    VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1380:    using a local ordering of the nodes. 

1382:    Not Collective

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

1393:    Level: intermediate

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

1399:    Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES 
1400:    options cannot be mixed without intervening calls to the assembly
1401:    routines.

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

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


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

1411: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(), 
1412:            VecSetLocalToGlobalMappingBlocked()
1413: @*/
1414: int VecSetValuesBlockedLocal(Vec x,int ni,const int ix[],const PetscScalar y[],InsertMode iora)
1415: {
1416:   int ierr,lixp[128],*lix = lixp;

1423:   if (!x->bmapping) {
1424:     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMappingBlocked()");
1425:   }
1426:   if (ni > 128) {
1427:     PetscMalloc(ni*sizeof(int),&lix);
1428:   }

1430:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1431:   ISLocalToGlobalMappingApply(x->bmapping,ni,(int*)ix,lix);
1432:   (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1433:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1434:   if (ni > 128) {
1435:     PetscFree(lix);
1436:   }
1437:   return(0);
1438: }

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

1444:    Collective on Vec

1446:    Input Parameter:
1447: .  vec - the vector

1449:    Level: beginner

1451:    Concepts: assembly^vectors

1453: .seealso: VecAssemblyEnd(), VecSetValues()
1454: @*/
1455: int VecAssemblyBegin(Vec vec)
1456: {
1457:   int        ierr;
1458:   PetscTruth flg;


1464:   PetscOptionsHasName(vec->prefix,"-vec_view_stash",&flg);
1465:   if (flg) {
1466:     VecStashView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1467:   }

1469:   PetscLogEventBegin(VEC_AssemblyBegin,vec,0,0,0);
1470:   if (vec->ops->assemblybegin) {
1471:     (*vec->ops->assemblybegin)(vec);
1472:   }
1473:   PetscLogEventEnd(VEC_AssemblyBegin,vec,0,0,0);
1474:   return(0);
1475: }

1477: /*@
1478:    VecAssemblyEnd - Completes assembling the vector.  This routine should
1479:    be called after VecAssemblyBegin().

1481:    Collective on Vec

1483:    Input Parameter:
1484: .  vec - the vector

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

1497: .seealso: VecAssemblyBegin(), VecSetValues()
1498: @*/
1499: int VecAssemblyEnd(Vec vec)
1500: {
1501:   int        ierr;
1502:   PetscTruth flg;

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


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

1557:    Collective on Vec

1559:    Input Parameters:
1560: +  nv - number of vectors
1561: .  x - one vector
1562: -  y - array of vectors.  Note that vectors are pointers

1564:    Output Parameter:
1565: .  val - array of the dot products

1567:    Notes for Users of Complex Numbers:
1568:    For complex vectors, VecMTDot() computes the indefinite form
1569: $      val = (x,y) = y^T x,
1570:    where y^T denotes the transpose of y.

1572:    Use VecMDot() for the inner product
1573: $      val = (x,y) = y^H x,
1574:    where y^H denotes the conjugate transpose of y.

1576:    Level: intermediate

1578:    Concepts: inner product^multiple
1579:    Concepts: vector^multiple inner products

1581: .seealso: VecMDot(), VecTDot()
1582: @*/
1583: int VecMTDot(int nv,Vec x,const Vec y[],PetscScalar *val)
1584: {

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

1598:   PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1599:   (*x->ops->mtdot)(nv,x,y,val);
1600:   PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1601:   return(0);
1602: }

1604: /*@C
1605:    VecMDot - Computes vector multiple dot products. 

1607:    Collective on Vec

1609:    Input Parameters:
1610: +  nv - number of vectors
1611: .  x - one vector
1612: -  y - array of vectors. 

1614:    Output Parameter:
1615: .  val - array of the dot products

1617:    Notes for Users of Complex Numbers:
1618:    For complex vectors, VecMDot() computes 
1619: $     val = (x,y) = y^H x,
1620:    where y^H denotes the conjugate transpose of y.

1622:    Use VecMTDot() for the indefinite form
1623: $     val = (x,y) = y^T x,
1624:    where y^T denotes the transpose of y.

1626:    Level: intermediate

1628:    Concepts: inner product^multiple
1629:    Concepts: vector^multiple inner products

1631: .seealso: VecMTDot(), VecDot()
1632: @*/
1633: int VecMDot(int nv,Vec x,const Vec y[],PetscScalar *val)
1634: {

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

1648:   PetscLogEventBarrierBegin(VEC_MDotBarrier,x,*y,0,0,x->comm);
1649:   (*x->ops->mdot)(nv,x,y,val);
1650:   PetscLogEventBarrierEnd(VEC_MDotBarrier,x,*y,0,0,x->comm);
1651:   return(0);
1652: }

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

1657:    Collective on Vec

1659:    Input Parameters:
1660: +  nv - number of scalars and x-vectors
1661: .  alpha - array of scalars
1662: .  y - one vector
1663: -  x - array of vectors

1665:    Level: intermediate

1667:    Concepts: BLAS

1669: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1670: @*/
1671: int  VecMAXPY(int nv,const PetscScalar *alpha,Vec y,Vec *x)
1672: {

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

1686:   PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1687:   (*y->ops->maxpy)(nv,alpha,y,x);
1688:   PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1689:   return(0);
1690: }

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

1701:    Not Collective

1703:    Input Parameter:
1704: .  x - the vector

1706:    Output Parameter:
1707: .  a - location to put pointer to the array

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

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

1727:    Level: beginner

1729:    Concepts: vector^accessing local values

1731: .seealso: VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(), VecGetArray2d()
1732: @*/
1733: int VecGetArray(Vec x,PetscScalar *a[])
1734: {

1741:   (*x->ops->getarray)(x,a);
1742:   return(0);
1743: }


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

1751:    Not Collective

1753:    Input Parameter:
1754: +  x - the vectors
1755: -  n - the number of vectors

1757:    Output Parameter:
1758: .  a - location to put pointer to the array

1760:    Fortran Note:
1761:    This routine is not supported in Fortran.

1763:    Level: intermediate

1765: .seealso: VecGetArray(), VecRestoreArrays()
1766: @*/
1767: int VecGetArrays(const Vec x[],int n,PetscScalar **a[])
1768: {
1769:   int         i,ierr;
1770:   PetscScalar **q;

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

1784: /*@C
1785:    VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1786:    has been called.

1788:    Not Collective

1790:    Input Parameters:
1791: +  x - the vector
1792: .  n - the number of vectors
1793: -  a - location of pointer to arrays obtained from VecGetArrays()

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

1801:    Fortran Note:
1802:    This routine is not supported in Fortran.

1804:    Level: intermediate

1806: .seealso: VecGetArrays(), VecRestoreArray()
1807: @*/
1808: int VecRestoreArrays(const Vec x[],int n,PetscScalar **a[])
1809: {
1810:   int         i,ierr;
1811:   PetscScalar **q = *a;


1817:   for(i=0;i<n;++i) {
1818:     VecRestoreArray(x[i],&q[i]);
1819:   }
1820:   PetscFree(q);
1821:   return(0);
1822: }

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

1827:    Not Collective

1829:    Input Parameters:
1830: +  x - the vector
1831: -  a - location of pointer to array obtained from VecGetArray()

1833:    Level: beginner

1835:    Notes:
1836:    For regular PETSc vectors this routine does not involve any copies. For
1837:    any special vectors that do not store local vector data in a contiguous
1838:    array, this routine will copy the data back into the underlying 
1839:    vector data structure from the array obtained with VecGetArray().

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

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

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

1863: .seealso: VecGetArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(), VecRestoreArray2d()
1864: @*/
1865: int VecRestoreArray(Vec x,PetscScalar *a[])
1866: {

1873:   if (x->ops->restorearray) {
1874:     (*x->ops->restorearray)(x,a);
1875:   }
1876:   return(0);
1877: }

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

1882:   Collective on Vec

1884:   Input Parameters:
1885: . vec   - The vector
1886: . title - The title

1888:   Level: intermediate

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

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

1944: /*@C
1945:    VecView - Views a vector object. 

1947:    Collective on Vec

1949:    Input Parameters:
1950: +  v - the vector
1951: -  viewer - an optional visualization context

1953:    Notes:
1954:    The available visualization contexts include
1955: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
1956: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1957:          output where only the first processor opens
1958:          the file.  All other processors send their 
1959:          data to the first processor to print. 

1961:    You can change the format the vector is printed using the 
1962:    option PetscViewerSetFormat().

1964:    The user can open alternative visualization contexts with
1965: +    PetscViewerASCIIOpen() - Outputs vector to a specified file
1966: .    PetscViewerBinaryOpen() - Outputs vector in binary to a
1967:          specified file; corresponding input uses VecLoad()
1968: .    PetscViewerDrawOpen() - Outputs vector to an X window display
1969: -    PetscViewerSocketOpen() - Outputs vector to Socket viewer

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

1980:    Level: beginner

1982:    Concepts: vector^printing
1983:    Concepts: vector^saving to disk

1985: .seealso: PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscDrawLGCreate(),
1986:           PetscViewerSocketOpen(), PetscViewerBinaryOpen(), VecLoad(), PetscViewerCreate(),
1987:           PetscRealView(), PetscScalarView(), PetscIntView()
1988: @*/
1989: int VecView(Vec vec,PetscViewer viewer)
1990: {
1991:   int               ierr;
1992:   PetscViewerFormat format;

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

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

2016: /*@
2017:    VecGetSize - Returns the global number of elements of the vector.

2019:    Not Collective

2021:    Input Parameter:
2022: .  x - the vector

2024:    Output Parameters:
2025: .  size - the global length of the vector

2027:    Level: beginner

2029:    Concepts: vector^local size

2031: .seealso: VecGetLocalSize()
2032: @*/
2033: int VecGetSize(Vec x,int *size)
2034: {

2041:   (*x->ops->getsize)(x,size);
2042:   return(0);
2043: }

2045: /*@
2046:    VecGetLocalSize - Returns the number of elements of the vector stored 
2047:    in local memory. This routine may be implementation dependent, so use 
2048:    with care.

2050:    Not Collective

2052:    Input Parameter:
2053: .  x - the vector

2055:    Output Parameter:
2056: .  size - the length of the local piece of the vector

2058:    Level: beginner

2060:    Concepts: vector^size

2062: .seealso: VecGetSize()
2063: @*/
2064: int VecGetLocalSize(Vec x,int *size)
2065: {

2072:   (*x->ops->getlocalsize)(x,size);
2073:   return(0);
2074: }

2076: /*@C
2077:    VecGetOwnershipRange - Returns the range of indices owned by 
2078:    this processor, assuming that the vectors are laid out with the
2079:    first n1 elements on the first processor, next n2 elements on the
2080:    second, etc.  For certain parallel layouts this range may not be 
2081:    well defined. 

2083:    Not Collective

2085:    Input Parameter:
2086: .  x - the vector

2088:    Output Parameters:
2089: +  low - the first local element, pass in PETSC_NULL if not interested
2090: -  high - one more than the last local element, pass in PETSC_NULL if not interested

2092:    Note:
2093:    The high argument is one more than the last element stored locally.

2095:    Fortran: PETSC_NULL_INTEGER should be used instead of PETSC_NULL

2097:    Level: beginner

2099:    Concepts: ownership^of vectors
2100:    Concepts: vector^ownership of elements

2102: @*/
2103: int VecGetOwnershipRange(Vec x,int *low,int *high)
2104: {
2105:   int      ierr;

2112:   PetscMapGetLocalRange(x->map,low,high);
2113:   return(0);
2114: }

2116: /*@C
2117:    VecGetPetscMap - Returns the map associated with the vector

2119:    Not Collective

2121:    Input Parameter:
2122: .  x - the vector

2124:    Output Parameters:
2125: .  map - the map

2127:    Level: developer

2129: @*/
2130: int VecGetPetscMap(Vec x,PetscMap *map)
2131: {
2135:   *map = x->map;
2136:   return(0);
2137: }

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

2142:    Collective on Vec

2144:    Input Parameter:
2145: +  x - the vector
2146: -  op - the option

2148:    Supported Options:
2149: +     VEC_IGNORE_OFF_PROC_ENTRIES, which causes VecSetValues() to ignore 
2150:       entries destined to be stored on a seperate processor. This can be used
2151:       to eliminate the global reduction in the VecAssemblyXXXX() if you know 
2152:       that you have only used VecSetValues() to set local elements
2153: -   VEC_TREAT_OFF_PROC_ENTRIES restores the treatment of off processor entries.

2155:    Level: intermediate

2157: @*/
2158: int VecSetOption(Vec x,VecOption op)
2159: {

2165:   if (x->ops->setoption) {
2166:     (*x->ops->setoption)(x,op);
2167:   }
2168:   return(0);
2169: }

2171: /* Default routines for obtaining and releasing; */
2172: /* may be used by any implementation */
2173: int VecDuplicateVecs_Default(Vec w,int m,Vec *V[])
2174: {
2175:   int  i,ierr;

2180:   if (m <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %d",m);
2181:   PetscMalloc(m*sizeof(Vec*),V);
2182:   for (i=0; i<m; i++) {VecDuplicate(w,*V+i);}
2183:   return(0);
2184: }

2186: int VecDestroyVecs_Default(const Vec v[], int m)
2187: {
2188:   int i,ierr;

2192:   if (m <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %d",m);
2193:   for (i=0; i<m; i++) {VecDestroy(v[i]);}
2194:   PetscFree((Vec*)v);
2195:   return(0);
2196: }

2198: /*@
2199:    VecPlaceArray - Allows one to replace the array in a vector with an
2200:    array provided by the user. This is useful to avoid copying an array
2201:    into a vector.

2203:    Not Collective

2205:    Input Parameters:
2206: +  vec - the vector
2207: -  array - the array

2209:    Notes:
2210:    You can return to the original array with a call to VecResetArray()

2212:    Level: developer

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

2216: @*/
2217: int VecPlaceArray(Vec vec,const PetscScalar array[])
2218: {

2224:   if (vec->ops->placearray) {
2225:     (*vec->ops->placearray)(vec,array);
2226:   } else {
2227:     SETERRQ(1,"Cannot place array in this type of vector");
2228:   }
2229:   return(0);
2230: }

2232: /*@
2233:    VecResetArray - Resets a vector to use its default memory. Call this 
2234:    after the use of VecPlaceArray().

2236:    Not Collective

2238:    Input Parameters:
2239: .  vec - the vector

2241:    Level: developer

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

2245: @*/
2246: int VecResetArray(Vec vec)
2247: {

2253:   if (vec->ops->resetarray) {
2254:     (*vec->ops->resetarray)(vec);
2255:   } else {
2256:     SETERRQ(1,"Cannot reset array in this type of vector");
2257:   }
2258:   return(0);
2259: }

2261: /*@C
2262:    VecReplaceArray - Allows one to replace the array in a vector with an
2263:    array provided by the user. This is useful to avoid copying an array
2264:    into a vector.

2266:    Not Collective

2268:    Input Parameters:
2269: +  vec - the vector
2270: -  array - the array

2272:    Notes:
2273:    This permanently replaces the array and frees the memory associated
2274:    with the old array.

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

2279:    Not supported from Fortran

2281:    Level: developer

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

2285: @*/
2286: int VecReplaceArray(Vec vec,const PetscScalar array[])
2287: {

2293:   if (vec->ops->replacearray) {
2294:     (*vec->ops->replacearray)(vec,array);
2295:   } else {
2296:     SETERRQ(1,"Cannot replace array in this type of vector");
2297:   }
2298:   return(0);
2299: }

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

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

2308:     Collective on Vec

2310:     Input Parameters:
2311: +   x - a vector to mimic
2312: -   n - the number of vectors to obtain

2314:     Output Parameters:
2315: +   y - Fortran90 pointer to the array of vectors
2316: -   ierr - error code

2318:     Example of Usage: 
2319: .vb
2320:     Vec x
2321:     Vec, pointer :: y(:)
2322:     ....
2323:     call VecDuplicateVecsF90(x,2,y,ierr)
2324:     call VecSet(alpha,y(2),ierr)
2325:     call VecSet(alpha,y(2),ierr)
2326:     ....
2327:     call VecDestroyVecsF90(y,2,ierr)
2328: .ve

2330:     Notes:
2331:     Not yet supported for all F90 compilers

2333:     Use VecDestroyVecsF90() to free the space.

2335:     Level: beginner

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

2339: M*/

2341: /*MC
2342:     VecRestoreArrayF90 - Restores a vector to a usable state after a call to
2343:     VecGetArrayF90().

2345:     Synopsis:
2346:     VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2348:     Not collective

2350:     Input Parameters:
2351: +   x - vector
2352: -   xx_v - the Fortran90 pointer to the array

2354:     Output Parameter:
2355: .   ierr - error code

2357:     Example of Usage: 
2358: .vb
2359:     PetscScalar, pointer :: xx_v(:)
2360:     ....
2361:     call VecGetArrayF90(x,xx_v,ierr)
2362:     a = xx_v(3)
2363:     call VecRestoreArrayF90(x,xx_v,ierr)
2364: .ve
2365:    
2366:     Notes:
2367:     Not yet supported for all F90 compilers

2369:     Level: beginner

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

2373: M*/

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

2378:     Synopsis:
2379:     VecDestroyVecsF90({Vec, pointer :: x(:)},integer n,integer ierr)

2381:     Input Parameters:
2382: +   x - pointer to array of vector pointers
2383: -   n - the number of vectors previously obtained

2385:     Output Parameter:
2386: .   ierr - error code

2388:     Notes:
2389:     Not yet supported for all F90 compilers

2391:     Level: beginner

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

2395: M*/

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

2403:     Synopsis:
2404:     VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2406:     Not Collective 

2408:     Input Parameter:
2409: .   x - vector

2411:     Output Parameters:
2412: +   xx_v - the Fortran90 pointer to the array
2413: -   ierr - error code

2415:     Example of Usage: 
2416: .vb
2417:     PetscScalar, pointer :: xx_v(:)
2418:     ....
2419:     call VecGetArrayF90(x,xx_v,ierr)
2420:     a = xx_v(3)
2421:     call VecRestoreArrayF90(x,xx_v,ierr)
2422: .ve

2424:     Notes:
2425:     Not yet supported for all F90 compilers

2427:     Level: beginner

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

2431: M*/

2433: /*@C 
2434:   VecLoadIntoVector - Loads a vector that has been stored in binary format
2435:   with VecView().

2437:   Collective on PetscViewer 

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

2443:   Level: intermediate

2445:   Notes:
2446:   The input file must contain the full global vector, as
2447:   written by the routine VecView().

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

2451:   Notes for advanced users:
2452:   Most users should not need to know the details of the binary storage
2453:   format, since VecLoad() and VecView() completely hide these details.
2454:   But for anyone who's interested, the standard binary matrix storage
2455:   format is
2456: .vb
2457:      int    VEC_FILE_COOKIE
2458:      int    number of rows
2459:      PetscScalar *values of all nonzeros
2460: .ve

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

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

2475:    Concepts: vector^loading from file

2477: .seealso: PetscViewerBinaryOpen(), VecView(), MatLoad(), VecLoad() 
2478: @*/
2479: int VecLoadIntoVector(PetscViewer viewer,Vec vec)
2480: {

2487:   if (!vec->ops->loadintovector) {
2488:     SETERRQ(1,"Vector does not support load");
2489:   }
2490:   (*vec->ops->loadintovector)(viewer,vec);
2491:   return(0);
2492: }

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

2497:    Collective on Vec

2499:    Input Parameter:
2500: .  v - the vector 

2502:    Output Parameter:
2503: .  v - the vector reciprocal

2505:    Level: intermediate

2507:    Concepts: vector^reciprocal

2509: @*/
2510: int VecReciprocal(Vec vec)
2511: {
2512:   int    ierr;

2517:   if (!vec->ops->reciprocal) {
2518:     SETERRQ(1,"Vector does not support reciprocal operation");
2519:   }
2520:   (*vec->ops->reciprocal)(vec);
2521:   return(0);
2522: }

2524: int VecSetOperation(Vec vec,VecOperation op, void (*f)(void))
2525: {
2528:   /* save the native version of the viewer */
2529:   if (op == VECOP_VIEW && !vec->ops->viewnative) {
2530:     vec->ops->viewnative = vec->ops->view;
2531:   }
2532:   (((void(**)(void))vec->ops)[(int)op]) = f;
2533:   return(0);
2534: }

2536: /*@
2537:    VecSetStashInitialSize - sets the sizes of the vec-stash, that is
2538:    used during the assembly process to store values that belong to 
2539:    other processors.

2541:    Collective on Vec

2543:    Input Parameters:
2544: +  vec   - the vector
2545: .  size  - the initial size of the stash.
2546: -  bsize - the initial size of the block-stash(if used).

2548:    Options Database Keys:
2549: +   -vecstash_initial_size <size> or <size0,size1,...sizep-1>
2550: -   -vecstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1>

2552:    Level: intermediate

2554:    Notes: 
2555:      The block-stash is used for values set with VecSetValuesBlocked() while
2556:      the stash is used for values set with VecSetValues()

2558:      Run with the option -log_info and look for output of the form
2559:      VecAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
2560:      to determine the appropriate value, MM, to use for size and 
2561:      VecAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
2562:      to determine the value, BMM to use for bsize

2564:    Concepts: vector^stash
2565:    Concepts: stash^vector

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

2569: @*/
2570: int VecSetStashInitialSize(Vec vec,int size,int bsize)
2571: {

2576:   VecStashSetInitialSize_Private(&vec->stash,size);
2577:   VecStashSetInitialSize_Private(&vec->bstash,bsize);
2578:   return(0);
2579: }

2581: /*@
2582:    VecStashView - Prints the entries in the vector stash and block stash.

2584:    Collective on Vec

2586:    Input Parameters:
2587: +  vec   - the vector
2588: -  viewer - the viewer

2590:    Level: advanced

2592:    Concepts: vector^stash
2593:    Concepts: stash^vector

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

2597: @*/
2598: int VecStashView(Vec v,PetscViewer viewer)
2599: {
2600:   int          ierr,rank,i,j;
2601:   PetscTruth   match;
2602:   VecStash     *s;
2603:   PetscScalar  val;


2610:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&match);
2611:   if (!match) SETERRQ1(1,"Stash viewer only works with ASCII viewer not %sn",((PetscObject)v)->type_name);
2612:   PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
2613:   MPI_Comm_rank(v->comm,&rank);
2614:   s = &v->bstash;

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

2632:   s = &v->stash;

2634:   /* print basic stash */
2635:   PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector stash size %dn",rank,s->n);
2636:   for (i=0; i<s->n; i++) {
2637:     val = s->array[i];
2638: #if defined(PETSC_USE_COMPLEX)
2639:       PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d (%18.16e %18.16e) ",rank,s->idx[i],PetscRealPart(val),PetscImaginaryPart(val));
2640: #else
2641:     PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d %18.16en",rank,s->idx[i],val);
2642: #endif
2643:   }
2644:   PetscViewerFlush(viewer);

2646:   PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
2647:   return(0);
2648: }

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

2655:    Not Collective

2657:    Input Parameter:
2658: +  x - the vector
2659: .  m - first dimension of two dimensional array
2660: .  n - second dimension of two dimensional array
2661: .  mstart - first index you will use in first coordinate direction (often 0)
2662: -  nstart - first index in the second coordinate direction (often 0)

2664:    Output Parameter:
2665: .  a - location to put pointer to the array

2667:    Level: beginner

2669:   Notes:
2670:    For a vector obtained from DACreateLocalVector() mstart and nstart are likely
2671:    obtained from the corner indices obtained from DAGetGhostCorners() while for
2672:    DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
2673:    the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray2d().
2674:    
2675:    For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

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

2679: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2680:           VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2681:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2682: @*/
2683: int VecGetArray2d(Vec x,int m,int n,int mstart,int nstart,PetscScalar **a[])
2684: {
2685:   int         i,ierr,N;
2686:   PetscScalar *aa;

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

2696:   PetscMalloc(m*sizeof(PetscScalar*),a);
2697:   for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
2698:   *a -= mstart;
2699:   return(0);
2700: }

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

2705:    Not Collective

2707:    Input Parameters:
2708: +  x - the vector
2709: .  m - first dimension of two dimensional array
2710: .  n - second dimension of the two dimensional array
2711: .  mstart - first index you will use in first coordinate direction (often 0)
2712: .  nstart - first index in the second coordinate direction (often 0)
2713: -  a - location of pointer to array obtained from VecGetArray2d()

2715:    Level: beginner

2717:    Notes:
2718:    For regular PETSc vectors this routine does not involve any copies. For
2719:    any special vectors that do not store local vector data in a contiguous
2720:    array, this routine will copy the data back into the underlying 
2721:    vector data structure from the array obtained with VecGetArray().

2723:    This routine actually zeros out the a pointer. 

2725: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2726:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2727:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2728: @*/
2729: int VecRestoreArray2d(Vec x,int m,int n,int mstart,int nstart,PetscScalar **a[])
2730: {

2736:   PetscFree(*a + mstart);
2737:   VecRestoreArray(x,PETSC_NULL);
2738:   return(0);
2739: }

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

2746:    Not Collective

2748:    Input Parameter:
2749: +  x - the vector
2750: .  m - first dimension of two dimensional array
2751: -  mstart - first index you will use in first coordinate direction (often 0)

2753:    Output Parameter:
2754: .  a - location to put pointer to the array

2756:    Level: beginner

2758:   Notes:
2759:    For a vector obtained from DACreateLocalVector() mstart are likely
2760:    obtained from the corner indices obtained from DAGetGhostCorners() while for
2761:    DACreateGlobalVector() they are the corner indices from DAGetCorners(). 
2762:    
2763:    For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2765: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2766:           VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2767:           VecGetarray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2768: @*/
2769: int VecGetArray1d(Vec x,int m,int mstart,PetscScalar *a[])
2770: {
2771:   int ierr,N;

2777:   VecGetLocalSize(x,&N);
2778:   if (m != N) SETERRQ2(1,"Local array size %d does not match 1d array dimensions %d",N,m);
2779:   VecGetArray(x,a);
2780:   *a  -= mstart;
2781:   return(0);
2782: }

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

2787:    Not Collective

2789:    Input Parameters:
2790: +  x - the vector
2791: .  m - first dimension of two dimensional array
2792: .  mstart - first index you will use in first coordinate direction (often 0)
2793: -  a - location of pointer to array obtained from VecGetArray21()

2795:    Level: beginner

2797:    Notes:
2798:    For regular PETSc vectors this routine does not involve any copies. For
2799:    any special vectors that do not store local vector data in a contiguous
2800:    array, this routine will copy the data back into the underlying 
2801:    vector data structure from the array obtained with VecGetArray1d().

2803:    This routine actually zeros out the a pointer. 

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

2807: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2808:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2809:           VecGetarray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2810: @*/
2811: int VecRestoreArray1d(Vec x,int m,int mstart,PetscScalar *a[])
2812: {

2818:   VecRestoreArray(x,PETSC_NULL);
2819:   return(0);
2820: }

2822: /*@C
2823:    VecConjugate - Conjugates a vector.

2825:    Collective on Vec

2827:    Input Parameters:
2828: .  x - the vector

2830:    Level: intermediate

2832:    Concepts: vector^conjugate

2834: @*/
2835: int VecConjugate(Vec x)
2836: {
2837: #ifdef PETSC_USE_COMPLEX

2843:   (*x->ops->conjugate)(x);
2844:   return(0);
2845: #else
2846:   return(0);
2847: #endif
2848: }

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

2855:    Not Collective

2857:    Input Parameter:
2858: +  x - the vector
2859: .  m - first dimension of three dimensional array
2860: .  n - second dimension of three dimensional array
2861: .  p - third dimension of three dimensional array
2862: .  mstart - first index you will use in first coordinate direction (often 0)
2863: .  nstart - first index in the second coordinate direction (often 0)
2864: -  pstart - first index in the third coordinate direction (often 0)

2866:    Output Parameter:
2867: .  a - location to put pointer to the array

2869:    Level: beginner

2871:   Notes:
2872:    For a vector obtained from DACreateLocalVector() mstart, nstart, and pstart are likely
2873:    obtained from the corner indices obtained from DAGetGhostCorners() while for
2874:    DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
2875:    the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray3d().
2876:    
2877:    For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

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

2881: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2882:           VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2883:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2884: @*/
2885: int VecGetArray3d(Vec x,int m,int n,int p,int mstart,int nstart,int pstart,PetscScalar ***a[])
2886: {
2887:   int         i,ierr,N,j;
2888:   PetscScalar *aa,**b;

2894:   VecGetLocalSize(x,&N);
2895:   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);
2896:   VecGetArray(x,&aa);

2898:   PetscMalloc(m*sizeof(PetscScalar**)+m*n*sizeof(PetscScalar*),a);
2899:   b    = (PetscScalar **)((*a) + m);
2900:   for (i=0; i<m; i++)   (*a)[i] = b + i*n - nstart;
2901:   for (i=0; i<m; i++) {
2902:     for (j=0; j<n; j++) {
2903:       b[i*n+j] = aa + i*n*p + j*p - pstart;
2904:     }
2905:   }
2906:   *a -= mstart;
2907:   return(0);
2908: }

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

2913:    Not Collective

2915:    Input Parameters:
2916: +  x - the vector
2917: .  m - first dimension of three dimensional array
2918: .  n - second dimension of the three dimensional array
2919: .  p - third dimension of the three dimensional array
2920: .  mstart - first index you will use in first coordinate direction (often 0)
2921: .  nstart - first index in the second coordinate direction (often 0)
2922: .  pstart - first index in the third coordinate direction (often 0)
2923: -  a - location of pointer to array obtained from VecGetArray3d()

2925:    Level: beginner

2927:    Notes:
2928:    For regular PETSc vectors this routine does not involve any copies. For
2929:    any special vectors that do not store local vector data in a contiguous
2930:    array, this routine will copy the data back into the underlying 
2931:    vector data structure from the array obtained with VecGetArray().

2933:    This routine actually zeros out the a pointer. 

2935: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2936:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2937:           VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2938: @*/
2939: int VecRestoreArray3d(Vec x,int m,int n,int p,int mstart,int nstart,int pstart,PetscScalar ***a[])
2940: {

2946:   PetscFree(*a + mstart);
2947:   VecRestoreArray(x,PETSC_NULL);
2948:   return(0);
2949: }