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-2009, 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: typedef struct {
 33:   lobpcg_Tolerance           tol;
 34:   lobpcg_BLASLAPACKFunctions blap_fn;
 35:   mv_MultiVectorPtr          eigenvectors;
 36:   mv_MultiVectorPtr          Y;
 37:   mv_InterfaceInterpreter    ii;
 38:   KSP                        ksp;
 39: } EPS_BLOPEX;


 44: static void Precond_FnSingleVector(void *data,void *x,void *y)
 45: {
 46:   PetscErrorCode  ierr;
 47:   EPS             eps = (EPS)data;
 48:   EPS_BLOPEX      *blopex = (EPS_BLOPEX*)eps->data;
 49:   PetscInt        lits;
 50: 
 52:   KSPSolve(blopex->ksp,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
 53:   KSPGetIterationNumber(blopex->ksp, &lits); CHKERRABORT(PETSC_COMM_WORLD,ierr);
 54:   eps->OP->lineariterations+= lits;

 56:   PetscFunctionReturnVoid();
 57: }

 61: static void Precond_FnMultiVector(void *data,void *x,void *y)
 62: {
 63:   EPS             eps = (EPS)data;
 64:   EPS_BLOPEX      *blopex = (EPS_BLOPEX*)eps->data;

 67:   blopex->ii.Eval(Precond_FnSingleVector,data,x,y);
 68:   PetscFunctionReturnVoid();
 69: }

 73: static void OperatorASingleVector(void *data,void *x,void *y)
 74: {
 75:   PetscErrorCode  ierr;
 76:   EPS                  eps = (EPS)data;
 77:   Mat             A;
 78: 
 80:   STGetOperators(eps->OP,&A,PETSC_NULL); CHKERRABORT(PETSC_COMM_WORLD,ierr);
 81:   MatMult(A,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
 82:   PetscFunctionReturnVoid();
 83: }

 87: static void OperatorAMultiVector(void *data,void *x,void *y)
 88: {
 89:   EPS                  eps = (EPS)data;
 90:   EPS_BLOPEX          *blopex = (EPS_BLOPEX*)eps->data;

 93:   blopex->ii.Eval(OperatorASingleVector,data,x,y);
 94:   PetscFunctionReturnVoid();
 95: }

 99: static void OperatorBSingleVector(void *data,void *x,void *y)
100: {
101:   PetscErrorCode  ierr;
102:   EPS                  eps = (EPS)data;
103:   Mat             B;
104: 
106:   STGetOperators(eps->OP,PETSC_NULL,&B); CHKERRABORT(PETSC_COMM_WORLD,ierr);
107:   MatMult(B,(Vec)x,(Vec)y); CHKERRABORT(PETSC_COMM_WORLD,ierr);
108:   PetscFunctionReturnVoid();
109: }

113: static void OperatorBMultiVector(void *data,void *x,void *y)
114: {
115:   EPS                  eps = (EPS)data;
116:   EPS_BLOPEX          *blopex = (EPS_BLOPEX*)eps->data;

119:   blopex->ii.Eval(OperatorBSingleVector,data,x,y);
120:   PetscFunctionReturnVoid();
121: }


126: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
127: {
128:   PetscErrorCode  ierr;
129:   EPS_BLOPEX      *blopex = (EPS_BLOPEX *)eps->data;
130:   Mat             A,B;
131:   PetscInt        N;
132:   PetscTruth      isShift;

135:   if (!eps->ishermitian) {
136:     SETERRQ(PETSC_ERR_SUP,"blopex only works for hermitian problems");
137:   }
138:   PetscTypeCompare((PetscObject)eps->OP,STSHIFT,&isShift);
139:   if (!isShift) {
140:     SETERRQ(PETSC_ERR_SUP,"blopex only works with shift spectral transformation");
141:   }
142:   if (eps->which!=EPS_SMALLEST_REAL) {
143:     SETERRQ(1,"Wrong value of eps->which");
144:   }

146:   KSPSetType(blopex->ksp, KSPPREONLY);
147:   KSPSetFromOptions(blopex->ksp);
148:   STGetOperators(eps->OP,&A,&B);
149:   KSPSetOperators(blopex->ksp,A,A,DIFFERENT_NONZERO_PATTERN);

151:   VecGetSize(eps->vec_initial,&N);
152:   eps->ncv = eps->nev = PetscMin(eps->nev,N);
153:   if (eps->mpd) PetscInfo(eps,"Warning: parameter mpd ignored\n");
154:   if (!eps->max_it) eps->max_it = PetscMax(100,2*N/eps->ncv);

156:   EPSAllocateSolution(eps);
157:   EPSDefaultGetWork(eps,1);
158: 
159:   blopex->tol.absolute = eps->tol;
160:   blopex->tol.relative = 1e-50;
161: 
162:   LOBPCG_InitRandomContext();
163:   SLEPCSetupInterpreter(&blopex->ii);
164:   blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->V);
165:   mv_MultiVectorSetRandom(blopex->eigenvectors,1234);

167:   if (eps->nds > 0) {
168:     blopex->Y = mv_MultiVectorCreateFromSampleVector(&blopex->ii, eps->nds, eps->DS);
169:   } else
170:     blopex->Y = PETSC_NULL;
171: 
172:   blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
173:   blopex->blap_fn.dsygv = PETSC_dsygv_interface;

175:   if (eps->extraction) {
176:      PetscInfo(eps,"Warning: extraction type ignored\n");
177:   }

179:   return(0);
180: }

184: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
185: {
186:   EPS_BLOPEX      *blopex = (EPS_BLOPEX *)eps->data;
187:   int             info,its;
188: 
190: 
191:   info = lobpcg_solve(blopex->eigenvectors,eps,OperatorAMultiVector,
192:           eps->isgeneralized?eps:PETSC_NULL,eps->isgeneralized?OperatorBMultiVector:PETSC_NULL,
193:         eps,Precond_FnMultiVector,blopex->Y,
194:         blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
195:         eps->eigr,PETSC_NULL,0,eps->errest,PETSC_NULL,0);
196:   if (info>0) SETERRQ1(PETSC_ERR_LIB,"Error in blopex (code=%d)",info);

198:   eps->its = its;
199:   eps->nconv = eps->ncv;
200:   if (info==-1) eps->reason = EPS_DIVERGED_ITS;
201:   else eps->reason = EPS_CONVERGED_TOL;
202: 
203:   return(0);
204: }

208: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
209: {
211:   EPS_BLOPEX     *blopex = (EPS_BLOPEX *)eps->data;

214:   KSPDestroy(blopex->ksp);
215:   LOBPCG_DestroyRandomContext();
216:   SLEPCSetupInterpreterForDignifiedDeath(&blopex->ii);
217:   mv_MultiVectorDestroy(blopex->eigenvectors);
218:   mv_MultiVectorDestroy(blopex->Y);
219:   PetscFree(eps->data);
220:   EPSDestroy_Default(eps);
221:   return(0);
222: }

227: PetscErrorCode EPSCreate_BLOPEX(EPS eps)
228: {
230:   EPS_BLOPEX     *blopex;
231:   const char*    prefix;

234:   PetscNew(EPS_BLOPEX,&blopex);
235:   PetscLogObjectMemory(eps,sizeof(EPS_BLOPEX));
236:   KSPCreate(((PetscObject)eps)->comm,&blopex->ksp);
237:   EPSGetOptionsPrefix(eps,&prefix);
238:   KSPSetOptionsPrefix(blopex->ksp,prefix);
239:   KSPAppendOptionsPrefix(blopex->ksp,"eps_blopex_");
240:   eps->data                      = (void *) blopex;
241:   eps->ops->solve                = EPSSolve_BLOPEX;
242:   eps->ops->setup                = EPSSetUp_BLOPEX;
243:   eps->ops->setfromoptions       = PETSC_NULL;
244:   eps->ops->destroy              = EPSDestroy_BLOPEX;
245:   eps->ops->backtransform        = EPSBackTransform_Default;
246:   eps->ops->computevectors       = EPSComputeVectors_Default;
247:   eps->which = EPS_SMALLEST_REAL;
248:   return(0);
249: }