Actual source code: slepcbv.h

slepc-main 2024-11-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: */

 14: #pragma once

 16: #include <slepcsys.h>

 18: /* SUBMANSEC = BV */

 20: SLEPC_EXTERN PetscErrorCode BVInitializePackage(void);
 21: SLEPC_EXTERN PetscErrorCode BVFinalizePackage(void);

 23: /*S
 24:     BV - Basis vectors, SLEPc object representing a collection of vectors
 25:     that typically constitute a basis of a subspace.

 27:     Level: beginner

 29: .seealso:  BVCreate()
 30: S*/
 31: typedef struct _p_BV* BV;

 33: /*J
 34:     BVType - String with the name of the type of BV. Each type differs in
 35:     the way data is stored internally.

 37:     Level: beginner

 39: .seealso: BVSetType(), BV
 40: J*/
 41: typedef const char* BVType;
 42: #define BVMAT        "mat"
 43: #define BVSVEC       "svec"
 44: #define BVVECS       "vecs"
 45: #define BVCONTIGUOUS "contiguous"
 46: #define BVTENSOR     "tensor"

 48: /* Logging support */
 49: SLEPC_EXTERN PetscClassId BV_CLASSID;

 51: /*E
 52:     BVOrthogType - Determines the method used in the orthogonalization
 53:     of vectors

 55:     Level: advanced

 57: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn(), BVOrthogRefineType
 58: E*/
 59: typedef enum { BV_ORTHOG_CGS,
 60:                BV_ORTHOG_MGS } BVOrthogType;
 61: SLEPC_EXTERN const char *BVOrthogTypes[];

 63: /*E
 64:     BVOrthogRefineType - Determines what type of refinement to use
 65:     during orthogonalization of vectors

 67:     Level: advanced

 69: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
 70: E*/
 71: typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
 72:                BV_ORTHOG_REFINE_NEVER,
 73:                BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
 74: SLEPC_EXTERN const char *BVOrthogRefineTypes[];

 76: /*E
 77:     BVOrthogBlockType - Determines the method used in block
 78:     orthogonalization (simultaneous orthogonalization of a set of vectors)

 80:     Level: advanced

 82: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalize()
 83: E*/
 84: typedef enum { BV_ORTHOG_BLOCK_GS,
 85:                BV_ORTHOG_BLOCK_CHOL,
 86:                BV_ORTHOG_BLOCK_TSQR,
 87:                BV_ORTHOG_BLOCK_TSQRCHOL,
 88:                BV_ORTHOG_BLOCK_SVQB     } BVOrthogBlockType;
 89: SLEPC_EXTERN const char *BVOrthogBlockTypes[];

 91: /*E
 92:    BVMatMultType - Different ways of performing the BVMatMult() operation

 94:    Notes:
 95:    Allowed values are
 96: +  BV_MATMULT_VECS - perform a matrix-vector multiply per each column
 97: .  BV_MATMULT_MAT - carry out a Mat-Mat product with a dense matrix
 98: -  BV_MATMULT_MAT_SAVE - this case is deprecated

100:    The default is BV_MATMULT_MAT except in the case of BVVECS.

102:    Level: advanced

104: .seealso: BVSetMatMultMethod(), BVMatMult()
105: E*/
106: typedef enum { BV_MATMULT_VECS,
107:                BV_MATMULT_MAT,
108:                BV_MATMULT_MAT_SAVE } BVMatMultType;
109: SLEPC_EXTERN const char *BVMatMultTypes[];

111: SLEPC_EXTERN PetscErrorCode BVCreate(MPI_Comm,BV*);
112: SLEPC_EXTERN PetscErrorCode BVDestroy(BV*);
113: SLEPC_EXTERN PetscErrorCode BVSetType(BV,BVType);
114: SLEPC_EXTERN PetscErrorCode BVGetType(BV,BVType*);
115: SLEPC_EXTERN PetscErrorCode BVSetSizes(BV,PetscInt,PetscInt,PetscInt);
116: SLEPC_EXTERN PetscErrorCode BVSetSizesFromVec(BV,Vec,PetscInt);
117: SLEPC_EXTERN PetscErrorCode BVSetLeadingDimension(BV,PetscInt);
118: SLEPC_EXTERN PetscErrorCode BVGetLeadingDimension(BV,PetscInt*);
119: SLEPC_EXTERN PetscErrorCode BVGetSizes(BV,PetscInt*,PetscInt*,PetscInt*);
120: SLEPC_EXTERN PetscErrorCode BVResize(BV,PetscInt,PetscBool);
121: SLEPC_EXTERN PetscErrorCode BVSetFromOptions(BV);
122: SLEPC_EXTERN PetscErrorCode BVView(BV,PetscViewer);
123: SLEPC_EXTERN PetscErrorCode BVViewFromOptions(BV,PetscObject,const char[]);

125: SLEPC_EXTERN PetscErrorCode BVGetColumn(BV,PetscInt,Vec*);
126: SLEPC_EXTERN PetscErrorCode BVRestoreColumn(BV,PetscInt,Vec*);
127: SLEPC_EXTERN PetscErrorCode BVGetSplit(BV,BV*,BV*);
128: SLEPC_EXTERN PetscErrorCode BVRestoreSplit(BV,BV*,BV*);
129: SLEPC_EXTERN PetscErrorCode BVGetSplitRows(BV,IS,IS,BV*,BV*);
130: SLEPC_EXTERN PetscErrorCode BVRestoreSplitRows(BV,IS,IS,BV*,BV*);
131: SLEPC_EXTERN PetscErrorCode BVGetArray(BV,PetscScalar**);
132: SLEPC_EXTERN PetscErrorCode BVRestoreArray(BV,PetscScalar**);
133: SLEPC_EXTERN PetscErrorCode BVGetArrayRead(BV,const PetscScalar**);
134: SLEPC_EXTERN PetscErrorCode BVRestoreArrayRead(BV,const PetscScalar**);
135: SLEPC_EXTERN PetscErrorCode BVCreateVec(BV,Vec*);
136: SLEPC_EXTERN PetscErrorCode BVCreateVecEmpty(BV,Vec*);
137: SLEPC_EXTERN PetscErrorCode BVSetVecType(BV,VecType);
138: SLEPC_EXTERN PetscErrorCode BVGetVecType(BV,VecType*);
139: SLEPC_EXTERN PetscErrorCode BVSetActiveColumns(BV,PetscInt,PetscInt);
140: SLEPC_EXTERN PetscErrorCode BVGetActiveColumns(BV,PetscInt*,PetscInt*);
141: SLEPC_EXTERN PetscErrorCode BVInsertVec(BV,PetscInt,Vec);
142: SLEPC_EXTERN PetscErrorCode BVInsertVecs(BV,PetscInt,PetscInt*,Vec*,PetscBool);
143: SLEPC_EXTERN PetscErrorCode BVInsertConstraints(BV,PetscInt*,Vec*);
144: SLEPC_EXTERN PetscErrorCode BVSetNumConstraints(BV,PetscInt);
145: SLEPC_EXTERN PetscErrorCode BVGetNumConstraints(BV,PetscInt*);
146: SLEPC_EXTERN PetscErrorCode BVSetDefiniteTolerance(BV,PetscReal);
147: SLEPC_EXTERN PetscErrorCode BVGetDefiniteTolerance(BV,PetscReal*);
148: SLEPC_EXTERN PetscErrorCode BVDuplicate(BV,BV*);
149: SLEPC_EXTERN PetscErrorCode BVDuplicateResize(BV,PetscInt,BV*);
150: SLEPC_EXTERN PetscErrorCode BVCopy(BV,BV);
151: SLEPC_EXTERN PetscErrorCode BVCopyVec(BV,PetscInt,Vec);
152: SLEPC_EXTERN PetscErrorCode BVCopyColumn(BV,PetscInt,PetscInt);
153: SLEPC_EXTERN PetscErrorCode BVSetMatrix(BV,Mat,PetscBool);
154: SLEPC_EXTERN PetscErrorCode BVGetMatrix(BV,Mat*,PetscBool*);
155: SLEPC_EXTERN PetscErrorCode BVApplyMatrix(BV,Vec,Vec);
156: SLEPC_EXTERN PetscErrorCode BVApplyMatrixBV(BV,BV);
157: SLEPC_EXTERN PetscErrorCode BVGetCachedBV(BV,BV*);
158: SLEPC_EXTERN PetscErrorCode BVSetSignature(BV,Vec);
159: SLEPC_EXTERN PetscErrorCode BVGetSignature(BV,Vec);
160: SLEPC_EXTERN PetscErrorCode BVSetBufferVec(BV,Vec);
161: SLEPC_EXTERN PetscErrorCode BVGetBufferVec(BV,Vec*);

163: SLEPC_EXTERN PetscErrorCode BVMult(BV,PetscScalar,PetscScalar,BV,Mat);
164: SLEPC_EXTERN PetscErrorCode BVMultVec(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
165: SLEPC_EXTERN PetscErrorCode BVMultColumn(BV,PetscScalar,PetscScalar,PetscInt,PetscScalar*);
166: SLEPC_EXTERN PetscErrorCode BVMultInPlace(BV,Mat,PetscInt,PetscInt);
167: SLEPC_EXTERN PetscErrorCode BVMultInPlaceHermitianTranspose(BV,Mat,PetscInt,PetscInt);
168: PETSC_DEPRECATED_FUNCTION(3, 16, 0, "BVMultInPlaceHermitianTranspose()", ) static inline PetscErrorCode BVMultInPlaceTranspose(BV bv,Mat A,PetscInt s,PetscInt e) {return BVMultInPlaceHermitianTranspose(bv,A,s,e);}
169: SLEPC_EXTERN PetscErrorCode BVMatMult(BV,Mat,BV);
170: SLEPC_EXTERN PetscErrorCode BVMatMultTranspose(BV,Mat,BV);
171: SLEPC_EXTERN PetscErrorCode BVMatMultHermitianTranspose(BV,Mat,BV);
172: SLEPC_EXTERN PetscErrorCode BVMatMultColumn(BV,Mat,PetscInt);
173: SLEPC_EXTERN PetscErrorCode BVMatMultTransposeColumn(BV,Mat,PetscInt);
174: SLEPC_EXTERN PetscErrorCode BVMatMultHermitianTransposeColumn(BV,Mat,PetscInt);
175: SLEPC_EXTERN PetscErrorCode BVMatProject(BV,Mat,BV,Mat);
176: SLEPC_EXTERN PetscErrorCode BVMatArnoldi(BV,Mat,Mat,PetscInt,PetscInt*,PetscReal*,PetscBool*);
177: SLEPC_EXTERN PetscErrorCode BVMatLanczos(BV,Mat,Mat,PetscInt,PetscInt*,PetscReal*,PetscBool*);

179: SLEPC_EXTERN PetscErrorCode BVDot(BV,BV,Mat);
180: SLEPC_EXTERN PetscErrorCode BVDotVec(BV,Vec,PetscScalar*);
181: SLEPC_EXTERN PetscErrorCode BVDotVecBegin(BV,Vec,PetscScalar*);
182: SLEPC_EXTERN PetscErrorCode BVDotVecEnd(BV,Vec,PetscScalar*);
183: SLEPC_EXTERN PetscErrorCode BVDotColumn(BV,PetscInt,PetscScalar*);
184: SLEPC_EXTERN PetscErrorCode BVDotColumnBegin(BV,PetscInt,PetscScalar*);
185: SLEPC_EXTERN PetscErrorCode BVDotColumnEnd(BV,PetscInt,PetscScalar*);
186: SLEPC_EXTERN PetscErrorCode BVScale(BV,PetscScalar);
187: SLEPC_EXTERN PetscErrorCode BVScaleColumn(BV,PetscInt,PetscScalar);
188: SLEPC_EXTERN PetscErrorCode BVNorm(BV,NormType,PetscReal*);
189: SLEPC_EXTERN PetscErrorCode BVNormVec(BV,Vec,NormType,PetscReal*);
190: SLEPC_EXTERN PetscErrorCode BVNormVecBegin(BV,Vec,NormType,PetscReal*);
191: SLEPC_EXTERN PetscErrorCode BVNormVecEnd(BV,Vec,NormType,PetscReal*);
192: SLEPC_EXTERN PetscErrorCode BVNormColumn(BV,PetscInt,NormType,PetscReal*);
193: SLEPC_EXTERN PetscErrorCode BVNormColumnBegin(BV,PetscInt,NormType,PetscReal*);
194: SLEPC_EXTERN PetscErrorCode BVNormColumnEnd(BV,PetscInt,NormType,PetscReal*);
195: SLEPC_EXTERN PetscErrorCode BVNormalize(BV,PetscScalar*);
196: SLEPC_EXTERN PetscErrorCode BVSetRandom(BV);
197: SLEPC_EXTERN PetscErrorCode BVSetRandomNormal(BV);
198: SLEPC_EXTERN PetscErrorCode BVSetRandomSign(BV);
199: SLEPC_EXTERN PetscErrorCode BVSetRandomColumn(BV,PetscInt);
200: SLEPC_EXTERN PetscErrorCode BVSetRandomCond(BV,PetscReal);
201: SLEPC_EXTERN PetscErrorCode BVSetRandomContext(BV,PetscRandom);
202: SLEPC_EXTERN PetscErrorCode BVGetRandomContext(BV,PetscRandom*);

204: SLEPC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal,BVOrthogBlockType);
205: SLEPC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*,BVOrthogBlockType*);
206: SLEPC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
207: SLEPC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar*,PetscReal*,PetscBool*);
208: SLEPC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar*,PetscReal*,PetscBool*);
209: SLEPC_EXTERN PetscErrorCode BVOrthonormalizeColumn(BV,PetscInt,PetscBool,PetscReal*,PetscBool*);
210: SLEPC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar*,PetscReal*,PetscBool*);
211: SLEPC_EXTERN PetscErrorCode BVBiorthogonalizeColumn(BV,BV,PetscInt);
212: SLEPC_EXTERN PetscErrorCode BVBiorthonormalizeColumn(BV,BV,PetscInt,PetscReal*);
213: SLEPC_EXTERN PetscErrorCode BVSetMatMultMethod(BV,BVMatMultType);
214: SLEPC_EXTERN PetscErrorCode BVGetMatMultMethod(BV,BVMatMultType*);

216: SLEPC_EXTERN PetscErrorCode BVCreateFromMat(Mat,BV*);
217: SLEPC_EXTERN PetscErrorCode BVCreateMat(BV,Mat*);
218: SLEPC_EXTERN PetscErrorCode BVGetMat(BV,Mat*);
219: SLEPC_EXTERN PetscErrorCode BVRestoreMat(BV,Mat*);

221: SLEPC_EXTERN PetscErrorCode BVScatter(BV,BV,VecScatter,Vec);
222: SLEPC_EXTERN PetscErrorCode BVSumQuadrature(BV,BV,PetscInt,PetscInt,PetscInt,PetscScalar*,PetscScalar*,VecScatter,PetscSubcomm,PetscInt,PetscBool);
223: SLEPC_EXTERN PetscErrorCode BVDotQuadrature(BV,BV,PetscScalar*,PetscInt,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscSubcomm,PetscInt,PetscBool);
224: SLEPC_EXTERN PetscErrorCode BVTraceQuadrature(BV,BV,PetscInt,PetscInt,PetscScalar*,VecScatter,PetscSubcomm,PetscInt,PetscBool,PetscReal*);

226: /*E
227:    BVSVDMethod - Different methods for computing the SVD of a BV

229:    Notes:
230:    Allowed values are
231: +  BV_SVD_METHOD_REFINE - based on the SVD of the cross product matrix S'*S, with refinement
232: .  BV_SVD_METHOD_QR     - based on the SVD of the triangular factor of qr(S)
233: -  BV_SVD_METHOD_QR_CAA - variant of QR intended for use in cammunication-avoiding Arnoldi

235:    Level: developer

237: .seealso: BVSVDAndRank()
238: E*/
239: typedef enum { BV_SVD_METHOD_REFINE,
240:                BV_SVD_METHOD_QR,
241:                BV_SVD_METHOD_QR_CAA } BVSVDMethod;
242: SLEPC_EXTERN const char *BVSVDMethods[];

244: SLEPC_EXTERN PetscErrorCode BVSVDAndRank(BV,PetscInt,PetscInt,PetscReal,BVSVDMethod,PetscScalar*,PetscReal*,PetscInt*);
245: SLEPC_EXTERN PetscErrorCode BVCISSResizeBases(BV,BV,BV,PetscInt,PetscInt,PetscInt,PetscInt);

247: SLEPC_EXTERN PetscErrorCode BVCreateTensor(BV,PetscInt,BV*);
248: SLEPC_EXTERN PetscErrorCode BVTensorBuildFirstColumn(BV,PetscInt);
249: SLEPC_EXTERN PetscErrorCode BVTensorCompress(BV,PetscInt);
250: SLEPC_EXTERN PetscErrorCode BVTensorGetDegree(BV,PetscInt*);
251: SLEPC_EXTERN PetscErrorCode BVTensorGetFactors(BV,BV*,Mat*);
252: SLEPC_EXTERN PetscErrorCode BVTensorRestoreFactors(BV,BV*,Mat*);

254: SLEPC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char*);
255: SLEPC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char*);
256: SLEPC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);

258: SLEPC_EXTERN PetscFunctionList BVList;
259: SLEPC_EXTERN PetscErrorCode BVRegister(const char[],PetscErrorCode(*)(BV));