Actual source code: bvec3.c

petsc-master 2019-08-15
Report Typos and Errors

  2: /*
  3:    Implements the sequential vectors.
  4: */

  6:  #include <../src/vec/vec/impls/dvecimpl.h>
  7: /*MC
  8:    VECSEQ - VECSEQ = "seq" - The basic sequential vector

 10:    Options Database Keys:
 11: . -vec_type seq - sets the vector type to VECSEQ during a call to VecSetFromOptions()

 13:   Level: beginner

 15: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateSeqWithArray(), VECMPI, VecType, VecCreateMPI(), VecCreateSeq()
 16: M*/

 18: #if defined(PETSC_USE_MIXED_PRECISION)
 19: extern PetscErrorCode VecCreate_Seq_Private(Vec,const float*);
 20: extern PetscErrorCode VecCreate_Seq_Private(Vec,const double*);
 21: #endif

 23: PETSC_EXTERN PetscErrorCode VecCreate_Seq(Vec V)
 24: {
 25:   Vec_Seq        *s;
 26:   PetscScalar    *array;
 28:   PetscInt       n = PetscMax(V->map->n,V->map->N);
 29:   PetscMPIInt    size;

 32:   MPI_Comm_size(PetscObjectComm((PetscObject)V),&size);
 33:   if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
 34: #if !defined(PETSC_USE_MIXED_PRECISION)
 35:   PetscMalloc1(n,&array);
 36:   PetscLogObjectMemory((PetscObject)V, n*sizeof(PetscScalar));
 37:   VecCreate_Seq_Private(V,array);

 39:   s                  = (Vec_Seq*)V->data;
 40:   s->array_allocated = array;

 42:   VecSet(V,0.0);
 43: #else
 44:   switch (((PetscObject)V)->precision) {
 45:   case PETSC_PRECISION_SINGLE: {
 46:     float *aarray;

 48:     PetscCalloc1(n,&aarray);
 49:     PetscLogObjectMemory((PetscObject)V, n*sizeof(float));
 50:     VecCreate_Seq_Private(V,aarray);

 52:     s                  = (Vec_Seq*)V->data;
 53:     s->array_allocated = (PetscScalar*)aarray;
 54:   } break;
 55:   case PETSC_PRECISION_DOUBLE: {
 56:     double *aarray;

 58:     PetscCalloc1(n,&aarray);
 59:     PetscLogObjectMemory((PetscObject)V, n*sizeof(double));
 60:     VecCreate_Seq_Private(V,aarray);

 62:     s                  = (Vec_Seq*)V->data;
 63:     s->array_allocated = (PetscScalar*)aarray;
 64:   } break;
 65:   default: SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"No support for mixed precision %d",(int)(((PetscObject)V)->precision));
 66:   }
 67: #endif
 68:   return(0);
 69: }