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-2010, Universidad 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
25: #include private/epsimpl.h
26: #include "slepc-interface.h"
27: #include "lobpcg.h"
28: #include "interpreter.h"
29: #include "multivector.h"
30: #include "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: } EPS_BLOPEX;
46: static void Precond_FnSingleVector(void *data,void *x,void *y)
47: {
48: PetscErrorCode ierr;
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;
58: PetscFunctionReturnVoid();
59: }
63: static void Precond_FnMultiVector(void *data,void *x,void *y)
64: {
65: EPS eps = (EPS)data;
66: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
69: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
70: PetscFunctionReturnVoid();
71: }
75: static void OperatorASingleVector(void *data,void *x,void *y)
76: {
77: PetscErrorCode ierr;
78: EPS eps = (EPS)data;
79: Mat A;
80:
82: STGetOperators(eps->OP,&A,PETSC_NULL); CHKERRABORT(PETSC_COMM_WORLD,ierr);
83: MatMult(A,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
84: PetscFunctionReturnVoid();
85: }
89: static void OperatorAMultiVector(void *data,void *x,void *y)
90: {
91: EPS eps = (EPS)data;
92: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
95: blopex->ii.Eval(OperatorASingleVector,data,x,y);
96: PetscFunctionReturnVoid();
97: }
101: static void OperatorBSingleVector(void *data,void *x,void *y)
102: {
103: PetscErrorCode ierr;
104: EPS eps = (EPS)data;
105: Mat B;
106:
108: STGetOperators(eps->OP,PETSC_NULL,&B); CHKERRABORT(PETSC_COMM_WORLD,ierr);
109: MatMult(B,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
110: PetscFunctionReturnVoid();
111: }
115: static void OperatorBMultiVector(void *data,void *x,void *y)
116: {
117: EPS eps = (EPS)data;
118: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
121: blopex->ii.Eval(OperatorBSingleVector,data,x,y);
122: PetscFunctionReturnVoid();
123: }
128: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
129: {
130: PetscErrorCode ierr;
131: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
132: PetscTruth isPrecond, isPreonly;
135: if (!eps->ishermitian) {
136: SETERRQ(PETSC_ERR_SUP,"blopex only works for hermitian problems");
137: }
138: if (!eps->which) eps->which = EPS_SMALLEST_REAL;
139: if (eps->which!=EPS_SMALLEST_REAL) {
140: SETERRQ(1,"Wrong value of eps->which");
141: }
143: /* Change the default sigma to inf if necessary */
144: if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL ||
145: eps->which == EPS_LARGEST_IMAGINARY) {
146: STSetDefaultShift(eps->OP, 3e300);
147: }
149: STSetUp(eps->OP);
150: PetscTypeCompare((PetscObject)eps->OP, STPRECOND, &isPrecond);
151: if (!isPrecond) SETERRQ(PETSC_ERR_SUP, "blopex only works with precond spectral transformation");
152: STGetKSP(eps->OP, &blopex->ksp);
153: PetscTypeCompare((PetscObject)blopex->ksp, KSPPREONLY, &isPreonly);
154: if (!isPreonly)
155: SETERRQ(PETSC_ERR_SUP, "blopex only works with preonly ksp of the spectral transformation");
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);
161: EPSAllocateSolution(eps);
162: EPSDefaultGetWork(eps,1);
163:
164: blopex->tol.absolute = eps->tol;
165: blopex->tol.relative = 1e-50;
166:
167: LOBPCG_InitRandomContext();
168: SLEPCSetupInterpreter(&blopex->ii);
169: blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->V);
170: mv_MultiVectorSetRandom(blopex->eigenvectors,1234);
172: if (eps->nds > 0) {
173: blopex->Y = mv_MultiVectorCreateFromSampleVector(&blopex->ii, eps->nds, eps->DS);
174: } else
175: blopex->Y = PETSC_NULL;
177: #ifdef PETSC_USE_COMPLEX
178: blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
179: blopex->blap_fn.zhegv = PETSC_zsygv_interface;
180: #else
181: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
182: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
183: #endif
185: if (eps->extraction) {
186: PetscInfo(eps,"Warning: extraction type ignored\n");
187: }
189: /* dispatch solve method */
190: if (eps->leftvecs) SETERRQ(PETSC_ERR_SUP,"Left vectors not supported in this solver");
191: eps->ops->solve = EPSSolve_BLOPEX;
192: return(0);
193: }
197: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
198: {
199: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
200: int info,its;
201:
203:
204: #ifdef PETSC_USE_COMPLEX
205: info = lobpcg_solve_complex(blopex->eigenvectors,eps,OperatorAMultiVector,
206: eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
207: eps,Precond_FnMultiVector,blopex->Y,
208: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
209: (komplex*)eps->eigr,PETSC_NULL,0,eps->errest,PETSC_NULL,0);
210: #else
211: info = lobpcg_solve_double(blopex->eigenvectors,eps,OperatorAMultiVector,
212: eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
213: eps,Precond_FnMultiVector,blopex->Y,
214: blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
215: eps->eigr,PETSC_NULL,0,eps->errest,PETSC_NULL,0);
216: #endif
217: if (info>0) SETERRQ1(PETSC_ERR_LIB,"Error in blopex (code=%d)",info);
219: eps->its = its;
220: eps->nconv = eps->ncv;
221: if (info==-1) eps->reason = EPS_DIVERGED_ITS;
222: else eps->reason = EPS_CONVERGED_TOL;
223:
224: return(0);
225: }
229: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
230: {
232: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
235: LOBPCG_DestroyRandomContext();
236: SLEPCSetupInterpreterForDignifiedDeath(&blopex->ii);
237: mv_MultiVectorDestroy(blopex->eigenvectors);
238: mv_MultiVectorDestroy(blopex->Y);
239: PetscFree(eps->data);
240: EPSDestroy_Default(eps);
241: return(0);
242: }
247: PetscErrorCode EPSCreate_BLOPEX(EPS eps)
248: {
250: EPS_BLOPEX *blopex;
251: const char* prefix;
255: STSetType(eps->OP, STPRECOND);
256: STPrecondSetKSPHasMat(eps->OP, PETSC_TRUE);
258: PetscNew(EPS_BLOPEX,&blopex);
259: PetscLogObjectMemory(eps,sizeof(EPS_BLOPEX));
260: KSPCreate(((PetscObject)eps)->comm,&blopex->ksp);
261: EPSGetOptionsPrefix(eps,&prefix);
262: KSPSetOptionsPrefix(blopex->ksp,prefix);
263: KSPAppendOptionsPrefix(blopex->ksp,"eps_blopex_");
264: eps->data = (void *) blopex;
265: eps->ops->setup = EPSSetUp_BLOPEX;
266: eps->ops->setfromoptions = PETSC_NULL;
267: eps->ops->destroy = EPSDestroy_BLOPEX;
268: eps->ops->backtransform = EPSBackTransform_Default;
269: eps->ops->computevectors = EPSComputeVectors_Default;
270: return(0);
271: }