Actual source code: blopex.c
slepc-3.7.1 2016-05-27
1: /*
2: This file implements a wrapper to the BLOPEX package
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
25: #include "slepc-interface.h"
26: #include <blopex_lobpcg.h>
27: #include <blopex_interpreter.h>
28: #include <blopex_multivector.h>
29: #include <blopex_temp_multivector.h>
31: PetscInt slepc_blopex_useconstr = -1;
33: PetscErrorCode EPSSolve_BLOPEX(EPS);
35: typedef struct {
36: lobpcg_Tolerance tol;
37: lobpcg_BLASLAPACKFunctions blap_fn;
38: mv_InterfaceInterpreter ii;
39: ST st;
40: Vec w;
41: PetscInt bs; /* block size */
42: } EPS_BLOPEX;
46: static void Precond_FnSingleVector(void *data,void *x,void *y)
47: {
49: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
50: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
51: KSP ksp;
54: STGetKSP(blopex->st,&ksp);CHKERRABORT(comm,ierr);
55: KSPSolve(ksp,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
56: PetscFunctionReturnVoid();
57: }
61: static void Precond_FnMultiVector(void *data,void *x,void *y)
62: {
63: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
66: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
67: PetscFunctionReturnVoid();
68: }
72: static void OperatorASingleVector(void *data,void *x,void *y)
73: {
75: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
76: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
77: Mat A,B;
78: PetscScalar sigma;
79: PetscInt nmat;
82: STGetNumMatrices(blopex->st,&nmat);CHKERRABORT(comm,ierr);
83: STGetOperators(blopex->st,0,&A);CHKERRABORT(comm,ierr);
84: if (nmat>1) { STGetOperators(blopex->st,1,&B);CHKERRABORT(comm,ierr); }
85: MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
86: STGetShift(blopex->st,&sigma);CHKERRABORT(comm,ierr);
87: if (sigma != 0.0) {
88: if (nmat>1) {
89: MatMult(B,(Vec)x,blopex->w);CHKERRABORT(comm,ierr);
90: } else {
91: VecCopy((Vec)x,blopex->w);CHKERRABORT(comm,ierr);
92: }
93: VecAXPY((Vec)y,-sigma,blopex->w);CHKERRABORT(comm,ierr);
94: }
95: PetscFunctionReturnVoid();
96: }
100: static void OperatorAMultiVector(void *data,void *x,void *y)
101: {
102: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
105: blopex->ii.Eval(OperatorASingleVector,data,x,y);
106: PetscFunctionReturnVoid();
107: }
111: static void OperatorBSingleVector(void *data,void *x,void *y)
112: {
114: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
115: MPI_Comm comm = PetscObjectComm((PetscObject)blopex->st);
116: Mat B;
119: STGetOperators(blopex->st,1,&B);CHKERRABORT(comm,ierr);
120: MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(comm,ierr);
121: PetscFunctionReturnVoid();
122: }
126: static void OperatorBMultiVector(void *data,void *x,void *y)
127: {
128: EPS_BLOPEX *blopex = (EPS_BLOPEX*)data;
131: blopex->ii.Eval(OperatorBSingleVector,data,x,y);
132: PetscFunctionReturnVoid();
133: }
137: PetscErrorCode EPSSetDimensions_BLOPEX(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
138: {
140: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
141: PetscInt k;
144: k = ((eps->nev-1)/ctx->bs+1)*ctx->bs;
145: if (*ncv) { /* ncv set */
146: if (*ncv<k) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv is not sufficiently large");
147: } else *ncv = k;
148: if (!*mpd) *mpd = *ncv;
149: else { PetscInfo(eps,"Warning: given value of mpd ignored\n"); }
150: return(0);
151: }
155: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
156: {
157: #if defined(PETSC_MISSING_LAPACK_POTRF) || defined(PETSC_MISSING_LAPACK_SYGV)
159: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF/SYGV - Lapack routine is unavailable");
160: #else
162: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
163: PetscBool isPrecond,istrivial,flg;
166: if (!eps->ishermitian || (eps->isgeneralized && !eps->ispositive)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"blopex only works for Hermitian problems");
167: if (!blopex->bs) blopex->bs = PetscMin(16,eps->nev);
168: EPSSetDimensions_BLOPEX(eps,eps->nev,&eps->ncv,&eps->mpd);
169: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
170: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
171: if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");
172: if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
173: if (eps->stopping!=EPSStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"External packages do not support user-defined stopping test");
174: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
175: RGIsTrivial(eps->rg,&istrivial);
176: if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");
178: STSetUp(eps->st);
179: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&isPrecond);
180: if (!isPrecond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"blopex only works with STPRECOND");
181: blopex->st = eps->st;
183: if (eps->converged == EPSConvergedRelative) {
184: blopex->tol.absolute = 0.0;
185: blopex->tol.relative = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
186: } else if (eps->converged == EPSConvergedAbsolute) {
187: blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
188: blopex->tol.relative = 0.0;
189: } else SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Convergence test not supported in this solver");
191: SLEPCSetupInterpreter(&blopex->ii);
193: /* allocate memory */
194: if (!eps->V) { EPSGetBV(eps,&eps->V); }
195: PetscObjectTypeCompareAny((PetscObject)eps->V,&flg,BVVECS,BVCONTIGUOUS,"");
196: if (!flg) { /* blopex only works with BVVECS or BVCONTIGUOUS */
197: BVSetType(eps->V,BVCONTIGUOUS);
198: }
199: EPSAllocateSolution(eps,0);
200: BVCreateVec(eps->V,&blopex->w);
201: PetscLogObjectParent((PetscObject)eps,(PetscObject)blopex->w);
203: #if defined(PETSC_USE_COMPLEX)
204: blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
205: blopex->blap_fn.zhegv = PETSC_zsygv_interface;
206: #else
207: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
208: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
209: #endif
211: /* dispatch solve method */
212: eps->ops->solve = EPSSolve_BLOPEX;
213: return(0);
214: #endif
215: }
219: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
220: {
221: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
222: PetscScalar sigma,*eigr=NULL;
223: PetscReal *errest=NULL;
224: int i,j,info,its,nconv;
225: double *residhist=NULL;
226: PetscErrorCode ierr;
227: mv_MultiVectorPtr eigenvectors,constraints;
228: #if defined(PETSC_USE_COMPLEX)
229: komplex *lambda=NULL,*lambdahist=NULL;
230: #else
231: double *lambda=NULL,*lambdahist=NULL;
232: #endif
235: STGetShift(eps->st,&sigma);
236: PetscMalloc1(blopex->bs,&lambda);
237: if (eps->numbermonitors>0) {
238: PetscMalloc4(blopex->bs*(eps->max_it+1),&lambdahist,eps->ncv,&eigr,blopex->bs*(eps->max_it+1),&residhist,eps->ncv,&errest);
239: }
241: /* Complete the initial basis with random vectors */
242: for (i=eps->nini;i<eps->ncv;i++) {
243: BVSetRandomColumn(eps->V,i);
244: }
246: while (eps->reason == EPS_CONVERGED_ITERATING) {
248: /* Create multivector of constraints from leading columns of V */
249: PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,1);
250: BVSetActiveColumns(eps->V,0,eps->nconv);
251: constraints = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds+eps->nconv,eps->V);
253: /* Create multivector where eigenvectors of this run will be stored */
254: PetscObjectComposedDataSetInt((PetscObject)eps->V,slepc_blopex_useconstr,0);
255: BVSetActiveColumns(eps->V,eps->nconv,eps->nconv+blopex->bs);
256: eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,blopex->bs,eps->V);
258: #if defined(PETSC_USE_COMPLEX)
259: info = lobpcg_solve_complex(eigenvectors,blopex,OperatorAMultiVector,
260: eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
261: blopex,Precond_FnMultiVector,constraints,
262: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
263: lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs);
264: #else
265: info = lobpcg_solve_double(eigenvectors,blopex,OperatorAMultiVector,
266: eps->isgeneralized?blopex:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
267: blopex,Precond_FnMultiVector,constraints,
268: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
269: lambda,lambdahist,blopex->bs,eps->errest+eps->nconv,residhist,blopex->bs);
270: #endif
271: if (info>0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"BLOPEX failed with exit code=%d",info);
272: mv_MultiVectorDestroy(constraints);
273: mv_MultiVectorDestroy(eigenvectors);
275: for (j=0;j<blopex->bs;j++) {
276: #if defined(PETSC_USE_COMPLEX)
277: eps->eigr[eps->nconv+j] = lambda[j].real+PETSC_i*lambda[j].imag;
278: #else
279: eps->eigr[eps->nconv+j] = lambda[j];
280: #endif
281: }
283: if (eps->numbermonitors>0) {
284: for (i=0;i<its;i++) {
285: nconv = 0;
286: for (j=0;j<blopex->bs;j++) {
287: #if defined(PETSC_USE_COMPLEX)
288: eigr[eps->nconv+j] = lambdahist[j+i*blopex->bs].real+PETSC_i*lambdahist[j+i*blopex->bs].imag;
289: #else
290: eigr[eps->nconv+j] = lambdahist[j+i*blopex->bs];
291: #endif
292: errest[eps->nconv+j] = residhist[j+i*blopex->bs];
293: if (residhist[j+i*blopex->bs]<=eps->tol) nconv++;
294: }
295: EPSMonitor(eps,eps->its+i,eps->nconv+nconv,eigr,eps->eigi,errest,eps->nconv+blopex->bs);
296: }
297: }
299: eps->its += its;
300: if (info==-1) {
301: eps->reason = EPS_DIVERGED_ITS;
302: break;
303: } else {
304: for (i=0;i<blopex->bs;i++) {
305: if (sigma != 0.0) eps->eigr[eps->nconv+i] += sigma;
306: }
307: eps->nconv += blopex->bs;
308: if (eps->nconv>=eps->nev) eps->reason = EPS_CONVERGED_TOL;
309: }
310: }
312: PetscFree(lambda);
313: if (eps->numbermonitors>0) {
314: PetscFree4(lambdahist,eigr,residhist,errest);
315: }
316: return(0);
317: }
321: static PetscErrorCode EPSBLOPEXSetBlockSize_BLOPEX(EPS eps,PetscInt bs)
322: {
323: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
326: ctx->bs = bs;
327: return(0);
328: }
332: /*@
333: EPSBLOPEXSetBlockSize - Sets the block size of the BLOPEX solver.
335: Logically Collective on EPS
337: Input Parameters:
338: + eps - the eigenproblem solver context
339: - bs - the block size
341: Options Database Key:
342: . -eps_blopex_blocksize - Sets the block size
344: Level: advanced
346: .seealso: EPSBLOPEXGetBlockSize()
347: @*/
348: PetscErrorCode EPSBLOPEXSetBlockSize(EPS eps,PetscInt bs)
349: {
355: PetscTryMethod(eps,"EPSBLOPEXSetBlockSize_C",(EPS,PetscInt),(eps,bs));
356: return(0);
357: }
361: static PetscErrorCode EPSBLOPEXGetBlockSize_BLOPEX(EPS eps,PetscInt *bs)
362: {
363: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
366: *bs = ctx->bs;
367: return(0);
368: }
372: /*@
373: EPSBLOPEXGetBlockSize - Gets the block size used in the BLOPEX solver.
375: Not Collective
377: Input Parameter:
378: . eps - the eigenproblem solver context
380: Output Parameter:
381: . bs - the block size
383: Level: advanced
385: .seealso: EPSBLOPEXSetBlockSize()
386: @*/
387: PetscErrorCode EPSBLOPEXGetBlockSize(EPS eps,PetscInt *bs)
388: {
394: PetscUseMethod(eps,"EPSBLOPEXGetBlockSize_C",(EPS,PetscInt*),(eps,bs));
395: return(0);
396: }
400: PetscErrorCode EPSReset_BLOPEX(EPS eps)
401: {
403: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
406: VecDestroy(&blopex->w);
407: return(0);
408: }
412: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
413: {
417: LOBPCG_DestroyRandomContext();
418: PetscFree(eps->data);
419: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",NULL);
420: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",NULL);
421: return(0);
422: }
426: PetscErrorCode EPSView_BLOPEX(EPS eps,PetscViewer viewer)
427: {
429: EPS_BLOPEX *ctx = (EPS_BLOPEX*)eps->data;
430: PetscBool isascii;
433: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
434: if (isascii) {
435: PetscViewerASCIIPrintf(viewer," BLOPEX: block size %D\n",ctx->bs);
436: }
437: return(0);
438: }
442: PetscErrorCode EPSSetFromOptions_BLOPEX(PetscOptionItems *PetscOptionsObject,EPS eps)
443: {
445: KSP ksp;
446: PetscBool flg;
447: PetscInt bs;
450: PetscOptionsHead(PetscOptionsObject,"EPS BLOPEX Options");
451: PetscOptionsInt("-eps_blopex_blocksize","BLOPEX block size","EPSBLOPEXSetBlockSize",20,&bs,&flg);
452: if (flg) {
453: EPSBLOPEXSetBlockSize(eps,bs);
454: }
455: LOBPCG_SetFromOptionsRandomContext();
457: /* Set STPrecond as the default ST */
458: if (!((PetscObject)eps->st)->type_name) {
459: STSetType(eps->st,STPRECOND);
460: }
461: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
463: /* Set the default options of the KSP */
464: STGetKSP(eps->st,&ksp);
465: if (!((PetscObject)ksp)->type_name) {
466: KSPSetType(ksp,KSPPREONLY);
467: }
468: PetscOptionsTail();
469: return(0);
470: }
474: PETSC_EXTERN PetscErrorCode EPSCreate_BLOPEX(EPS eps)
475: {
476: EPS_BLOPEX *ctx;
480: PetscNewLog(eps,&ctx);
481: eps->data = (void*)ctx;
483: eps->ops->setup = EPSSetUp_BLOPEX;
484: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
485: eps->ops->destroy = EPSDestroy_BLOPEX;
486: eps->ops->reset = EPSReset_BLOPEX;
487: eps->ops->view = EPSView_BLOPEX;
488: eps->ops->backtransform = EPSBackTransform_Default;
489: LOBPCG_InitRandomContext(PetscObjectComm((PetscObject)eps),NULL);
490: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXSetBlockSize_C",EPSBLOPEXSetBlockSize_BLOPEX);
491: PetscObjectComposeFunction((PetscObject)eps,"EPSBLOPEXGetBlockSize_C",EPSBLOPEXGetBlockSize_BLOPEX);
492: if (slepc_blopex_useconstr < 0) { PetscObjectComposedDataRegister(&slepc_blopex_useconstr); }
493: return(0);
494: }