Actual source code: slepcbv.h

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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;
 65: PETSC_EXTERN const char *BVOrthogTypes[];

 67: /*E
 68:     BVOrthogRefineType - Determines what type of refinement to use
 69:     during orthogonalization of vectors

 71:     Level: advanced

 73: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalizeColumn()
 74: E*/
 75: typedef enum { BV_ORTHOG_REFINE_IFNEEDED,
 76:                BV_ORTHOG_REFINE_NEVER,
 77:                BV_ORTHOG_REFINE_ALWAYS } BVOrthogRefineType;
 78: PETSC_EXTERN const char *BVOrthogRefineTypes[];

 80: /*E
 81:     BVOrthogBlockType - Determines the method used in block
 82:     orthogonalization (simultaneous orthogonalization of a set of vectors)

 84:     Level: advanced

 86: .seealso: BVSetOrthogonalization(), BVGetOrthogonalization(), BVOrthogonalize()
 87: E*/
 88: typedef enum { BV_ORTHOG_BLOCK_GS,
 89:                BV_ORTHOG_BLOCK_CHOL } BVOrthogBlockType;
 90: PETSC_EXTERN const char *BVOrthogBlockTypes[];

 92: /*E
 93:     BVMatMultType - Determines how to perform the BVMatMult() operation:
 94:        BV_MATMULT_VECS: perform a matrix-vector multiply per each column;
 95:        BV_MATMULT_MAT: carry out a MatMatMult() product with a dense matrix (default);
 96:        BV_MATMULT_MAT_SAVE: call MatMatMult() and keep auxiliary matrices
 97:          (more efficient but needs more memory)

 99:     Level: advanced

101: .seealso: 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);

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

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

176: PETSC_EXTERN PetscErrorCode BVSetOrthogonalization(BV,BVOrthogType,BVOrthogRefineType,PetscReal,BVOrthogBlockType);
177: PETSC_EXTERN PetscErrorCode BVGetOrthogonalization(BV,BVOrthogType*,BVOrthogRefineType*,PetscReal*,BVOrthogBlockType*);
178: PETSC_EXTERN PetscErrorCode BVOrthogonalize(BV,Mat);
179: PETSC_EXTERN PetscErrorCode BVOrthogonalizeVec(BV,Vec,PetscScalar*,PetscReal*,PetscBool*);
180: PETSC_EXTERN PetscErrorCode BVOrthogonalizeColumn(BV,PetscInt,PetscScalar*,PetscReal*,PetscBool*);
181: PETSC_EXTERN PetscErrorCode BVOrthogonalizeSomeColumn(BV,PetscInt,PetscBool*,PetscScalar*,PetscReal*,PetscBool*);
182: PETSC_EXTERN PetscErrorCode BVSetMatMultMethod(BV,BVMatMultType);
183: PETSC_EXTERN PetscErrorCode BVGetMatMultMethod(BV,BVMatMultType*);

185: PETSC_EXTERN PetscErrorCode BVSetOptionsPrefix(BV,const char*);
186: PETSC_EXTERN PetscErrorCode BVAppendOptionsPrefix(BV,const char*);
187: PETSC_EXTERN PetscErrorCode BVGetOptionsPrefix(BV,const char*[]);

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

192: #endif