Actual source code: svdprimme.c

slepc-3.17.0 2022-03-31
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: {
 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 = 0;
 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:   for (i=0;i<*blockSize;i++) {
131:     if (*transpose) {
132:       PetscObjectComm((PetscObject)svd),VecPlaceArray(y,(PetscScalar*)xa+(*ldx)*i);
133:       PetscObjectComm((PetscObject)svd),VecPlaceArray(x,(PetscScalar*)ya+(*ldy)*i);
134:       PetscObjectComm((PetscObject)svd),MatMult(svd->AT,y,x);
135:     } else {
136:       PetscObjectComm((PetscObject)svd),VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);
137:       PetscObjectComm((PetscObject)svd),VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);
138:       PetscObjectComm((PetscObject)svd),MatMult(svd->A,x,y);
139:     }
140:     PetscObjectComm((PetscObject)svd),VecResetArray(x);
141:     PetscObjectComm((PetscObject)svd),VecResetArray(y);
142:   }
143:   return;
144: }

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

153:   SVDCheckStandard(svd);
154:   MPI_Comm_size(PetscObjectComm((PetscObject)svd),&numProcs);
155:   MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&procID);

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

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

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

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


213:   svd->mpd = primme->maxBasisSize;
214:   svd->ncv = (primme->locking?svd->nsv:0)+primme->maxBasisSize;
215:   ops->bs  = primme->maxBlockSize;

217:   /* Set workspace */
218:   SVDAllocateSolution(svd,0);

220:   /* Prepare auxiliary vectors */
221:   if (!ops->x) {
222:     MatCreateVecsEmpty(svd->A,&ops->x,&ops->y);
223:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ops->x);
224:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ops->y);
225:   }
226:   PetscFunctionReturn(0);
227: }

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

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:   PetscCalloc1(ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal),&svecs);
246:   PetscLogObjectMemory((PetscObject)svd,sizeof(PetscReal)*ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal));

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

255:   /* Copy left and right singular vectors into svd */
256:   BVGetArray(svd->U,&a);
257:   PetscArraycpy(a,svecs,ops->primme.mLocal*ops->primme.initSize);
258:   BVRestoreArray(svd->U,&a);

260:   BVGetArray(svd->V,&a);
261:   PetscArraycpy(a,svecs+ops->primme.mLocal*ops->primme.initSize,ops->primme.nLocal*ops->primme.initSize);
262:   BVRestoreArray(svd->V,&a);

264:   PetscFree(svecs);

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

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

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

291:   primme_svds_free(&ops->primme);
292:   VecDestroy(&ops->x);
293:   VecDestroy(&ops->y);
294:   PetscFunctionReturn(0);
295: }

297: PetscErrorCode SVDDestroy_PRIMME(SVD svd)
298: {
299:   PetscFree(svd->data);
300:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL);
301:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL);
302:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL);
303:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL);
304:   PetscFunctionReturn(0);
305: }

307: PetscErrorCode SVDView_PRIMME(SVD svd,PetscViewer viewer)
308: {
309:   PetscBool      isascii;
310:   SVD_PRIMME     *ctx = (SVD_PRIMME*)svd->data;
311:   PetscMPIInt    rank;

313:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
314:   if (isascii) {
315:     PetscViewerASCIIPrintf(viewer,"  block size=%" PetscInt_FMT "\n",ctx->bs);
316:     PetscViewerASCIIPrintf(viewer,"  solver method: %s\n",SVDPRIMMEMethods[(SVDPRIMMEMethod)ctx->method]);

318:     /* Display PRIMME params */
319:     MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank);
320:     if (!rank) primme_svds_display_params(ctx->primme);
321:   }
322:   PetscFunctionReturn(0);
323: }

325: PetscErrorCode SVDSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,SVD svd)
326: {
327:   SVD_PRIMME      *ctx = (SVD_PRIMME*)svd->data;
328:   PetscInt        bs;
329:   SVDPRIMMEMethod meth;
330:   PetscBool       flg;

332:   PetscOptionsHead(PetscOptionsObject,"SVD PRIMME Options");

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

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

340:   PetscOptionsTail();
341:   PetscFunctionReturn(0);
342: }

344: static PetscErrorCode SVDPRIMMESetBlockSize_PRIMME(SVD svd,PetscInt bs)
345: {
346:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

348:   if (bs == PETSC_DEFAULT) ops->bs = 0;
349:   else {
351:     ops->bs = bs;
352:   }
353:   PetscFunctionReturn(0);
354: }

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

359:    Logically Collective on svd

361:    Input Parameters:
362: +  svd - the singular value solver context
363: -  bs - block size

365:    Options Database Key:
366: .  -svd_primme_blocksize - Sets the max allowed block size value

368:    Notes:
369:    If the block size is not set, the value established by primme_svds_initialize
370:    is used.

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

377:    Level: advanced

379: .seealso: SVDPRIMMEGetBlockSize()
380: @*/
381: PetscErrorCode SVDPRIMMESetBlockSize(SVD svd,PetscInt bs)
382: {
385:   PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs));
386:   PetscFunctionReturn(0);
387: }

389: static PetscErrorCode SVDPRIMMEGetBlockSize_PRIMME(SVD svd,PetscInt *bs)
390: {
391:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

393:   *bs = ops->bs;
394:   PetscFunctionReturn(0);
395: }

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

400:    Not Collective

402:    Input Parameter:
403: .  svd - the singular value solver context

405:    Output Parameter:
406: .  bs - returned block size

408:    Level: advanced

410: .seealso: SVDPRIMMESetBlockSize()
411: @*/
412: PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs)
413: {
416:   PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs));
417:   PetscFunctionReturn(0);
418: }

420: static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method)
421: {
422:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

424:   ops->method = (primme_svds_preset_method)method;
425:   PetscFunctionReturn(0);
426: }

428: /*@
429:    SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver.

431:    Logically Collective on svd

433:    Input Parameters:
434: +  svd - the singular value solver context
435: -  method - method that will be used by PRIMME

437:    Options Database Key:
438: .  -svd_primme_method - Sets the method for the PRIMME SVD solver

440:    Note:
441:    If not set, the method defaults to SVD_PRIMME_HYBRID.

443:    Level: advanced

445: .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod
446: @*/
447: PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method)
448: {
451:   PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method));
452:   PetscFunctionReturn(0);
453: }

455: static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method)
456: {
457:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

459:   *method = (SVDPRIMMEMethod)ops->method;
460:   PetscFunctionReturn(0);
461: }

463: /*@
464:    SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver.

466:    Not Collective

468:    Input Parameter:
469: .  svd - the singular value solver context

471:    Output Parameter:
472: .  method - method that will be used by PRIMME

474:    Level: advanced

476: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod
477: @*/
478: PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method)
479: {
482:   PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method));
483:   PetscFunctionReturn(0);
484: }

486: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd)
487: {
488:   SVD_PRIMME     *primme;

490:   PetscNewLog(svd,&primme);
491:   svd->data = (void*)primme;

493:   primme_svds_initialize(&primme->primme);
494:   primme->bs = 0;
495:   primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID;
496:   primme->svd = svd;

498:   svd->ops->solve          = SVDSolve_PRIMME;
499:   svd->ops->setup          = SVDSetUp_PRIMME;
500:   svd->ops->setfromoptions = SVDSetFromOptions_PRIMME;
501:   svd->ops->destroy        = SVDDestroy_PRIMME;
502:   svd->ops->reset          = SVDReset_PRIMME;
503:   svd->ops->view           = SVDView_PRIMME;

505:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME);
506:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME);
507:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME);
508:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME);
509:   PetscFunctionReturn(0);
510: }