Actual source code: lobpcg.c
2: /*
3: This file implements a wrapper to the LOBPCG solver from HYPRE package
4: */
5: #include src/eps/epsimpl.h
7: #include "HYPRE_lobpcg.h"
8: #include "HYPRE_parcsr_int.h"
10: #include "multivector.h"
11: /* #include "temp_multivector.h" */
12: typedef struct
13: {
14: long numVectors;
15: int* mask; void** vector;
16: int ownsVectors;
17: int ownsMask;
18:
19: HYPRE_InterfaceInterpreter* interpreter;
20:
21: } hypre_TempMultiVector;
25: EXTERN PetscErrorCode VecHYPRE_IJVectorCreate(Vec,HYPRE_IJVector*);
26: EXTERN PetscErrorCode VecHYPRE_IJVectorCopy(Vec v,HYPRE_IJVector ij);
27: EXTERN PetscErrorCode VecHYPRE_IJVectorCopyFrom(HYPRE_IJVector ij,Vec v);
29: typedef struct {
30: HYPRE_InterfaceInterpreter interpreter;
31: HYPRE_Solver solver;
32: int bsize;
33: hypre_MultiVectorPtr eigenvectors;
34: double *eigval;
35: KSP ksp;
36: Vec x,b;
37: HYPRE_IJVector ijx,ijb;
38: } EPS_LOBPCG;
40: /* Nasty global variable to access EPS data from FunctSolver and FunctA */
41: static EPS globaleps;
45: int EPSFunctSolver(HYPRE_Solver s, HYPRE_Matrix m, HYPRE_Vector x, HYPRE_Vector y)
46: { /* solves A*x=b */
48: PetscErrorCode ierr;
49: EPS_LOBPCG *lobpcg = (EPS_LOBPCG *)globaleps->data;
50: HYPRE_ParVector px,pb;
54: HYPRE_IJVectorGetObject(lobpcg->ijb,(void**)&pb);
55: HYPRE_ParVectorCopy((HYPRE_ParVector)x,pb);
56: VecHYPRE_IJVectorCopyFrom(lobpcg->ijb,lobpcg->b);
57:
58: KSPSolve(lobpcg->ksp,lobpcg->b,lobpcg->x);
59:
60: VecHYPRE_IJVectorCopy(lobpcg->x,lobpcg->ijx);
61: HYPRE_IJVectorGetObject(lobpcg->ijx,(void**)&px);
62: HYPRE_ParVectorCopy(px,(HYPRE_ParVector)y);
63:
64: return(0);
65: }
69: int EPSFunctA(void *matvec_data, double alpha, void *A, void *x, double beta, void *y)
70: {
71: /*
72: computes y <- alpha * A * x + beta * y
73: only works if alpha=1 and beta=0
74: */
75: PetscErrorCode ierr;
76: EPS_LOBPCG *lobpcg = (EPS_LOBPCG *)globaleps->data;
77: HYPRE_ParVector px,py;
80:
81: HYPRE_IJVectorGetObject(lobpcg->ijx,(void**)&px);
82: HYPRE_ParVectorCopy((HYPRE_ParVector)x,px);
83: VecHYPRE_IJVectorCopyFrom(lobpcg->ijx,lobpcg->x);
84:
85: STApply(globaleps->OP,lobpcg->x,lobpcg->b);
86:
87: VecHYPRE_IJVectorCopy(lobpcg->b,lobpcg->ijb);
88: HYPRE_IJVectorGetObject(lobpcg->ijb,(void**)&py);
89: HYPRE_ParVectorCopy(py,(HYPRE_ParVector)y);
90:
91: return(0);
92: }
96: PetscErrorCode EPSSetUp_LOBPCG(EPS eps)
97: {
98: PetscErrorCode ierr;
99: EPS_LOBPCG *lobpcg = (EPS_LOBPCG *)eps->data;
100: Mat A;
101: HYPRE_ParVector px;
102: int i,N;
103: PetscTruth isShift;
106: if (!eps->ishermitian || eps->isgeneralized) {
107: SETERRQ(PETSC_ERR_SUP,"LOBPCG only works for standard symmetric problems");
108: }
109: PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
110: if (!isShift) {
111: SETERRQ(PETSC_ERR_SUP,"LOBPCG only works with shift spectral transformation");
112: }
113: if (eps->which!=EPS_SMALLEST_REAL) {
114: SETERRQ(1,"Wrong value of eps->which");
115: }
116: STGetOperators(eps->OP,&A,PETSC_NULL);
117: MatGetVecs(A,&lobpcg->x,&lobpcg->b);
118: VecHYPRE_IJVectorCreate(lobpcg->x,&lobpcg->ijx);
119: VecHYPRE_IJVectorCreate(lobpcg->b,&lobpcg->ijb);
120: KSPSetOperators(lobpcg->ksp,A,A,DIFFERENT_NONZERO_PATTERN);
121: KSPSetUp(lobpcg->ksp);
122: VecGetSize(lobpcg->x,&N);
124: if (!eps->max_it) eps->max_it = PetscMax(100,N);
125: HYPRE_LOBPCGSetMaxIter(lobpcg->solver,eps->max_it);
126: if (!eps->tol) eps->tol = 1.e-7;
127: HYPRE_LOBPCGSetTol(lobpcg->solver,eps->tol);
128: HYPRE_LOBPCGSetPrintLevel(lobpcg->solver,0);
130: HYPRE_LOBPCGSetPrecond(lobpcg->solver,EPSFunctSolver,PETSC_NULL,PETSC_NULL);
131: HYPRE_LOBPCGSetup(lobpcg->solver,PETSC_NULL,PETSC_NULL,PETSC_NULL);
133: lobpcg->bsize = eps->ncv = eps->nev;
134: PetscMalloc(lobpcg->bsize*sizeof(PetscScalar),&lobpcg->eigval);
135: PetscMemzero(lobpcg->eigval,lobpcg->bsize*sizeof(PetscScalar));
137: HYPRE_IJVectorGetObject(lobpcg->ijx,(void**)&px);
138: lobpcg->eigenvectors = hypre_MultiVectorCreateFromSampleVector(&lobpcg->interpreter,lobpcg->bsize,px);
139: for (i=0;i<lobpcg->bsize;i++) {
140: if (i==0) {
141: VecHYPRE_IJVectorCopy(eps->vec_initial,lobpcg->ijx);
142: } else {
143: SlepcVecSetRandom(lobpcg->x);
144: VecHYPRE_IJVectorCopy(lobpcg->x,lobpcg->ijx);
145: }
146: HYPRE_ParVectorCopy(px,((hypre_TempMultiVector *) ((hypre_MultiVector *) lobpcg->eigenvectors)->data)->vector[i]);
147: }
148:
149: EPSAllocateSolution(eps);
150: return(0);
151: }
155: PetscErrorCode EPSSolve_LOBPCG(EPS eps)
156: {
157: PetscErrorCode ierr;
158: EPS_LOBPCG *lobpcg = (EPS_LOBPCG *)eps->data;
159: int i;
160: HYPRE_ParVector px;
163: globaleps = eps;
164:
165: HYPRE_LOBPCGSolve(lobpcg->solver,PETSC_NULL,lobpcg->eigenvectors,lobpcg->eigval);
166: HYPRE_IJVectorGetObject(lobpcg->ijx,(void**)&px);
167: for (i=0;i<lobpcg->bsize;i++) {
168: eps->eigr[i] = lobpcg->eigval[i];
169: eps->eigi[i] = 0;
170: HYPRE_ParVectorCopy(((hypre_TempMultiVector *) ((hypre_MultiVector *) lobpcg->eigenvectors)->data)->vector[i],px);
171: VecHYPRE_IJVectorCopyFrom(lobpcg->ijx,eps->V[i]);
172: }
173: eps->nconv = lobpcg->bsize;
174: eps->reason = EPS_CONVERGED_TOL;
175: eps->its = HYPRE_LOBPCGIterations(lobpcg->solver);
176:
177: return(0);
178: }
182: PetscErrorCode EPSSetFromOptions_LOBPCG(EPS eps)
183: {
185: EPS_LOBPCG *lobpcg = (EPS_LOBPCG *)eps->data;
188: KSPSetFromOptions(lobpcg->ksp);
189: return(0);
190: }
194: PetscErrorCode EPSDestroy_LOBPCG(EPS eps)
195: {
197: EPS_LOBPCG *lobpcg = (EPS_LOBPCG *)eps->data;
200: HYPRE_LOBPCGDestroy(lobpcg->solver);
201: hypre_MultiVectorDestroy(lobpcg->eigenvectors);
202: PetscFree(lobpcg->eigval);
203: HYPRE_IJVectorDestroy(lobpcg->ijx);
204: HYPRE_IJVectorDestroy(lobpcg->ijb);
205:
206: KSPDestroy(lobpcg->ksp);
207: VecDestroy(lobpcg->x);
208: VecDestroy(lobpcg->b);
209: PetscFree(eps->data);
210: EPSFreeSolution(eps);
211: return(0);
212: }
217: PetscErrorCode EPSCreate_LOBPCG(EPS eps)
218: {
220: EPS_LOBPCG *lobpcg;
221: const char* prefix;
224: PetscNew(EPS_LOBPCG,&lobpcg);
225: PetscLogObjectMemory(eps,sizeof(EPS_LOBPCG));
226: HYPRE_ParCSRSetupInterpreter(&lobpcg->interpreter);
227: lobpcg->interpreter.Matvec = EPSFunctA;
228: HYPRE_LOBPCGCreate(&lobpcg->interpreter,&lobpcg->solver);
229: KSPCreate(eps->comm,&lobpcg->ksp);
230: EPSGetOptionsPrefix(eps,&prefix);
231: KSPSetOptionsPrefix(lobpcg->ksp,prefix);
232: KSPAppendOptionsPrefix(lobpcg->ksp,"eps_lobpcg_");
233: eps->data = (void *) lobpcg;
234: eps->ops->solve = EPSSolve_LOBPCG;
235: eps->ops->setup = EPSSetUp_LOBPCG;
236: eps->ops->setfromoptions = EPSSetFromOptions_LOBPCG;
237: eps->ops->destroy = EPSDestroy_LOBPCG;
238: eps->ops->backtransform = EPSBackTransform_Default;
239: eps->ops->computevectors = EPSComputeVectors_Default;
240: eps->which = EPS_SMALLEST_REAL;
241: return(0);
242: }