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-2011, 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 <private/stimpl.h> /*I "slepcst.h" I*/
25: #include <private/epsimpl.h> /*I "slepceps.h" I*/
26: #include "slepc-interface.h"
27: #include <blopex_lobpcg.h>
28: #include <blopex_interpreter.h>
29: #include <blopex_multivector.h>
30: #include <blopex_temp_multivector.h>
32: PetscErrorCode EPSSolve_BLOPEX(EPS);
34: typedef struct {
35: lobpcg_Tolerance tol;
36: lobpcg_BLASLAPACKFunctions blap_fn;
37: mv_MultiVectorPtr eigenvectors;
38: mv_MultiVectorPtr Y;
39: mv_InterfaceInterpreter ii;
40: KSP ksp;
41: Vec w;
42: } EPS_BLOPEX;
46: static void Precond_FnSingleVector(void *data,void *x,void *y)
47: {
49: EPS eps = (EPS)data;
50: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
51: PetscInt lits;
52:
54: KSPSolve(blopex->ksp,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
55: KSPGetIterationNumber(blopex->ksp,&lits); CHKERRABORT(PETSC_COMM_WORLD,ierr);
56: eps->OP->lineariterations+= lits;
57: PetscFunctionReturnVoid();
58: }
62: static void Precond_FnMultiVector(void *data,void *x,void *y)
63: {
64: EPS eps = (EPS)data;
65: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
68: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
69: PetscFunctionReturnVoid();
70: }
74: static void OperatorASingleVector(void *data,void *x,void *y)
75: {
77: EPS eps = (EPS)data;
78: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
79: Mat A,B;
80:
82: STGetOperators(eps->OP,&A,&B);CHKERRABORT(PETSC_COMM_WORLD,ierr);
83: MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
84: if (eps->OP->sigma != 0.0) {
85: if (B) { MatMult(B,(Vec)x,blopex->w);CHKERRABORT(PETSC_COMM_WORLD,ierr); }
86: else { VecCopy((Vec)x,blopex->w);CHKERRABORT(PETSC_COMM_WORLD,ierr); }
87: VecAXPY((Vec)y,-eps->OP->sigma,blopex->w);CHKERRABORT(PETSC_COMM_WORLD,ierr);
88: }
89: PetscFunctionReturnVoid();
90: }
94: static void OperatorAMultiVector(void *data,void *x,void *y)
95: {
96: EPS eps = (EPS)data;
97: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
100: blopex->ii.Eval(OperatorASingleVector,data,x,y);
101: PetscFunctionReturnVoid();
102: }
106: static void OperatorBSingleVector(void *data,void *x,void *y)
107: {
109: EPS eps = (EPS)data;
110: Mat B;
111:
113: STGetOperators(eps->OP,PETSC_NULL,&B);CHKERRABORT(PETSC_COMM_WORLD,ierr);
114: MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(PETSC_COMM_WORLD,ierr);
115: PetscFunctionReturnVoid();
116: }
120: static void OperatorBMultiVector(void *data,void *x,void *y)
121: {
122: EPS eps = (EPS)data;
123: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
126: blopex->ii.Eval(OperatorBSingleVector,data,x,y);
127: PetscFunctionReturnVoid();
128: }
132: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
133: {
135: PetscInt i;
136: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
137: PetscBool isPrecond;
140: if (!eps->ishermitian) {
141: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"blopex only works for hermitian problems");
142: }
143: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
144: if (eps->which!=EPS_SMALLEST_REAL) {
145: SETERRQ(((PetscObject)eps)->comm,1,"Wrong value of eps->which");
146: }
148: /* Change the default sigma to inf if necessary */
149: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
150: eps->which == EPS_LARGEST_IMAGINARY) {
151: STSetDefaultShift(eps->OP,3e300);
152: }
154: STSetUp(eps->OP);
155: PetscTypeCompare((PetscObject)eps->OP,STPRECOND,&isPrecond);
156: if (!isPrecond) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"blopex only works with STPRECOND");
157: STGetKSP(eps->OP,&blopex->ksp);
159: eps->ncv = eps->nev = PetscMin(eps->nev,eps->n);
160: if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
161: if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
163: EPSAllocateSolution(eps);
164: EPSDefaultGetWork(eps,1);
165:
166: blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
167: blopex->tol.relative = 1e-50;
168:
169: SLEPCSetupInterpreter(&blopex->ii);
170: blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->V);
171: for (i=0;i<eps->ncv;i++) { PetscObjectReference((PetscObject)eps->V[i]); }
172: mv_MultiVectorSetRandom(blopex->eigenvectors,1234);
173: VecDuplicate(eps->V[0],&blopex->w);
175: if (eps->nds > 0) {
176: blopex->Y = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->nds,eps->DS);
177: for (i=0;i<eps->nds;i++) { PetscObjectReference((PetscObject)eps->DS[i]); }
178: } else
179: blopex->Y = PETSC_NULL;
181: #if defined(PETSC_USE_COMPLEX)
182: blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
183: blopex->blap_fn.zhegv = PETSC_zsygv_interface;
184: #else
185: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
186: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
187: #endif
189: if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
191: /* dispatch solve method */
192: if (eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Left vectors not supported in this solver");
193: eps->ops->solve = EPSSolve_BLOPEX;
194: return(0);
195: }
199: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
200: {
201: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
202: int i,j,info,its,nconv;
203: double *lambdahist=PETSC_NULL,*residhist=PETSC_NULL;
205:
207: if (eps->numbermonitors>0) {
208: PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(double),&lambdahist);
209: PetscMalloc(eps->ncv*(eps->max_it+1)*sizeof(double),&residhist);
210: }
212: #if defined(PETSC_USE_COMPLEX)
213: info = lobpcg_solve_complex(blopex->eigenvectors,eps,OperatorAMultiVector,
214: eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
215: eps,Precond_FnMultiVector,blopex->Y,
216: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
217: (komplex*)eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
218: #else
219: info = lobpcg_solve_double(blopex->eigenvectors,eps,OperatorAMultiVector,
220: eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
221: eps,Precond_FnMultiVector,blopex->Y,
222: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
223: eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
224: #endif
225: if (info>0) SETERRQ1(((PetscObject)eps)->comm,PETSC_ERR_LIB,"Error in blopex (code=%d)",info);
227: if (eps->numbermonitors>0) {
228: for (i=0;i<its;i++) {
229: nconv = 0;
230: for (j=0;j<eps->ncv;j++) { if (residhist[j+i*eps->ncv]>eps->tol) break; else nconv++; }
231: EPSMonitor(eps,i,nconv,lambdahist+i*eps->ncv,eps->eigi,residhist+i*eps->ncv,eps->ncv);
232: }
233: PetscFree(lambdahist);
234: PetscFree(residhist);
235: }
237: eps->its = its;
238: eps->nconv = eps->ncv;
239: if (eps->OP->sigma != 0.0) {
240: for (i=0;i<eps->nconv;i++) eps->eigr[i]+=eps->OP->sigma;
241: }
242: if (info==-1) eps->reason = EPS_DIVERGED_ITS;
243: else eps->reason = EPS_CONVERGED_TOL;
244: return(0);
245: }
249: PetscErrorCode EPSReset_BLOPEX(EPS eps)
250: {
252: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
255: mv_MultiVectorDestroy(blopex->eigenvectors);
256: mv_MultiVectorDestroy(blopex->Y);
257: VecDestroy(&blopex->w);
258: EPSReset_Default(eps);
259: return(0);
260: }
264: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
265: {
269: LOBPCG_DestroyRandomContext();
270: PetscFree(eps->data);
271: return(0);
272: }
276: PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps)
277: {
278: PetscErrorCode ierr;
281: PetscOptionsHead("EPS BLOPEX Options");
282: LOBPCG_SetFromOptionsRandomContext();
283: PetscOptionsTail();
284: return(0);
285: }
290: PetscErrorCode EPSCreate_BLOPEX(EPS eps)
291: {
295: PetscNewLog(eps,EPS_BLOPEX,&eps->data);
296: eps->ops->setup = EPSSetUp_BLOPEX;
297: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
298: eps->ops->destroy = EPSDestroy_BLOPEX;
299: eps->ops->reset = EPSReset_BLOPEX;
300: eps->ops->backtransform = EPSBackTransform_Default;
301: eps->ops->computevectors = EPSComputeVectors_Default;
302: STSetType(eps->OP,STPRECOND);
303: STPrecondSetKSPHasMat(eps->OP,PETSC_TRUE);
304: LOBPCG_InitRandomContext(((PetscObject)eps)->comm);
305: return(0);
306: }