Actual source code: slepcbv.h
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
24: #include <slepcsys.h>
26: PETSC_EXTERN PetscErrorCode BVInitializePackage(void);
28: /*S
29: BV - Basis vectors, SLEPc object representing a collection of vectors
30: that typically constitute a basis of a subspace.
32: Level: beginner
34: .seealso: BVCreate()
35: S*/
36: typedef struct _p_BV* BV;
38: /*J
39: BVType - String with the name of the type of BV. Each type differs in
40: the way data is stored internally.
42: Level: beginner
44: .seealso: BVSetType(), BV
45: J*/
46: typedef const char* BVType;
47: #define BVMAT "mat"
48: #define BVSVEC "svec"
49: #define BVVECS "vecs"
50: #define BVCONTIGUOUS "contiguous"
52: /* Logging support */
53: PETSC_EXTERN PetscClassId BV_CLASSID;
55: /*E
56: BVOrthogType - Determines what type of orthogonalization to use
58: Level: advanced
60: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
61: E*/
62: typedef enum { BV_ORTHOG_CGS,
63: BV_ORTHOG_MGS } BVOrthogType;
65: /*E
66: BVOrthogRefineType - Determines what type of refinement
67: to use during orthogonalization
69: Level: advanced
71: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
72: E*/
73: typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
74: BV_ORTHOG_REFINE_NEVER,
75: BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
77: PETSC_EXTERN PetscErrorCode BVCreate(MPI_Comm,BV*);
78: PETSC_EXTERN PetscErrorCode BVDestroy(BV*);
79: PETSC_EXTERN PetscErrorCode BVSetType(BV,BVType);
80: PETSC_EXTERN PetscErrorCode BVGetType(BV,BVType*);
81: PETSC_EXTERN PetscErrorCode BVSetSizes(BV,PetscInt,PetscInt,PetscInt);
82: PETSC_EXTERN PetscErrorCode BVSetSizesFromVec(BV,Vec,PetscInt);
83: PETSC_EXTERN PetscErrorCode BVGetSizes(BV,PetscInt*,PetscInt*,PetscInt*);
84: PETSC_EXTERN PetscErrorCode BVResize(BV,PetscInt,PetscBool);
85: PETSC_EXTERN PetscErrorCode BVSetFromOptions(BV);
86: PETSC_EXTERN PetscErrorCode BVView(BV,PetscViewer);
88: PETSC_EXTERN PetscErrorCode BVGetColumn(BV,PetscInt,Vec*);
89: PETSC_EXTERN PetscErrorCode BVRestoreColumn(BV,PetscInt,Vec*);
90: PETSC_EXTERN PetscErrorCode BVGetArray(BV,PetscScalar**);
91: PETSC_EXTERN PetscErrorCode BVRestoreArray(BV,PetscScalar**);
92: PETSC_EXTERN PetscErrorCode BVGetVec(BV,Vec*);
93: PETSC_EXTERN PetscErrorCode BVSetActiveColumns(BV,PetscInt,PetscInt);
94: PETSC_EXTERN PetscErrorCode BVGetActiveColumns(BV,PetscInt*,PetscInt*);
95: PETSC_EXTERN PetscErrorCode BVInsertVec(BV,PetscInt,Vec);
96: PETSC_EXTERN PetscErrorCode BVInsertVecs(BV,PetscInt,PetscInt*,Vec*,PetscBool);
97: PETSC_EXTERN PetscErrorCode BVInsertConstraints(BV,PetscInt*,Vec*);
98: PETSC_EXTERN PetscErrorCode BVSetNumConstraints(BV,PetscInt);
99: PETSC_EXTERN PetscErrorCode BVGetNumConstraints(BV,PetscInt*);
100: PETSC_EXTERN PetscErrorCode BVDuplicate(BV,BV*);
101: PETSC_EXTERN PetscErrorCode BVDuplicateResize(BV,PetscInt,BV*);
102: PETSC_EXTERN PetscErrorCode BVCopy(BV,BV);
103: PETSC_EXTERN PetscErrorCode BVCopyVec(BV,PetscInt,Vec);
104: PETSC_EXTERN PetscErrorCode BVCopyColumn(BV,PetscInt,PetscInt);
105: PETSC_EXTERN PetscErrorCode BVSetMatrix(BV,Mat,PetscBool);
106: PETSC_EXTERN PetscErrorCode BVGetMatrix(BV,Mat*,PetscBool*);
107: PETSC_EXTERN PetscErrorCode BVApplyMatrix(BV,Vec,Vec);
108: PETSC_EXTERN PetscErrorCode BVSetSignature(BV,Vec);
109: PETSC_EXTERN PetscErrorCode BVGetSignature(BV,Vec);
111: PETSC_EXTERN PetscErrorCode BVMult(BV,PetscScalar,PetscScalar,BV,Mat);
112: PETSC_EXTERN PetscErrorCode BVMultVec(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
113: PETSC_EXTERN PetscErrorCode BVMultColumn(BV,PetscScalar,PetscScalar,PetscInt,PetscScalar*);
114: PETSC_EXTERN PetscErrorCode BVMultInPlace(BV,Mat,PetscInt,PetscInt);
115: PETSC_EXTERN PetscErrorCode BVMultInPlaceTranspose(BV,Mat,PetscInt,PetscInt);
116: PETSC_EXTERN PetscErrorCode BVMatMult(BV,Mat,BV);
117: PETSC_EXTERN PetscErrorCode BVMatMultColumn(BV,Mat,PetscInt);
118: PETSC_EXTERN PetscErrorCode BVMatProject(BV,Mat,BV,Mat);
119: PETSC_EXTERN PetscErrorCode BVAXPY(BV,PetscScalar,BV);
120: PETSC_EXTERN PetscErrorCode BVDot(BV,BV,Mat);
121: PETSC_EXTERN PetscErrorCode BVDotVec(BV,Vec,PetscScalar*);
122: PETSC_EXTERN PetscErrorCode BVDotColumn(BV,PetscInt,PetscScalar*);
123: PETSC_EXTERN PetscErrorCode BVScale(BV,PetscScalar);
124: PETSC_EXTERN PetscErrorCode BVScaleColumn(BV,PetscInt,PetscScalar);
125: PETSC_EXTERN PetscErrorCode BVNorm(BV,NormType,PetscReal*);
126: PETSC_EXTERN PetscErrorCode BVNormVec(BV,Vec,NormType,PetscReal*);
127: PETSC_EXTERN PetscErrorCode BVNormColumn(BV,PetscInt,NormType,PetscReal*);
128: PETSC_EXTERN PetscErrorCode BVSetRandom(BV,PetscRandom);
129: PETSC_EXTERN PetscErrorCode BVSetRandomColumn(BV,PetscInt,PetscRandom);
131: PETSC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal);
132: PETSC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*);
133: PETSC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
134: PETSC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar*,PetscReal*,PetscBool*);
135: PETSC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar*,PetscReal*,PetscBool*);
136: PETSC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar*,PetscReal*,PetscBool*);
138: PETSC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char*);
139: PETSC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char*);
140: PETSC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);
142: PETSC_EXTERN PetscFunctionList BVList;
143: PETSC_EXTERN PetscBool BVRegisterAllCalled;
144: PETSC_EXTERN PetscErrorCode BVRegisterAll(void);
145: PETSC_EXTERN PetscErrorCode BVRegister(const char[],PetscErrorCode(*)(BV));
147: #endif