Actual source code: blopex.c
2: /*
3: This file implements a wrapper to the BLOPEX solver
4: */
5: #include src/eps/epsimpl.h
6: #include "src/contrib/blopex/petsc-interface/petsc-interface.h"
7: #include "lobpcg.h"
8: #include "interpreter.h"
9: #include "multivector.h"
10: #include "temp_multivector.h"
12: typedef struct {
13: lobpcg_Tolerance tol;
14: lobpcg_BLASLAPACKFunctions blap_fn;
15: mv_MultiVectorPtr eigenvectors;
16: mv_InterfaceInterpreter ii;
17: KSP ksp;
18: } EPS_BLOPEX;
22: static void Precond_FnSingleVector(void *data,void *x,void *y)
23: {
24: PetscErrorCode ierr;
25: EPS eps = (EPS)data;
26: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
27:
29: KSPSolve(blopex->ksp,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
30: PetscFunctionReturnVoid();
31: }
35: static void Precond_FnMultiVector(void *data,void *x,void *y)
36: {
37: EPS eps = (EPS)data;
38: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
41: blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
42: PetscFunctionReturnVoid();
43: }
47: static void OperatorASingleVector(void *data,void *x,void *y)
48: {
49: PetscErrorCode ierr;
50: EPS eps = (EPS)data;
51:
53: STApply(eps->OP,(Vec)x,(Vec)y);
54: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,(Vec)y,PETSC_NULL,PETSC_NULL,PETSC_NULL); CHKERRABORT(PETSC_COMM_WORLD,ierr);
55: PetscFunctionReturnVoid();
56: }
60: static void OperatorAMultiVector(void *data,void *x,void *y)
61: {
62: EPS eps = (EPS)data;
63: EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;
66: blopex->ii.Eval(OperatorASingleVector,data,x,y);
67: PetscFunctionReturnVoid();
68: }
72: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
73: {
74: PetscErrorCode ierr;
75: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
76: Mat A;
77: int N;
78: PetscTruth isShift;
81: if (!eps->ishermitian || eps->isgeneralized) {
82: SETERRQ(PETSC_ERR_SUP,"blopex only works for standard symmetric problems");
83: }
84: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
85: if (!isShift) {
86: SETERRQ(PETSC_ERR_SUP,"blopex only works with shift spectral transformation");
87: }
88: if (eps->which!=EPS_SMALLEST_REAL) {
89: SETERRQ(1,"Wrong value of eps->which");
90: }
91: STGetOperators(eps->OP,&A,PETSC_NULL);
92: KSPSetOperators(blopex->ksp,A,A,DIFFERENT_NONZERO_PATTERN);
93: KSPSetUp(blopex->ksp);
95: VecGetSize(eps->vec_initial,&N);
96: eps->ncv = eps->nev = PetscMin(eps->nev,N);
97: if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);
98:
99: blopex->tol.absolute = eps->tol;
100: blopex->tol.relative = 1e-50;
101:
102: LOBPCG_InitRandomContext();
103: PETSCSetupInterpreter(&blopex->ii);
104: blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->vec_initial);
105: mv_MultiVectorSetRandom(blopex->eigenvectors,1234);
106:
107: blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
108: blopex->blap_fn.dsygv = PETSC_dsygv_interface;
110: EPSAllocateSolution(eps);
111: return(0);
112: }
116: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
117: {
118: PetscErrorCode ierr;
119: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
120: int info,i;
121: mv_TempMultiVector *mv;
122:
124: info = lobpcg_solve(blopex->eigenvectors,eps,OperatorAMultiVector,
125: PETSC_NULL,PETSC_NULL,eps,Precond_FnMultiVector,NULL,
126: blopex->blap_fn,blopex->tol,eps->max_it,0,&eps->its,
127: eps->eigr,PETSC_NULL,0,eps->errest,PETSC_NULL,0);
128: if (info>0) SETERRQ1(PETSC_ERR_LIB,"Error in blopex (code=%d)",info);
130: eps->nconv = eps->ncv;
131: if (info==-1) eps->reason = EPS_DIVERGED_ITS;
132: else eps->reason = EPS_CONVERGED_TOL;
133:
134: mv = (mv_TempMultiVector*)mv_MultiVectorGetData(blopex->eigenvectors);
135: for (i=0;i<eps->nconv;i++) {
136: VecCopy((Vec)mv->vector[i],eps->V[i]);
137: }
138: return(0);
139: }
143: PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps)
144: {
146: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
149: KSPSetFromOptions(blopex->ksp);
150: return(0);
151: }
155: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
156: {
158: EPS_BLOPEX *blopex = (EPS_BLOPEX *)eps->data;
161: KSPDestroy(blopex->ksp);
162: LOBPCG_DestroyRandomContext();
163: mv_MultiVectorDestroy(blopex->eigenvectors);
164: PetscFree(eps->data);
165: EPSFreeSolution(eps);
166: return(0);
167: }
172: PetscErrorCode EPSCreate_BLOPEX(EPS eps)
173: {
175: EPS_BLOPEX *blopex;
176: const char* prefix;
179: PetscNew(EPS_BLOPEX,&blopex);
180: PetscLogObjectMemory(eps,sizeof(EPS_BLOPEX));
181: KSPCreate(eps->comm,&blopex->ksp);
182: EPSGetOptionsPrefix(eps,&prefix);
183: KSPSetOptionsPrefix(blopex->ksp,prefix);
184: KSPAppendOptionsPrefix(blopex->ksp,"eps_blopex_");
185: eps->data = (void *) blopex;
186: eps->ops->solve = EPSSolve_BLOPEX;
187: eps->ops->setup = EPSSetUp_BLOPEX;
188: eps->ops->setfromoptions = EPSSetFromOptions_BLOPEX;
189: eps->ops->destroy = EPSDestroy_BLOPEX;
190: eps->ops->backtransform = EPSBackTransform_Default;
191: eps->ops->computevectors = EPSComputeVectors_Default;
192: eps->which = EPS_SMALLEST_REAL;
193: return(0);
194: }