Actual source code: blopex.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    This file implements a wrapper to the BLOPEX package

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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 <slepc-private/epsimpl.h>
 25: #include <slepc-private/stimpl.h>
 26: #include "slepc-interface.h"
 27: #include <blopex_lobpcg.h>
 28: #include <blopex_interpreter.h>
 29: #include <blopex_multivector.h>
 30: #include <blopex_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:   ST                         st;
 41:   Vec                        w;
 42: } EPS_BLOPEX;

 46: static void Precond_FnSingleVector(void *data,void *x,void *y)
 47: {
 49:   EPS            eps = (EPS)data;
 50:   EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;

 53:   KSPSolve(blopex->st->ksp,(Vec)x,(Vec)y);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 54:   PetscFunctionReturnVoid();
 55: }

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

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

 71: static void OperatorASingleVector(void *data,void *x,void *y)
 72: {
 74:   EPS            eps = (EPS)data;
 75:   EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;
 76:   Mat            A,B;
 77:   PetscScalar    sigma;
 78:   PetscInt       nmat;

 81:   STGetNumMatrices(eps->st,&nmat);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 82:   STGetOperators(eps->st,0,&A);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 83:   if (nmat>1) { STGetOperators(eps->st,1,&B);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr); }
 84:   MatMult(A,(Vec)x,(Vec)y);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 85:   STGetShift(eps->st,&sigma);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 86:   if (sigma != 0.0) {
 87:     if (nmat>1) {
 88:       MatMult(B,(Vec)x,blopex->w);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 89:     } else {
 90:       VecCopy((Vec)x,blopex->w);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 91:     }
 92:     VecAXPY((Vec)y,-sigma,blopex->w);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
 93:   }
 94:   PetscFunctionReturnVoid();
 95: }

 99: static void OperatorAMultiVector(void *data,void *x,void *y)
100: {
101:   EPS        eps = (EPS)data;
102:   EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;

105:   blopex->ii.Eval(OperatorASingleVector,data,x,y);
106:   PetscFunctionReturnVoid();
107: }

111: static void OperatorBSingleVector(void *data,void *x,void *y)
112: {
114:   EPS            eps = (EPS)data;
115:   Mat            B;

118:   STGetOperators(eps->st,1,&B);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
119:   MatMult(B,(Vec)x,(Vec)y);CHKERRABORT(PetscObjectComm((PetscObject)eps),ierr);
120:   PetscFunctionReturnVoid();
121: }

125: static void OperatorBMultiVector(void *data,void *x,void *y)
126: {
127:   EPS        eps = (EPS)data;
128:   EPS_BLOPEX *blopex = (EPS_BLOPEX*)eps->data;

131:   blopex->ii.Eval(OperatorBSingleVector,data,x,y);
132:   PetscFunctionReturnVoid();
133: }

137: PetscErrorCode EPSSetUp_BLOPEX(EPS eps)
138: {
139: #if defined(PETSC_MISSING_LAPACK_POTRF) || defined(PETSC_MISSING_LAPACK_SYGV)
141:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF/SYGV - Lapack routine is unavailable");
142: #else
144:   EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;
145:   PetscBool      isPrecond,istrivial,flg;
146:   BV             Y;
147:   PetscInt       k;

150:   if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"blopex only works for hermitian problems");
151:   if (!eps->which) eps->which = EPS_SMALLEST_REAL;
152:   if (eps->which!=EPS_SMALLEST_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");

154:   STSetUp(eps->st);
155:   PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&isPrecond);
156:   if (!isPrecond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"blopex only works with STPRECOND");
157:   blopex->st = eps->st;

159:   eps->ncv = eps->nev = PetscMin(eps->nev,eps->n);
160:   if (eps->mpd) { PetscInfo(eps,"Warning: parameter mpd ignored\n"); }
161:   if (!eps->max_it) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
162:   if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not supported in this solver");

164:   /* blopex only works with BVVECS or BVCONTIGUOUS, if different set to CONTIGUOUS */
165:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
166:   PetscObjectTypeCompareAny((PetscObject)eps->V,&flg,BVVECS,BVCONTIGUOUS,"");
167:   if (!flg) {
168:     BVSetType(eps->V,BVCONTIGUOUS);
169:   }

171:   EPSAllocateSolution(eps,0);
172:   EPSSetWorkVecs(eps,1);

174:   if (eps->converged == EPSConvergedEigRelative) {
175:     blopex->tol.absolute = 0.0;
176:     blopex->tol.relative = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
177:   } else if (eps->converged == EPSConvergedAbsolute) {
178:     blopex->tol.absolute = eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol;
179:     blopex->tol.relative = 0.0;
180:   } else {
181:     SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Convergence test not supported in this solver");
182:   }

184:   SLEPCSetupInterpreter(&blopex->ii);
185:   blopex->eigenvectors = mv_MultiVectorCreateFromSampleVector(&blopex->ii,eps->ncv,eps->V);

187:   BVGetVec(eps->V,&blopex->w);
188:   PetscLogObjectParent((PetscObject)eps,(PetscObject)blopex->w);
189:   if (eps->nds<0) {
190:     k = -eps->nds;
191:     BVCreate(PetscObjectComm((PetscObject)eps),&Y);
192:     BVSetSizesFromVec(Y,blopex->w,k);
193:     BVSetType(Y,BVVECS);
194:     BVInsertVecs(Y,0,&k,eps->defl,PETSC_FALSE);
195:     SlepcBasisDestroy_Private(&eps->nds,&eps->defl);
196:     blopex->Y = mv_MultiVectorCreateFromSampleVector(&blopex->ii,k,Y);
197:     BVDestroy(&Y);
198:   } else blopex->Y = NULL;

200: #if defined(PETSC_USE_COMPLEX)
201:   blopex->blap_fn.zpotrf = PETSC_zpotrf_interface;
202:   blopex->blap_fn.zhegv = PETSC_zsygv_interface;
203: #else
204:   blopex->blap_fn.dpotrf = PETSC_dpotrf_interface;
205:   blopex->blap_fn.dsygv = PETSC_dsygv_interface;
206: #endif

208:   if (eps->extraction) { PetscInfo(eps,"Warning: extraction type ignored\n"); }
209:   RGIsTrivial(eps->rg,&istrivial);
210:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support region filtering");

212:   /* dispatch solve method */
213:   eps->ops->solve = EPSSolve_BLOPEX;
214:   return(0);
215: #endif
216: }

220: PetscErrorCode EPSSolve_BLOPEX(EPS eps)
221: {
222:   EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;
223:   PetscScalar    sigma;
224:   int            i,j,info,its,nconv;
225:   double         *residhist=NULL;
227: #if defined(PETSC_USE_COMPLEX)
228:   komplex        *lambdahist=NULL;
229: #else
230:   double         *lambdahist=NULL;
231: #endif

234:   /* Complete the initial basis with random vectors */
235:   for (i=eps->nini;i<eps->ncv;i++) {
236:     BVSetRandomColumn(eps->V,i,eps->rand);
237:   }

239:   if (eps->numbermonitors>0) {
240:     PetscMalloc2(eps->ncv*(eps->max_it+1),&lambdahist,eps->ncv*(eps->max_it+1),&residhist);
241:   }

243: #if defined(PETSC_USE_COMPLEX)
244:   info = lobpcg_solve_complex(blopex->eigenvectors,eps,OperatorAMultiVector,
245:         eps->isgeneralized?eps:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
246:         eps,Precond_FnMultiVector,blopex->Y,
247:         blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
248:         (komplex*)eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
249: #else
250:   info = lobpcg_solve_double(blopex->eigenvectors,eps,OperatorAMultiVector,
251:         eps->isgeneralized?eps:NULL,eps->isgeneralized?OperatorBMultiVector:NULL,
252:         eps,Precond_FnMultiVector,blopex->Y,
253:         blopex->blap_fn,blopex->tol,eps->max_it,0,&its,
254:         eps->eigr,lambdahist,eps->ncv,eps->errest,residhist,eps->ncv);
255: #endif
256:   if (info>0) SETERRQ1(PetscObjectComm((PetscObject)eps),PETSC_ERR_LIB,"Error in blopex (code=%d)",info);

258:   if (eps->numbermonitors>0) {
259:     for (i=0;i<its;i++) {
260:       nconv = 0;
261:       for (j=0;j<eps->ncv;j++) {
262:         if (residhist[j+i*eps->ncv]>eps->tol) break;
263:         else nconv++;
264:       }
265:       EPSMonitor(eps,i,nconv,(PetscScalar*)lambdahist+i*eps->ncv,eps->eigi,residhist+i*eps->ncv,eps->ncv);
266:     }
267:     PetscFree2(lambdahist,residhist);
268:   }

270:   eps->its = its;
271:   eps->nconv = eps->ncv;
272:   STGetShift(eps->st,&sigma);
273:   if (sigma != 0.0) {
274:     for (i=0;i<eps->nconv;i++) eps->eigr[i]+=sigma;
275:   }
276:   if (info==-1) eps->reason = EPS_DIVERGED_ITS;
277:   else eps->reason = EPS_CONVERGED_TOL;
278:   return(0);
279: }

283: PetscErrorCode EPSReset_BLOPEX(EPS eps)
284: {
286:   EPS_BLOPEX     *blopex = (EPS_BLOPEX*)eps->data;

289:   mv_MultiVectorDestroy(blopex->eigenvectors);
290:   mv_MultiVectorDestroy(blopex->Y);
291:   VecDestroy(&blopex->w);
292:   return(0);
293: }

297: PetscErrorCode EPSDestroy_BLOPEX(EPS eps)
298: {

302:   LOBPCG_DestroyRandomContext();
303:   PetscFree(eps->data);
304:   return(0);
305: }

309: PetscErrorCode EPSSetFromOptions_BLOPEX(EPS eps)
310: {
311:   PetscErrorCode  ierr;
312:   KSP             ksp;

315:   PetscOptionsHead("EPS BLOPEX Options");
316:   LOBPCG_SetFromOptionsRandomContext();

318:   /* Set STPrecond as the default ST */
319:   if (!((PetscObject)eps->st)->type_name) {
320:     STSetType(eps->st,STPRECOND);
321:   }
322:   STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);

324:   /* Set the default options of the KSP */
325:   STGetKSP(eps->st,&ksp);
326:   if (!((PetscObject)ksp)->type_name) {
327:     KSPSetType(ksp,KSPPREONLY);
328:   }
329:   PetscOptionsTail();
330:   return(0);
331: }

335: PETSC_EXTERN PetscErrorCode EPSCreate_BLOPEX(EPS eps)
336: {
337:   EPS_BLOPEX     *ctx;

341:   PetscNewLog(eps,&ctx);
342:   eps->data = (void*)ctx;

344:   eps->ops->setup                = EPSSetUp_BLOPEX;
345:   eps->ops->setfromoptions       = EPSSetFromOptions_BLOPEX;
346:   eps->ops->destroy              = EPSDestroy_BLOPEX;
347:   eps->ops->reset                = EPSReset_BLOPEX;
348:   eps->ops->backtransform        = EPSBackTransform_Default;
349:   LOBPCG_InitRandomContext(PetscObjectComm((PetscObject)eps),eps->rand);
350:   return(0);
351: }