Actual source code: svdprimme.c
slepc-3.18.3 2023-03-24
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: SVDCheckDefinite(svd);
155: MPI_Comm_size(PetscObjectComm((PetscObject)svd),&numProcs);
156: MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&procID);
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 = SlepcDefaultTol(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;
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) MatCreateVecsEmpty(svd->A,&ops->x,&ops->y);
223: return 0;
224: }
226: PetscErrorCode SVDSolve_PRIMME(SVD svd)
227: {
228: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
229: PetscScalar *svecs, *a;
230: PetscInt i,ierrprimme;
231: PetscReal *svals,*rnorms;
233: /* Reset some parameters left from previous runs */
234: ops->primme.aNorm = 0.0;
235: ops->primme.initSize = svd->nini;
236: ops->primme.iseed[0] = -1;
237: ops->primme.iseed[1] = -1;
238: ops->primme.iseed[2] = -1;
239: ops->primme.iseed[3] = -1;
241: /* Allocating left and right singular vectors contiguously */
242: PetscCalloc1(ops->primme.numSvals*(ops->primme.mLocal+ops->primme.nLocal),&svecs);
244: /* Call PRIMME solver */
245: PetscMalloc2(svd->ncv,&svals,svd->ncv,&rnorms);
246: ierrprimme = PRIMME_DRIVER(svals,svecs,rnorms,&ops->primme);
247: for (i=0;i<svd->ncv;i++) svd->sigma[i] = svals[i];
248: for (i=0;i<svd->ncv;i++) svd->errest[i] = rnorms[i];
249: PetscFree2(svals,rnorms);
251: /* Copy left and right singular vectors into svd */
252: BVGetArray(svd->U,&a);
253: PetscArraycpy(a,svecs,ops->primme.mLocal*ops->primme.initSize);
254: BVRestoreArray(svd->U,&a);
256: BVGetArray(svd->V,&a);
257: PetscArraycpy(a,svecs+ops->primme.mLocal*ops->primme.initSize,ops->primme.nLocal*ops->primme.initSize);
258: BVRestoreArray(svd->V,&a);
260: PetscFree(svecs);
262: svd->nconv = ops->primme.initSize >= 0 ? ops->primme.initSize : 0;
263: svd->reason = svd->nconv >= svd->nsv ? SVD_CONVERGED_TOL: SVD_DIVERGED_ITS;
264: svd->its = ops->primme.stats.numOuterIterations;
266: /* Process PRIMME error code */
267: if (ierrprimme != 0) {
268: switch (ierrprimme%100) {
269: case -1:
270: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": unexpected error",ierrprimme);
271: case -2:
272: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_LIB,"PRIMME library failed with error code=%" PetscInt_FMT ": allocation error",ierrprimme);
273: case -3: /* stop due to maximum number of iterations or matvecs */
274: break;
275: default:
278: }
279: }
280: return 0;
281: }
283: PetscErrorCode SVDReset_PRIMME(SVD svd)
284: {
285: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
287: primme_svds_free(&ops->primme);
288: VecDestroy(&ops->x);
289: VecDestroy(&ops->y);
290: return 0;
291: }
293: PetscErrorCode SVDDestroy_PRIMME(SVD svd)
294: {
295: PetscFree(svd->data);
296: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",NULL);
297: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",NULL);
298: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",NULL);
299: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",NULL);
300: return 0;
301: }
303: PetscErrorCode SVDView_PRIMME(SVD svd,PetscViewer viewer)
304: {
305: PetscBool isascii;
306: SVD_PRIMME *ctx = (SVD_PRIMME*)svd->data;
307: PetscMPIInt rank;
309: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
310: if (isascii) {
311: PetscViewerASCIIPrintf(viewer," block size=%" PetscInt_FMT "\n",ctx->bs);
312: PetscViewerASCIIPrintf(viewer," solver method: %s\n",SVDPRIMMEMethods[(SVDPRIMMEMethod)ctx->method]);
314: /* Display PRIMME params */
315: MPI_Comm_rank(PetscObjectComm((PetscObject)svd),&rank);
316: if (!rank) primme_svds_display_params(ctx->primme);
317: }
318: return 0;
319: }
321: PetscErrorCode SVDSetFromOptions_PRIMME(SVD svd,PetscOptionItems *PetscOptionsObject)
322: {
323: SVD_PRIMME *ctx = (SVD_PRIMME*)svd->data;
324: PetscInt bs;
325: SVDPRIMMEMethod meth;
326: PetscBool flg;
328: PetscOptionsHeadBegin(PetscOptionsObject,"SVD PRIMME Options");
330: PetscOptionsInt("-svd_primme_blocksize","Maximum block size","SVDPRIMMESetBlockSize",ctx->bs,&bs,&flg);
331: if (flg) SVDPRIMMESetBlockSize(svd,bs);
333: PetscOptionsEnum("-svd_primme_method","Method for solving the singular value problem","SVDPRIMMESetMethod",SVDPRIMMEMethods,(PetscEnum)ctx->method,(PetscEnum*)&meth,&flg);
334: if (flg) SVDPRIMMESetMethod(svd,meth);
336: PetscOptionsHeadEnd();
337: return 0;
338: }
340: static PetscErrorCode SVDPRIMMESetBlockSize_PRIMME(SVD svd,PetscInt bs)
341: {
342: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
344: if (bs == PETSC_DEFAULT) ops->bs = 0;
345: else {
347: ops->bs = bs;
348: }
349: return 0;
350: }
352: /*@
353: SVDPRIMMESetBlockSize - The maximum block size that PRIMME will try to use.
355: Logically Collective on svd
357: Input Parameters:
358: + svd - the singular value solver context
359: - bs - block size
361: Options Database Key:
362: . -svd_primme_blocksize - Sets the max allowed block size value
364: Notes:
365: If the block size is not set, the value established by primme_svds_initialize
366: is used.
368: The user should set the block size based on the architecture specifics
369: of the target computer, as well as any a priori knowledge of multiplicities.
370: The code does NOT require bs > 1 to find multiple eigenvalues. For some
371: methods, keeping bs = 1 yields the best overall performance.
373: Level: advanced
375: .seealso: SVDPRIMMEGetBlockSize()
376: @*/
377: PetscErrorCode SVDPRIMMESetBlockSize(SVD svd,PetscInt bs)
378: {
381: PetscTryMethod(svd,"SVDPRIMMESetBlockSize_C",(SVD,PetscInt),(svd,bs));
382: return 0;
383: }
385: static PetscErrorCode SVDPRIMMEGetBlockSize_PRIMME(SVD svd,PetscInt *bs)
386: {
387: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
389: *bs = ops->bs;
390: return 0;
391: }
393: /*@
394: SVDPRIMMEGetBlockSize - Get the maximum block size the code will try to use.
396: Not Collective
398: Input Parameter:
399: . svd - the singular value solver context
401: Output Parameter:
402: . bs - returned block size
404: Level: advanced
406: .seealso: SVDPRIMMESetBlockSize()
407: @*/
408: PetscErrorCode SVDPRIMMEGetBlockSize(SVD svd,PetscInt *bs)
409: {
412: PetscUseMethod(svd,"SVDPRIMMEGetBlockSize_C",(SVD,PetscInt*),(svd,bs));
413: return 0;
414: }
416: static PetscErrorCode SVDPRIMMESetMethod_PRIMME(SVD svd,SVDPRIMMEMethod method)
417: {
418: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
420: ops->method = (primme_svds_preset_method)method;
421: return 0;
422: }
424: /*@
425: SVDPRIMMESetMethod - Sets the method for the PRIMME SVD solver.
427: Logically Collective on svd
429: Input Parameters:
430: + svd - the singular value solver context
431: - method - method that will be used by PRIMME
433: Options Database Key:
434: . -svd_primme_method - Sets the method for the PRIMME SVD solver
436: Note:
437: If not set, the method defaults to SVD_PRIMME_HYBRID.
439: Level: advanced
441: .seealso: SVDPRIMMEGetMethod(), SVDPRIMMEMethod
442: @*/
443: PetscErrorCode SVDPRIMMESetMethod(SVD svd,SVDPRIMMEMethod method)
444: {
447: PetscTryMethod(svd,"SVDPRIMMESetMethod_C",(SVD,SVDPRIMMEMethod),(svd,method));
448: return 0;
449: }
451: static PetscErrorCode SVDPRIMMEGetMethod_PRIMME(SVD svd,SVDPRIMMEMethod *method)
452: {
453: SVD_PRIMME *ops = (SVD_PRIMME*)svd->data;
455: *method = (SVDPRIMMEMethod)ops->method;
456: return 0;
457: }
459: /*@
460: SVDPRIMMEGetMethod - Gets the method for the PRIMME SVD solver.
462: Not Collective
464: Input Parameter:
465: . svd - the singular value solver context
467: Output Parameter:
468: . method - method that will be used by PRIMME
470: Level: advanced
472: .seealso: SVDPRIMMESetMethod(), SVDPRIMMEMethod
473: @*/
474: PetscErrorCode SVDPRIMMEGetMethod(SVD svd,SVDPRIMMEMethod *method)
475: {
478: PetscUseMethod(svd,"SVDPRIMMEGetMethod_C",(SVD,SVDPRIMMEMethod*),(svd,method));
479: return 0;
480: }
482: SLEPC_EXTERN PetscErrorCode SVDCreate_PRIMME(SVD svd)
483: {
484: SVD_PRIMME *primme;
486: PetscNew(&primme);
487: svd->data = (void*)primme;
489: primme_svds_initialize(&primme->primme);
490: primme->bs = 0;
491: primme->method = (primme_svds_preset_method)SVD_PRIMME_HYBRID;
492: primme->svd = svd;
494: svd->ops->solve = SVDSolve_PRIMME;
495: svd->ops->setup = SVDSetUp_PRIMME;
496: svd->ops->setfromoptions = SVDSetFromOptions_PRIMME;
497: svd->ops->destroy = SVDDestroy_PRIMME;
498: svd->ops->reset = SVDReset_PRIMME;
499: svd->ops->view = SVDView_PRIMME;
501: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetBlockSize_C",SVDPRIMMESetBlockSize_PRIMME);
502: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetBlockSize_C",SVDPRIMMEGetBlockSize_PRIMME);
503: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMESetMethod_C",SVDPRIMMESetMethod_PRIMME);
504: PetscObjectComposeFunction((PetscObject)svd,"SVDPRIMMEGetMethod_C",SVDPRIMMEGetMethod_PRIMME);
505: return 0;
506: }