Actual source code: bvimpl.h
slepc-3.9.2 2018-07-02
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: */
11: #if !defined(_BVIMPL)
12: #define _BVIMPL
14: #include <slepcbv.h>
15: #include <slepc/private/slepcimpl.h>
17: PETSC_EXTERN PetscBool BVRegisterAllCalled;
18: PETSC_EXTERN PetscErrorCode BVRegisterAll(void);
20: 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;
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 (*matmult)(BV,Mat,BV);
40: PetscErrorCode (*copy)(BV,BV);
41: PetscErrorCode (*copycolumn)(BV,PetscInt,PetscInt);
42: PetscErrorCode (*resize)(BV,PetscInt,PetscBool);
43: PetscErrorCode (*getcolumn)(BV,PetscInt,Vec*);
44: PetscErrorCode (*restorecolumn)(BV,PetscInt,Vec*);
45: PetscErrorCode (*getarray)(BV,PetscScalar**);
46: PetscErrorCode (*restorearray)(BV,PetscScalar**);
47: PetscErrorCode (*getarrayread)(BV,const PetscScalar**);
48: PetscErrorCode (*restorearrayread)(BV,const PetscScalar**);
49: PetscErrorCode (*restoresplit)(BV,BV*,BV*);
50: PetscErrorCode (*gramschmidt)(BV,PetscInt,Vec,PetscBool*,PetscScalar*,PetscScalar*,PetscReal*,PetscReal*);
51: PetscErrorCode (*duplicate)(BV,BV);
52: PetscErrorCode (*create)(BV);
53: PetscErrorCode (*setfromoptions)(PetscOptionItems*,BV);
54: PetscErrorCode (*view)(BV,PetscViewer);
55: PetscErrorCode (*destroy)(BV);
56: };
58: struct _p_BV {
59: PETSCHEADER(struct _BVOps);
60: /*------------------------- User parameters --------------------------*/
61: Vec t; /* template vector */
62: PetscInt n,N; /* dimensions of vectors (local, global) */
63: PetscInt m; /* number of vectors */
64: PetscInt l; /* number of leading columns */
65: PetscInt k; /* number of active columns */
66: PetscInt nc; /* number of constraints */
67: BVOrthogType orthog_type; /* the method of vector orthogonalization */
68: BVOrthogRefineType orthog_ref; /* refinement method */
69: PetscReal orthog_eta; /* refinement threshold */
70: BVOrthogBlockType orthog_block; /* the method of block orthogonalization */
71: Mat matrix; /* inner product matrix */
72: PetscBool indef; /* matrix is indefinite */
73: BVMatMultType vmm; /* version of matmult operation */
74: PetscBool rrandom; /* reproducible random vectors */
76: /*---------------------- Cached data and workspace -------------------*/
77: Vec buffer; /* buffer vector used in orthogonalization */
78: Mat Abuffer; /* auxiliary seqdense matrix that wraps the buffer */
79: Vec Bx; /* result of matrix times a vector x */
80: PetscObjectId xid; /* object id of vector x */
81: PetscObjectState xstate; /* state of vector x */
82: Vec cv[2]; /* column vectors obtained with BVGetColumn() */
83: PetscInt ci[2]; /* column indices of obtained vectors */
84: PetscObjectState st[2]; /* state of obtained vectors */
85: PetscObjectId id[2]; /* object id of obtained vectors */
86: PetscScalar *h,*c; /* orthogonalization coefficients */
87: Vec omega; /* signature matrix values for indefinite case */
88: PetscBool defersfo; /* deferred call to setfromoptions */
89: BV cached; /* cached BV to store result of matrix times BV */
90: PetscObjectState bvstate; /* state of BV when BVApplyMatrixBV() was called */
91: BV L,R; /* BV objects obtained with BVGetSplit() */
92: PetscObjectState lstate,rstate;/* state of L and R when BVGetSplit() was called */
93: PetscInt lsplit; /* the value of l when BVGetSplit() was called */
94: PetscInt issplit; /* >0 if this BV has been created by splitting (1=left, 2=right) */
95: BV splitparent; /* my parent if I am a split BV */
96: PetscRandom rand; /* random number generator */
97: Mat Acreate; /* matrix given at BVCreateFromMat() */
98: Mat Aget; /* matrix returned for BVGetMat() */
99: PetscBool cuda; /* true if GPU must be used in SVEC */
100: PetscScalar *work;
101: PetscInt lwork;
102: void *data;
103: };
105: /*
106: BV_SafeSqrt - Computes the square root of a scalar value alpha, which is
107: assumed to be z'*B*z. The result is
108: if definite inner product: res = sqrt(alpha)
109: if indefinite inner product: res = sgn(alpha)*sqrt(abs(alpha))
110: */
111: PETSC_STATIC_INLINE PetscErrorCode BV_SafeSqrt(BV bv,PetscScalar alpha,PetscReal *res)
112: {
114: PetscReal absal,realp;
117: absal = PetscAbsScalar(alpha);
118: realp = PetscRealPart(alpha);
119: if (absal<PETSC_MACHINE_EPSILON) {
120: PetscInfo(bv,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
121: }
122: #if defined(PETSC_USE_COMPLEX)
123: if (PetscAbsReal(PetscImaginaryPart(alpha))>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));
124: #endif
125: if (bv->indef) {
126: *res = (realp<0.0)? -PetscSqrtReal(-realp): PetscSqrtReal(realp);
127: } else {
128: if (realp<-10*PETSC_MACHINE_EPSILON) SETERRQ(PetscObjectComm((PetscObject)bv),1,"The inner product is not well defined: indefinite matrix");
129: *res = (realp<0.0)? 0.0: PetscSqrtReal(realp);
130: }
131: return(0);
132: }
134: /*
135: BV_IPMatMult - Multiply a vector x by the inner-product matrix, cache the
136: result in Bx.
137: */
138: PETSC_STATIC_INLINE PetscErrorCode BV_IPMatMult(BV bv,Vec x)
139: {
143: if (((PetscObject)x)->id != bv->xid || ((PetscObject)x)->state != bv->xstate) {
144: if (!bv->Bx) {
145: MatCreateVecs(bv->matrix,&bv->Bx,NULL);
146: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
147: }
148: MatMult(bv->matrix,x,bv->Bx);
149: bv->xid = ((PetscObject)x)->id;
150: bv->xstate = ((PetscObject)x)->state;
151: }
152: return(0);
153: }
155: /*
156: BV_IPMatMultBV - Multiply BV by the inner-product matrix, cache the
157: result internally in bv->cached.
158: */
159: PETSC_STATIC_INLINE PetscErrorCode BV_IPMatMultBV(BV bv)
160: {
164: BVGetCachedBV(bv,&bv->cached);
165: if (((PetscObject)bv)->state != bv->bvstate || bv->l != bv->cached->l || bv->k != bv->cached->k) {
166: BVSetActiveColumns(bv->cached,bv->l,bv->k);
167: if (bv->matrix) {
168: BVMatMult(bv,bv->matrix,bv->cached);
169: } else {
170: BVCopy(bv,bv->cached);
171: }
172: bv->bvstate = ((PetscObject)bv)->state;
173: }
174: return(0);
175: }
177: /*
178: BV_AllocateCoeffs - Allocate orthogonalization coefficients if not done already.
179: */
180: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateCoeffs(BV bv)
181: {
185: if (!bv->h) {
186: PetscMalloc2(bv->nc+bv->m,&bv->h,bv->nc+bv->m,&bv->c);
187: PetscLogObjectMemory((PetscObject)bv,2*bv->m*sizeof(PetscScalar));
188: }
189: return(0);
190: }
192: /*
193: BV_AllocateSignature - Allocate signature coefficients if not done already.
194: */
195: PETSC_STATIC_INLINE PetscErrorCode BV_AllocateSignature(BV bv)
196: {
200: if (bv->indef && !bv->omega) {
201: if (bv->cuda) {
202: #if defined(PETSC_HAVE_VECCUDA)
203: VecCreateSeqCUDA(PETSC_COMM_SELF,bv->nc+bv->m,&bv->omega);
204: #else
205: SETERRQ(PetscObjectComm((PetscObject)bv),1,"Something wrong happened");
206: #endif
207: } else {
208: VecCreateSeq(PETSC_COMM_SELF,bv->nc+bv->m,&bv->omega);
209: }
210: PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->omega);
211: VecSet(bv->omega,1.0);
212: }
213: return(0);
214: }
216: /*
217: BVAvailableVec: First (0) or second (1) vector available for
218: getcolumn operation (or -1 if both vectors already fetched).
219: */
220: #define BVAvailableVec (((bv->ci[0]==-bv->nc-1)? 0: (bv->ci[1]==-bv->nc-1)? 1: -1))
222: /*
223: Macros to test valid BV arguments
224: */
225: #if !defined(PETSC_USE_DEBUG)
227: #define BVCheckSizes(h,arg) do {} while (0)
229: #else
231: #define BVCheckSizes(h,arg) \
232: do { \
233: if (!h->m) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"BV sizes have not been defined: Parameter #%d",arg); \
234: } while (0)
236: #endif
238: PETSC_INTERN PetscErrorCode BVView_Vecs(BV,PetscViewer);
240: PETSC_INTERN PetscErrorCode BVAllocateWork_Private(BV,PetscInt);
242: PETSC_INTERN PetscErrorCode BVMult_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
243: PETSC_INTERN PetscErrorCode BVMultVec_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
244: PETSC_INTERN PetscErrorCode BVMultInPlace_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar*,const PetscScalar*,PetscBool);
245: PETSC_INTERN PetscErrorCode BVMultInPlace_Vecs_Private(BV,PetscInt,PetscInt,PetscInt,Vec*,const PetscScalar*,PetscBool);
246: PETSC_INTERN PetscErrorCode BVAXPY_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,PetscScalar,PetscScalar*);
247: PETSC_INTERN PetscErrorCode BVDot_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
248: PETSC_INTERN PetscErrorCode BVDotVec_BLAS_Private(BV,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
249: PETSC_INTERN PetscErrorCode BVScale_BLAS_Private(BV,PetscInt,PetscScalar*,PetscScalar);
250: PETSC_INTERN PetscErrorCode BVNorm_LAPACK_Private(BV,PetscInt,PetscInt,const PetscScalar*,NormType,PetscReal*,PetscBool);
251: PETSC_INTERN PetscErrorCode BVMatCholInv_LAPACK_Private(BV,Mat,Mat);
252: PETSC_INTERN PetscErrorCode BVMatTriInv_LAPACK_Private(BV,Mat,Mat);
253: PETSC_INTERN PetscErrorCode BVMatSVQB_LAPACK_Private(BV,Mat,Mat);
254: PETSC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_TSQR(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscInt);
255: PETSC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_TSQR_OnlyR(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscInt);
257: /* reduction operation used in BVOrthogonalize */
258: PETSC_EXTERN MPI_Op MPIU_TSQR;
259: PETSC_EXTERN void SlepcGivensPacked(void*,void*,PetscMPIInt*,MPI_Datatype*);
261: #if defined(PETSC_HAVE_VECCUDA)
262: #include <petsccuda.h>
263: #include <cublas_v2.h>
265: #define WaitForGPU() PetscCUDASynchronize ? cudaThreadSynchronize() : 0
267: /* complex single */
268: #if defined(PETSC_USE_COMPLEX)
269: #if defined(PETSC_USE_REAL_SINGLE)
270: #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))
271: #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))
272: #define cublasXscal(a,b,c,d,e) cublasCscal(a,b,(const cuComplex*)(c),(cuComplex*)(d),e)
273: #define cublasXnrm2(a,b,c,d,e) cublasScnrm2(a,b,(const cuComplex*)(c),d,e)
274: #define cublasXaxpy(a,b,c,d,e,f,g) cublasCaxpy((a),(b),(cuComplex*)(c),(cuComplex*)(d),(e),(cuComplex*)(f),(g))
275: #define cublasXdotc(a,b,c,d,e,f,g) cublasCdotc((a),(b),(const cuComplex *)(c),(d),(const cuComplex *)(e),(f),(cuComplex *)(g))
276: #else /* complex double */
277: #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))
278: #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))
279: #define cublasXscal(a,b,c,d,e) cublasZscal(a,b,(const cuDoubleComplex*)(c),(cuDoubleComplex*)(d),e)
280: #define cublasXnrm2(a,b,c,d,e) cublasDznrm2(a,b,(const cuDoubleComplex*)(c),d,e)
281: #define cublasXaxpy(a,b,c,d,e,f,g) cublasZaxpy((a),(b),(cuDoubleComplex*)(c),(cuDoubleComplex*)(d),(e),(cuDoubleComplex*)(f),(g))
282: #define cublasXdotc(a,b,c,d,e,f,g) cublasZdotc((a),(b),(const cuDoubleComplex *)(c),(d),(const cuDoubleComplex *)(e),(f),(cuDoubleComplex *)(g))
283: #endif
284: #else /* real single */
285: #if defined(PETSC_USE_REAL_SINGLE)
286: #define cublasXgemm cublasSgemm
287: #define cublasXgemv cublasSgemv
288: #define cublasXscal cublasSscal
289: #define cublasXnrm2 cublasSnrm2
290: #define cublasXaxpy cublasSaxpy
291: #define cublasXdotc cublasSdot
292: #else /* real double */
293: #define cublasXgemm cublasDgemm
294: #define cublasXgemv cublasDgemv
295: #define cublasXscal cublasDscal
296: #define cublasXnrm2 cublasDnrm2
297: #define cublasXaxpy cublasDaxpy
298: #define cublasXdotc cublasDdot
299: #endif
300: #endif
302: PETSC_INTERN PetscErrorCode BV_CleanCoefficients_CUDA(BV,PetscInt,PetscScalar*);
303: PETSC_INTERN PetscErrorCode BV_AddCoefficients_CUDA(BV,PetscInt,PetscScalar*,PetscScalar*);
304: PETSC_INTERN PetscErrorCode BV_SetValue_CUDA(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar);
305: PETSC_INTERN PetscErrorCode BV_SquareSum_CUDA(BV,PetscInt,PetscScalar*,PetscReal*);
306: PETSC_INTERN PetscErrorCode BV_ApplySignature_CUDA(BV,PetscInt,PetscScalar*,PetscBool);
307: PETSC_INTERN PetscErrorCode BV_SquareRoot_CUDA(BV,PetscInt,PetscScalar*,PetscReal*);
308: PETSC_INTERN PetscErrorCode BV_StoreCoefficients_CUDA(BV,PetscInt,PetscScalar*,PetscScalar*);
310: #endif /* PETSC_HAVE_VECCUDA */
312: /*
313: BV_CleanCoefficients_Default - Sets to zero all entries of column j of the bv buffer
314: */
315: PETSC_STATIC_INLINE PetscErrorCode BV_CleanCoefficients_Default(BV bv,PetscInt j,PetscScalar *h)
316: {
318: PetscScalar *hh=h,*a;
319: PetscInt i;
322: if (!h) {
323: VecGetArray(bv->buffer,&a);
324: hh = a + j*(bv->nc+bv->m);
325: }
326: for (i=0;i<bv->nc+j;i++) hh[i] = 0.0;
327: if (!h) { VecRestoreArray(bv->buffer,&a); }
328: return(0);
329: }
331: /*
332: BV_AddCoefficients_Default - Add the contents of the scratch (0-th column) of the bv buffer
333: into column j of the bv buffer
334: */
335: PETSC_STATIC_INLINE PetscErrorCode BV_AddCoefficients_Default(BV bv,PetscInt j,PetscScalar *h,PetscScalar *c)
336: {
338: PetscScalar *hh=h,*cc=c;
339: PetscInt i;
342: if (!h) {
343: VecGetArray(bv->buffer,&cc);
344: hh = cc + j*(bv->nc+bv->m);
345: }
346: for (i=0;i<bv->nc+j;i++) hh[i] += cc[i];
347: if (!h) { VecRestoreArray(bv->buffer,&cc); }
348: return(0);
349: }
351: /*
352: BV_SetValue_Default - Sets value in row j (counted after the constraints) of column k
353: of the coefficients array
354: */
355: PETSC_STATIC_INLINE PetscErrorCode BV_SetValue_Default(BV bv,PetscInt j,PetscInt k,PetscScalar *h,PetscScalar value)
356: {
358: PetscScalar *hh=h,*a;
361: if (!h) {
362: VecGetArray(bv->buffer,&a);
363: hh = a + k*(bv->nc+bv->m);
364: }
365: hh[bv->nc+j] = value;
366: if (!h) { VecRestoreArray(bv->buffer,&a); }
367: return(0);
368: }
370: /*
371: BV_SquareSum_Default - Returns the value h'*h, where h represents the contents of the
372: coefficients array (up to position j)
373: */
374: PETSC_STATIC_INLINE PetscErrorCode BV_SquareSum_Default(BV bv,PetscInt j,PetscScalar *h,PetscReal *sum)
375: {
377: PetscScalar *hh=h;
378: PetscInt i;
381: *sum = 0.0;
382: if (!h) { VecGetArray(bv->buffer,&hh); }
383: for (i=0;i<bv->nc+j;i++) *sum += PetscRealPart(hh[i]*PetscConj(hh[i]));
384: if (!h) { VecRestoreArray(bv->buffer,&hh); }
385: return(0);
386: }
388: /*
389: BV_ApplySignature_Default - Computes the pointwise product h*omega, where h represents
390: the contents of the coefficients array (up to position j) and omega is the signature;
391: if inverse=TRUE then the operation is h/omega
392: */
393: PETSC_STATIC_INLINE PetscErrorCode BV_ApplySignature_Default(BV bv,PetscInt j,PetscScalar *h,PetscBool inverse)
394: {
395: PetscErrorCode ierr;
396: PetscScalar *hh=h;
397: PetscInt i;
398: const PetscScalar *omega;
401: if (!(bv->nc+j)) return(0);
402: if (!h) { VecGetArray(bv->buffer,&hh); }
403: VecGetArrayRead(bv->omega,&omega);
404: if (inverse) for (i=0;i<bv->nc+j;i++) hh[i] /= PetscRealPart(omega[i]);
405: else for (i=0;i<bv->nc+j;i++) hh[i] *= PetscRealPart(omega[i]);
406: VecRestoreArrayRead(bv->omega,&omega);
407: if (!h) { VecRestoreArray(bv->buffer,&hh); }
408: return(0);
409: }
411: /*
412: BV_SquareRoot_Default - Returns the square root of position j (counted after the constraints)
413: of the coefficients array
414: */
415: PETSC_STATIC_INLINE PetscErrorCode BV_SquareRoot_Default(BV bv,PetscInt j,PetscScalar *h,PetscReal *beta)
416: {
418: PetscScalar *hh=h;
421: if (!h) { VecGetArray(bv->buffer,&hh); }
422: BV_SafeSqrt(bv,hh[bv->nc+j],beta);
423: if (!h) { VecRestoreArray(bv->buffer,&hh); }
424: return(0);
425: }
427: /*
428: BV_StoreCoefficients_Default - Copy the contents of the coefficients array to an array dest
429: provided by the caller (only values from l to j are copied)
430: */
431: PETSC_STATIC_INLINE PetscErrorCode BV_StoreCoefficients_Default(BV bv,PetscInt j,PetscScalar *h,PetscScalar *dest)
432: {
434: PetscScalar *hh=h,*a;
435: PetscInt i;
438: if (!h) {
439: VecGetArray(bv->buffer,&a);
440: hh = a + j*(bv->nc+bv->m);
441: }
442: for (i=bv->l;i<j;i++) dest[i-bv->l] = hh[bv->nc+i];
443: if (!h) { VecRestoreArray(bv->buffer,&a); }
444: return(0);
445: }
447: #if defined(PETSC_HAVE_VECCUDA)
448: #define BV_CleanCoefficients(a,b,c) ((a)->cuda?BV_CleanCoefficients_CUDA:BV_CleanCoefficients_Default)((a),(b),(c))
449: #define BV_AddCoefficients(a,b,c,d) ((a)->cuda?BV_AddCoefficients_CUDA:BV_AddCoefficients_Default)((a),(b),(c),(d))
450: #define BV_SetValue(a,b,c,d,e) ((a)->cuda?BV_SetValue_CUDA:BV_SetValue_Default)((a),(b),(c),(d),(e))
451: #define BV_SquareSum(a,b,c,d) ((a)->cuda?BV_SquareSum_CUDA:BV_SquareSum_Default)((a),(b),(c),(d))
452: #define BV_ApplySignature(a,b,c,d) ((a)->cuda?BV_ApplySignature_CUDA:BV_ApplySignature_Default)((a),(b),(c),(d))
453: #define BV_SquareRoot(a,b,c,d) ((a)->cuda?BV_SquareRoot_CUDA:BV_SquareRoot_Default)((a),(b),(c),(d))
454: #define BV_StoreCoefficients(a,b,c,d) ((a)->cuda?BV_StoreCoefficients_CUDA:BV_StoreCoefficients_Default)((a),(b),(c),(d))
455: #else
456: #define BV_CleanCoefficients(a,b,c) BV_CleanCoefficients_Default((a),(b),(c))
457: #define BV_AddCoefficients(a,b,c,d) BV_AddCoefficients_Default((a),(b),(c),(d))
458: #define BV_SetValue(a,b,c,d,e) BV_SetValue_Default((a),(b),(c),(d),(e))
459: #define BV_SquareSum(a,b,c,d) BV_SquareSum_Default((a),(b),(c),(d))
460: #define BV_ApplySignature(a,b,c,d) BV_ApplySignature_Default((a),(b),(c),(d))
461: #define BV_SquareRoot(a,b,c,d) BV_SquareRoot_Default((a),(b),(c),(d))
462: #define BV_StoreCoefficients(a,b,c,d) BV_StoreCoefficients_Default((a),(b),(c),(d))
463: #endif /* PETSC_HAVE_VECCUDA */
465: #endif