Actual source code: svdprimme.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: #if defined(PETSC_HAVE_MPI_IN_PLACE)
 50:     *MPI_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
 51: #else
 52: #error Cannot use SLEPc-PRIMME interface if your MPI does not support MPI_IN_PLACE
 53: #endif
 54:   } else {
 55:     *MPI_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));
 56:   }
 57: }

 59: #if defined(SLEPC_HAVE_PRIMME3)
 60: static void par_broadcastReal(void *buf,int *count,primme_svds_params *primme,int *ierr)
 61: {
 62:   *MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));
 63: }
 64: #endif

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

 78:   (*svd->converged)(svd,sigma,r,&errest,svd->convergedctx);
 79:   if (ierr) *err = 1;
 80:   else {
 81:     *isconv = (errest<=svd->tol?1:0);
 82:     *err = 0;
 83:   }
 84: }

 86: 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,
 87: #if defined(SLEPC_HAVE_PRIMME3)
 88:                        const char *msg,double *time,
 89: #endif
 90:                        primme_event *event,int *stage,struct primme_svds_params *primme,int *err)
 91: {

 93:   PetscErrorCode 0;
 94:   SVD            svd = (SVD)primme->commInfo;
 95:   PetscInt       i,k,nerrest;

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

127: 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)
128: {
129:   PetscInt   i;
130:   SVD_PRIMME *ops = (SVD_PRIMME*)primme->matrix;
131:   Vec        x = ops->x,y = ops->y;
132:   SVD        svd = ops->svd;

135:   for (i=0;i<*blockSize;i++) {
136:     if (*transpose) {
137:       *VecPlaceArray(y,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
138:       *VecPlaceArray(x,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
139:       *MatMult(svd->AT,y,x);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
140:     } else {
141:       *VecPlaceArray(x,(PetscScalar*)xa+(*ldx)*i);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
142:       *VecPlaceArray(y,(PetscScalar*)ya+(*ldy)*i);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
143:       *MatMult(svd->A,x,y);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
144:     }
145:     *VecResetArray(x);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
146:     *VecResetArray(y);CHKERRABORT(PetscObjectComm((PetscObject)svd),*ierr);
147:   }
148:   PetscFunctionReturnVoid();
149: }

151: PetscErrorCode SVDSetUp_PRIMME(SVD svd)
152: {
153:   PetscErrorCode     ierr;
154:   PetscMPIInt        numProcs,procID;
155:   PetscInt           n,m,nloc,mloc;
156:   SVD_PRIMME         *ops = (SVD_PRIMME*)svd->data;
157:   primme_svds_params *primme = &ops->primme;

160:   SVDCheckStandard(svd);
161:   MPI_Comm_size(PetscObjectComm((PetscObject)svd),&numProcs);
162:   MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&procID);

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

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

203:   switch (svd->which) {
204:     case SVD_LARGEST:
205:       primme->target = primme_svds_largest;
206:       break;
207:     case SVD_SMALLEST:
208:       primme->target = primme_svds_smallest;
209:       break;
210:   }

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

218:   if (primme_svds_set_method(ops->method,(primme_preset_method)EPS_PRIMME_DEFAULT_MIN_TIME,PRIMME_DEFAULT_METHOD,primme) < 0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"PRIMME method not valid");

220:   svd->mpd = primme->maxBasisSize;
221:   svd->ncv = (primme->locking?svd->nsv:0)+primme->maxBasisSize;
222:   ops->bs  = primme->maxBlockSize;

224:   /* Set workspace */
225:   SVDAllocateSolution(svd,0);

227:   /* Prepare auxiliary vectors */
228:   if (!ops->x) {
229:     MatCreateVecsEmpty(svd->A,&ops->x,&ops->y);
230:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ops->x);
231:     PetscLogObjectParent((PetscObject)svd,(PetscObject)ops->y);
232:   }
233:   return(0);
234: }

236: PetscErrorCode SVDSolve_PRIMME(SVD svd)
237: {
239:   SVD_PRIMME     *ops = (SVD_PRIMME*)svd->data;
240:   PetscScalar    *svecs, *a;
241:   PetscInt       i,ierrprimme;
242:   PetscReal      *svals,*rnorms;

245:   /* Reset some parameters left from previous runs */
246:   ops->primme.aNorm    = 0.0;
247:   ops->primme.initSize = svd->nini;
248:   ops->primme.iseed[0] = -1;
249:   ops->primme.iseed[1] = -1;
250:   ops->primme.iseed[2] = -1;
251:   ops->primme.iseed[3] = -1;

253:   /* Allocating left and right singular vectors contiguously */
254:   PetscCalloc1(ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal),&svecs);
255:   PetscLogObjectMemory((PetscObject)svd,sizeof(PetscReal)*ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal));

257:   /* Call PRIMME solver */
258:   PetscMalloc2(svd->ncv,&svals,svd->ncv,&rnorms);
259:   ierrprimme = PRIMME_DRIVER(svals,svecs,rnorms,&ops->primme);
260:   for (i=0;i<svd->ncv;i++) svd->sigma[i] = svals[i];
261:   for (i=0;i<svd->ncv;i++) svd->errest[i] = rnorms[i];
262:   PetscFree2(svals,rnorms);

264:   /* Copy left and right singular vectors into svd */
265:   BVGetArray(svd->U,&a);
266:   PetscArraycpy(a,svecs,ops->primme.mLocal*ops->primme.initSize);
267:   BVRestoreArray(svd->U,&a);

269:   BVGetArray(svd->V,&a);
270:   PetscArraycpy(a,svecs+ops->primme.mLocal*ops->primme.initSize,ops->primme.nLocal*ops->primme.initSize);
271:   BVRestoreArray(svd->V,&a);

273:   PetscFree(svecs);

275:   svd->nconv  = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
276:   svd->reason = svd->nconv >= svd->nsv ? SVD_CONVERGED_TOL: SVD_DIVERGED_ITS;
277:   svd->its    = ops->primme.stats.numOuterIterations;

279:   /* Process PRIMME error code */
280:   if (ierrprimme == 0) {
281:     /* no error */
282:   } else if (ierrprimme%100 == -1) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: unexpected error",ierrprimme);
283:   else if (ierrprimme%100 == -2) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: allocation error",ierrprimme);
284:   else if (ierrprimme%100 == -3) {
285:     /* stop by maximum number of iteration or matvecs */
286:   } else if (ierrprimme%100 >= -39) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: configuration error; check PRIMME's manual",ierrprimme);
287:   else SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: runtime error; check PRIMME's manual",ierrprimme);
288:   return(0);
289: }

291: PetscErrorCode SVDReset_PRIMME(SVD svd)
292: {
294:   SVD_PRIMME     *ops = (SVD_PRIMME*)svd->data;

297:   primme_svds_free(&ops->primme);
298:   VecDestroy(&ops->x);
299:   VecDestroy(&ops->y);
300:   return(0);
301: }

303: PetscErrorCode SVDDestroy_PRIMME(SVD svd)
304: {

308:   PetscFree(svd->data);
309:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL);
310:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL);
311:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL);
312:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL);
313:   return(0);
314: }

316: PetscErrorCode SVDView_PRIMME(SVD svd,PetscViewer viewer)
317: {
319:   PetscBool      isascii;
320:   SVD_PRIMME     *ctx = (SVD_PRIMME*)svd->data;
321:   PetscMPIInt    rank;

324:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
325:   if (isascii) {
326:     PetscViewerASCIIPrintf(viewer,"  block size=%D\n",ctx->bs);
327:     PetscViewerASCIIPrintf(viewer,"  solver method: %s\n",SVDPRIMMEMethods[(SVDPRIMMEMethod)ctx->method]);

329:     /* Display PRIMME params */
330:     MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank);
331:     if (!rank) primme_svds_display_params(ctx->primme);
332:   }
333:   return(0);
334: }

336: PetscErrorCode SVDSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,SVD svd)
337: {
338:   PetscErrorCode  ierr;
339:   SVD_PRIMME      *ctx = (SVD_PRIMME*)svd->data;
340:   PetscInt        bs;
341:   SVDPRIMMEMethod meth;
342:   PetscBool       flg;

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

347:     PetscOptionsInt("-svd_primme_blocksize","Maximum block size","SVDPRIMMESetBlockSize",ctx->bs,&bs,&flg);
348:     if (flg) { SVDPRIMMESetBlockSize(svd,bs); }

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

353:   PetscOptionsTail();
354:   return(0);
355: }

357: static PetscErrorCode SVDPRIMMESetBlockSize_PRIMME(SVD svd,PetscInt bs)
358: {
359:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

362:   if (bs == PETSC_DEFAULT) ops->bs = 0;
363:   else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
364:   else ops->bs = bs;
365:   return(0);
366: }

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

371:    Logically Collective on svd

373:    Input Parameters:
374: +  svd - the singular value solver context
375: -  bs - block size

377:    Options Database Key:
378: .  -svd_primme_blocksize - Sets the max allowed block size value

380:    Notes:
381:    If the block size is not set, the value established by primme_svds_initialize
382:    is used.

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

389:    Level: advanced

391: .seealso: SVDPRIMMEGetBlockSize()
392: @*/
393: PetscErrorCode SVDPRIMMESetBlockSize(SVD svd,PetscInt bs)
394: {

400:   PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs));
401:   return(0);
402: }

404: static PetscErrorCode SVDPRIMMEGetBlockSize_PRIMME(SVD svd,PetscInt *bs)
405: {
406:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

409:   *bs = ops->bs;
410:   return(0);
411: }

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

416:    Not Collective

418:    Input Parameter:
419: .  svd - the singular value solver context

421:    Output Parameter:
422: .  bs - returned block size

424:    Level: advanced

426: .seealso: SVDPRIMMESetBlockSize()
427: @*/
428: PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs)
429: {

435:   PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs));
436:   return(0);
437: }

439: static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method)
440: {
441:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

444:   ops->method = (primme_svds_preset_method)method;
445:   return(0);
446: }

448: /*@
449:    SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver.

451:    Logically Collective on svd

453:    Input Parameters:
454: +  svd - the singular value solver context
455: -  method - method that will be used by PRIMME

457:    Options Database Key:
458: .  -svd_primme_method - Sets the method for the PRIMME SVD solver

460:    Note:
461:    If not set, the method defaults to SVD_PRIMME_HYBRID.

463:    Level: advanced

465: .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod
466: @*/
467: PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method)
468: {

474:   PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method));
475:   return(0);
476: }

478: static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method)
479: {
480:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

483:   *method = (SVDPRIMMEMethod)ops->method;
484:   return(0);
485: }

487: /*@
488:    SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver.

490:    Not Collective

492:    Input Parameter:
493: .  svd - the singular value solver context

495:    Output Parameter:
496: .  method - method that will be used by PRIMME

498:    Level: advanced

500: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod
501: @*/
502: PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method)
503: {

509:   PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method));
510:   return(0);
511: }

513: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd)
514: {
516:   SVD_PRIMME     *primme;

519:   PetscNewLog(svd,&primme);
520:   svd->data = (void*)primme;

522:   primme_svds_initialize(&primme->primme);
523:   primme->bs = 0;
524:   primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID;
525:   primme->svd = svd;

527:   svd->ops->solve          = SVDSolve_PRIMME;
528:   svd->ops->setup          = SVDSetUp_PRIMME;
529:   svd->ops->setfromoptions = SVDSetFromOptions_PRIMME;
530:   svd->ops->destroy        = SVDDestroy_PRIMME;
531:   svd->ops->reset          = SVDReset_PRIMME;
532:   svd->ops->view           = SVDView_PRIMME;

534:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME);
535:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME);
536:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME);
537:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME);
538:   return(0);
539: }