Actual source code: dvecimpl.h
1: /*
2: This should not be included in users code.
4: Includes definition of structure for seqential double precision vectors
6: These are shared by dvec1.c dvec2.c dvec3.c bvec1.c bvec2.c
7: pvectors/pvec.c pvectors/pbvec.c
8: */
10: #ifndef __DVECIMPL
13: #include vecimpl.h
15: typedef struct {
16: VECHEADER
17: } Vec_Seq;
19: EXTERN PetscErrorCode VecMDot_Seq(PetscInt,Vec,const Vec[],PetscScalar *);
20: EXTERN PetscErrorCode VecMTDot_Seq(PetscInt,Vec,const Vec[],PetscScalar *);
21: EXTERN PetscErrorCode VecMin_Seq(Vec,PetscInt*,PetscReal *);
22: EXTERN PetscErrorCode VecSet_Seq(const PetscScalar*,Vec);
23: EXTERN PetscErrorCode VecSetRandom_Seq(PetscRandom,Vec);
24: EXTERN PetscErrorCode VecMAXPY_Seq(PetscInt,const PetscScalar *,Vec,Vec *);
25: EXTERN PetscErrorCode VecAYPX_Seq(const PetscScalar *,Vec,Vec);
26: EXTERN PetscErrorCode VecWAXPY_Seq(const PetscScalar*,Vec,Vec,Vec);
27: EXTERN PetscErrorCode VecPointwiseMult_Seq(Vec,Vec,Vec);
28: EXTERN PetscErrorCode VecPointwiseDivide_Seq(Vec,Vec,Vec);
29: EXTERN PetscErrorCode VecMaxPointwiseDivide_Seq(Vec,Vec,PetscReal*);
30: EXTERN PetscErrorCode VecGetArray_Seq(Vec,PetscScalar *[]);
31: EXTERN PetscErrorCode VecRestoreArray_Seq(Vec,PetscScalar *[]);
32: EXTERN PetscErrorCode VecPlaceArray_Seq(Vec,const PetscScalar *);
33: EXTERN PetscErrorCode VecResetArray_Seq(Vec);
34: EXTERN PetscErrorCode VecReplaceArray_Seq(Vec,const PetscScalar *);
35: EXTERN PetscErrorCode VecGetSize_Seq(Vec,PetscInt *);
36: EXTERN PetscErrorCode VecDot_Seq(Vec,Vec,PetscScalar *);
37: EXTERN PetscErrorCode VecTDot_Seq(Vec,Vec,PetscScalar *);
38: EXTERN PetscErrorCode VecScale_Seq(const PetscScalar *,Vec);
39: EXTERN PetscErrorCode VecCopy_Seq(Vec,Vec);
40: EXTERN PetscErrorCode VecSwap_Seq(Vec,Vec);
41: EXTERN PetscErrorCode VecAXPY_Seq(const PetscScalar *,Vec,Vec);
42: EXTERN PetscErrorCode VecAXPBY_Seq(const PetscScalar *,const PetscScalar *,Vec,Vec);
43: EXTERN PetscErrorCode VecMax_Seq(Vec,PetscInt*,PetscReal *);
44: EXTERN PetscErrorCode VecDuplicate_Seq(Vec,Vec *);
45: EXTERN PetscErrorCode VecConjugate_Seq(Vec);
46: EXTERN PetscErrorCode VecNorm_Seq(Vec,NormType,PetscReal*);
47: #endif