1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407
|
/*
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
SLEPc - Scalable Library for Eigenvalue Problem Computations
Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
This file is part of SLEPc.
SLEPc is distributed under a 2-clause BSD license (see LICENSE).
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
*/
/*
User interface for the basis vectors object in SLEPc
*/
#pragma once
#include <slepcsys.h>
/* SUBMANSEC = BV */
SLEPC_EXTERN PetscErrorCode BVInitializePackage(void);
SLEPC_EXTERN PetscErrorCode BVFinalizePackage(void);
/*S
BV - Basis vectors, SLEPc object representing a collection of vectors
that typically constitute a basis of a subspace.
Level: beginner
.seealso: [](sec:bv), `BVCreate()`
S*/
typedef struct _p_BV* BV;
/*J
BVType - String with the name of the type of `BV`. Each type differs in
the way data is stored internally.
Level: beginner
.seealso: [](sec:bv), `BVSetType()`, `BV`
J*/
typedef const char *BVType;
#define BVMAT "mat"
#define BVSVEC "svec"
#define BVVECS "vecs"
#define BVCONTIGUOUS "contiguous"
#define BVTENSOR "tensor"
/* Logging support */
SLEPC_EXTERN PetscClassId BV_CLASSID;
/*E
BVOrthogType - Determines the method used in the orthogonalization of vectors.
Values:
+ `BV_ORTHOG_CGS` - Classical Gram-Schmidt
- `BV_ORTHOG_MGS` - Modified Gram-Schmidt
Note:
The default is to use CGS. Both CGS and MGS provide similar numerical accuracy when
combined with iterative refinement (the default). MGS results in poorer performance,
both sequential and in parallel, so it is not recommended unless iterative
refinement is disabled.
Level: advanced
.seealso: [](sec:bv), `BVSetOrthogonalization()`, `BVGetOrthogonalization()`, `BVOrthogonalizeColumn()`, `BVOrthogRefineType`
E*/
typedef enum { BV_ORTHOG_CGS,
BV_ORTHOG_MGS } BVOrthogType;
SLEPC_EXTERN const char *BVOrthogTypes[];
/*MC
BV_ORTHOG_CGS - Use Classical Gram-Schmidt for orthogonalization of vectors.
Level: advanced
Note:
Although CGS is not numerically stable, when combined with iterative refinement
it yields a sufficiently accurate result.
.seealso: [](sec:bv), `BVOrthogType`, `BVSetOrthogonalization()`, `BV_ORTHOG_MGS`
M*/
/*MC
BV_ORTHOG_MGS - Use Modified Gram-Schmidt for orthogonalization of vectors.
Level: advanced
Note:
MGS results in poorer performance, both sequential and in parallel, so it is
not recommended unless iterative refinement is disabled.
.seealso: [](sec:bv), `BVOrthogType`, `BVSetOrthogonalization()`, `BV_ORTHOG_CGS`
M*/
/*E
BVOrthogRefineType - Determines what type of iterative refinement to use
during orthogonalization of vectors.
Values:
+ `BV_ORTHOG_REFINE_IFNEEDED` - Refine only if a certain criterion is satisfied
. `BV_ORTHOG_REFINE_NEVER` - Never refine, do plain CGS or MGS
- `BV_ORTHOG_REFINE_ALWAYS` - Always refine, i.e., run CGS2 or MGS2
Notes:
The default is to perform one step of iterative refinement if a certain
numerical criterion holds. In ill-conditioned cases, a third orthogonalization
can also be done.
Never refining is not recommended because it will generally make the eigensolver
numerically unstable.
Always refining is numerically stable, but it often performs unnecessary
computation.
Level: advanced
.seealso: [](sec:bv), `BVSetOrthogonalization()`, `BVGetOrthogonalization()`, `BVOrthogonalizeColumn()`
E*/
typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
BV_ORTHOG_REFINE_NEVER,
BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
SLEPC_EXTERN const char *BVOrthogRefineTypes[];
/*MC
BV_ORTHOG_REFINE_IFNEEDED - In the orthogonalization of vectors, do iterative
refinement only if a certain criterion is satisfied.
Note:
This is the default. One step of iterative refinement will be done if a certain
numerical criterion holds. In ill-conditioned cases, a third orthogonalization
can also be done. The criterion is similar to the one proposed in {cite:p}`Dan76`.
The parameter $\eta$ for the criterion can be set in `BVSetOrthogonalization()`.
Level: advanced
.seealso: [](sec:bv), `BVOrthogRefineType`, `BVSetOrthogonalization()`, `BV_ORTHOG_REFINE_NEVER`, `BV_ORTHOG_REFINE_ALWAYS`
M*/
/*MC
BV_ORTHOG_REFINE_NEVER - In the orthogonalization of vectors, never do an iterative
refinement step, i.e., do plain CGS or MGS.
Note:
Never refining is not recommended because it will generally make the eigensolver
numerically unstable.
Level: advanced
.seealso: [](sec:bv), `BVOrthogRefineType`, `BVSetOrthogonalization()`, `BV_ORTHOG_REFINE_IFNEEDED`, `BV_ORTHOG_REFINE_ALWAYS`
M*/
/*MC
BV_ORTHOG_REFINE_ALWAYS - In the orthogonalization of vectors, always do an iterative
refinement step, i.e., run CGS2 or MGS2.
Note:
Always refining is numerically stable, but it often performs unnecessary
computation.
Level: advanced
.seealso: [](sec:bv), `BVOrthogRefineType`, `BVSetOrthogonalization()`, `BV_ORTHOG_REFINE_IFNEEDED`, `BV_ORTHOG_REFINE_NEVER`
M*/
/*E
BVOrthogBlockType - Determines the method used in block
orthogonalization (simultaneous orthogonalization of a set of vectors).
Values:
+ `BV_ORTHOG_BLOCK_GS` - Gram-Schmidt, column by column
. `BV_ORTHOG_BLOCK_CHOL` - Cholesky QR method
. `BV_ORTHOG_BLOCK_TSQR` - Tall-Skinny QR method
. `BV_ORTHOG_BLOCK_TSQRCHOL` - TSQR, but computing the triangular factor only
- `BV_ORTHOG_BLOCK_SVQB` - SVQB method
Level: advanced
.seealso: [](sec:bv), `BVSetOrthogonalization()`, `BVGetOrthogonalization()`, `BVOrthogonalize()`
E*/
typedef enum { BV_ORTHOG_BLOCK_GS,
BV_ORTHOG_BLOCK_CHOL,
BV_ORTHOG_BLOCK_TSQR,
BV_ORTHOG_BLOCK_TSQRCHOL,
BV_ORTHOG_BLOCK_SVQB } BVOrthogBlockType;
SLEPC_EXTERN const char *BVOrthogBlockTypes[];
/*E
BVMatMultType - Different ways of performing the `BVMatMult()` operation.
Values:
Allowed values are
+ `BV_MATMULT_VECS` - perform a matrix-vector multiply per each column
. `BV_MATMULT_MAT` - carry out a Mat-Mat product with a dense matrix
- `BV_MATMULT_MAT_SAVE` - this case is deprecated
Note:
The default is `BV_MATMULT_MAT` except in the case of `BVVECS`.
Level: advanced
.seealso: [](sec:bv), `BVSetMatMultMethod()`, `BVMatMult()`
E*/
typedef enum { BV_MATMULT_VECS,
BV_MATMULT_MAT,
BV_MATMULT_MAT_SAVE } BVMatMultType;
SLEPC_EXTERN const char *BVMatMultTypes[];
/*MC
BV_MATMULT_VECS - Perform a matrix-vector multiply per each column.
Level: advanced
.seealso: [](sec:bv), `BVMatMultType`, `BVSetMatMultMethod()`, `BVMatMult()`, `BV_MATMULT_MAT`, `BV_MATMULT_MAT_SAVE`
M*/
/*MC
BV_MATMULT_MAT - Perform a `Mat`-`Mat` product with a dense matrix.
Level: advanced
.seealso: [](sec:bv), `BVMatMultType`, `BVSetMatMultMethod()`, `BVMatMult()`, `BV_MATMULT_VECS`, `BV_MATMULT_MAT_SAVE`
M*/
/*MC
BV_MATMULT_MAT_SAVE - This case is deprecated, do not use.
Level: advanced
.seealso: [](sec:bv), `BVMatMultType`, `BVSetMatMultMethod()`, `BVMatMult()`, `BV_MATMULT_VECS`, `BV_MATMULT_MAT`
M*/
SLEPC_EXTERN PetscErrorCode BVCreate(MPI_Comm,BV*);
SLEPC_EXTERN PetscErrorCode BVDestroy(BV*);
SLEPC_EXTERN PetscErrorCode BVSetType(BV,BVType);
SLEPC_EXTERN PetscErrorCode BVGetType(BV,BVType*);
SLEPC_EXTERN PetscErrorCode BVSetSizes(BV,PetscInt,PetscInt,PetscInt);
SLEPC_EXTERN PetscErrorCode BVSetSizesFromVec(BV,Vec,PetscInt);
SLEPC_EXTERN PetscErrorCode BVSetLeadingDimension(BV,PetscInt);
SLEPC_EXTERN PetscErrorCode BVGetLeadingDimension(BV,PetscInt*);
SLEPC_EXTERN PetscErrorCode BVGetSizes(BV,PetscInt*,PetscInt*,PetscInt*);
SLEPC_EXTERN PetscErrorCode BVResize(BV,PetscInt,PetscBool);
SLEPC_EXTERN PetscErrorCode BVSetFromOptions(BV);
SLEPC_EXTERN PetscErrorCode BVView(BV,PetscViewer);
SLEPC_EXTERN PetscErrorCode BVViewFromOptions(BV,PetscObject,const char[]);
SLEPC_EXTERN PetscErrorCode BVGetColumn(BV,PetscInt,Vec*);
SLEPC_EXTERN PetscErrorCode BVRestoreColumn(BV,PetscInt,Vec*);
SLEPC_EXTERN PetscErrorCode BVGetSplit(BV,BV*,BV*);
SLEPC_EXTERN PetscErrorCode BVRestoreSplit(BV,BV*,BV*);
SLEPC_EXTERN PetscErrorCode BVGetSplitRows(BV,IS,IS,BV*,BV*);
SLEPC_EXTERN PetscErrorCode BVRestoreSplitRows(BV,IS,IS,BV*,BV*);
SLEPC_EXTERN PetscErrorCode BVGetArray(BV,PetscScalar*[]);
SLEPC_EXTERN PetscErrorCode BVRestoreArray(BV,PetscScalar*[]);
SLEPC_EXTERN PetscErrorCode BVGetArrayRead(BV,const PetscScalar*[]);
SLEPC_EXTERN PetscErrorCode BVRestoreArrayRead(BV,const PetscScalar*[]);
SLEPC_EXTERN PetscErrorCode BVCreateVec(BV,Vec*);
SLEPC_EXTERN PetscErrorCode BVCreateVecEmpty(BV,Vec*);
SLEPC_EXTERN PetscErrorCode BVSetVecType(BV,VecType);
SLEPC_EXTERN PetscErrorCode BVGetVecType(BV,VecType*);
SLEPC_EXTERN PetscErrorCode BVSetActiveColumns(BV,PetscInt,PetscInt);
SLEPC_EXTERN PetscErrorCode BVGetActiveColumns(BV,PetscInt*,PetscInt*);
SLEPC_EXTERN PetscErrorCode BVInsertVec(BV,PetscInt,Vec);
SLEPC_EXTERN PetscErrorCode BVInsertVecs(BV,PetscInt,PetscInt*,Vec[],PetscBool);
SLEPC_EXTERN PetscErrorCode BVInsertConstraints(BV,PetscInt*,Vec[]);
SLEPC_EXTERN PetscErrorCode BVSetNumConstraints(BV,PetscInt);
SLEPC_EXTERN PetscErrorCode BVGetNumConstraints(BV,PetscInt*);
SLEPC_EXTERN PetscErrorCode BVSetDefiniteTolerance(BV,PetscReal);
SLEPC_EXTERN PetscErrorCode BVGetDefiniteTolerance(BV,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVDuplicate(BV,BV*);
SLEPC_EXTERN PetscErrorCode BVDuplicateResize(BV,PetscInt,BV*);
SLEPC_EXTERN PetscErrorCode BVCopy(BV,BV);
SLEPC_EXTERN PetscErrorCode BVCopyVec(BV,PetscInt,Vec);
SLEPC_EXTERN PetscErrorCode BVCopyColumn(BV,PetscInt,PetscInt);
SLEPC_EXTERN PetscErrorCode BVSetMatrix(BV,Mat,PetscBool);
SLEPC_EXTERN PetscErrorCode BVGetMatrix(BV,Mat*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVApplyMatrix(BV,Vec,Vec);
SLEPC_EXTERN PetscErrorCode BVApplyMatrixBV(BV,BV);
SLEPC_EXTERN PetscErrorCode BVGetCachedBV(BV,BV*);
SLEPC_EXTERN PetscErrorCode BVSetSignature(BV,Vec);
SLEPC_EXTERN PetscErrorCode BVGetSignature(BV,Vec);
SLEPC_EXTERN PetscErrorCode BVSetBufferVec(BV,Vec);
SLEPC_EXTERN PetscErrorCode BVGetBufferVec(BV,Vec*);
SLEPC_EXTERN PetscErrorCode BVMult(BV,PetscScalar,PetscScalar,BV,Mat);
SLEPC_EXTERN PetscErrorCode BVMultVec(BV,PetscScalar,PetscScalar,Vec,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVMultColumn(BV,PetscScalar,PetscScalar,PetscInt,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVMultInPlace(BV,Mat,PetscInt,PetscInt);
SLEPC_EXTERN PetscErrorCode BVMultInPlaceHermitianTranspose(BV,Mat,PetscInt,PetscInt);
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);}
SLEPC_EXTERN PetscErrorCode BVMatMult(BV,Mat,BV);
SLEPC_EXTERN PetscErrorCode BVMatMultTranspose(BV,Mat,BV);
SLEPC_EXTERN PetscErrorCode BVMatMultHermitianTranspose(BV,Mat,BV);
SLEPC_EXTERN PetscErrorCode BVMatMultColumn(BV,Mat,PetscInt);
SLEPC_EXTERN PetscErrorCode BVMatMultTransposeColumn(BV,Mat,PetscInt);
SLEPC_EXTERN PetscErrorCode BVMatMultHermitianTransposeColumn(BV,Mat,PetscInt);
SLEPC_EXTERN PetscErrorCode BVMatProject(BV,Mat,BV,Mat);
SLEPC_EXTERN PetscErrorCode BVMatArnoldi(BV,Mat,Mat,PetscInt,PetscInt*,PetscReal*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVMatLanczos(BV,Mat,Mat,PetscInt,PetscInt*,PetscReal*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVDot(BV,BV,Mat);
SLEPC_EXTERN PetscErrorCode BVDotVec(BV,Vec,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVDotVecBegin(BV,Vec,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVDotVecEnd(BV,Vec,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVDotColumn(BV,PetscInt,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVDotColumnBegin(BV,PetscInt,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVDotColumnEnd(BV,PetscInt,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVScale(BV,PetscScalar);
SLEPC_EXTERN PetscErrorCode BVScaleColumn(BV,PetscInt,PetscScalar);
SLEPC_EXTERN PetscErrorCode BVNorm(BV,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormVec(BV,Vec,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormVecBegin(BV,Vec,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormVecEnd(BV,Vec,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormColumn(BV,PetscInt,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormColumnBegin(BV,PetscInt,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormColumnEnd(BV,PetscInt,NormType,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVNormalize(BV,PetscScalar[]);
SLEPC_EXTERN PetscErrorCode BVSetRandom(BV);
SLEPC_EXTERN PetscErrorCode BVSetRandomNormal(BV);
SLEPC_EXTERN PetscErrorCode BVSetRandomSign(BV);
SLEPC_EXTERN PetscErrorCode BVSetRandomColumn(BV,PetscInt);
SLEPC_EXTERN PetscErrorCode BVSetRandomCond(BV,PetscReal);
SLEPC_EXTERN PetscErrorCode BVSetRandomContext(BV,PetscRandom);
SLEPC_EXTERN PetscErrorCode BVGetRandomContext(BV,PetscRandom*);
SLEPC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal,BVOrthogBlockType);
SLEPC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*,BVOrthogBlockType*);
SLEPC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
SLEPC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar[],PetscReal*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar[],PetscReal*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVOrthonormalizeColumn(BV,PetscInt,PetscBool,PetscReal*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar[],PetscReal*,PetscBool*);
SLEPC_EXTERN PetscErrorCode BVBiorthogonalizeColumn(BV,BV,PetscInt);
SLEPC_EXTERN PetscErrorCode BVBiorthonormalizeColumn(BV,BV,PetscInt,PetscReal*);
SLEPC_EXTERN PetscErrorCode BVSetMatMultMethod(BV,BVMatMultType);
SLEPC_EXTERN PetscErrorCode BVGetMatMultMethod(BV,BVMatMultType*);
SLEPC_EXTERN PetscErrorCode BVCreateFromMat(Mat,BV*);
SLEPC_EXTERN PetscErrorCode BVCreateMat(BV,Mat*);
SLEPC_EXTERN PetscErrorCode BVGetMat(BV,Mat*);
SLEPC_EXTERN PetscErrorCode BVRestoreMat(BV,Mat*);
SLEPC_EXTERN PetscErrorCode BVScatter(BV,BV,VecScatter,Vec);
SLEPC_EXTERN PetscErrorCode BVSumQuadrature(BV,BV,PetscInt,PetscInt,PetscInt,PetscScalar[],PetscScalar[],VecScatter,PetscSubcomm,PetscInt,PetscBool);
SLEPC_EXTERN PetscErrorCode BVDotQuadrature(BV,BV,PetscScalar[],PetscInt,PetscInt,PetscInt,PetscScalar[],PetscScalar[],PetscSubcomm,PetscInt,PetscBool);
SLEPC_EXTERN PetscErrorCode BVTraceQuadrature(BV,BV,PetscInt,PetscInt,PetscScalar[],VecScatter,PetscSubcomm,PetscInt,PetscBool,PetscReal*);
/*E
BVSVDMethod - Different methods for computing the SVD of a `BV`.
Values:
+ `BV_SVD_METHOD_REFINE` - based on the SVD of the cross product matrix $S^*S$, with refinement
. `BV_SVD_METHOD_QR` - based on the SVD of the triangular factor of qr(S)
- `BV_SVD_METHOD_QR_CAA` - variant of QR intended for use in communication-avoiding Arnoldi
Level: developer
.seealso: [](sec:bv), `BVSVDAndRank()`
E*/
typedef enum { BV_SVD_METHOD_REFINE,
BV_SVD_METHOD_QR,
BV_SVD_METHOD_QR_CAA } BVSVDMethod;
SLEPC_EXTERN const char *BVSVDMethods[];
/*MC
BV_SVD_METHOD_REFINE - Compute the SVD of the `BV` via the cross product matrix $S^*S$,
with iterative refinement.
Level: advanced
.seealso: [](sec:bv), `BVSVDMethod`, `BVSVDAndRank()`, `BV_SVD_METHOD_QR`, `BV_SVD_METHOD_QR_CAA`
M*/
/*MC
BV_SVD_METHOD_QR - Compute the SVD of the `BV` via the SVD of the triangular factor of
the QR factorization of $S$.
Level: advanced
.seealso: [](sec:bv), `BVSVDMethod`, `BVSVDAndRank()`, `BV_SVD_METHOD_REFINE`, `BV_SVD_METHOD_QR_CAA`
M*/
/*MC
BV_SVD_METHOD_QR_CAA - Compute the SVD of the `BV` employing a variant of
`BV_SVD_METHOD_QR` intended for use in communication-avoiding Arnoldi.
Level: advanced
.seealso: [](sec:bv), `BVSVDMethod`, `BVSVDAndRank()`, `BV_SVD_METHOD_REFINE`, `BV_SVD_METHOD_QR`
M*/
SLEPC_EXTERN PetscErrorCode BVSVDAndRank(BV,PetscInt,PetscInt,PetscReal,BVSVDMethod,PetscScalar[],PetscReal[],PetscInt*);
SLEPC_EXTERN PetscErrorCode BVCISSResizeBases(BV,BV,BV,PetscInt,PetscInt,PetscInt,PetscInt);
SLEPC_EXTERN PetscErrorCode BVCreateTensor(BV,PetscInt,BV*);
SLEPC_EXTERN PetscErrorCode BVTensorBuildFirstColumn(BV,PetscInt);
SLEPC_EXTERN PetscErrorCode BVTensorCompress(BV,PetscInt);
SLEPC_EXTERN PetscErrorCode BVTensorGetDegree(BV,PetscInt*);
SLEPC_EXTERN PetscErrorCode BVTensorGetFactors(BV,BV*,Mat*);
SLEPC_EXTERN PetscErrorCode BVTensorRestoreFactors(BV,BV*,Mat*);
SLEPC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char[]);
SLEPC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char[]);
SLEPC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);
SLEPC_EXTERN PetscFunctionList BVList;
SLEPC_EXTERN PetscErrorCode BVRegister(const char[],PetscErrorCode(*)(BV));
|