Actual source code: blopex.c
slepc-3.13.0 2020-03-31
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 BLOPEX package
12: */
14: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
15: #include "blopex.h"
16: #include <lobpcg.h>
17: #include <interpreter.h>
18: #include <multivector.h>
19: #include <temp_multivector.h>
21: PetscInt slepc_blopex_useconstr = -1;
23: typedef struct {
24: lobpcg_Tolerance tol;
25: lobpcg_BLASLAPACKFunctions blap_fn;
26: mv_InterfaceInterpreter ii;
27: ST st;
28: Vec w;
29: PetscInt bs; /* block size */
30: } EPS_BLOPEX;
32: static void Precond_FnSingleVector(void *data,void *x,void *y)
33: {
35: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
36: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
37: KSP ksp;
40: STGetKSP(blopex->st,&ksp);CHKERRABORT(comm,ierr);
41: KSPSolve(ksp,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
42: PetscFunctionReturnVoid();
43: }
45: static void Precond_FnMultiVector(void *data,void *x,void *y)
46: {
47: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
50: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
51: PetscFunctionReturnVoid();
52: }
54: static void OperatorASingleVector(void *data,void *x,void *y)
55: {
57: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
58: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
59: Mat A,B;
60: PetscScalar sigma;
61: PetscInt nmat;
64: STGetNumMatrices(blopex->st,&nmat);CHKERRABORT(comm,ierr);
65: STGetMatrix(blopex->st,0,&A);CHKERRABORT(comm,ierr);
66: if (nmat>1) { STGetMatrix(blopex->st,1,&B);CHKERRABORT(comm,ierr); }
67: MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
68: STGetShift(blopex->st,&sigma);CHKERRABORT(comm,ierr);
69: if (sigma != 0.0) {
70: if (nmat>1) {
71: MatMult(B,(Vec)x,blopex->w);CHKERRABORT(comm,ierr);
72: } else {
73: VecCopy((Vec)x,blopex->w);CHKERRABORT(comm,ierr);
74: }
75: VecAXPY((Vec)y,-sigma,blopex->w);CHKERRABORT(comm,ierr);
76: }
77: PetscFunctionReturnVoid();
78: }
80: static void OperatorAMultiVector(void *data,void *x,void *y)
81: {
82: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
85: blopex->ii.Eval(OperatorASingleVector,data,x,y);
86: PetscFunctionReturnVoid();
87: }
89: static void OperatorBSingleVector(void *data,void *x,void *y)
90: {
92: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
93: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
94: Mat B;
97: STGetMatrix(blopex->st,1,&B);CHKERRABORT(comm,ierr);
98: MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
99: PetscFunctionReturnVoid();
100: }
102: static void OperatorBMultiVector(void *data,void *x,void *y)
103: {
104: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
107: blopex->ii.Eval(OperatorBSingleVector,data,x,y);
108: PetscFunctionReturnVoid();
109: }
111: PetscErrorCode EPSSetDimensions_BLOPEX(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
112: {
114: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
115: PetscInt k;
118: k = ((eps->nev-1)/ctx->bs+1)*ctx->bs;
119: if (*ncv) { /* ncv set */
120: if (*ncv<k) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv is not sufficiently large");
121: } else *ncv = k;
122: if (!*mpd) *mpd = *ncv;
123: else { PetscInfo(eps,"Warning: given value of mpd ignored\n"); }
124: return(0);
125: }
127: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
128: {
130: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
131: PetscBool istrivial,flg;
132: KSP ksp;
135: if (!eps->ishermitian || (eps->isgeneralized && !eps->ispositive)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"blopex only works for Hermitian problems");
136: if (!blopex->bs) blopex->bs = PetscMin(16,eps->nev);
137: EPSSetDimensions_BLOPEX(eps,eps->nev,&eps->ncv,&eps->mpd);
138: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
139: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
140: if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
141: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
142: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
143: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
144: RGIsTrivial(eps->rg,&istrivial);
145: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
147: blopex->st = eps->st;
149: if (eps->converged == EPSConvergedRelative) {
150: blopex->tol.absolute = 0.0;
151: blopex->tol.relative = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
152: } else if (eps->converged == EPSConvergedAbsolute) {
153: blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
154: blopex->tol.relative = 0.0;
155: } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Convergence test not supported in this solver");
157: SLEPCSetupInterpreter(&blopex->ii);
159: STGetKSP(eps->st,&ksp);
160: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
161: if (!flg) { PetscInfo(eps,"Warning: ignoring KSP, should use KSPPREONLY\n"); }
163: /* allocate memory */
164: if (!eps->V) { EPSGetBV(eps,&eps->V); }
165: PetscObjectTypeCompareAny((PetscObject)eps->V,&flg,BVVECS,BVCONTIGUOUS,"");
166: if (!flg) { /* blopex only works with BVVECS or BVCONTIGUOUS */
167: BVSetType(eps->V,BVCONTIGUOUS);
168: }
169: EPSAllocateSolution(eps,0);
170: if (!blopex->w) {
171: BVCreateVec(eps->V,&blopex->w);
172: PetscLogObjectParent((PetscObject)eps,(PetscObject)blopex->w);
173: }
175: #if defined(PETSC_USE_COMPLEX)
176: blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
177: blopex->blap_fn.zhegv = PETSC_zsygv_interface;
178: #else
179: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
180: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
181: #endif
182: return(0);
183: }
185: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
186: {
187: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
188: PetscScalar sigma,*eigr=NULL;
189: PetscReal *errest=NULL;
190: int i,j,info,its,nconv;
191: double *residhist=NULL;
192: PetscErrorCode ierr;
193: mv_MultiVectorPtr eigenvectors,constraints;
194: #if defined(PETSC_USE_COMPLEX)
195: komplex *lambda=NULL,*lambdahist=NULL;
196: #else
197: double *lambda=NULL,*lambdahist=NULL;
198: #endif
201: STGetShift(eps->st,&sigma);
202: PetscMalloc1(blopex->bs,&lambda);
203: if (eps->numbermonitors>0) {
204: PetscMalloc4(blopex->bs*(eps->max_it+1),&lambdahist,eps->ncv,&eigr,blopex->bs*(eps->max_it+1),&residhist,eps->ncv,&errest);
205: }
207: /* Complete the initial basis with random vectors */
208: for (i=0;i<eps->nini;i++) { /* in case the initial vectors were also set with VecSetRandom */
209: BVSetRandomColumn(eps->V,eps->nini);
210: }
211: for (i=eps->nini;i<eps->ncv;i++) {
212: BVSetRandomColumn(eps->V,i);
213: }
215: while (eps->reason == EPS_CONVERGED_ITERATING) {
217: /* Create multivector of constraints from leading columns of V */
218: PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,1);
219: BVSetActiveColumns(eps->V,0,eps->nconv);
220: constraints = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds+eps->nconv,eps->V);
222: /* Create multivector where eigenvectors of this run will be stored */
223: PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,0);
224: BVSetActiveColumns(eps->V,eps->nconv,eps->nconv+blopex->bs);
225: eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,blopex->bs,eps->V);
227: #if defined(PETSC_USE_COMPLEX)
228: info = lobpcg_solve_complex(eigenvectors,blopex,OperatorAMultiVector,
229: eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
230: blopex,Precond_FnMultiVector,constraints,
231: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
232: lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs);
233: #else
234: info = lobpcg_solve_double(eigenvectors,blopex,OperatorAMultiVector,
235: eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
236: blopex,Precond_FnMultiVector,constraints,
237: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
238: lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs);
239: #endif
240: if (info>0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"BLOPEX failed with exit code=%d",info);
241: mv_MultiVectorDestroy(constraints);
242: mv_MultiVectorDestroy(eigenvectors);
244: for (j=0;j<blopex->bs;j++) {
245: #if defined(PETSC_USE_COMPLEX)
246: eps->eigr[eps->nconv+j] = PetscCMPLX(lambda[j].real,lambda[j].imag);
247: #else
248: eps->eigr[eps->nconv+j] = lambda[j];
249: #endif
250: }
252: if (eps->numbermonitors>0) {
253: for (i=0;i<its;i++) {
254: nconv = 0;
255: for (j=0;j<blopex->bs;j++) {
256: #if defined(PETSC_USE_COMPLEX)
257: eigr[eps->nconv+j] = PetscCMPLX(lambdahist[j+i*blopex->bs].real,lambdahist[j+i*blopex->bs].imag);
258: #else
259: eigr[eps->nconv+j] = lambdahist[j+i*blopex->bs];
260: #endif
261: errest[eps->nconv+j] = residhist[j+i*blopex->bs];
262: if (residhist[j+i*blopex->bs]<=eps->tol) nconv++;
263: }
264: EPSMonitor(eps,eps->its+i,eps->nconv+nconv,eigr,eps->eigi,errest,eps->nconv+blopex->bs);
265: }
266: }
268: eps->its += its;
269: if (info==-1) {
270: eps->reason = EPS_DIVERGED_ITS;
271: break;
272: } else {
273: for (i=0;i<blopex->bs;i++) {
274: if (sigma != 0.0) eps->eigr[eps->nconv+i] += sigma;
275: }
276: eps->nconv += blopex->bs;
277: if (eps->nconv>=eps->nev) eps->reason = EPS_CONVERGED_TOL;
278: }
279: }
281: PetscFree(lambda);
282: if (eps->numbermonitors>0) {
283: PetscFree4(lambdahist,eigr,residhist,errest);
284: }
285: return(0);
286: }
288: static PetscErrorCode EPSBLOPEXSetBlockSize_BLOPEX(EPS eps,PetscInt bs)
289: {
290: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
293: if (bs==PETSC_DEFAULT) {
294: ctx->bs = 0;
295: eps->state = EPS_STATE_INITIAL;
296: } else {
297: if (bs<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Block size must be >0");
298: ctx->bs = bs;
299: }
300: return(0);
301: }
303: /*@
304: EPSBLOPEXSetBlockSize - Sets the block size of the BLOPEX solver.
306: Logically Collective on eps
308: Input Parameters:
309: + eps - the eigenproblem solver context
310: - bs - the block size
312: Options Database Key:
313: . -eps_blopex_blocksize - Sets the block size
315: Level: advanced
317: .seealso: EPSBLOPEXGetBlockSize()
318: @*/
319: PetscErrorCode EPSBLOPEXSetBlockSize(EPS eps,PetscInt bs)
320: {
326: PetscTryMethod(eps,"EPSBLOPEXSetBlockSize_C",(EPS,PetscInt),(eps,bs));
327: return(0);
328: }
330: static PetscErrorCode EPSBLOPEXGetBlockSize_BLOPEX(EPS eps,PetscInt *bs)
331: {
332: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
335: *bs = ctx->bs;
336: return(0);
337: }
339: /*@
340: EPSBLOPEXGetBlockSize - Gets the block size used in the BLOPEX solver.
342: Not Collective
344: Input Parameter:
345: . eps - the eigenproblem solver context
347: Output Parameter:
348: . bs - the block size
350: Level: advanced
352: .seealso: EPSBLOPEXSetBlockSize()
353: @*/
354: PetscErrorCode EPSBLOPEXGetBlockSize(EPS eps,PetscInt *bs)
355: {
361: PetscUseMethod(eps,"EPSBLOPEXGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
362: return(0);
363: }
365: PetscErrorCode EPSReset_BLOPEX(EPS eps)
366: {
368: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
371: VecDestroy(&blopex->w);
372: return(0);
373: }
375: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
376: {
380: LOBPCG_DestroyRandomContext();
381: PetscFree(eps->data);
382: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",NULL);
383: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",NULL);
384: return(0);
385: }
387: PetscErrorCode EPSView_BLOPEX(EPS eps,PetscViewer viewer)
388: {
390: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
391: PetscBool isascii;
394: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
395: if (isascii) {
396: PetscViewerASCIIPrintf(viewer," block size %D\n",ctx->bs);
397: }
398: return(0);
399: }
401: PetscErrorCode EPSSetFromOptions_BLOPEX(PetscOptionItems *PetscOptionsObject,EPS eps)
402: {
404: PetscBool flg;
405: PetscInt bs;
408: PetscOptionsHead(PetscOptionsObject,"EPS BLOPEX Options");
410: PetscOptionsInt("-eps_blopex_blocksize","Block size","EPSBLOPEXSetBlockSize",20,&bs,&flg);
411: if (flg) { EPSBLOPEXSetBlockSize(eps,bs); }
413: PetscOptionsTail();
415: LOBPCG_SetFromOptionsRandomContext();
416: return(0);
417: }
419: SLEPC_EXTERN PetscErrorCode EPSCreate_BLOPEX(EPS eps)
420: {
421: EPS_BLOPEX *ctx;
425: PetscNewLog(eps,&ctx);
426: eps->data = (void*)ctx;
428: eps->categ = EPS_CATEGORY_PRECOND;
430: eps->ops->solve = EPSSolve_BLOPEX;
431: eps->ops->setup = EPSSetUp_BLOPEX;
432: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
433: eps->ops->destroy = EPSDestroy_BLOPEX;
434: eps->ops->reset = EPSReset_BLOPEX;
435: eps->ops->view = EPSView_BLOPEX;
436: eps->ops->backtransform = EPSBackTransform_Default;
437: eps->ops->setdefaultst = EPSSetDefaultST_GMRES;
439: LOBPCG_InitRandomContext(PetscObjectComm((PetscObject)eps),NULL);
440: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",EPSBLOPEXSetBlockSize_BLOPEX);
441: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",EPSBLOPEXGetBlockSize_BLOPEX);
442: if (slepc_blopex_useconstr < 0) { PetscObjectComposedDataRegister(&slepc_blopex_useconstr); }
443: return(0);
444: }