Actual source code: blopex.c
1: /*
2: This file implements a wrapper to the BLOPEX solver
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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: PetscErrorCode EPSSolve_BLOPEX(EPS);
33: typedef struct {
34: lobpcg_Tolerance tol;
35: lobpcg_BLASLAPACKFunctions blap_fn;
36: mv_MultiVectorPtr eigenvectors;
37: mv_MultiVectorPtr Y;
38: mv_InterfaceInterpreter ii;
39: ST st;
40: Vec w;
41: } EPS_BLOPEX;
45: static void Precond_FnSingleVector(void *data,void *x,void *y)
46: {
48: EPS eps = (EPS)data;
49: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
50:
52: STAssociatedKSPSolve(blopex->st,(Vec)x,(Vec)y); CHKERRABORT(((PetscObject)eps)->comm,ierr);
53: PetscFunctionReturnVoid();
54: }
58: static void Precond_FnMultiVector(void *data,void *x,void *y)
59: {
60: EPS eps = (EPS)data;
61: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
64: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
65: PetscFunctionReturnVoid();
66: }
70: static void OperatorASingleVector(void *data,void *x,void *y)
71: {
73: EPS eps = (EPS)data;
74: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
75: Mat A,B;
76: PetscScalar sigma;
77:
79: STGetOperators(eps->OP,&A,&B);CHKERRABORT(((PetscObject)eps)->comm,ierr);
80: MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(((PetscObject)eps)->comm,ierr);
81: STGetShift(eps->OP,&sigma);CHKERRABORT(((PetscObject)eps)->comm,ierr);
82: if (sigma != 0.0) {
83: if (B) { MatMult(B,(Vec)x,blopex->w);CHKERRABORT(((PetscObject)eps)->comm,ierr); }
84: else { VecCopy((Vec)x,blopex->w);CHKERRABORT(((PetscObject)eps)->comm,ierr); }
85: VecAXPY((Vec)y,-sigma,blopex->w);CHKERRABORT(((PetscObject)eps)->comm,ierr);
86: }
87: PetscFunctionReturnVoid();
88: }
92: static void OperatorAMultiVector(void *data,void *x,void *y)
93: {
94: EPS eps = (EPS)data;
95: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
98: blopex->ii.Eval(OperatorASingleVector,data,x,y);
99: PetscFunctionReturnVoid();
100: }
104: static void OperatorBSingleVector(void *data,void *x,void *y)
105: {
107: EPS eps = (EPS)data;
108: Mat B;
109:
111: STGetOperators(eps->OP,PETSC_NULL,&B);CHKERRABORT(((PetscObject)eps)->comm,ierr);
112: MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(((PetscObject)eps)->comm,ierr);
113: PetscFunctionReturnVoid();
114: }
118: static void OperatorBMultiVector(void *data,void *x,void *y)
119: {
120: EPS eps = (EPS)data;
121: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
124: blopex->ii.Eval(OperatorBSingleVector,data,x,y);
125: PetscFunctionReturnVoid();
126: }
130: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
131: {
132: #if defined(PETSC_MISSING_LAPACK_POTRF) || defined(PETSC_MISSING_LAPACK_SYGV)
134: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF/SYGV - Lapack routine is unavailable");
135: #else
137: PetscInt i;
138: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
139: PetscBool isPrecond;
142: if (!eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"blopex only works for hermitian problems");
143: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
144: if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
146: /* Change the default sigma to inf if necessary */
147: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
148: eps->which == EPS_LARGEST_IMAGINARY) {
149: STSetDefaultShift(eps->OP,3e300);
150: }
152: STSetUp(eps->OP);
153: PetscObjectTypeCompare((PetscObject)eps->OP,STPRECOND,&isPrecond);
154: if (!isPrecond) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"blopex only works with STPRECOND");
155: blopex->st = eps->OP;
157: eps->ncv = eps->nev = PetscMin(eps->nev,eps->n);
158: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
159: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
160: if (eps->arbit_func) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");
162: EPSAllocateSolution(eps);
163: EPSDefaultGetWork(eps,1);
164:
165: blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
166: blopex->tol.relative = 1e-50;
167:
168: SLEPCSetupInterpreter(&blopex->ii);
169: blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->V);
170: for (i=0;i<eps->ncv;i++) { PetscObjectReference((PetscObject)eps->V[i]); }
171: mv_MultiVectorSetRandom(blopex->eigenvectors,1234);
172: VecDuplicate(eps->V[0],&blopex->w);
174: if (eps->nds > 0) {
175: blopex->Y = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds,eps->defl);
176: for (i=0;i<eps->nds;i++) { PetscObjectReference((PetscObject)eps->defl[i]); }
177: } else
178: blopex->Y = PETSC_NULL;
180: #if defined(PETSC_USE_COMPLEX)
181: blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
182: blopex->blap_fn.zhegv = PETSC_zsygv_interface;
183: #else
184: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
185: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
186: #endif
188: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
190: /* dispatch solve method */
191: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
192: eps->ops->solve = EPSSolve_BLOPEX;
193: return(0);
194: #endif
195: }
199: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
200: {
201: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
202: PetscScalar sigma;
203: int i,j,info,its,nconv;
204: double *residhist=PETSC_NULL;
206: #if defined(PETSC_USE_COMPLEX)
207: komplex *lambdahist=PETSC_NULL;
208: #else
209: double *lambdahist=PETSC_NULL;
210: #endif
211:
213: if (eps->numbermonitors>0) {
214: #if defined(PETSC_USE_COMPLEX)
215: PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(komplex),&lambdahist);
216: #else
217: PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(double),&lambdahist);
218: #endif
219: PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(double),&residhist);
220: }
222: #if defined(PETSC_USE_COMPLEX)
223: info = lobpcg_solve_complex(blopex->eigenvectors,eps,OperatorAMultiVector,
224: eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
225: eps,Precond_FnMultiVector,blopex->Y,
226: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
227: (komplex*)eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
228: #else
229: info = lobpcg_solve_double(blopex->eigenvectors,eps,OperatorAMultiVector,
230: eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
231: eps,Precond_FnMultiVector,blopex->Y,
232: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
233: eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
234: #endif
235: if (info>0) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error in blopex (code=%d)",info);
237: if (eps->numbermonitors>0) {
238: for (i=0;i<its;i++) {
239: nconv = 0;
240: for (j=0;j<eps->ncv;j++) { if (residhist[j+i*eps->ncv]>eps->tol) break; else nconv++; }
241: EPSMonitor(eps,i,nconv,(PetscScalar*)lambdahist+i*eps->ncv,eps->eigi,residhist+i*eps->ncv,eps->ncv);
242: }
243: PetscFree(lambdahist);
244: PetscFree(residhist);
245: }
247: eps->its = its;
248: eps->nconv = eps->ncv;
249: STGetShift(eps->OP,&sigma);
250: if (sigma != 0.0) {
251: for (i=0;i<eps->nconv;i++) eps->eigr[i]+=sigma;
252: }
253: if (info==-1) eps->reason = EPS_DIVERGED_ITS;
254: else eps->reason = EPS_CONVERGED_TOL;
255: return(0);
256: }
260: PetscErrorCode EPSReset_BLOPEX(EPS eps)
261: {
263: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
266: mv_MultiVectorDestroy(blopex->eigenvectors);
267: mv_MultiVectorDestroy(blopex->Y);
268: VecDestroy(&blopex->w);
269: EPSReset_Default(eps);
270: return(0);
271: }
275: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
276: {
280: LOBPCG_DestroyRandomContext();
281: PetscFree(eps->data);
282: return(0);
283: }
287: PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps)
288: {
289: PetscErrorCode ierr;
292: PetscOptionsHead("EPS BLOPEX Options");
293: LOBPCG_SetFromOptionsRandomContext();
294: PetscOptionsTail();
295: return(0);
296: }
298: EXTERN_C_BEGIN
301: PetscErrorCode EPSCreate_BLOPEX(EPS eps)
302: {
306: PetscNewLog(eps,EPS_BLOPEX,&eps->data);
307: eps->ops->setup = EPSSetUp_BLOPEX;
308: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
309: eps->ops->destroy = EPSDestroy_BLOPEX;
310: eps->ops->reset = EPSReset_BLOPEX;
311: eps->ops->backtransform = EPSBackTransform_Default;
312: eps->ops->computevectors = EPSComputeVectors_Default;
313: STSetType(eps->OP,STPRECOND);
314: STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);
315: LOBPCG_InitRandomContext(((PetscObject)eps)->comm,eps->rand);
316: return(0);
317: }
318: EXTERN_C_END