Actual source code: bvimpl.h
slepc-3.7.0 2016-05-16
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) , 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: */
22: #if !defined(_BVIMPL)
23: #define _BVIMPL
25: #include <slepcbv.h>
26: #include <slepc/private/slepcimpl.h>
28: PETSC_EXTERN PetscBool BVRegisterAllCalled;
29: PETSC_EXTERN PetscErrorCode BVRegisterAll(void);
31: PETSC_EXTERN PetscLogEvent BV_Create,BV_Copy,BV_Mult,BV_MultVec,BV_MultInPlace,BV_Dot,BV_DotVec,BV_Orthogonalize,BV_OrthogonalizeVec,BV_Scale,BV_Norm,BV_NormVec,BV_SetRandom,BV_MatMult,BV_MatMultVec,BV_MatProject;
33: typedef struct _BVOps *BVOps;
35: struct _BVOps {
36: PetscErrorCode (*mult)(BV,PetscScalar,PetscScalar,BV,Mat);
37: PetscErrorCode (*multvec)(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
38: PetscErrorCode (*multinplace)(BV,Mat,PetscInt,PetscInt);
39: PetscErrorCode (*multinplacetrans)(BV,Mat,PetscInt,PetscInt);
40: PetscErrorCode (*dot)(BV,BV,Mat);
41: PetscErrorCode (*dotvec)(BV,Vec,PetscScalar*);
42: PetscErrorCode (*dotvec_local)(BV,Vec,PetscScalar*);
43: PetscErrorCode (*dotvec_begin)(BV,Vec,PetscScalar*);
44: PetscErrorCode (*dotvec_end)(BV,Vec,PetscScalar*);
45: PetscErrorCode (*scale)(BV,PetscInt,PetscScalar);
46: PetscErrorCode (*norm)(BV,PetscInt,NormType,PetscReal*);
47: PetscErrorCode (*norm_local)(BV,PetscInt,NormType,PetscReal*);
48: PetscErrorCode (*norm_begin)(BV,PetscInt,NormType,PetscReal*);
49: PetscErrorCode (*norm_end)(BV,PetscInt,NormType,PetscReal*);
50: PetscErrorCode (*orthogonalize)(BV,Mat);
51: PetscErrorCode (*matmult)(BV,Mat,BV);
52: PetscErrorCode (*copy)(BV,BV);
53: PetscErrorCode (*resize)(BV,PetscInt,PetscBool);
54: PetscErrorCode (*getcolumn)(BV,PetscInt,Vec*);
55: PetscErrorCode (*restorecolumn)(BV,PetscInt,Vec*);
56: PetscErrorCode (*getarray)(BV,PetscScalar**);
57: PetscErrorCode (*restorearray)(BV,PetscScalar**);
58: PetscErrorCode (*getarrayread)(BV,const PetscScalar**);
59: PetscErrorCode (*restorearrayread)(BV,const PetscScalar**);
60: PetscErrorCode (*duplicate)(BV,BV*);
61: PetscErrorCode (*create)(BV);
62: PetscErrorCode (*setfromoptions)(PetscOptionItems*,BV);
63: PetscErrorCode (*view)(BV,PetscViewer);
64: PetscErrorCode (*destroy)(BV);
65: };
67: struct _p_BV {
68: PETSCHEADER(struct _BVOps);
69: /*------------------------- User parameters --------------------------*/
70: Vec t; /* template vector */
71: PetscInt n,N; /* dimensions of vectors (local, global) */
72: PetscInt m; /* number of vectors */
73: PetscInt l; /* number of leading columns */
74: PetscInt k; /* number of active columns */
75: PetscInt nc; /* number of constraints */
76: BVOrthogType orthog_type; /* the method of vector orthogonalization */
77: BVOrthogRefineType orthog_ref; /* refinement method */
78: PetscReal orthog_eta; /* refinement threshold */
79: BVOrthogBlockType orthog_block; /* the method of block orthogonalization */
80: Mat matrix; /* inner product matrix */
81: PetscBool indef; /* matrix is indefinite */
82: BVMatMultType vmm; /* version of matmult operation */
84: /*---------------------- Cached data and workspace -------------------*/
85: Vec Bx; /* result of matrix times a vector x */
86: PetscObjectId xid; /* object id of vector x */
87: PetscObjectState xstate; /* state of vector x */
88: Vec cv[2]; /* column vectors obtained with BVGetColumn() */
89: PetscInt ci[2]; /* column indices of obtained vectors */
90: PetscObjectState st[2]; /* state of obtained vectors */
91: PetscObjectId id[2]; /* object id of obtained vectors */
92: PetscScalar *h,*c; /* orthogonalization coefficients */
93: PetscReal *omega; /* signature matrix values for indefinite case */
94: Mat B,C; /* auxiliary dense matrices for matmult operation */
95: PetscObjectId Aid; /* object id of matrix A of matmult operation */
96: PetscBool defersfo; /* deferred call to setfromoptions */
97: BV cached; /* cached BV to store result of matrix times BV */
98: PetscObjectState bvstate; /* state of BV when BVApplyMatrixBV() was called */
99: PetscRandom rand; /* random number generator */
100: PetscBool rrandom; /* reproducible random vectors */
101: PetscScalar *work;
102: PetscInt lwork;
103: void *data;
104: };
108: /*
109: BV_SafeSqrt - Computes the square root of a scalar value alpha, which is
110: assumed to be z'*B*z. The result is
111: if definite inner product: res = sqrt(alpha)
112: if indefinite inner product: res = sgn(alpha)*sqrt(abs(alpha))
113: */
114: PETSC_STATIC_INLINE PetscErrorCode BV_SafeSqrt(BV bv,PetscScalar alpha,PetscReal *res)
115: {
117: PetscReal absal,realp;
120: absal = PetscAbsScalar(alpha);
121: realp = PetscRealPart(alpha);
122: if (absal<PETSC_MACHINE_EPSILON) {
123: PetscInfo(bv,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
124: }
125: #if defined(PETSC_USE_COMPLEX)
126: if (PetscAbsReal(PetscImaginaryPart(alpha))/absal>PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)bv),1,"The inner product is not well defined: nonzero imaginary part");
127: #endif
128: if (bv->indef) {
129: *res = (realp<0.0)? -PetscSqrtReal(-realp): PetscSqrtReal(realp);
130: } else {
131: if (realp<0.0) SETERRQ(PetscObjectComm((PetscObject)bv),1,"The inner product is not well defined: indefinite matrix");
132: *res = PetscSqrtReal(realp);
133: }
134: return(0);
135: }
139: /*
140: BV_IPMatMult - Multiply a vector x by the inner-product matrix, cache the
141: result in Bx.
142: */
143: PETSC_STATIC_INLINE PetscErrorCode BV_IPMatMult(BV bv,Vec x)
144: {
148: if (((PetscObject)x)->id != bv->xid || ((PetscObject)x)->state != bv->xstate) {
149: MatMult(bv->matrix,x,bv->Bx);
150: bv->xid = ((PetscObject)x)->id;
151: bv->xstate = ((PetscObject)x)->state;
152: }
153: return(0);
154: }
158: /*
159: BV_AllocateCachedBV - Allocate auxiliary BV required for BVApplyMatrixBV if not available.
160: */
161: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateCachedBV(BV V)
162: {
166: if (!V->cached) {
167: BVCreate(PetscObjectComm((PetscObject)V),&V->cached);
168: BVSetSizesFromVec(V->cached,V->t,V->m);
169: BVSetType(V->cached,((PetscObject)V)->type_name);
170: BVSetOrthogonalization(V->cached,V->orthog_type,V->orthog_ref,V->orthog_eta,V->orthog_block);
171: }
172: return(0);
173: }
177: /*
178: BV_IPMatMultBV - Multiply BV by the inner-product matrix, cache the
179: result internally in bv->cached.
180: */
181: PETSC_STATIC_INLINE PetscErrorCode BV_IPMatMultBV(BV bv)
182: {
186: BV_AllocateCachedBV(bv);
187: BVSetActiveColumns(bv->cached,bv->l,bv->k);
188: if (((PetscObject)bv)->state != bv->bvstate) {
189: if (bv->matrix) {
190: BVMatMult(bv,bv->matrix,bv->cached);
191: } else {
192: BVCopy(bv,bv->cached);
193: }
194: bv->bvstate = ((PetscObject)bv)->state;
195: }
196: return(0);
197: }
201: /*
202: BV_AllocateCoeffs - Allocate orthogonalization coefficients if not done already.
203: */
204: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateCoeffs(BV bv)
205: {
209: if (!bv->h) {
210: PetscMalloc2(bv->nc+bv->m,&bv->h,bv->nc+bv->m,&bv->c);
211: PetscLogObjectMemory((PetscObject)bv,2*bv->m*sizeof(PetscScalar));
212: }
213: return(0);
214: }
218: /*
219: BV_AllocateSignature - Allocate signature coefficients if not done already.
220: */
221: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateSignature(BV bv)
222: {
224: PetscInt i;
227: if (bv->indef && !bv->omega) {
228: PetscMalloc1(bv->nc+bv->m,&bv->omega);
229: PetscLogObjectMemory((PetscObject)bv,bv->m*sizeof(PetscReal));
230: for (i=-bv->nc;i<bv->m;i++) bv->omega[i] = 1.0;
231: }
232: return(0);
233: }
237: /*
238: BV_AllocateMatMult - Allocate auxiliary matrices required for BVMatMult if not available.
239: */
240: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateMatMult(BV bv,Mat A,PetscInt m)
241: {
243: PetscObjectId Aid;
244: PetscBool create=PETSC_FALSE;
245: PetscInt cols;
248: if (!bv->B) create=PETSC_TRUE;
249: else {
250: MatGetSize(bv->B,NULL,&cols);
251: PetscObjectGetId((PetscObject)A,&Aid);
252: if (cols!=m || bv->Aid!=Aid) {
253: MatDestroy(&bv->B);
254: MatDestroy(&bv->C);
255: create=PETSC_TRUE;
256: }
257: }
258: if (create) {
259: MatCreateDense(PetscObjectComm((PetscObject)bv),bv->n,PETSC_DECIDE,bv->N,m,NULL,&bv->B);
260: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->B);
261: MatAssemblyBegin(bv->B,MAT_FINAL_ASSEMBLY);
262: MatAssemblyEnd(bv->B,MAT_FINAL_ASSEMBLY);
263: }
264: return(0);
265: }
267: /*
268: BVAvailableVec: First (0) or second (1) vector available for
269: getcolumn operation (or -1 if both vectors already fetched).
270: */
271: #define BVAvailableVec (((bv->ci[0]==-bv->nc-1)? 0: (bv->ci[1]==-bv->nc-1)? 1: -1))
273: /*
274: Macros to test valid BV arguments
275: */
276: #if !defined(PETSC_USE_DEBUG)
278: #define BVCheckSizes(h,arg) do {} while (0)
280: #else
282: #define BVCheckSizes(h,arg) \
283: do { \
284: if (!h->m) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"BV sizes have not been defined: Parameter #%d",arg); \
285: } while (0)
287: #endif
289: PETSC_INTERN PetscErrorCode BVView_Vecs(BV,PetscViewer);
291: PETSC_INTERN PetscErrorCode BVAllocateWork_Private(BV,PetscInt);
293: PETSC_INTERN PetscErrorCode BVMult_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
294: PETSC_INTERN PetscErrorCode BVMultVec_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
295: PETSC_INTERN PetscErrorCode BVMultInPlace_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar*,const PetscScalar*,PetscBool);
296: PETSC_INTERN PetscErrorCode BVMultInPlace_Vecs_Private(BV,PetscInt,PetscInt,PetscInt,Vec*,const PetscScalar*,PetscBool);
297: PETSC_INTERN PetscErrorCode BVAXPY_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,PetscScalar,PetscScalar*);
298: PETSC_INTERN PetscErrorCode BVDot_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
299: PETSC_INTERN PetscErrorCode BVDotVec_BLAS_Private(BV,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
300: PETSC_INTERN PetscErrorCode BVScale_BLAS_Private(BV,PetscInt,PetscScalar*,PetscScalar);
301: PETSC_INTERN PetscErrorCode BVNorm_LAPACK_Private(BV,PetscInt,PetscInt,const PetscScalar*,NormType,PetscReal*,PetscBool);
302: PETSC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_Private(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscBool);
304: #endif