Actual source code: bvimpl.h

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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,BV_SVDAndRank;

 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 */
 78:   PetscReal          deftol;       /* tolerance for BV_SafeSqrt */

 80:   /*---------------------- Cached data and workspace -------------------*/
 81:   Vec                buffer;       /* buffer vector used in orthogonalization */
 82:   Mat                Abuffer;      /* auxiliary seqdense matrix that wraps the buffer */
 83:   Vec                Bx;           /* result of matrix times a vector x */
 84:   PetscObjectId      xid;          /* object id of vector x */
 85:   PetscObjectState   xstate;       /* state of vector x */
 86:   Vec                cv[2];        /* column vectors obtained with BVGetColumn() */
 87:   PetscInt           ci[2];        /* column indices of obtained vectors */
 88:   PetscObjectState   st[2];        /* state of obtained vectors */
 89:   PetscObjectId      id[2];        /* object id of obtained vectors */
 90:   PetscScalar        *h,*c;        /* orthogonalization coefficients */
 91:   Vec                omega;        /* signature matrix values for indefinite case */
 92:   PetscBool          defersfo;     /* deferred call to setfromoptions */
 93:   BV                 cached;       /* cached BV to store result of matrix times BV */
 94:   PetscObjectState   bvstate;      /* state of BV when BVApplyMatrixBV() was called */
 95:   BV                 L,R;          /* BV objects obtained with BVGetSplit() */
 96:   PetscObjectState   lstate,rstate;/* state of L and R when BVGetSplit() was called */
 97:   PetscInt           lsplit;       /* the value of l when BVGetSplit() was called */
 98:   PetscInt           issplit;      /* >0 if this BV has been created by splitting (1=left, 2=right) */
 99:   BV                 splitparent;  /* my parent if I am a split BV */
100:   PetscRandom        rand;         /* random number generator */
101:   Mat                Acreate;      /* matrix given at BVCreateFromMat() */
102:   Mat                Aget;         /* matrix returned for BVGetMat() */
103:   PetscBool          cuda;         /* true if GPU must be used in SVEC */
104:   PetscBool          sfocalled;    /* setfromoptions has been called */
105:   PetscScalar        *work;
106:   PetscInt           lwork;
107:   void               *data;
108: };

110: /*
111:   BV_SafeSqrt - Computes the square root of a scalar value alpha, which is
112:   assumed to be z'*B*z. The result is
113:     if definite inner product:     res = sqrt(alpha)
114:     if indefinite inner product:   res = sgn(alpha)*sqrt(abs(alpha))
115: */
116: static inline PetscErrorCode BV_SafeSqrt(BV bv,PetscScalar alpha,PetscReal *res)
117: {
118:   PetscReal      absal,realp;

120:   absal = PetscAbsScalar(alpha);
121:   realp = PetscRealPart(alpha);
122:   if (PetscUnlikely(absal<PETSC_MACHINE_EPSILON)) PetscInfo(bv,"Zero norm, either the vector is zero or a semi-inner product is being used\n");
123: #if defined(PETSC_USE_COMPLEX)
125: #endif
126:   if (PetscUnlikely(bv->indef)) {
127:     *res = (realp<0.0)? -PetscSqrtReal(-realp): PetscSqrtReal(realp);
128:   } else {
130:     *res = (realp<0.0)? 0.0: PetscSqrtReal(realp);
131:   }
132:   PetscFunctionReturn(0);
133: }

135: /*
136:   BV_IPMatMult - Multiply a vector x by the inner-product matrix, cache the
137:   result in Bx.
138: */
139: static inline PetscErrorCode BV_IPMatMult(BV bv,Vec x)
140: {
141:   if (((PetscObject)x)->id != bv->xid || ((PetscObject)x)->state != bv->xstate) {
142:     if (PetscUnlikely(!bv->Bx)) {
143:       MatCreateVecs(bv->matrix,&bv->Bx,NULL);
144:       PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->Bx);
145:     }
146:     MatMult(bv->matrix,x,bv->Bx);
147:     PetscObjectGetId((PetscObject)x,&bv->xid);
148:     PetscObjectStateGet((PetscObject)x,&bv->xstate);
149:   }
150:   PetscFunctionReturn(0);
151: }

153: /*
154:   BV_IPMatMultBV - Multiply BV by the inner-product matrix, cache the
155:   result internally in bv->cached.
156: */
157: static inline PetscErrorCode BV_IPMatMultBV(BV bv)
158: {
159:   BVGetCachedBV(bv,&bv->cached);
160:   if (((PetscObject)bv)->state != bv->bvstate || bv->l != bv->cached->l || bv->k != bv->cached->k) {
161:     BVSetActiveColumns(bv->cached,bv->l,bv->k);
162:     if (bv->matrix) BVMatMult(bv,bv->matrix,bv->cached);
163:     else BVCopy(bv,bv->cached);
164:     bv->bvstate = ((PetscObject)bv)->state;
165:   }
166:   PetscFunctionReturn(0);
167: }

169: /*
170:   BV_AllocateCoeffs - Allocate orthogonalization coefficients if not done already.
171: */
172: static inline PetscErrorCode BV_AllocateCoeffs(BV bv)
173: {
174:   if (!bv->h) {
175:     PetscMalloc2(bv->nc+bv->m,&bv->h,bv->nc+bv->m,&bv->c);
176:     PetscLogObjectMemory((PetscObject)bv,2*bv->m*sizeof(PetscScalar));
177:   }
178:   PetscFunctionReturn(0);
179: }

181: /*
182:   BV_AllocateSignature - Allocate signature coefficients if not done already.
183: */
184: static inline PetscErrorCode BV_AllocateSignature(BV bv)
185: {
186:   if (bv->indef && !bv->omega) {
187:     if (bv->cuda) {
188: #if defined(PETSC_HAVE_CUDA)
189:       VecCreateSeqCUDA(PETSC_COMM_SELF,bv->nc+bv->m,&bv->omega);
190: #else
191:       SETERRQ(PetscObjectComm((PetscObject)bv),PETSC_ERR_PLIB,"Something wrong happened");
192: #endif
193:     } else VecCreateSeq(PETSC_COMM_SELF,bv->nc+bv->m,&bv->omega);
194:     PetscLogObjectParent((PetscObject)bv,(PetscObject)bv->omega);
195:     VecSet(bv->omega,1.0);
196:   }
197:   PetscFunctionReturn(0);
198: }

200: /*
201:   BVAvailableVec: First (0) or second (1) vector available for
202:   getcolumn operation (or -1 if both vectors already fetched).
203: */
204: #define BVAvailableVec (((bv->ci[0]==-bv->nc-1)? 0: (bv->ci[1]==-bv->nc-1)? 1: -1))

206: /*
207:     Macros to test valid BV arguments
208: */
209: #if !defined(PETSC_USE_DEBUG)

211: #define BVCheckSizes(h,arg) do {(void)(h);} while (0)
212: #define BVCheckOp(h,arg,op) do {(void)(h);} while (0)

214: #else

216: #define BVCheckSizes(h,arg) \
217:   do { \
219:   } while (0)

221: #define BVCheckOp(h,arg,op) \
222:   do { \
224:   } while (0)

226: #endif

228: SLEPC_INTERN PetscErrorCode BVView_Vecs(BV,PetscViewer);

230: SLEPC_INTERN PetscErrorCode BVAllocateWork_Private(BV,PetscInt);

232: SLEPC_INTERN PetscErrorCode BVMult_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
233: SLEPC_INTERN PetscErrorCode BVMultVec_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,const PetscScalar*,PetscScalar,PetscScalar*);
234: SLEPC_INTERN PetscErrorCode BVMultInPlace_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar*,const PetscScalar*,PetscBool);
235: SLEPC_INTERN PetscErrorCode BVMultInPlace_Vecs_Private(BV,PetscInt,PetscInt,PetscInt,Vec*,const PetscScalar*,PetscBool);
236: SLEPC_INTERN PetscErrorCode BVAXPY_BLAS_Private(BV,PetscInt,PetscInt,PetscScalar,const PetscScalar*,PetscScalar,PetscScalar*);
237: SLEPC_INTERN PetscErrorCode BVDot_BLAS_Private(BV,PetscInt,PetscInt,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
238: SLEPC_INTERN PetscErrorCode BVDotVec_BLAS_Private(BV,PetscInt,PetscInt,const PetscScalar*,const PetscScalar*,PetscScalar*,PetscBool);
239: SLEPC_INTERN PetscErrorCode BVScale_BLAS_Private(BV,PetscInt,PetscScalar*,PetscScalar);
240: SLEPC_INTERN PetscErrorCode BVNorm_LAPACK_Private(BV,PetscInt,PetscInt,const PetscScalar*,NormType,PetscReal*,PetscBool);
241: SLEPC_INTERN PetscErrorCode BVNormalize_LAPACK_Private(BV,PetscInt,PetscInt,const PetscScalar*,PetscScalar*,PetscBool);
242: SLEPC_INTERN PetscErrorCode BVMatCholInv_LAPACK_Private(BV,Mat,Mat);
243: SLEPC_INTERN PetscErrorCode BVMatTriInv_LAPACK_Private(BV,Mat,Mat);
244: SLEPC_INTERN PetscErrorCode BVMatSVQB_LAPACK_Private(BV,Mat,Mat);
245: SLEPC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_TSQR(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscInt);
246: SLEPC_INTERN PetscErrorCode BVOrthogonalize_LAPACK_TSQR_OnlyR(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscInt);

248: /* reduction operations used in BVOrthogonalize and BVNormalize */
249: SLEPC_EXTERN MPI_Op MPIU_TSQR, MPIU_LAPY2;
250: SLEPC_EXTERN void MPIAPI SlepcGivensPacked(void*,void*,PetscMPIInt*,MPI_Datatype*);
251: SLEPC_EXTERN void MPIAPI SlepcPythag(void*,void*,PetscMPIInt*,MPI_Datatype*);

253: /*
254:    BV_CleanCoefficients_Default - Sets to zero all entries of column j of the bv buffer
255: */
256: static inline PetscErrorCode BV_CleanCoefficients_Default(BV bv,PetscInt j,PetscScalar *h)
257: {
258:   PetscScalar    *hh=h,*a;
259:   PetscInt       i;

261:   if (!h) {
262:     VecGetArray(bv->buffer,&a);
263:     hh = a + j*(bv->nc+bv->m);
264:   }
265:   for (i=0;i<bv->nc+j;i++) hh[i] = 0.0;
266:   if (!h) VecRestoreArray(bv->buffer,&a);
267:   PetscFunctionReturn(0);
268: }

270: /*
271:    BV_AddCoefficients_Default - Add the contents of the scratch (0-th column) of the bv buffer
272:    into column j of the bv buffer
273: */
274: static inline PetscErrorCode BV_AddCoefficients_Default(BV bv,PetscInt j,PetscScalar *h,PetscScalar *c)
275: {
276:   PetscScalar    *hh=h,*cc=c;
277:   PetscInt       i;

279:   if (!h) {
280:     VecGetArray(bv->buffer,&cc);
281:     hh = cc + j*(bv->nc+bv->m);
282:   }
283:   for (i=0;i<bv->nc+j;i++) hh[i] += cc[i];
284:   if (!h) VecRestoreArray(bv->buffer,&cc);
285:   PetscLogFlops(1.0*(bv->nc+j));
286:   PetscFunctionReturn(0);
287: }

289: /*
290:    BV_SetValue_Default - Sets value in row j (counted after the constraints) of column k
291:    of the coefficients array
292: */
293: static inline PetscErrorCode BV_SetValue_Default(BV bv,PetscInt j,PetscInt k,PetscScalar *h,PetscScalar value)
294: {
295:   PetscScalar    *hh=h,*a;

297:   if (!h) {
298:     VecGetArray(bv->buffer,&a);
299:     hh = a + k*(bv->nc+bv->m);
300:   }
301:   hh[bv->nc+j] = value;
302:   if (!h) VecRestoreArray(bv->buffer,&a);
303:   PetscFunctionReturn(0);
304: }

306: /*
307:    BV_SquareSum_Default - Returns the value h'*h, where h represents the contents of the
308:    coefficients array (up to position j)
309: */
310: static inline PetscErrorCode BV_SquareSum_Default(BV bv,PetscInt j,PetscScalar *h,PetscReal *sum)
311: {
312:   PetscScalar    *hh=h;
313:   PetscInt       i;

315:   *sum = 0.0;
316:   if (!h) VecGetArray(bv->buffer,&hh);
317:   for (i=0;i<bv->nc+j;i++) *sum += PetscRealPart(hh[i]*PetscConj(hh[i]));
318:   if (!h) VecRestoreArray(bv->buffer,&hh);
319:   PetscLogFlops(2.0*(bv->nc+j));
320:   PetscFunctionReturn(0);
321: }

323: /*
324:    BV_ApplySignature_Default - Computes the pointwise product h*omega, where h represents
325:    the contents of the coefficients array (up to position j) and omega is the signature;
326:    if inverse=TRUE then the operation is h/omega
327: */
328: static inline PetscErrorCode BV_ApplySignature_Default(BV bv,PetscInt j,PetscScalar *h,PetscBool inverse)
329: {
330:   PetscScalar       *hh=h;
331:   PetscInt          i;
332:   const PetscScalar *omega;

334:   if (PetscUnlikely(!(bv->nc+j))) PetscFunctionReturn(0);
335:   if (!h) VecGetArray(bv->buffer,&hh);
336:   VecGetArrayRead(bv->omega,&omega);
337:   if (inverse) for (i=0;i<bv->nc+j;i++) hh[i] /= PetscRealPart(omega[i]);
338:   else for (i=0;i<bv->nc+j;i++) hh[i] *= PetscRealPart(omega[i]);
339:   VecRestoreArrayRead(bv->omega,&omega);
340:   if (!h) VecRestoreArray(bv->buffer,&hh);
341:   PetscLogFlops(1.0*(bv->nc+j));
342:   PetscFunctionReturn(0);
343: }

345: /*
346:    BV_SquareRoot_Default - Returns the square root of position j (counted after the constraints)
347:    of the coefficients array
348: */
349: static inline PetscErrorCode BV_SquareRoot_Default(BV bv,PetscInt j,PetscScalar *h,PetscReal *beta)
350: {
351:   PetscScalar    *hh=h;

353:   if (!h) VecGetArray(bv->buffer,&hh);
354:   BV_SafeSqrt(bv,hh[bv->nc+j],beta);
355:   if (!h) VecRestoreArray(bv->buffer,&hh);
356:   PetscFunctionReturn(0);
357: }

359: /*
360:    BV_StoreCoefficients_Default - Copy the contents of the coefficients array to an array dest
361:    provided by the caller (only values from l to j are copied)
362: */
363: static inline PetscErrorCode BV_StoreCoefficients_Default(BV bv,PetscInt j,PetscScalar *h,PetscScalar *dest)
364: {
365:   PetscScalar    *hh=h,*a;
366:   PetscInt       i;

368:   if (!h) {
369:     VecGetArray(bv->buffer,&a);
370:     hh = a + j*(bv->nc+bv->m);
371:   }
372:   for (i=bv->l;i<j;i++) dest[i-bv->l] = hh[bv->nc+i];
373:   if (!h) VecRestoreArray(bv->buffer,&a);
374:   PetscFunctionReturn(0);
375: }

377: /*
378:   BV_GetEigenvector - retrieves k-th eigenvector from basis vectors V.
379:   The argument eigi is the imaginary part of the corresponding eigenvalue.
380: */
381: static inline PetscErrorCode BV_GetEigenvector(BV V,PetscInt k,PetscScalar eigi,Vec Vr,Vec Vi)
382: {
383: #if defined(PETSC_USE_COMPLEX)
384:   if (Vr) BVCopyVec(V,k,Vr);
385:   if (Vi) VecSet(Vi,0.0);
386: #else
387:   if (eigi > 0.0) { /* first value of conjugate pair */
388:     if (Vr) BVCopyVec(V,k,Vr);
389:     if (Vi) BVCopyVec(V,k+1,Vi);
390:   } else if (eigi < 0.0) { /* second value of conjugate pair */
391:     if (Vr) BVCopyVec(V,k-1,Vr);
392:     if (Vi) {
393:       BVCopyVec(V,k,Vi);
394:       VecScale(Vi,-1.0);
395:     }
396:   } else { /* real eigenvalue */
397:     if (Vr) BVCopyVec(V,k,Vr);
398:     if (Vi) VecSet(Vi,0.0);
399:   }
400: #endif
401:   PetscFunctionReturn(0);
402: }

404: /*
405:    BV_OrthogonalizeColumn_Safe - this is intended for cases where we know that
406:    the resulting vector is going to be numerically zero, so normalization or
407:    iterative refinement may cause problems in parallel (collective operation
408:    not being called by all processes)
409: */
410: static inline PetscErrorCode BV_OrthogonalizeColumn_Safe(BV bv,PetscInt j,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
411: {
412:   BVOrthogRefineType orthog_ref;

414:   PetscInfo(bv,"Orthogonalizing column %" PetscInt_FMT " without refinement\n",j);
415:   orthog_ref     = bv->orthog_ref;
416:   bv->orthog_ref = BV_ORTHOG_REFINE_NEVER;  /* avoid refinement */
417:   BVOrthogonalizeColumn(bv,j,H,NULL,NULL);
418:   bv->orthog_ref = orthog_ref;  /* restore refinement setting */
419:   if (norm)   *norm  = 0.0;
420:   if (lindep) *lindep = PETSC_TRUE;
421:   PetscFunctionReturn(0);
422: }

424: #if defined(PETSC_HAVE_CUDA)
425: SLEPC_INTERN PetscErrorCode BV_CleanCoefficients_CUDA(BV,PetscInt,PetscScalar*);
426: SLEPC_INTERN PetscErrorCode BV_AddCoefficients_CUDA(BV,PetscInt,PetscScalar*,PetscScalar*);
427: SLEPC_INTERN PetscErrorCode BV_SetValue_CUDA(BV,PetscInt,PetscInt,PetscScalar*,PetscScalar);
428: SLEPC_INTERN PetscErrorCode BV_SquareSum_CUDA(BV,PetscInt,PetscScalar*,PetscReal*);
429: SLEPC_INTERN PetscErrorCode BV_ApplySignature_CUDA(BV,PetscInt,PetscScalar*,PetscBool);
430: SLEPC_INTERN PetscErrorCode BV_SquareRoot_CUDA(BV,PetscInt,PetscScalar*,PetscReal*);
431: SLEPC_INTERN PetscErrorCode BV_StoreCoefficients_CUDA(BV,PetscInt,PetscScalar*,PetscScalar*);
432: #define BV_CleanCoefficients(a,b,c)   ((a)->cuda?BV_CleanCoefficients_CUDA:BV_CleanCoefficients_Default)((a),(b),(c))
433: #define BV_AddCoefficients(a,b,c,d)   ((a)->cuda?BV_AddCoefficients_CUDA:BV_AddCoefficients_Default)((a),(b),(c),(d))
434: #define BV_SetValue(a,b,c,d,e)        ((a)->cuda?BV_SetValue_CUDA:BV_SetValue_Default)((a),(b),(c),(d),(e))
435: #define BV_SquareSum(a,b,c,d)         ((a)->cuda?BV_SquareSum_CUDA:BV_SquareSum_Default)((a),(b),(c),(d))
436: #define BV_ApplySignature(a,b,c,d)    ((a)->cuda?BV_ApplySignature_CUDA:BV_ApplySignature_Default)((a),(b),(c),(d))
437: #define BV_SquareRoot(a,b,c,d)        ((a)->cuda?BV_SquareRoot_CUDA:BV_SquareRoot_Default)((a),(b),(c),(d))
438: #define BV_StoreCoefficients(a,b,c,d) ((a)->cuda?BV_StoreCoefficients_CUDA:BV_StoreCoefficients_Default)((a),(b),(c),(d))
439: #else
440: #define BV_CleanCoefficients(a,b,c)   BV_CleanCoefficients_Default((a),(b),(c))
441: #define BV_AddCoefficients(a,b,c,d)   BV_AddCoefficients_Default((a),(b),(c),(d))
442: #define BV_SetValue(a,b,c,d,e)        BV_SetValue_Default((a),(b),(c),(d),(e))
443: #define BV_SquareSum(a,b,c,d)         BV_SquareSum_Default((a),(b),(c),(d))
444: #define BV_ApplySignature(a,b,c,d)    BV_ApplySignature_Default((a),(b),(c),(d))
445: #define BV_SquareRoot(a,b,c,d)        BV_SquareRoot_Default((a),(b),(c),(d))
446: #define BV_StoreCoefficients(a,b,c,d) BV_StoreCoefficients_Default((a),(b),(c),(d))
447: #endif /* PETSC_HAVE_CUDA */

449: #endif