Actual source code: slepcbv.h
slepc-3.6.3 2016-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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 the method used in the orthogonalization
57: of vectors
59: Level: advanced
61: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn(), BVOrthogRefineType
62: E*/
63: typedef enum { BV_ORTHOG_CGS,
64: BV_ORTHOG_MGS } BVOrthogType;
66: /*E
67: BVOrthogRefineType - Determines what type of refinement to use
68: during orthogonalization of vectors
70: Level: advanced
72: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
73: E*/
74: typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
75: BV_ORTHOG_REFINE_NEVER,
76: BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
78: /*E
79: BVOrthogBlockType - Determines the method used in block
80: orthogonalization (simultaneous orthogonalization of a set of vectors)
82: Level: advanced
84: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalize()
85: E*/
86: typedef enum { BV_ORTHOG_BLOCK_GS,
87: BV_ORTHOG_BLOCK_CHOL } BVOrthogBlockType;
89: /*E
90: BVMatMultType - Determines how to perform the BVMatMult() operation:
91: BV_MATMULT_VECS: perform a matrix-vector multiply per each column;
92: BV_MATMULT_MAT: carry out a MatMatMult() product with a dense matrix (default);
93: BV_MATMULT_MAT_SAVE: call MatMatMult() and keep auxiliary matrices
94: (more efficient but needs more memory)
96: Level: advanced
98: .seealso: BVMatMult()
99: E*/
100: typedef enum { BV_MATMULT_VECS,
101: BV_MATMULT_MAT,
102: BV_MATMULT_MAT_SAVE } BVMatMultType;
104: PETSC_EXTERN PetscErrorCode BVCreate(MPI_Comm,BV*);
105: PETSC_EXTERN PetscErrorCode BVDestroy(BV*);
106: PETSC_EXTERN PetscErrorCode BVSetType(BV,BVType);
107: PETSC_EXTERN PetscErrorCode BVGetType(BV,BVType*);
108: PETSC_EXTERN PetscErrorCode BVSetSizes(BV,PetscInt,PetscInt,PetscInt);
109: PETSC_EXTERN PetscErrorCode BVSetSizesFromVec(BV,Vec,PetscInt);
110: PETSC_EXTERN PetscErrorCode BVGetSizes(BV,PetscInt*,PetscInt*,PetscInt*);
111: PETSC_EXTERN PetscErrorCode BVResize(BV,PetscInt,PetscBool);
112: PETSC_EXTERN PetscErrorCode BVSetFromOptions(BV);
113: PETSC_EXTERN PetscErrorCode BVView(BV,PetscViewer);
115: PETSC_EXTERN PetscErrorCode BVGetColumn(BV,PetscInt,Vec*);
116: PETSC_EXTERN PetscErrorCode BVRestoreColumn(BV,PetscInt,Vec*);
117: PETSC_EXTERN PetscErrorCode BVGetArray(BV,PetscScalar**);
118: PETSC_EXTERN PetscErrorCode BVRestoreArray(BV,PetscScalar**);
119: PETSC_EXTERN PetscErrorCode BVCreateVec(BV,Vec*);
120: PETSC_EXTERN PetscErrorCode BVSetActiveColumns(BV,PetscInt,PetscInt);
121: PETSC_EXTERN PetscErrorCode BVGetActiveColumns(BV,PetscInt*,PetscInt*);
122: PETSC_EXTERN PetscErrorCode BVInsertVec(BV,PetscInt,Vec);
123: PETSC_EXTERN PetscErrorCode BVInsertVecs(BV,PetscInt,PetscInt*,Vec*,PetscBool);
124: PETSC_EXTERN PetscErrorCode BVInsertConstraints(BV,PetscInt*,Vec*);
125: PETSC_EXTERN PetscErrorCode BVSetNumConstraints(BV,PetscInt);
126: PETSC_EXTERN PetscErrorCode BVGetNumConstraints(BV,PetscInt*);
127: PETSC_EXTERN PetscErrorCode BVDuplicate(BV,BV*);
128: PETSC_EXTERN PetscErrorCode BVDuplicateResize(BV,PetscInt,BV*);
129: PETSC_EXTERN PetscErrorCode BVCopy(BV,BV);
130: PETSC_EXTERN PetscErrorCode BVCopyVec(BV,PetscInt,Vec);
131: PETSC_EXTERN PetscErrorCode BVCopyColumn(BV,PetscInt,PetscInt);
132: PETSC_EXTERN PetscErrorCode BVSetMatrix(BV,Mat,PetscBool);
133: PETSC_EXTERN PetscErrorCode BVGetMatrix(BV,Mat*,PetscBool*);
134: PETSC_EXTERN PetscErrorCode BVApplyMatrix(BV,Vec,Vec);
135: PETSC_EXTERN PetscErrorCode BVApplyMatrixBV(BV,BV);
136: PETSC_EXTERN PetscErrorCode BVGetCachedBV(BV,BV*);
137: PETSC_EXTERN PetscErrorCode BVSetSignature(BV,Vec);
138: PETSC_EXTERN PetscErrorCode BVGetSignature(BV,Vec);
140: PETSC_EXTERN PetscErrorCode BVMult(BV,PetscScalar,PetscScalar,BV,Mat);
141: PETSC_EXTERN PetscErrorCode BVMultVec(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
142: PETSC_EXTERN PetscErrorCode BVMultColumn(BV,PetscScalar,PetscScalar,PetscInt,PetscScalar*);
143: PETSC_EXTERN PetscErrorCode BVMultInPlace(BV,Mat,PetscInt,PetscInt);
144: PETSC_EXTERN PetscErrorCode BVMultInPlaceTranspose(BV,Mat,PetscInt,PetscInt);
145: PETSC_EXTERN PetscErrorCode BVMatMult(BV,Mat,BV);
146: PETSC_EXTERN PetscErrorCode BVMatMultHermitianTranspose(BV,Mat,BV);
147: PETSC_EXTERN PetscErrorCode BVMatMultColumn(BV,Mat,PetscInt);
148: PETSC_EXTERN PetscErrorCode BVMatProject(BV,Mat,BV,Mat);
149: PETSC_EXTERN PetscErrorCode BVAXPY(BV,PetscScalar,BV);
150: PETSC_EXTERN PetscErrorCode BVDot(BV,BV,Mat);
151: PETSC_EXTERN PetscErrorCode BVDotVec(BV,Vec,PetscScalar*);
152: PETSC_EXTERN PetscErrorCode BVDotVecBegin(BV,Vec,PetscScalar*);
153: PETSC_EXTERN PetscErrorCode BVDotVecEnd(BV,Vec,PetscScalar*);
154: PETSC_EXTERN PetscErrorCode BVDotColumn(BV,PetscInt,PetscScalar*);
155: PETSC_EXTERN PetscErrorCode BVDotColumnBegin(BV,PetscInt,PetscScalar*);
156: PETSC_EXTERN PetscErrorCode BVDotColumnEnd(BV,PetscInt,PetscScalar*);
157: PETSC_EXTERN PetscErrorCode BVScale(BV,PetscScalar);
158: PETSC_EXTERN PetscErrorCode BVScaleColumn(BV,PetscInt,PetscScalar);
159: PETSC_EXTERN PetscErrorCode BVNorm(BV,NormType,PetscReal*);
160: PETSC_EXTERN PetscErrorCode BVNormVec(BV,Vec,NormType,PetscReal*);
161: PETSC_EXTERN PetscErrorCode BVNormVecBegin(BV,Vec,NormType,PetscReal*);
162: PETSC_EXTERN PetscErrorCode BVNormVecEnd(BV,Vec,NormType,PetscReal*);
163: PETSC_EXTERN PetscErrorCode BVNormColumn(BV,PetscInt,NormType,PetscReal*);
164: PETSC_EXTERN PetscErrorCode BVNormColumnBegin(BV,PetscInt,NormType,PetscReal*);
165: PETSC_EXTERN PetscErrorCode BVNormColumnEnd(BV,PetscInt,NormType,PetscReal*);
166: PETSC_EXTERN PetscErrorCode BVSetRandom(BV,PetscRandom);
167: PETSC_EXTERN PetscErrorCode BVSetRandomColumn(BV,PetscInt,PetscRandom);
169: PETSC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal,BVOrthogBlockType);
170: PETSC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*,BVOrthogBlockType*);
171: PETSC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
172: PETSC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar*,PetscReal*,PetscBool*);
173: PETSC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar*,PetscReal*,PetscBool*);
174: PETSC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar*,PetscReal*,PetscBool*);
175: PETSC_EXTERN PetscErrorCode BVSetMatMultMethod(BV,BVMatMultType);
176: PETSC_EXTERN PetscErrorCode BVGetMatMultMethod(BV,BVMatMultType*);
178: PETSC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char*);
179: PETSC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char*);
180: PETSC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);
182: PETSC_EXTERN PetscFunctionList BVList;
183: PETSC_EXTERN PetscErrorCode BVRegister(const char[],PetscErrorCode(*)(BV));
185: #endif