Actual source code: svdprimme.c
slepc-3.14.1 2020-12-08
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: *SVDMatMult(svd,PETSC_TRUE,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: *SVDMatMult(svd,PETSC_FALSE,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: 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: SVDMatGetSize(svd,&m,&n);
159: SVDMatGetLocalSize(svd,&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 = svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL: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: default:
204: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"'which' value not supported by PRIMME");
205: }
207: /* If user sets mpd or ncv, maxBasisSize is modified */
208: if (svd->mpd!=PETSC_DEFAULT) {
209: primme->maxBasisSize = svd->mpd;
210: if (svd->ncv!=PETSC_DEFAULT) { PetscInfo(svd,"Warning: 'ncv' is ignored by PRIMME\n"); }
211: } else if (svd->ncv!=PETSC_DEFAULT) primme->maxBasisSize = svd->ncv;
213: 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");
215: svd->mpd = primme->maxBasisSize;
216: svd->ncv = (primme->locking?svd->nsv:0)+primme->maxBasisSize;
217: ops->bs = primme->maxBlockSize;
219: /* Set workspace */
220: SVDAllocateSolution(svd,0);
222: /* Prepare auxiliary vectors */
223: if (!ops->x) {
224: MatCreateVecsEmpty(svd->A,&ops->x,&ops->y);
225: PetscLogObjectParent((PetscObject)svd,(PetscObject)ops->x);
226: PetscLogObjectParent((PetscObject)svd,(PetscObject)ops->y);
227: }
228: return(0);
229: }
231: PetscErrorCode SVDSolve_PRIMME(SVD svd)
232: {
234: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
235: PetscScalar *svecs, *a;
236: PetscInt i,ierrprimme;
237: PetscReal *svals,*rnorms;
240: /* Reset some parameters left from previous runs */
241: ops->primme.aNorm = 0.0;
242: ops->primme.initSize = svd->nini;
243: ops->primme.iseed[0] = -1;
244: ops->primme.iseed[1] = -1;
245: ops->primme.iseed[2] = -1;
246: ops->primme.iseed[3] = -1;
248: /* Allocating left and right singular vectors contiguously */
249: PetscCalloc1(ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal),&svecs);
250: PetscLogObjectMemory((PetscObject)svd,sizeof(PetscReal)*ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal));
252: /* Call PRIMME solver */
253: PetscMalloc2(svd->ncv,&svals,svd->ncv,&rnorms);
254: ierrprimme = PRIMME_DRIVER(svals,svecs,rnorms,&ops->primme);
255: for (i=0;i<svd->ncv;i++) svd->sigma[i] = svals[i];
256: for (i=0;i<svd->ncv;i++) svd->errest[i] = rnorms[i];
257: PetscFree2(svals,rnorms);
259: /* Copy left and right singular vectors into svd */
260: BVGetArray(svd->U,&a);
261: PetscArraycpy(a,svecs,ops->primme.mLocal*ops->primme.initSize);
262: BVRestoreArray(svd->U,&a);
264: BVGetArray(svd->V,&a);
265: PetscArraycpy(a,svecs+ops->primme.mLocal*ops->primme.initSize,ops->primme.nLocal*ops->primme.initSize);
266: BVRestoreArray(svd->V,&a);
268: PetscFree(svecs);
270: svd->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
271: svd->reason = svd->nconv >= svd->nsv ? SVD_CONVERGED_TOL: SVD_DIVERGED_ITS;
272: svd->its = ops->primme.stats.numOuterIterations;
274: /* Process PRIMME error code */
275: if (ierrprimme == 0) {
276: /* no error */
277: } else if (ierrprimme%100 == -1) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: unexpected error",ierrprimme);
278: else if (ierrprimme%100 == -2) SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: allocation error",ierrprimme);
279: else if (ierrprimme%100 == -3) {
280: /* stop by maximum number of iteration or matvecs */
281: } 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);
282: else SETERRQ1(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%d: runtime error; check PRIMME's manual",ierrprimme);
283: return(0);
284: }
286: PetscErrorCode SVDReset_PRIMME(SVD svd)
287: {
289: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
292: primme_svds_free(&ops->primme);
293: VecDestroy(&ops->x);
294: VecDestroy(&ops->y);
295: return(0);
296: }
298: PetscErrorCode SVDDestroy_PRIMME(SVD svd)
299: {
303: PetscFree(svd->data);
304: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL);
305: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL);
306: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL);
307: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL);
308: return(0);
309: }
311: PetscErrorCode SVDView_PRIMME(SVD svd,PetscViewer viewer)
312: {
314: PetscBool isascii;
315: SVD_PRIMME *ctx = (SVD_PRIMME*)svd->data;
316: PetscMPIInt rank;
319: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
320: if (isascii) {
321: PetscViewerASCIIPrintf(viewer," block size=%D\n",ctx->bs);
322: PetscViewerASCIIPrintf(viewer," solver method: %s\n",SVDPRIMMEMethods[(SVDPRIMMEMethod)ctx->method]);
324: /* Display PRIMME params */
325: MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank);
326: if (!rank) primme_svds_display_params(ctx->primme);
327: }
328: return(0);
329: }
331: PetscErrorCode SVDSetFromOptions_PRIMME(PetscOptionItems *PetscOptionsObject,SVD svd)
332: {
333: PetscErrorCode ierr;
334: SVD_PRIMME *ctx = (SVD_PRIMME*)svd->data;
335: PetscInt bs;
336: SVDPRIMMEMethod meth;
337: PetscBool flg;
340: PetscOptionsHead(PetscOptionsObject,"SVD PRIMME Options");
342: PetscOptionsInt("-svd_primme_blocksize","Maximum block size","SVDPRIMMESetBlockSize",ctx->bs,&bs,&flg);
343: if (flg) { SVDPRIMMESetBlockSize(svd,bs); }
345: PetscOptionsEnum("-svd_primme_method","Method for solving the singular value problem","SVDPRIMMESetMethod",SVDPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
346: if (flg) { SVDPRIMMESetMethod(svd,meth); }
348: PetscOptionsTail();
349: return(0);
350: }
352: static PetscErrorCode SVDPRIMMESetBlockSize_PRIMME(SVD svd,PetscInt bs)
353: {
354: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
357: if (bs == PETSC_DEFAULT) ops->bs = 0;
358: else if (bs <= 0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"PRIMME: block size must be positive");
359: else ops->bs = bs;
360: return(0);
361: }
363: /*@
364: SVDPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
366: Logically Collective on svd
368: Input Parameters:
369: + svd - the singular value solver context
370: - bs - block size
372: Options Database Key:
373: . -svd_primme_blocksize - Sets the max allowed block size value
375: Notes:
376: If the block size is not set, the value established by primme_svds_initialize
377: is used.
379: The user should set the block size based on the architecture specifics
380: of the target computer, as well as any a priori knowledge of multiplicities.
381: The code does NOT require bs > 1 to find multiple eigenvalues. For some
382: methods, keeping bs = 1 yields the best overall performance.
384: Level: advanced
386: .seealso: SVDPRIMMEGetBlockSize()
387: @*/
388: PetscErrorCode SVDPRIMMESetBlockSize(SVD svd,PetscInt bs)
389: {
395: PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs));
396: return(0);
397: }
399: static PetscErrorCode SVDPRIMMEGetBlockSize_PRIMME(SVD svd,PetscInt *bs)
400: {
401: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
404: *bs = ops->bs;
405: return(0);
406: }
408: /*@
409: SVDPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
411: Not Collective
413: Input Parameter:
414: . svd - the singular value solver context
416: Output Parameter:
417: . bs - returned block size
419: Level: advanced
421: .seealso: SVDPRIMMESetBlockSize()
422: @*/
423: PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs)
424: {
430: PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs));
431: return(0);
432: }
434: static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method)
435: {
436: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
439: ops->method = (primme_svds_preset_method)method;
440: return(0);
441: }
443: /*@
444: SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver.
446: Logically Collective on svd
448: Input Parameters:
449: + svd - the singular value solver context
450: - method - method that will be used by PRIMME
452: Options Database Key:
453: . -svd_primme_method - Sets the method for the PRIMME SVD solver
455: Note:
456: If not set, the method defaults to SVD_PRIMME_HYBRID.
458: Level: advanced
460: .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod
461: @*/
462: PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method)
463: {
469: PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method));
470: return(0);
471: }
473: static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method)
474: {
475: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
478: *method = (SVDPRIMMEMethod)ops->method;
479: return(0);
480: }
482: /*@
483: SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver.
485: Not Collective
487: Input Parameter:
488: . svd - the singular value solver context
490: Output Parameter:
491: . method - method that will be used by PRIMME
493: Level: advanced
495: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod
496: @*/
497: PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method)
498: {
504: PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method));
505: return(0);
506: }
508: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd)
509: {
511: SVD_PRIMME *primme;
514: PetscNewLog(svd,&primme);
515: svd->data = (void*)primme;
517: primme_svds_initialize(&primme->primme);
518: primme->bs = 0;
519: primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID;
520: primme->svd = svd;
522: svd->ops->solve = SVDSolve_PRIMME;
523: svd->ops->setup = SVDSetUp_PRIMME;
524: svd->ops->setfromoptions = SVDSetFromOptions_PRIMME;
525: svd->ops->destroy = SVDDestroy_PRIMME;
526: svd->ops->reset = SVDReset_PRIMME;
527: svd->ops->view = SVDView_PRIMME;
529: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME);
530: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME);
531: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME);
532: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME);
533: return(0);
534: }