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-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/eps/epsimpl.h
14: #include "src/contrib/blopex/petsc-interface/petsc-interface.h"
15: #include "lobpcg.h"
16: #include "interpreter.h"
17: #include "multivector.h"
18: #include "temp_multivector.h"
20: typedef struct {
21: lobpcg_Tolerance tol;
22: lobpcg_BLASLAPACKFunctions blap_fn;
23: mv_MultiVectorPtr eigenvectors;
24: mv_InterfaceInterpreter ii;
25: KSP ksp;
26: } EPS_BLOPEX;
30: static void Precond_FnSingleVector(void *data,void *x,void *y)
31: {
32: PetscErrorCode ierr;
33: EPS eps = (EPS)data;
34: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
35:
37: KSPSolve(blopex->ksp,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
38: PetscFunctionReturnVoid();
39: }
43: static void Precond_FnMultiVector(void *data,void *x,void *y)
44: {
45: EPS eps = (EPS)data;
46: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
49: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
50: PetscFunctionReturnVoid();
51: }
55: static void OperatorASingleVector(void *data,void *x,void *y)
56: {
57: PetscErrorCode ierr;
58: EPS eps = (EPS)data;
59:
61: STApply(eps->OP,(Vec)x,(Vec)y);
62: IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,(Vec)y,PETSC_NULL,PETSC_NULL,PETSC_NULL,eps->work[0]); CHKERRABORT(PETSC_COMM_WORLD,ierr);
63: PetscFunctionReturnVoid();
64: }
68: static void OperatorAMultiVector(void *data,void *x,void *y)
69: {
70: EPS eps = (EPS)data;
71: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
74: blopex->ii.Eval(OperatorASingleVector,data,x,y);
75: PetscFunctionReturnVoid();
76: }
80: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
81: {
82: PetscErrorCode ierr;
83: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
84: Mat A;
85: int N;
86: PetscTruth isShift;
89: if (!eps->ishermitian || eps->isgeneralized) {
90: SETERRQ(PETSC_ERR_SUP,"blopex only works for standard symmetric problems");
91: }
92: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
93: if (!isShift) {
94: SETERRQ(PETSC_ERR_SUP,"blopex only works with shift spectral transformation");
95: }
96: if (eps->which!=EPS_SMALLEST_REAL) {
97: SETERRQ(1,"Wrong value of eps->which");
98: }
99: STGetOperators(eps->OP,&A,PETSC_NULL);
100: KSPSetOperators(blopex->ksp,A,A,DIFFERENT_NONZERO_PATTERN);
101: KSPSetUp(blopex->ksp);
103: VecGetSize(eps->vec_initial,&N);
104: eps->ncv = eps->nev = PetscMin(eps->nev,N);
105: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
106:
107: blopex->tol.absolute = eps->tol;
108: blopex->tol.relative = 1e-50;
109:
110: LOBPCG_InitRandomContext();
111: PETSCSetupInterpreter(&blopex->ii);
112: blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->vec_initial);
113: mv_MultiVectorSetRandom(blopex->eigenvectors,1234);
114:
115: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
116: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
118: EPSAllocateSolution(eps);
119: EPSDefaultGetWork(eps,1);
120: return(0);
121: }
125: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
126: {
127: PetscErrorCode ierr;
128: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
129: int info,i;
130: mv_TempMultiVector *mv;
131:
133: info = lobpcg_solve(blopex->eigenvectors,eps,OperatorAMultiVector,
134: PETSC_NULL,PETSC_NULL,eps,Precond_FnMultiVector,NULL,
135: blopex->blap_fn,blopex->tol,eps->max_it,0,&eps->its,
136: eps->eigr,PETSC_NULL,0,eps->errest,PETSC_NULL,0);
137: if (info>0) SETERRQ1(PETSC_ERR_LIB,"Error in blopex (code=%d)",info);
139: eps->nconv = eps->ncv;
140: if (info==-1) eps->reason = EPS_DIVERGED_ITS;
141: else eps->reason = EPS_CONVERGED_TOL;
142:
143: mv = (mv_TempMultiVector*)mv_MultiVectorGetData(blopex->eigenvectors);
144: for (i=0;i<eps->nconv;i++) {
145: VecCopy((Vec)mv->vector[i],eps->V[i]);
146: }
147: return(0);
148: }
152: PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps)
153: {
155: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
158: KSPSetFromOptions(blopex->ksp);
159: return(0);
160: }
164: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
165: {
167: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
170: KSPDestroy(blopex->ksp);
171: LOBPCG_DestroyRandomContext();
172: mv_MultiVectorDestroy(blopex->eigenvectors);
173: PetscFree(eps->data);
174: EPSFreeSolution(eps);
175: return(0);
176: }
181: PetscErrorCode EPSCreate_BLOPEX(EPS eps)
182: {
184: EPS_BLOPEX *blopex;
185: const char* prefix;
188: PetscNew(EPS_BLOPEX,&blopex);
189: PetscLogObjectMemory(eps,sizeof(EPS_BLOPEX));
190: KSPCreate(eps->comm,&blopex->ksp);
191: EPSGetOptionsPrefix(eps,&prefix);
192: KSPSetOptionsPrefix(blopex->ksp,prefix);
193: KSPAppendOptionsPrefix(blopex->ksp,"eps_blopex_");
194: eps->data = (void *) blopex;
195: eps->ops->solve = EPSSolve_BLOPEX;
196: eps->ops->setup = EPSSetUp_BLOPEX;
197: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
198: eps->ops->destroy = EPSDestroy_BLOPEX;
199: eps->ops->backtransform = EPSBackTransform_Default;
200: eps->ops->computevectors = EPSComputeVectors_Default;
201: eps->which = EPS_SMALLEST_REAL;
202: return(0);
203: }