Actual source code: bvimpl.h
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: */
11: #if !defined(SLEPCBVIMPL_H)
12: #define SLEPCBVIMPL_H
14: #include <slepcbv.h>
15: #include <slepc/private/slepcimpl.h>
17: SLEPC_EXTERN PetscBool BVRegisterAllCalled;
18: SLEPC_EXTERN PetscErrorCode BVRegisterAll(void);
20: SLEPC_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_Normalize,BV_SetRandom,BV_MatMult,BV_MatMultVec,BV_MatProject;
22: typedef struct _BVOps *BVOps;
24: struct _BVOps {
25: PetscErrorCode (*mult)(BV,PetscScalar,PetscScalar,BV,Mat);
26: PetscErrorCode (*multvec)(BV,PetscScalar,PetscScalar,Vec,PetscScalar*);
27: PetscErrorCode (*multinplace)(BV,Mat,PetscInt,PetscInt);
28: PetscErrorCode (*multinplacetrans)(BV,Mat,PetscInt,PetscInt);
29: PetscErrorCode (*dot)(BV,BV,Mat);
30: PetscErrorCode (*dotvec)(BV,Vec,PetscScalar*);
31: PetscErrorCode (*dotvec_local)(BV,Vec,PetscScalar*);
32: PetscErrorCode (*dotvec_begin)(BV,Vec,PetscScalar*);
33: PetscErrorCode (*dotvec_end)(BV,Vec,PetscScalar*);
34: PetscErrorCode (*scale)(BV,PetscInt,PetscScalar);
35: PetscErrorCode (*norm)(BV,PetscInt,NormType,PetscReal*);
36: PetscErrorCode (*norm_local)(BV,PetscInt,NormType,PetscReal*);
37: PetscErrorCode (*norm_begin)(BV,PetscInt,NormType,PetscReal*);
38: PetscErrorCode (*norm_end)(BV,PetscInt,NormType,PetscReal*);
39: PetscErrorCode (*normalize)(BV,PetscScalar*);
40: PetscErrorCode (*matmult)(BV,Mat,BV);
41: PetscErrorCode (*copy)(BV,BV);
42: PetscErrorCode (*copycolumn)(BV,PetscInt,PetscInt);
43: PetscErrorCode (*resize)(BV,PetscInt,PetscBool);
44: PetscErrorCode (*getcolumn)(BV,PetscInt,Vec*);
45: PetscErrorCode (*restorecolumn)(BV,PetscInt,Vec*);
46: PetscErrorCode (*getarray)(BV,PetscScalar**);
47: PetscErrorCode (*restorearray)(BV,PetscScalar**);
48: PetscErrorCode (*getarrayread)(BV,const PetscScalar**);
49: PetscErrorCode (*restorearrayread)(BV,const PetscScalar**);
50: PetscErrorCode (*restoresplit)(BV,BV*,BV*);
51: PetscErrorCode (*gramschmidt)(BV,PetscInt,Vec,PetscBool*,PetscScalar*,PetscScalar*,PetscReal*,PetscReal*);
52: PetscErrorCode (*getmat)(BV,Mat*);
53: PetscErrorCode (*restoremat)(BV,Mat*);
54: PetscErrorCode (*duplicate)(BV,BV);
55: PetscErrorCode (*create)(BV);
56: PetscErrorCode (*setfromoptions)(PetscOptionItems*,BV);
57: PetscErrorCode (*view)(BV,PetscViewer);
58: PetscErrorCode (*destroy)(BV);
59: };
61: struct _p_BV {
62: PETSCHEADER(struct _BVOps);
63: /*------------------------- User parameters --------------------------*/
64: Vec t; /* template vector */
65: PetscInt n,N; /* dimensions of vectors (local, global) */
66: PetscInt m; /* number of vectors */
67: PetscInt l; /* number of leading columns */
68: PetscInt k; /* number of active columns */
69: PetscInt nc; /* number of constraints */
70: BVOrthogType orthog_type; /* the method of vector orthogonalization */
71: BVOrthogRefineType orthog_ref; /* refinement method */
72: PetscReal orthog_eta; /* refinement threshold */
73: BVOrthogBlockType orthog_block; /* the method of block orthogonalization */
74: Mat matrix; /* inner product matrix */
75: PetscBool indef; /* matrix is indefinite */
76: BVMatMultType vmm; /* version of matmult operation */
77: PetscBool rrandom; /* reproducible random vectors */
79: /*---------------------- Cached data and workspace -------------------*/
80: Vec buffer; /* buffer vector used in orthogonalization */
81: Mat Abuffer; /* auxiliary seqdense matrix that wraps the buffer */
82: Vec Bx; /* result of matrix times a vector x */
83: PetscObjectId xid; /* object id of vector x */
84: PetscObjectState xstate; /* state of vector x */
85: Vec cv[2]; /* column vectors obtained with BVGetColumn() */
86: PetscInt ci[2]; /* column indices of obtained vectors */
87: PetscObjectState st[2]; /* state of obtained vectors */
88: PetscObjectId id[2]; /* object id of obtained vectors */
89: PetscScalar *h,*c; /* orthogonalization coefficients */
90: Vec omega; /* signature matrix values for indefinite case */
91: PetscBool defersfo; /* deferred call to setfromoptions */
92: BV cached; /* cached BV to store result of matrix times BV */
93: PetscObjectState bvstate; /* state of BV when BVApplyMatrixBV() was called */
94: BV L,R; /* BV objects obtained with BVGetSplit() */
95: PetscObjectState lstate,rstate;/* state of L and R when BVGetSplit() was called */
96: PetscInt lsplit; /* the value of l when BVGetSplit() was called */
97: PetscInt issplit; /* >0 if this BV has been created by splitting (1=left, 2=right) */
98: BV splitparent; /* my parent if I am a split BV */
99: PetscRandom rand; /* random number generator */
100: Mat Acreate; /* matrix given at BVCreateFromMat() */
101: Mat Aget; /* matrix returned for BVGetMat() */
102: PetscBool cuda; /* true if GPU must be used in SVEC */
103: PetscScalar *work;
104: PetscInt lwork;
105: void *data;
106: };
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))>10*PETSC_MACHINE_EPSILON && PetscAbsReal(PetscImaginaryPart(alpha))/absal>100*PETSC_MACHINE_EPSILON) SETERRQ1(PetscObjectComm((PetscObject)bv),1,"The inner product is not well defined: nonzero imaginary part %g",PetscImaginaryPart(alpha));
127: #endif
128: if (bv->indef) {
129: *res = (realp<0.0)? -PetscSqrtReal(-realp): PetscSqrtReal(realp);
130: } else {
131: if (realp<-10*PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)bv),1,"The inner product is not well defined: indefinite matrix");
132: *res = (realp<0.0)? 0.0: PetscSqrtReal(realp);
133: }
134: return(0);
135: }
137: /*
138: BV_IPMatMult - Multiply a vector x by the inner-product matrix, cache the
139: result in Bx.
140: */
141: PETSC_STATIC_INLINE PetscErrorCode BV_IPMatMult(BV bv,Vec x)
142: {
146: if (((PetscObject)x)->id != bv->xid || ((PetscObject)x)->state != bv->xstate) {
147: if (!bv->Bx) {
148: MatCreateVecs(bv->matrix,&bv->Bx,NULL);
149: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
150: }
151: MatMult(bv->matrix,x,bv->Bx);
152: PetscObjectGetId((PetscObject)x,&bv->xid);
153: PetscObjectStateGet((PetscObject)x,&bv->xstate);
154: }
155: return(0);
156: }
158: /*
159: BV_IPMatMultBV - Multiply BV by the inner-product matrix, cache the
160: result internally in bv->cached.
161: */
162: PETSC_STATIC_INLINE PetscErrorCode BV_IPMatMultBV(BV bv)
163: {
167: BVGetCachedBV(bv,&bv->cached);
168: if (((PetscObject)bv)->state != bv->bvstate || bv->l != bv->cached->l || bv->k != bv->cached->k) {
169: BVSetActiveColumns(bv->cached,bv->l,bv->k);
170: if (bv->matrix) {
171: BVMatMult(bv,bv->matrix,bv->cached);
172: } else {
173: BVCopy(bv,bv->cached);
174: }
175: bv->bvstate = ((PetscObject)bv)->state;
176: }
177: return(0);
178: }
180: /*
181: BV_AllocateCoeffs - Allocate orthogonalization coefficients if not done already.
182: */
183: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateCoeffs(BV bv)
184: {
188: if (!bv->h) {
189: PetscMalloc2(bv->nc+bv->m,&bv->h,bv->nc+bv->m,&bv->c);
190: PetscLogObjectMemory((PetscObject)bv,2*bv->m*sizeof(PetscScalar));
191: }
192: return(0);
193: }
195: /*
196: BV_AllocateSignature - Allocate signature coefficients if not done already.
197: */
198: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateSignature(BV bv)
199: {
203: if (bv->indef && !bv->omega) {
204: if (bv->cuda) {
205: #if defined(PETSC_HAVE_CUDA)
206: VecCreateSeqCUDA(PETSC_COMM_SELF,bv->nc+bv->m,&bv->omega);
207: #else
208: SETERRQ(PetscObjectComm((PetscObject)bv),1,"Something wrong happened");
209: #endif
210: } else {
211: VecCreateSeq(PETSC_COMM_SELF,bv->nc+bv->m,&bv->omega);
212: }
213: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->omega);
214: VecSet(bv->omega,1.0);
215: }
216: return(0);
217: }
219: /*
220: BVAvailableVec: First (0) or second (1) vector available for
221: getcolumn operation (or -1 if both vectors already fetched).
222: */
223: #define BVAvailableVec (((bv->ci[0]==-bv->nc-1)? 0: (bv->ci[1]==-bv->nc-1)? 1: -1))
225: /*
226: Macros to test valid BV arguments
227: */
228: #if !defined(PETSC_USE_DEBUG)
230: #define BVCheckSizes(h,arg) do {} while (0)
231: #define BVCheckOp(h,arg,op) do {} while (0)
233: #else
235: #define BVCheckSizes(h,arg) \
236: do { \
237: if (!(h)->m) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"BV sizes have not been defined: Parameter #%d",arg); \
238: } while (0)
240: #define BVCheckOp(h,arg,op) \
241: do { \
242: if (!(h)->ops->op) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_SUP,"Operation not implemented in this BV type: Parameter #%d",arg); \
243: } while (0)
245: #endif
247: SLEPC_INTERN PetscErrorCode BVView_Vecs(BV,PetscViewer);
249: SLEPC_INTERN PetscErrorCode BVAllocateWork_Private(BV,PetscInt);
251: SLEPC_INTERN PetscErrorCode BVMult_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
252: SLEPC_INTERN PetscErrorCode BVMultVec_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
253: SLEPC_INTERN PetscErrorCode BVMultInPlace_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar*,const PetscScalar*,PetscBool);
254: SLEPC_INTERN PetscErrorCode BVMultInPlace_Vecs_Private(BV,PetscInt,PetscInt,PetscInt,Vec*,const PetscScalar*,PetscBool);
255: SLEPC_INTERN PetscErrorCode BVAXPY_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,PetscScalar,PetscScalar*);
256: SLEPC_INTERN PetscErrorCode BVDot_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
257: SLEPC_INTERN PetscErrorCode BVDotVec_BLAS_Private(BV,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
258: SLEPC_INTERN PetscErrorCode BVScale_BLAS_Private(BV,PetscInt,PetscScalar*,PetscScalar);
259: SLEPC_INTERN PetscErrorCode BVNorm_LAPACK_Private(BV,PetscInt,PetscInt,const PetscScalar*,NormType,PetscReal*,PetscBool);
260: SLEPC_INTERN PetscErrorCode BVNormalize_LAPACK_Private(BV,PetscInt,PetscInt,const PetscScalar*,PetscScalar*,PetscBool);
261: SLEPC_INTERN PetscErrorCode BVMatCholInv_LAPACK_Private(BV,Mat,Mat);
262: SLEPC_INTERN PetscErrorCode BVMatTriInv_LAPACK_Private(BV,Mat,Mat);
263: SLEPC_INTERN PetscErrorCode BVMatSVQB_LAPACK_Private(BV,Mat,Mat);
264: SLEPC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_TSQR(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscInt);
265: SLEPC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_TSQR_OnlyR(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscInt);
267: /* reduction operations used in BVOrthogonalize and BVNormalize */
268: SLEPC_EXTERN MPI_Op MPIU_TSQR, MPIU_LAPY2;
269: SLEPC_EXTERN void MPIAPI SlepcGivensPacked(void*,void*,PetscMPIInt*,MPI_Datatype*);
270: SLEPC_EXTERN void MPIAPI SlepcPythag(void*,void*,PetscMPIInt*,MPI_Datatype*);
272: #if defined(PETSC_HAVE_CUDA)
274: /* complex single */
275: #if defined(PETSC_USE_COMPLEX)
276: #if defined(PETSC_USE_REAL_SINGLE)
277: #define cublasXgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n) cublasCgemm((a),(b),(c),(d),(e),(f),(cuComplex*)(g),(cuComplex*)(h),(i),(cuComplex*)(j),(k),(cuComplex*)(l),(cuComplex*)(m),(n))
278: #define cublasXgemv(a,b,c,d,e,f,g,h,i,j,k,l) cublasCgemv((a),(b),(c),(d),(cuComplex*)(e),(cuComplex*)(f),(g),(cuComplex*)(h),(i),(cuComplex*)(j),(cuComplex*)(k),(l))
279: #define cublasXscal(a,b,c,d,e) cublasCscal(a,b,(const cuComplex*)(c),(cuComplex*)(d),e)
280: #define cublasXnrm2(a,b,c,d,e) cublasScnrm2(a,b,(const cuComplex*)(c),d,e)
281: #define cublasXaxpy(a,b,c,d,e,f,g) cublasCaxpy((a),(b),(cuComplex*)(c),(cuComplex*)(d),(e),(cuComplex*)(f),(g))
282: #define cublasXdotc(a,b,c,d,e,f,g) cublasCdotc((a),(b),(const cuComplex *)(c),(d),(const cuComplex *)(e),(f),(cuComplex *)(g))
283: #else /* complex double */
284: #define cublasXgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n) cublasZgemm((a),(b),(c),(d),(e),(f),(cuDoubleComplex*)(g),(cuDoubleComplex*)(h),(i),(cuDoubleComplex*)(j),(k),(cuDoubleComplex*)(l),(cuDoubleComplex*)(m),(n))
285: #define cublasXgemv(a,b,c,d,e,f,g,h,i,j,k,l) cublasZgemv((a),(b),(c),(d),(cuDoubleComplex*)(e),(cuDoubleComplex*)(f),(g),(cuDoubleComplex*)(h),(i),(cuDoubleComplex*)(j),(cuDoubleComplex*)(k),(l))
286: #define cublasXscal(a,b,c,d,e) cublasZscal(a,b,(const cuDoubleComplex*)(c),(cuDoubleComplex*)(d),e)
287: #define cublasXnrm2(a,b,c,d,e) cublasDznrm2(a,b,(const cuDoubleComplex*)(c),d,e)
288: #define cublasXaxpy(a,b,c,d,e,f,g) cublasZaxpy((a),(b),(cuDoubleComplex*)(c),(cuDoubleComplex*)(d),(e),(cuDoubleComplex*)(f),(g))
289: #define cublasXdotc(a,b,c,d,e,f,g) cublasZdotc((a),(b),(const cuDoubleComplex *)(c),(d),(const cuDoubleComplex *)(e),(f),(cuDoubleComplex *)(g))
290: #endif
291: #else /* real single */
292: #if defined(PETSC_USE_REAL_SINGLE)
293: #define cublasXgemm cublasSgemm
294: #define cublasXgemv cublasSgemv
295: #define cublasXscal cublasSscal
296: #define cublasXnrm2 cublasSnrm2
297: #define cublasXaxpy cublasSaxpy
298: #define cublasXdotc cublasSdot
299: #else /* real double */
300: #define cublasXgemm cublasDgemm
301: #define cublasXgemv cublasDgemv
302: #define cublasXscal cublasDscal
303: #define cublasXnrm2 cublasDnrm2
304: #define cublasXaxpy cublasDaxpy
305: #define cublasXdotc cublasDdot
306: #endif
307: #endif
309: SLEPC_INTERN PetscErrorCode BV_CleanCoefficients_CUDA(BV,PetscInt,PetscScalar*);
310: SLEPC_INTERN PetscErrorCode BV_AddCoefficients_CUDA(BV,PetscInt,PetscScalar*,PetscScalar*);
311: SLEPC_INTERN PetscErrorCode BV_SetValue_CUDA(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar);
312: SLEPC_INTERN PetscErrorCode BV_SquareSum_CUDA(BV,PetscInt,PetscScalar*,PetscReal*);
313: SLEPC_INTERN PetscErrorCode BV_ApplySignature_CUDA(BV,PetscInt,PetscScalar*,PetscBool);
314: SLEPC_INTERN PetscErrorCode BV_SquareRoot_CUDA(BV,PetscInt,PetscScalar*,PetscReal*);
315: SLEPC_INTERN PetscErrorCode BV_StoreCoefficients_CUDA(BV,PetscInt,PetscScalar*,PetscScalar*);
317: #endif /* PETSC_HAVE_CUDA */
319: /*
320: BV_CleanCoefficients_Default - Sets to zero all entries of column j of the bv buffer
321: */
322: PETSC_STATIC_INLINE PetscErrorCode BV_CleanCoefficients_Default(BV bv,PetscInt j,PetscScalar *h)
323: {
325: PetscScalar *hh=h,*a;
326: PetscInt i;
329: if (!h) {
330: VecGetArray(bv->buffer,&a);
331: hh = a + j*(bv->nc+bv->m);
332: }
333: for (i=0;i<bv->nc+j;i++) hh[i] = 0.0;
334: if (!h) { VecRestoreArray(bv->buffer,&a); }
335: return(0);
336: }
338: /*
339: BV_AddCoefficients_Default - Add the contents of the scratch (0-th column) of the bv buffer
340: into column j of the bv buffer
341: */
342: PETSC_STATIC_INLINE PetscErrorCode BV_AddCoefficients_Default(BV bv,PetscInt j,PetscScalar *h,PetscScalar *c)
343: {
345: PetscScalar *hh=h,*cc=c;
346: PetscInt i;
349: if (!h) {
350: VecGetArray(bv->buffer,&cc);
351: hh = cc + j*(bv->nc+bv->m);
352: }
353: for (i=0;i<bv->nc+j;i++) hh[i] += cc[i];
354: if (!h) { VecRestoreArray(bv->buffer,&cc); }
355: return(0);
356: }
358: /*
359: BV_SetValue_Default - Sets value in row j (counted after the constraints) of column k
360: of the coefficients array
361: */
362: PETSC_STATIC_INLINE PetscErrorCode BV_SetValue_Default(BV bv,PetscInt j,PetscInt k,PetscScalar *h,PetscScalar value)
363: {
365: PetscScalar *hh=h,*a;
368: if (!h) {
369: VecGetArray(bv->buffer,&a);
370: hh = a + k*(bv->nc+bv->m);
371: }
372: hh[bv->nc+j] = value;
373: if (!h) { VecRestoreArray(bv->buffer,&a); }
374: return(0);
375: }
377: /*
378: BV_SquareSum_Default - Returns the value h'*h, where h represents the contents of the
379: coefficients array (up to position j)
380: */
381: PETSC_STATIC_INLINE PetscErrorCode BV_SquareSum_Default(BV bv,PetscInt j,PetscScalar *h,PetscReal *sum)
382: {
384: PetscScalar *hh=h;
385: PetscInt i;
388: *sum = 0.0;
389: if (!h) { VecGetArray(bv->buffer,&hh); }
390: for (i=0;i<bv->nc+j;i++) *sum += PetscRealPart(hh[i]*PetscConj(hh[i]));
391: if (!h) { VecRestoreArray(bv->buffer,&hh); }
392: return(0);
393: }
395: /*
396: BV_ApplySignature_Default - Computes the pointwise product h*omega, where h represents
397: the contents of the coefficients array (up to position j) and omega is the signature;
398: if inverse=TRUE then the operation is h/omega
399: */
400: PETSC_STATIC_INLINE PetscErrorCode BV_ApplySignature_Default(BV bv,PetscInt j,PetscScalar *h,PetscBool inverse)
401: {
402: PetscErrorCode ierr;
403: PetscScalar *hh=h;
404: PetscInt i;
405: const PetscScalar *omega;
408: if (!(bv->nc+j)) return(0);
409: if (!h) { VecGetArray(bv->buffer,&hh); }
410: VecGetArrayRead(bv->omega,&omega);
411: if (inverse) for (i=0;i<bv->nc+j;i++) hh[i] /= PetscRealPart(omega[i]);
412: else for (i=0;i<bv->nc+j;i++) hh[i] *= PetscRealPart(omega[i]);
413: VecRestoreArrayRead(bv->omega,&omega);
414: if (!h) { VecRestoreArray(bv->buffer,&hh); }
415: return(0);
416: }
418: /*
419: BV_SquareRoot_Default - Returns the square root of position j (counted after the constraints)
420: of the coefficients array
421: */
422: PETSC_STATIC_INLINE PetscErrorCode BV_SquareRoot_Default(BV bv,PetscInt j,PetscScalar *h,PetscReal *beta)
423: {
425: PetscScalar *hh=h;
428: if (!h) { VecGetArray(bv->buffer,&hh); }
429: BV_SafeSqrt(bv,hh[bv->nc+j],beta);
430: if (!h) { VecRestoreArray(bv->buffer,&hh); }
431: return(0);
432: }
434: /*
435: BV_StoreCoefficients_Default - Copy the contents of the coefficients array to an array dest
436: provided by the caller (only values from l to j are copied)
437: */
438: PETSC_STATIC_INLINE PetscErrorCode BV_StoreCoefficients_Default(BV bv,PetscInt j,PetscScalar *h,PetscScalar *dest)
439: {
441: PetscScalar *hh=h,*a;
442: PetscInt i;
445: if (!h) {
446: VecGetArray(bv->buffer,&a);
447: hh = a + j*(bv->nc+bv->m);
448: }
449: for (i=bv->l;i<j;i++) dest[i-bv->l] = hh[bv->nc+i];
450: if (!h) { VecRestoreArray(bv->buffer,&a); }
451: return(0);
452: }
454: /*
455: BV_GetEigenvector - retrieves k-th eigenvector from basis vectors V.
456: The argument eigi is the imaginary part of the corresponding eigenvalue.
457: */
458: PETSC_STATIC_INLINE PetscErrorCode BV_GetEigenvector(BV V,PetscInt k,PetscScalar eigi,Vec Vr,Vec Vi)
459: {
463: #if defined(PETSC_USE_COMPLEX)
464: if (Vr) { BVCopyVec(V,k,Vr); }
465: if (Vi) { VecSet(Vi,0.0); }
466: #else
467: if (eigi > 0.0) { /* first value of conjugate pair */
468: if (Vr) { BVCopyVec(V,k,Vr); }
469: if (Vi) { BVCopyVec(V,k+1,Vi); }
470: } else if (eigi < 0.0) { /* second value of conjugate pair */
471: if (Vr) { BVCopyVec(V,k-1,Vr); }
472: if (Vi) {
473: BVCopyVec(V,k,Vi);
474: VecScale(Vi,-1.0);
475: }
476: } else { /* real eigenvalue */
477: if (Vr) { BVCopyVec(V,k,Vr); }
478: if (Vi) { VecSet(Vi,0.0); }
479: }
480: #endif
481: return(0);
482: }
484: /*
485: BV_OrthogonalizeColumn_Safe - this is intended for cases where we know that
486: the resulting vector is going to be numerically zero, so normalization or
487: iterative refinement may cause problems in parallel (collective operation
488: not being called by all processes)
489: */
490: PETSC_STATIC_INLINE PetscErrorCode BV_OrthogonalizeColumn_Safe(BV bv,PetscInt j,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
491: {
492: PetscErrorCode ierr;
493: BVOrthogRefineType orthog_ref;
496: PetscInfo1(bv,"Orthogonalizing column %D without refinement\n",j);
497: orthog_ref = bv->orthog_ref;
498: bv->orthog_ref = BV_ORTHOG_REFINE_NEVER; /* avoid refinement */
499: BVOrthogonalizeColumn(bv,j,H,NULL,NULL);
500: bv->orthog_ref = orthog_ref; /* restore refinement setting */
501: if (norm) *norm = 0.0;
502: if (lindep) *lindep = PETSC_TRUE;
503: return(0);
504: }
506: #if defined(PETSC_HAVE_CUDA)
507: #define BV_CleanCoefficients(a,b,c) ((a)->cuda?BV_CleanCoefficients_CUDA:BV_CleanCoefficients_Default)((a),(b),(c))
508: #define BV_AddCoefficients(a,b,c,d) ((a)->cuda?BV_AddCoefficients_CUDA:BV_AddCoefficients_Default)((a),(b),(c),(d))
509: #define BV_SetValue(a,b,c,d,e) ((a)->cuda?BV_SetValue_CUDA:BV_SetValue_Default)((a),(b),(c),(d),(e))
510: #define BV_SquareSum(a,b,c,d) ((a)->cuda?BV_SquareSum_CUDA:BV_SquareSum_Default)((a),(b),(c),(d))
511: #define BV_ApplySignature(a,b,c,d) ((a)->cuda?BV_ApplySignature_CUDA:BV_ApplySignature_Default)((a),(b),(c),(d))
512: #define BV_SquareRoot(a,b,c,d) ((a)->cuda?BV_SquareRoot_CUDA:BV_SquareRoot_Default)((a),(b),(c),(d))
513: #define BV_StoreCoefficients(a,b,c,d) ((a)->cuda?BV_StoreCoefficients_CUDA:BV_StoreCoefficients_Default)((a),(b),(c),(d))
514: #else
515: #define BV_CleanCoefficients(a,b,c) BV_CleanCoefficients_Default((a),(b),(c))
516: #define BV_AddCoefficients(a,b,c,d) BV_AddCoefficients_Default((a),(b),(c),(d))
517: #define BV_SetValue(a,b,c,d,e) BV_SetValue_Default((a),(b),(c),(d),(e))
518: #define BV_SquareSum(a,b,c,d) BV_SquareSum_Default((a),(b),(c),(d))
519: #define BV_ApplySignature(a,b,c,d) BV_ApplySignature_Default((a),(b),(c),(d))
520: #define BV_SquareRoot(a,b,c,d) BV_SquareRoot_Default((a),(b),(c),(d))
521: #define BV_StoreCoefficients(a,b,c,d) BV_StoreCoefficients_Default((a),(b),(c),(d))
522: #endif /* PETSC_HAVE_CUDA */
524: #endif