Actual source code: svdprimme.c

slepc-3.15.0 2021-03-31
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:     *MPIU_Allreduce(MPI_IN_PLACE,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
 50:   } else {
 51:     *MPIU_Allreduce(sendBuf,recvBuf,*count,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
 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:   *MPI_Bcast(buf,*count,MPIU_REAL,0/*root*/,PetscObjectComm((PetscObject)primme->commInfo));CHKERRABORT(PetscObjectComm((PetscObject)primme->commInfo),*ierr);
 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:   *err = 1;
 75:   (*svd->converged)(svd,sigma,r,&errest,svd->convergedctx);CHKERRABORT(PetscObjectComm((PetscObject)svd),ierr);
 76:   *isconv = (errest<=svd->tol?1:0);
 77:   *err = 0;
 78: }

 80: 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,
 81: #if defined(SLEPC_HAVE_PRIMME3)
 82:                        const char *msg,double *time,
 83: #endif
 84:                        primme_event *event,int *stage,struct primme_svds_params *primme,int *err)
 85: {

 88:   SVD            svd = (SVD)primme->commInfo;
 89:   PetscInt       i,k,nerrest;

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

121: 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)
122: {
123:   PetscInt   i;
124:   SVD_PRIMME *ops = (SVD_PRIMME*)primme->matrix;
125:   Vec        x = ops->x,y = ops->y;
126:   SVD        svd = ops->svd;

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

145: PetscErrorCode SVDSetUp_PRIMME(SVD svd)
146: {
147:   PetscErrorCode     ierr;
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;

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

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

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

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

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

212:   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");

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

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

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

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

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

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

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

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

263:   BVGetArray(svd->V,&a);
264:   PetscArraycpy(a,svecs+ops->primme.mLocal*ops->primme.initSize,ops->primme.nLocal*ops->primme.initSize);
265:   BVRestoreArray(svd->V,&a);

267:   PetscFree(svecs);

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

273:   /* Process PRIMME error code */
274:   if (ierrprimme == 0) {
275:     /* no error */
276:   } else if (ierrprimme%100 == -1) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: unexpected error",ierrprimme);
277:   else if (ierrprimme%100 == -2) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: allocation error",ierrprimme);
278:   else if (ierrprimme%100 == -3) {
279:     /* stop by maximum number of iteration or matvecs */
280:   } 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);
281:   else SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: runtime error; check PRIMME's manual",ierrprimme);
282:   return(0);
283: }

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

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

297: PetscErrorCode SVDDestroy_PRIMME(SVD svd)
298: {

302:   PetscFree(svd->data);
303:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL);
304:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL);
305:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL);
306:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL);
307:   return(0);
308: }

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

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

323:     /* Display PRIMME params */
324:     MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank);CHKERRMPI(ierr);
325:     if (!rank) primme_svds_display_params(ctx->primme);
326:   }
327:   return(0);
328: }

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

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

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

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

347:   PetscOptionsTail();
348:   return(0);
349: }

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

356:   if (bs == PETSC_DEFAULT) ops->bs = 0;
357:   else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
358:   else ops->bs = bs;
359:   return(0);
360: }

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

365:    Logically Collective on svd

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: {

394:   PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs));
395:   return(0);
396: }

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

403:   *bs = ops->bs;
404:   return(0);
405: }

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

410:    Not Collective

412:    Input Parameter:
413: .  svd - the singular value solver context

415:    Output Parameter:
416: .  bs - returned block size

418:    Level: advanced

420: .seealso: SVDPRIMMESetBlockSize()
421: @*/
422: PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs)
423: {

429:   PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs));
430:   return(0);
431: }

433: static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method)
434: {
435:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

438:   ops->method = (primme_svds_preset_method)method;
439:   return(0);
440: }

442: /*@
443:    SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver.

445:    Logically Collective on svd

447:    Input Parameters:
448: +  svd - the singular value solver context
449: -  method - method that will be used by PRIMME

451:    Options Database Key:
452: .  -svd_primme_method - Sets the method for the PRIMME SVD solver

454:    Note:
455:    If not set, the method defaults to SVD_PRIMME_HYBRID.

457:    Level: advanced

459: .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod
460: @*/
461: PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method)
462: {

468:   PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method));
469:   return(0);
470: }

472: static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method)
473: {
474:   SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;

477:   *method = (SVDPRIMMEMethod)ops->method;
478:   return(0);
479: }

481: /*@
482:    SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver.

484:    Not Collective

486:    Input Parameter:
487: .  svd - the singular value solver context

489:    Output Parameter:
490: .  method - method that will be used by PRIMME

492:    Level: advanced

494: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod
495: @*/
496: PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method)
497: {

503:   PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method));
504:   return(0);
505: }

507: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd)
508: {
510:   SVD_PRIMME     *primme;

513:   PetscNewLog(svd,&primme);
514:   svd->data = (void*)primme;

516:   primme_svds_initialize(&primme->primme);
517:   primme->bs = 0;
518:   primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID;
519:   primme->svd = svd;

521:   svd->ops->solve          = SVDSolve_PRIMME;
522:   svd->ops->setup          = SVDSetUp_PRIMME;
523:   svd->ops->setfromoptions = SVDSetFromOptions_PRIMME;
524:   svd->ops->destroy        = SVDDestroy_PRIMME;
525:   svd->ops->reset          = SVDReset_PRIMME;
526:   svd->ops->view           = SVDView_PRIMME;

528:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME);
529:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME);
530:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME);
531:   PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME);
532:   return(0);
533: }