Actual source code: slepcbv.h
slepc-3.8.2 2017-12-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: User interface for the basis vectors object in SLEPc
12: */
16: #include <slepcsys.h>
18: PETSC_EXTERN PetscErrorCode BVInitializePackage(void);
20: /*S
21: BV - Basis vectors, SLEPc object representing a collection of vectors
22: that typically constitute a basis of a subspace.
24: Level: beginner
26: .seealso: BVCreate()
27: S*/
28: typedef struct _p_BV* BV;
30: /*J
31: BVType - String with the name of the type of BV. Each type differs in
32: the way data is stored internally.
34: Level: beginner
36: .seealso: BVSetType(), BV
37: J*/
38: typedef const char* BVType;
39: #define BVMAT "mat"
40: #define BVSVEC "svec"
41: #define BVVECS "vecs"
42: #define BVCONTIGUOUS "contiguous"
44: /* Logging support */
45: PETSC_EXTERN PetscClassId BV_CLASSID;
47: /*E
48: BVOrthogType - Determines the method used in the orthogonalization
49: of vectors
51: Level: advanced
53: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn(), BVOrthogRefineType
54: E*/
55: typedef enum { BV_ORTHOG_CGS,
56: BV_ORTHOG_MGS } BVOrthogType;
57: PETSC_EXTERN const char *BVOrthogTypes[];
59: /*E
60: BVOrthogRefineType - Determines what type of refinement to use
61: during orthogonalization of vectors
63: Level: advanced
65: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
66: E*/
67: typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
68: BV_ORTHOG_REFINE_NEVER,
69: BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
70: PETSC_EXTERN const char *BVOrthogRefineTypes[];
72: /*E
73: BVOrthogBlockType - Determines the method used in block
74: orthogonalization (simultaneous orthogonalization of a set of vectors)
76: Level: advanced
78: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalize()
79: E*/
80: typedef enum { BV_ORTHOG_BLOCK_GS,
81: BV_ORTHOG_BLOCK_CHOL,
82: BV_ORTHOG_BLOCK_TSQR } BVOrthogBlockType;
83: PETSC_EXTERN const char *BVOrthogBlockTypes[];
85: /*E
86: BVMatMultType - Determines how to perform the BVMatMult() operation:
87: BV_MATMULT_VECS: perform a matrix-vector multiply per each column;
88: BV_MATMULT_MAT: carry out a MatMatMult() product with a dense matrix (default);
89: BV_MATMULT_MAT_SAVE: call MatMatMult() and keep auxiliary matrices
90: (more efficient but needs more memory)
92: Level: advanced
94: .seealso: BVMatMult()
95: E*/
96: typedef enum { BV_MATMULT_VECS,
97: BV_MATMULT_MAT,
98: BV_MATMULT_MAT_SAVE } BVMatMultType;
99: PETSC_EXTERN const char *BVMatMultTypes[];
101: PETSC_EXTERN PetscErrorCode BVCreate(MPI_Comm,BV*);
102: PETSC_EXTERN PetscErrorCode BVDestroy(BV*);
103: PETSC_EXTERN PetscErrorCode BVSetType(BV,BVType);
104: PETSC_EXTERN PetscErrorCode BVGetType(BV,BVType*);
105: PETSC_EXTERN PetscErrorCode BVSetSizes(BV,PetscInt,PetscInt,PetscInt);
106: PETSC_EXTERN PetscErrorCode BVSetSizesFromVec(BV,Vec,PetscInt);
107: PETSC_EXTERN PetscErrorCode BVGetSizes(BV,PetscInt*,PetscInt*,PetscInt*);
108: PETSC_EXTERN PetscErrorCode BVResize(BV,PetscInt,PetscBool);
109: PETSC_EXTERN PetscErrorCode BVSetFromOptions(BV);
110: PETSC_EXTERN PetscErrorCode BVView(BV,PetscViewer);
111: PETSC_STATIC_INLINE PetscErrorCode BVViewFromOptions(BV bv,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)bv,obj,name);}
113: PETSC_EXTERN PetscErrorCode BVGetColumn(BV,PetscInt,Vec*);
114: PETSC_EXTERN PetscErrorCode BVRestoreColumn(BV,PetscInt,Vec*);
115: PETSC_EXTERN PetscErrorCode BVGetArray(BV,PetscScalar**);
116: PETSC_EXTERN PetscErrorCode BVRestoreArray(BV,PetscScalar**);
117: PETSC_EXTERN PetscErrorCode BVGetArrayRead(BV,const PetscScalar**);
118: PETSC_EXTERN PetscErrorCode BVRestoreArrayRead(BV,const 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);
139: PETSC_EXTERN PetscErrorCode BVSetBufferVec(BV,Vec);
140: PETSC_EXTERN PetscErrorCode BVGetBufferVec(BV,Vec*);
142: PETSC_EXTERN PetscErrorCode BVMult(BV,PetscScalar,PetscScalar,BV,Mat);
143: PETSC_EXTERN PetscErrorCode BVMultVec(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
144: PETSC_EXTERN PetscErrorCode BVMultColumn(BV,PetscScalar,PetscScalar,PetscInt,PetscScalar*);
145: PETSC_EXTERN PetscErrorCode BVMultInPlace(BV,Mat,PetscInt,PetscInt);
146: PETSC_EXTERN PetscErrorCode BVMultInPlaceTranspose(BV,Mat,PetscInt,PetscInt);
147: PETSC_EXTERN PetscErrorCode BVMatMult(BV,Mat,BV);
148: PETSC_EXTERN PetscErrorCode BVMatMultHermitianTranspose(BV,Mat,BV);
149: PETSC_EXTERN PetscErrorCode BVMatMultColumn(BV,Mat,PetscInt);
150: PETSC_EXTERN PetscErrorCode BVMatProject(BV,Mat,BV,Mat);
151: PETSC_EXTERN PetscErrorCode BVDot(BV,BV,Mat);
152: PETSC_EXTERN PetscErrorCode BVDotVec(BV,Vec,PetscScalar*);
153: PETSC_EXTERN PetscErrorCode BVDotVecBegin(BV,Vec,PetscScalar*);
154: PETSC_EXTERN PetscErrorCode BVDotVecEnd(BV,Vec,PetscScalar*);
155: PETSC_EXTERN PetscErrorCode BVDotColumn(BV,PetscInt,PetscScalar*);
156: PETSC_EXTERN PetscErrorCode BVDotColumnBegin(BV,PetscInt,PetscScalar*);
157: PETSC_EXTERN PetscErrorCode BVDotColumnEnd(BV,PetscInt,PetscScalar*);
158: PETSC_EXTERN PetscErrorCode BVScale(BV,PetscScalar);
159: PETSC_EXTERN PetscErrorCode BVScaleColumn(BV,PetscInt,PetscScalar);
160: PETSC_EXTERN PetscErrorCode BVNorm(BV,NormType,PetscReal*);
161: PETSC_EXTERN PetscErrorCode BVNormVec(BV,Vec,NormType,PetscReal*);
162: PETSC_EXTERN PetscErrorCode BVNormVecBegin(BV,Vec,NormType,PetscReal*);
163: PETSC_EXTERN PetscErrorCode BVNormVecEnd(BV,Vec,NormType,PetscReal*);
164: PETSC_EXTERN PetscErrorCode BVNormColumn(BV,PetscInt,NormType,PetscReal*);
165: PETSC_EXTERN PetscErrorCode BVNormColumnBegin(BV,PetscInt,NormType,PetscReal*);
166: PETSC_EXTERN PetscErrorCode BVNormColumnEnd(BV,PetscInt,NormType,PetscReal*);
167: PETSC_EXTERN PetscErrorCode BVSetRandom(BV);
168: PETSC_EXTERN PetscErrorCode BVSetRandomColumn(BV,PetscInt);
169: PETSC_EXTERN PetscErrorCode BVSetRandomCond(BV,PetscReal);
170: PETSC_EXTERN PetscErrorCode BVSetRandomContext(BV,PetscRandom);
171: PETSC_EXTERN PetscErrorCode BVGetRandomContext(BV,PetscRandom*);
173: PETSC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal,BVOrthogBlockType);
174: PETSC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*,BVOrthogBlockType*);
175: PETSC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
176: PETSC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar*,PetscReal*,PetscBool*);
177: PETSC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar*,PetscReal*,PetscBool*);
178: PETSC_EXTERN PetscErrorCode BVOrthonormalizeColumn(BV,PetscInt,PetscBool,PetscReal*,PetscBool*);
179: PETSC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar*,PetscReal*,PetscBool*);
180: PETSC_EXTERN PetscErrorCode BVSetMatMultMethod(BV,BVMatMultType);
181: PETSC_EXTERN PetscErrorCode BVGetMatMultMethod(BV,BVMatMultType*);
183: PETSC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char*);
184: PETSC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char*);
185: PETSC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);
187: PETSC_EXTERN PetscFunctionList BVList;
188: PETSC_EXTERN PetscErrorCode BVRegister(const char[],PetscErrorCode(*)(BV));
190: PETSC_EXTERN PetscErrorCode BVCreateFromMat(Mat,BV*);
191: PETSC_EXTERN PetscErrorCode BVCreateMat(BV,Mat*);
193: #endif