Actual source code: bvec3.c
petsc-3.7.7 2017-09-25
2: /*
3: Implements the sequential vectors.
4: */
6: #include <../src/vec/vec/impls/dvecimpl.h> /*I "petscvec.h" I*/
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
25: PETSC_EXTERN PetscErrorCode VecCreate_Seq(Vec V)
26: {
27: Vec_Seq *s;
28: PetscScalar *array;
30: PetscInt n = PetscMax(V->map->n,V->map->N);
31: PetscMPIInt size;
34: MPI_Comm_size(PetscObjectComm((PetscObject)V),&size);
35: if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot create VECSEQ on more than one process");
36: #if !defined(PETSC_USE_MIXED_PRECISION)
37: PetscMalloc1(n,&array);
38: PetscLogObjectMemory((PetscObject)V, n*sizeof(PetscScalar));
39: VecCreate_Seq_Private(V,array);
41: s = (Vec_Seq*)V->data;
42: s->array_allocated = array;
44: VecSet(V,0.0);
45: #else
46: switch (((PetscObject)V)->precision) {
47: case PETSC_PRECISION_SINGLE: {
48: float *aarray;
50: PetscMalloc1(n,&aarray);
51: PetscLogObjectMemory((PetscObject)V, n*sizeof(float));
52: PetscMemzero(aarray,n*sizeof(float));
53: VecCreate_Seq_Private(V,aarray);
55: s = (Vec_Seq*)V->data;
56: s->array_allocated = (PetscScalar*)aarray;
57: } break;
58: case PETSC_PRECISION_DOUBLE: {
59: double *aarray;
61: PetscMalloc1(n,&aarray);
62: PetscLogObjectMemory((PetscObject)V, n*sizeof(double));
63: PetscMemzero(aarray,n*sizeof(double));
64: VecCreate_Seq_Private(V,aarray);
66: s = (Vec_Seq*)V->data;
67: s->array_allocated = (PetscScalar*)aarray;
68: } break;
69: default: SETERRQ1(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"No support for mixed precision %d",(int)(((PetscObject)V)->precision));
70: }
71: #endif
72: return(0);
73: }