Actual source code: svdprimme.c

slepc-3.21.0 2024-03-30
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: */
 10: /*
 11:    This file implements a wrapper to the PRIMME SVD solver
 12: */

 14: #include <slepc/private/svdimpl.h>

 16: #include <primme.h>

 18: #if defined(PETSC_USE_COMPLEX)
 19: #if defined(PETSC_USE_REAL_SINGLE)
 20: #define PRIMME_DRIVER cprimme_svds
 21: #else
 22: #define PRIMME_DRIVER zprimme_svds
 23: #endif
 24: #else
 25: #if defined(PETSC_USE_REAL_SINGLE)
 26: #define PRIMME_DRIVER sprimme_svds
 27: #else
 28: #define PRIMME_DRIVER dprimme_svds
 29: #endif
 30: #endif

 32: #if defined(PRIMME_VERSION_MAJOR) && PRIMME_VERSION_MAJOR*100+PRIMME_VERSION_MINOR >= 202
 33: #define SLEPC_HAVE_PRIMME2p2
 34: #endif

 36: typedef struct {
 37:   primme_svds_params        primme;   /* param struct */
 38:   PetscInt                  bs;       /* block size */
 39:   primme_svds_preset_method method;   /* primme method */
 40:   SVD                       svd;      /* reference to the solver */
 41:   Vec                       x,y;      /* auxiliary vectors */
 42: } SVD_PRIMME;

 44: static void multMatvec_PRIMME(void*,PRIMME_INT*,void*,PRIMME_INT*,int*,int*,struct primme_svds_params*,int*);

 46: static void par_GlobalSumReal(void *sendBuf,void *recvBuf,int *count,primme_svds_params *primme,int *ierr)
 47: {
 48:   if (sendBuf == recvBuf) {
 49:     *ierr = MPI_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
 50:   } else {
 51:     *ierr = MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
 52:   }
 53: }

 55: #if defined(SLEPC_HAVE_PRIMME3)
 56: static void par_broadcastReal(void *buf,int *count,primme_svds_params *primme,int *ierr)
 57: {
 58:   *ierr = MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));
 59: }
 60: #endif

 62: #if defined(SLEPC_HAVE_PRIMME2p2)
 63: static void convTestFun(double *sval,void *leftsvec,void *rightsvec,double *resNorm,
 64: #if defined(SLEPC_HAVE_PRIMME3)
 65:                         int *method,
 66: #endif
 67:                         int *isconv,struct primme_svds_params *primme,int *err)
 68: {
 69:   PetscErrorCode ierr;
 70:   SVD            svd = (SVD)primme->commInfo;
 71:   PetscReal      sigma = sval?*sval:0.0;
 72:   PetscReal      r = resNorm?*resNorm:HUGE_VAL,errest;

 74:   ierr = (*svd->converged)(svd,sigma,r,&errest,svd->convergedctx);
 75:   if (ierr) *err = 1;
 76:   else {
 77:     *isconv = (errest<=svd->tol?1:0);
 78:     *err = 0;
 79:   }
 80: }

 82: static void monitorFun(void *basisSvals,int *basisSize,int *basisFlags,int *iblock,int *blockSize,void *basisNorms,int *numConverged,void *lockedSvals,int *numLocked,int *lockedFlags,void *lockedNorms,int *inner_its,void *LSRes,
 83: #if defined(SLEPC_HAVE_PRIMME3)
 84:                        const char *msg,double *time,
 85: #endif
 86:                        primme_event *event,int *stage,struct primme_svds_params *primme,int *err)
 87: {

 89:   PetscErrorCode ierr = PETSC_SUCCESS;
 90:   SVD            svd = (SVD)primme->commInfo;
 91:   PetscInt       i,k,nerrest;

 93:   *err = 1;
 94:   switch (*event) {
 95:     case primme_event_outer_iteration:
 96:       /* Update SVD */
 97:       svd->its = primme->stats.numOuterIterations;
 98:       if (numConverged) svd->nconv = *numConverged;
 99:       k = 0;
100:       if (lockedSvals && numLocked) for (i=0; i<*numLocked && k<svd->ncv; i++) svd->sigma[k++] = ((PetscReal*)lockedSvals)[i];
101:       nerrest = k;
102:       if (iblock && blockSize) {
103:         for (i=0; i<*blockSize && k+iblock[i]<svd->ncv; i++) svd->errest[k+iblock[i]] = ((PetscReal*)basisNorms)[i];
104:         nerrest = k+(*blockSize>0?1+iblock[*blockSize-1]:0);
105:       }
106:       if (basisSvals && basisSize) for (i=0; i<*basisSize && k<svd->ncv; i++) svd->sigma[k++] = ((PetscReal*)basisSvals)[i];
107:       /* Show progress */
108:       ierr = SVDMonitor(svd,svd->its,numConverged?*numConverged:0,svd->sigma,svd->errest,nerrest);
109:       break;
110: #if defined(SLEPC_HAVE_PRIMME3)
111:     case primme_event_message:
112:       /* Print PRIMME information messages */
113:       ierr = PetscInfo(svd,"%s\n",msg);
114:       break;
115: #endif
116:     default:
117:       break;
118:   }
119:   *err = (ierr!=0)? 1: 0;
120: }
121: #endif /* SLEPC_HAVE_PRIMME2p2 */

123: static void multMatvec_PRIMME(void *xa,PRIMME_INT *ldx,void *ya,PRIMME_INT *ldy,int *blockSize,int *transpose,struct primme_svds_params *primme,int *ierr)
124: {
125:   PetscInt   i;
126:   SVD_PRIMME *ops = (SVD_PRIMME*)primme->matrix;
127:   Vec        x = ops->x,y = ops->y;
128:   SVD        svd = ops->svd;

130:   PetscFunctionBegin;
131:   for (i=0;i<*blockSize;i++) {
132:     if (*transpose) {
133:       PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(y,(PetscScalar*)xa+(*ldx)*i));
134:       PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(x,(PetscScalar*)ya+(*ldy)*i));
135:       PetscCallAbort(PetscObjectComm((PetscObject)svd),MatMult(svd->AT,y,x));
136:     } else {
137:       PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i));
138:       PetscCallAbort(PetscObjectComm((PetscObject)svd),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i));
139:       PetscCallAbort(PetscObjectComm((PetscObject)svd),MatMult(svd->A,x,y));
140:     }
141:     PetscCallAbort(PetscObjectComm((PetscObject)svd),VecResetArray(x));
142:     PetscCallAbort(PetscObjectComm((PetscObject)svd),VecResetArray(y));
143:   }
144:   PetscFunctionReturnVoid();
145: }

147: static PetscErrorCode SVDSetUp_PRIMME(SVD svd)
148: {
149:   PetscMPIInt        numProcs,procID;
150:   PetscInt           n,m,nloc,mloc;
151:   SVD_PRIMME         *ops = (SVD_PRIMME*)svd->data;
152:   primme_svds_params *primme = &ops->primme;

154:   PetscFunctionBegin;
155:   SVDCheckStandard(svd);
156:   SVDCheckDefinite(svd);
157:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)svd),&numProcs));
158:   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&procID));

160:   /* Check some constraints and set some default values */
161:   PetscCall(MatGetSize(svd->A,&m,&n));
162:   PetscCall(MatGetLocalSize(svd->A,&mloc,&nloc));
163:   PetscCall(SVDSetDimensions_Default(svd));
164:   if (svd->max_it==PETSC_DEFAULT) svd->max_it = PETSC_MAX_INT;
165:   svd->leftbasis = PETSC_TRUE;
166:   SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
167: #if !defined(SLEPC_HAVE_PRIMME2p2)
168:   if (svd->converged != SVDConvergedAbsolute) PetscCall(PetscInfo(svd,"Warning: using absolute convergence test\n"));
169: #endif

171:   /* Transfer SLEPc options to PRIMME options */
172:   primme_svds_free(primme);
173:   primme_svds_initialize(primme);
174:   primme->m             = m;
175:   primme->n             = n;
176:   primme->mLocal        = mloc;
177:   primme->nLocal        = nloc;
178:   primme->numSvals      = svd->nsv;
179:   primme->matrix        = ops;
180:   primme->commInfo      = svd;
181:   primme->maxMatvecs    = svd->max_it;
182: #if !defined(SLEPC_HAVE_PRIMME2p2)
183:   primme->eps           = SlepcDefaultTol(svd->tol);
184: #endif
185:   primme->numProcs      = numProcs;
186:   primme->procID        = procID;
187:   primme->printLevel    = 1;
188:   primme->matrixMatvec  = multMatvec_PRIMME;
189:   primme->globalSumReal = par_GlobalSumReal;
190: #if defined(SLEPC_HAVE_PRIMME3)
191:   primme->broadcastReal = par_broadcastReal;
192: #endif
193: #if defined(SLEPC_HAVE_PRIMME2p2)
194:   primme->convTestFun   = convTestFun;
195:   primme->monitorFun    = monitorFun;
196: #endif
197:   if (ops->bs > 0) primme->maxBlockSize = ops->bs;

199:   switch (svd->which) {
200:     case SVD_LARGEST:
201:       primme->target = primme_svds_largest;
202:       break;
203:     case SVD_SMALLEST:
204:       primme->target = primme_svds_smallest;
205:       break;
206:   }

208:   /* If user sets mpd or ncv, maxBasisSize is modified */
209:   if (svd->mpd!=PETSC_DEFAULT) {
210:     primme->maxBasisSize = svd->mpd;
211:     if (svd->ncv!=PETSC_DEFAULT) PetscCall(PetscInfo(svd,"Warning: 'ncv' is ignored by PRIMME\n"));
212:   } else if (svd->ncv!=PETSC_DEFAULT) primme->maxBasisSize = svd->ncv;

214:   PetscCheck(primme_svds_set_method(ops->method,(primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME,PRIMME_DEFAULT_METHOD,primme)>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"PRIMME method not valid");

216:   svd->mpd = primme->maxBasisSize;
217:   svd->ncv = (primme->locking?svd->nsv:0)+primme->maxBasisSize;
218:   ops->bs  = primme->maxBlockSize;

220:   /* Set workspace */
221:   PetscCall(SVDAllocateSolution(svd,0));

223:   /* Prepare auxiliary vectors */
224:   if (!ops->x) PetscCall(MatCreateVecsEmpty(svd->A,&ops->x,&ops->y));
225:   PetscFunctionReturn(PETSC_SUCCESS);
226: }

228: static PetscErrorCode SVDSolve_PRIMME(SVD svd)
229: {
230:   SVD_PRIMME     *ops = (SVD_PRIMME*)svd->data;
231:   PetscScalar    *svecs, *a;
232:   PetscInt       i,ierrprimme,ld;
233:   PetscReal      *svals,*rnorms;

235:   PetscFunctionBegin;
236:   /* Reset some parameters left from previous runs */
237:   ops->primme.aNorm    = 0.0;
238:   ops->primme.initSize = svd->nini;
239:   ops->primme.iseed[0] = -1;
240:   ops->primme.iseed[1] = -1;
241:   ops->primme.iseed[2] = -1;
242:   ops->primme.iseed[3] = -1;

244:   /* Allocating left and right singular vectors contiguously */
245:   PetscCall(PetscCalloc1(ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal),&svecs));

247:   /* Call PRIMME solver */
248:   PetscCall(PetscMalloc2(svd->ncv,&svals,svd->ncv,&rnorms));
249:   ierrprimme = PRIMME_DRIVER(svals,svecs,rnorms,&ops->primme);
250:   for (i=0;i<svd->ncv;i++) svd->sigma[i] = svals[i];
251:   for (i=0;i<svd->ncv;i++) svd->errest[i] = rnorms[i];
252:   PetscCall(PetscFree2(svals,rnorms));

254:   /* Copy left and right singular vectors into svd */
255:   PetscCall(BVGetLeadingDimension(svd->U,&ld));
256:   PetscCall(BVGetArray(svd->U,&a));
257:   for (i=0;i<ops->primme.initSize;i++) PetscCall(PetscArraycpy(a+i*ld,svecs+i*ops->primme.mLocal,ops->primme.mLocal));
258:   PetscCall(BVRestoreArray(svd->U,&a));

260:   PetscCall(BVGetLeadingDimension(svd->V,&ld));
261:   PetscCall(BVGetArray(svd->V,&a));
262:   for (i=0;i<ops->primme.initSize;i++) PetscCall(PetscArraycpy(a+i*ld,svecs+ops->primme.mLocal*ops->primme.initSize+i*ops->primme.nLocal,ops->primme.nLocal));
263:   PetscCall(BVRestoreArray(svd->V,&a));

265:   PetscCall(PetscFree(svecs));

267:   svd->nconv  = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
268:   svd->reason = svd->nconv >= svd->nsv ? SVD_CONVERGED_TOL: SVD_DIVERGED_ITS;
269:   svd->its    = ops->primme.stats.numOuterIterations;

271:   /* Process PRIMME error code */
272:   if (ierrprimme != 0) {
273:     switch (ierrprimme%100) {
274:       case -1:
275:         SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
276:       case -2:
277:         SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
278:       case -3: /* stop due to maximum number of iterations or matvecs */
279:         break;
280:       default:
281:         PetscCheck(ierrprimme<-39,PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": configuration error; check PRIMME's manual",ierrprimme);
282:         PetscCheck(ierrprimme>=-39,PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": runtime error; check PRIMME's manual",ierrprimme);
283:     }
284:   }
285:   PetscFunctionReturn(PETSC_SUCCESS);
286: }

288: static PetscErrorCode SVDReset_PRIMME(SVD svd)
289: {
290:   SVD_PRIMME     *ops = (SVD_PRIMME*)svd->data;

292:   PetscFunctionBegin;
293:   primme_svds_free(&ops->primme);
294:   PetscCall(VecDestroy(&ops->x));
295:   PetscCall(VecDestroy(&ops->y));
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: static PetscErrorCode SVDDestroy_PRIMME(SVD svd)
300: {
301:   PetscFunctionBegin;
302:   PetscCall(PetscFree(svd->data));
303:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL));
304:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL));
305:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL));
306:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL));
307:   PetscFunctionReturn(PETSC_SUCCESS);
308: }

310: static PetscErrorCode SVDView_PRIMME(SVD svd,PetscViewer viewer)
311: {
312:   PetscBool      isascii;
313:   SVD_PRIMME     *ctx = (SVD_PRIMME*)svd->data;
314:   PetscMPIInt    rank;

316:   PetscFunctionBegin;
317:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
318:   if (isascii) {
319:     PetscCall(PetscViewerASCIIPrintf(viewer,"  block size=%" PetscInt_FMT "\n",ctx->bs));
320:     PetscCall(PetscViewerASCIIPrintf(viewer,"  solver method: %s\n",SVDPRIMMEMethods[(SVDPRIMMEMethod)ctx->method]));

322:     /* Display PRIMME params */
323:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank));
324:     if (!rank) primme_svds_display_params(ctx->primme);
325:   }
326:   PetscFunctionReturn(PETSC_SUCCESS);
327: }

329: static PetscErrorCode SVDSetFromOptions_PRIMME(SVD svd,PetscOptionItems *PetscOptionsObject)
330: {
331:   SVD_PRIMME      *ctx = (SVD_PRIMME*)svd->data;
332:   PetscInt        bs;
333:   SVDPRIMMEMethod meth;
334:   PetscBool       flg;

336:   PetscFunctionBegin;
337:   PetscOptionsHeadBegin(PetscOptionsObject,"SVD PRIMME Options");

339:     PetscCall(PetscOptionsInt("-svd_primme_blocksize","Maximum block size","SVDPRIMMESetBlockSize",ctx->bs,&bs,&flg));
340:     if (flg) PetscCall(SVDPRIMMESetBlockSize(svd,bs));

342:     PetscCall(PetscOptionsEnum("-svd_primme_method","Method for solving the singular value problem","SVDPRIMMESetMethod",SVDPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg));
343:     if (flg) PetscCall(SVDPRIMMESetMethod(svd,meth));

345:   PetscOptionsHeadEnd();
346:   PetscFunctionReturn(PETSC_SUCCESS);
347: }

349: static PetscErrorCode SVDPRIMMESetBlockSize_PRIMME(SVD svd,PetscInt bs)
350: {
351:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

353:   PetscFunctionBegin;
354:   if (bs == PETSC_DEFAULT) ops->bs = 0;
355:   else {
356:     PetscCheck(bs>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
357:     ops->bs = bs;
358:   }
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: /*@
363:    SVDPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.

365:    Logically Collective

367:    Input Parameters:
368: +  svd - the singular value solver context
369: -  bs - block size

371:    Options Database Key:
372: .  -svd_primme_blocksize - Sets the max allowed block size value

374:    Notes:
375:    If the block size is not set, the value established by primme_svds_initialize
376:    is used.

378:    The user should set the block size based on the architecture specifics
379:    of the target computer, as well as any a priori knowledge of multiplicities.
380:    The code does NOT require bs > 1 to find multiple eigenvalues. For some
381:    methods, keeping bs = 1 yields the best overall performance.

383:    Level: advanced

385: .seealso: SVDPRIMMEGetBlockSize()
386: @*/
387: PetscErrorCode SVDPRIMMESetBlockSize(SVD svd,PetscInt bs)
388: {
389:   PetscFunctionBegin;
392:   PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs));
393:   PetscFunctionReturn(PETSC_SUCCESS);
394: }

396: static PetscErrorCode SVDPRIMMEGetBlockSize_PRIMME(SVD svd,PetscInt *bs)
397: {
398:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

400:   PetscFunctionBegin;
401:   *bs = ops->bs;
402:   PetscFunctionReturn(PETSC_SUCCESS);
403: }

405: /*@
406:    SVDPRIMMEGetBlockSize - Get the maximum block size the code will try to use.

408:    Not Collective

410:    Input Parameter:
411: .  svd - the singular value solver context

413:    Output Parameter:
414: .  bs - returned block size

416:    Level: advanced

418: .seealso: SVDPRIMMESetBlockSize()
419: @*/
420: PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs)
421: {
422:   PetscFunctionBegin;
424:   PetscAssertPointer(bs,2);
425:   PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs));
426:   PetscFunctionReturn(PETSC_SUCCESS);
427: }

429: static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method)
430: {
431:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

433:   PetscFunctionBegin;
434:   ops->method = (primme_svds_preset_method)method;
435:   PetscFunctionReturn(PETSC_SUCCESS);
436: }

438: /*@
439:    SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver.

441:    Logically Collective

443:    Input Parameters:
444: +  svd - the singular value solver context
445: -  method - method that will be used by PRIMME

447:    Options Database Key:
448: .  -svd_primme_method - Sets the method for the PRIMME SVD solver

450:    Note:
451:    If not set, the method defaults to SVD_PRIMME_HYBRID.

453:    Level: advanced

455: .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod
456: @*/
457: PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method)
458: {
459:   PetscFunctionBegin;
462:   PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method));
463:   PetscFunctionReturn(PETSC_SUCCESS);
464: }

466: static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method)
467: {
468:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

470:   PetscFunctionBegin;
471:   *method = (SVDPRIMMEMethod)ops->method;
472:   PetscFunctionReturn(PETSC_SUCCESS);
473: }

475: /*@
476:    SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver.

478:    Not Collective

480:    Input Parameter:
481: .  svd - the singular value solver context

483:    Output Parameter:
484: .  method - method that will be used by PRIMME

486:    Level: advanced

488: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod
489: @*/
490: PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method)
491: {
492:   PetscFunctionBegin;
494:   PetscAssertPointer(method,2);
495:   PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method));
496:   PetscFunctionReturn(PETSC_SUCCESS);
497: }

499: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd)
500: {
501:   SVD_PRIMME     *primme;

503:   PetscFunctionBegin;
504:   PetscCall(PetscNew(&primme));
505:   svd->data = (void*)primme;

507:   primme_svds_initialize(&primme->primme);
508:   primme->bs = 0;
509:   primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID;
510:   primme->svd = svd;

512:   svd->ops->solve          = SVDSolve_PRIMME;
513:   svd->ops->setup          = SVDSetUp_PRIMME;
514:   svd->ops->setfromoptions = SVDSetFromOptions_PRIMME;
515:   svd->ops->destroy        = SVDDestroy_PRIMME;
516:   svd->ops->reset          = SVDReset_PRIMME;
517:   svd->ops->view           = SVDView_PRIMME;

519:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME));
520:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME));
521:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME));
522:   PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME));
523:   PetscFunctionReturn(PETSC_SUCCESS);
524: }