Actual source code: slepcbv.h

slepc-3.9.2 2018-07-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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"
 43: #define BVTENSOR     "tensor"

 45: /* Logging support */
 46: PETSC_EXTERN PetscClassId BV_CLASSID;

 48: /*E
 49:     BVOrthogType - Determines the method used in the orthogonalization
 50:     of vectors

 52:     Level: advanced

 54: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn(), BVOrthogRefineType
 55: E*/
 56: typedef enum { BV_ORTHOG_CGS,
 57:                BV_ORTHOG_MGS } BVOrthogType;
 58: PETSC_EXTERN const char *BVOrthogTypes[];

 60: /*E
 61:     BVOrthogRefineType - Determines what type of refinement to use
 62:     during orthogonalization of vectors

 64:     Level: advanced

 66: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
 67: E*/
 68: typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
 69:                BV_ORTHOG_REFINE_NEVER,
 70:                BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
 71: PETSC_EXTERN const char *BVOrthogRefineTypes[];

 73: /*E
 74:     BVOrthogBlockType - Determines the method used in block
 75:     orthogonalization (simultaneous orthogonalization of a set of vectors)

 77:     Level: advanced

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

 88: /*E
 89:    BVMatMultType - Different ways of performing the BVMatMult() operation

 91:    Notes:
 92:    Allowed values are
 93: +  BV_MATMULT_VECS - perform a matrix-vector multiply per each column
 94: .  BV_MATMULT_MAT - carry out a MatMatMult() product with a dense matrix
 95: -  BV_MATMULT_MAT_SAVE - this case is deprecated

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

 99:    Level: advanced

101: .seealso: BVSetMatMultMethod(), BVMatMult()
102: E*/
103: typedef enum { BV_MATMULT_VECS,
104:                BV_MATMULT_MAT,
105:                BV_MATMULT_MAT_SAVE } BVMatMultType;
106: PETSC_EXTERN const char *BVMatMultTypes[];

108: PETSC_EXTERN PetscErrorCode BVCreate(MPI_Comm,BV*);
109: PETSC_EXTERN PetscErrorCode BVDestroy(BV*);
110: PETSC_EXTERN PetscErrorCode BVSetType(BV,BVType);
111: PETSC_EXTERN PetscErrorCode BVGetType(BV,BVType*);
112: PETSC_EXTERN PetscErrorCode BVSetSizes(BV,PetscInt,PetscInt,PetscInt);
113: PETSC_EXTERN PetscErrorCode BVSetSizesFromVec(BV,Vec,PetscInt);
114: PETSC_EXTERN PetscErrorCode BVGetSizes(BV,PetscInt*,PetscInt*,PetscInt*);
115: PETSC_EXTERN PetscErrorCode BVResize(BV,PetscInt,PetscBool);
116: PETSC_EXTERN PetscErrorCode BVSetFromOptions(BV);
117: PETSC_EXTERN PetscErrorCode BVView(BV,PetscViewer);
118: PETSC_STATIC_INLINE PetscErrorCode BVViewFromOptions(BV bv,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)bv,obj,name);}

120: PETSC_EXTERN PetscErrorCode BVGetColumn(BV,PetscInt,Vec*);
121: PETSC_EXTERN PetscErrorCode BVRestoreColumn(BV,PetscInt,Vec*);
122: PETSC_EXTERN PetscErrorCode BVGetSplit(BV,BV*,BV*);
123: PETSC_EXTERN PetscErrorCode BVRestoreSplit(BV,BV*,BV*);
124: PETSC_EXTERN PetscErrorCode BVGetArray(BV,PetscScalar**);
125: PETSC_EXTERN PetscErrorCode BVRestoreArray(BV,PetscScalar**);
126: PETSC_EXTERN PetscErrorCode BVGetArrayRead(BV,const PetscScalar**);
127: PETSC_EXTERN PetscErrorCode BVRestoreArrayRead(BV,const PetscScalar**);
128: PETSC_EXTERN PetscErrorCode BVCreateVec(BV,Vec*);
129: PETSC_EXTERN PetscErrorCode BVSetActiveColumns(BV,PetscInt,PetscInt);
130: PETSC_EXTERN PetscErrorCode BVGetActiveColumns(BV,PetscInt*,PetscInt*);
131: PETSC_EXTERN PetscErrorCode BVInsertVec(BV,PetscInt,Vec);
132: PETSC_EXTERN PetscErrorCode BVInsertVecs(BV,PetscInt,PetscInt*,Vec*,PetscBool);
133: PETSC_EXTERN PetscErrorCode BVInsertConstraints(BV,PetscInt*,Vec*);
134: PETSC_EXTERN PetscErrorCode BVSetNumConstraints(BV,PetscInt);
135: PETSC_EXTERN PetscErrorCode BVGetNumConstraints(BV,PetscInt*);
136: PETSC_EXTERN PetscErrorCode BVDuplicate(BV,BV*);
137: PETSC_EXTERN PetscErrorCode BVDuplicateResize(BV,PetscInt,BV*);
138: PETSC_EXTERN PetscErrorCode BVCopy(BV,BV);
139: PETSC_EXTERN PetscErrorCode BVCopyVec(BV,PetscInt,Vec);
140: PETSC_EXTERN PetscErrorCode BVCopyColumn(BV,PetscInt,PetscInt);
141: PETSC_EXTERN PetscErrorCode BVSetMatrix(BV,Mat,PetscBool);
142: PETSC_EXTERN PetscErrorCode BVGetMatrix(BV,Mat*,PetscBool*);
143: PETSC_EXTERN PetscErrorCode BVApplyMatrix(BV,Vec,Vec);
144: PETSC_EXTERN PetscErrorCode BVApplyMatrixBV(BV,BV);
145: PETSC_EXTERN PetscErrorCode BVGetCachedBV(BV,BV*);
146: PETSC_EXTERN PetscErrorCode BVSetSignature(BV,Vec);
147: PETSC_EXTERN PetscErrorCode BVGetSignature(BV,Vec);
148: PETSC_EXTERN PetscErrorCode BVSetBufferVec(BV,Vec);
149: PETSC_EXTERN PetscErrorCode BVGetBufferVec(BV,Vec*);

151: PETSC_EXTERN PetscErrorCode BVMult(BV,PetscScalar,PetscScalar,BV,Mat);
152: PETSC_EXTERN PetscErrorCode BVMultVec(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
153: PETSC_EXTERN PetscErrorCode BVMultColumn(BV,PetscScalar,PetscScalar,PetscInt,PetscScalar*);
154: PETSC_EXTERN PetscErrorCode BVMultInPlace(BV,Mat,PetscInt,PetscInt);
155: PETSC_EXTERN PetscErrorCode BVMultInPlaceTranspose(BV,Mat,PetscInt,PetscInt);
156: PETSC_EXTERN PetscErrorCode BVMatMult(BV,Mat,BV);
157: PETSC_EXTERN PetscErrorCode BVMatMultHermitianTranspose(BV,Mat,BV);
158: PETSC_EXTERN PetscErrorCode BVMatMultColumn(BV,Mat,PetscInt);
159: PETSC_EXTERN PetscErrorCode BVMatProject(BV,Mat,BV,Mat);
160: PETSC_EXTERN PetscErrorCode BVDot(BV,BV,Mat);
161: PETSC_EXTERN PetscErrorCode BVDotVec(BV,Vec,PetscScalar*);
162: PETSC_EXTERN PetscErrorCode BVDotVecBegin(BV,Vec,PetscScalar*);
163: PETSC_EXTERN PetscErrorCode BVDotVecEnd(BV,Vec,PetscScalar*);
164: PETSC_EXTERN PetscErrorCode BVDotColumn(BV,PetscInt,PetscScalar*);
165: PETSC_EXTERN PetscErrorCode BVDotColumnBegin(BV,PetscInt,PetscScalar*);
166: PETSC_EXTERN PetscErrorCode BVDotColumnEnd(BV,PetscInt,PetscScalar*);
167: PETSC_EXTERN PetscErrorCode BVScale(BV,PetscScalar);
168: PETSC_EXTERN PetscErrorCode BVScaleColumn(BV,PetscInt,PetscScalar);
169: PETSC_EXTERN PetscErrorCode BVNorm(BV,NormType,PetscReal*);
170: PETSC_EXTERN PetscErrorCode BVNormVec(BV,Vec,NormType,PetscReal*);
171: PETSC_EXTERN PetscErrorCode BVNormVecBegin(BV,Vec,NormType,PetscReal*);
172: PETSC_EXTERN PetscErrorCode BVNormVecEnd(BV,Vec,NormType,PetscReal*);
173: PETSC_EXTERN PetscErrorCode BVNormColumn(BV,PetscInt,NormType,PetscReal*);
174: PETSC_EXTERN PetscErrorCode BVNormColumnBegin(BV,PetscInt,NormType,PetscReal*);
175: PETSC_EXTERN PetscErrorCode BVNormColumnEnd(BV,PetscInt,NormType,PetscReal*);
176: PETSC_EXTERN PetscErrorCode BVSetRandom(BV);
177: PETSC_EXTERN PetscErrorCode BVSetRandomColumn(BV,PetscInt);
178: PETSC_EXTERN PetscErrorCode BVSetRandomCond(BV,PetscReal);
179: PETSC_EXTERN PetscErrorCode BVSetRandomContext(BV,PetscRandom);
180: PETSC_EXTERN PetscErrorCode BVGetRandomContext(BV,PetscRandom*);

182: PETSC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal,BVOrthogBlockType);
183: PETSC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*,BVOrthogBlockType*);
184: PETSC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
185: PETSC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar*,PetscReal*,PetscBool*);
186: PETSC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar*,PetscReal*,PetscBool*);
187: PETSC_EXTERN PetscErrorCode BVOrthonormalizeColumn(BV,PetscInt,PetscBool,PetscReal*,PetscBool*);
188: PETSC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar*,PetscReal*,PetscBool*);
189: PETSC_EXTERN PetscErrorCode BVSetMatMultMethod(BV,BVMatMultType);
190: PETSC_EXTERN PetscErrorCode BVGetMatMultMethod(BV,BVMatMultType*);

192: PETSC_EXTERN PetscErrorCode BVCreateFromMat(Mat,BV*);
193: PETSC_EXTERN PetscErrorCode BVCreateMat(BV,Mat*);
194: PETSC_EXTERN PetscErrorCode BVGetMat(BV,Mat*);
195: PETSC_EXTERN PetscErrorCode BVRestoreMat(BV,Mat*);

197: PETSC_EXTERN PetscErrorCode BVCreateTensor(BV,PetscInt,BV*);
198: PETSC_EXTERN PetscErrorCode BVTensorBuildFirstColumn(BV,PetscInt);
199: PETSC_EXTERN PetscErrorCode BVTensorCompress(BV,PetscInt);
200: PETSC_EXTERN PetscErrorCode BVTensorGetDegree(BV,PetscInt*);
201: PETSC_EXTERN PetscErrorCode BVTensorGetFactors(BV,BV*,Mat*);
202: PETSC_EXTERN PetscErrorCode BVTensorRestoreFactors(BV,BV*,Mat*);

204: PETSC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char*);
205: PETSC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char*);
206: PETSC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);

208: PETSC_EXTERN PetscFunctionList BVList;
209: PETSC_EXTERN PetscErrorCode BVRegister(const char[],PetscErrorCode(*)(BV));

211: #endif