Use LOBPCG with Epetra, for a generalized eigenvalue problem.This example computes the eigenvalues of largest magnitude of an generalized eigenvalue problem, using Anasazi's implementation of the LOBPCG method, with Epetra linear algebra.
#include "Epetra_CrsMatrix.h"
#include "Teuchos_StandardCatchMacros.hpp"
#ifdef HAVE_MPI
#include "Epetra_MpiComm.h"
#include <mpi.h>
#else
#include "Epetra_SerialComm.h"
#endif
#include "Epetra_Map.h"
#include "ModeLaplace2DQ2.h"
using namespace Anasazi;
int main(int argc, char *argv[]) {
#ifdef HAVE_MPI
  
  
  MPI_Init(&argc,&argv);
#endif
  bool success = false;
  try {
    
    
#ifdef HAVE_MPI
#else
#endif
    
    
    printer.
stream(
Errors) << Anasazi_Version() << std::endl << std::endl;
 
    
    
    std::string which("SM");
    cmdp.
setOption(
"sort",&which,
"Targetted eigenvalues (SM or LM).");
 
      throw -1;
    }
    
    const int space_dim = 2;
    
    std::vector<double> brick_dim( space_dim );
    brick_dim[0] = 1.0;
    brick_dim[1] = 1.0;
    
    std::vector<int> elements( space_dim );
    elements[0] = 10;
    elements[1] = 10;
    
      Teuchos::rcp( 
new ModeLaplace2DQ2(Comm, brick_dim[0], elements[0], brick_dim[1], elements[1]) );
 
    
    
    int nev = 10;
    int blockSize = 5;
    int maxIters = 500;
    double tol = 1.0e-8;
    ivec->Random();
    
    
    MyProblem->setHermitian(true);
    
    MyProblem->setNEV( nev );
    
    bool boolret = MyProblem->setProblem();
    if (boolret != true) {
      printer.
print(
Errors,
"Anasazi::BasicEigenproblem::setProblem() returned an error.\n");
 
      throw -1;
    }
    
    
    MyPL.
set( 
"Which", which );
 
    MyPL.
set( 
"Block Size", blockSize );
 
    MyPL.
set( 
"Maximum Iterations", maxIters );
 
    MyPL.
set( 
"Convergence Tolerance", tol );
 
    MyPL.
set( 
"Full Ortho", 
true );
 
    MyPL.
set( 
"Use Locking", 
true );
 
    MyPL.
set( 
"Verbosity", verbosity );
 
    
    
    
    
    
    
    std::vector<Value<double> > evals = sol.
Evals;
 
    
    
    std::vector<double> normR(sol.
numVecs);
 
      for (
int i=0; i<sol.
numVecs; i++) {
 
        T(i,i) = evals[i].realpart;
      }
      K->
Apply( *evecs, Kvec );
 
      M->
Apply( *evecs, Mvec );
 
      MVT::MvTimesMatAddMv( -1.0, Mvec, T, 1.0, Kvec );
      MVT::MvNorm( Kvec, normR );
    }
    
    
    std::ostringstream os;
    os.setf(std::ios_base::right, std::ios_base::adjustfield);
    os<<
"Solver manager returned " << (returnCode == 
Converged ? 
"converged." : 
"unconverged.") << std::endl;
 
    os<<std::endl;
    os<<"------------------------------------------------------"<<std::endl;
    os<<std::setw(16)<<"Eigenvalue"
      <<std::setw(18)<<"Direct Residual"
      <<std::endl;
    os<<"------------------------------------------------------"<<std::endl;
    for (
int i=0; i<sol.
numVecs; i++) {
 
      os<<std::setw(16)<<evals[i].realpart
        <<std::setw(18)<<normR[i]/evals[i].realpart
        <<std::endl;
    }
    os<<"------------------------------------------------------"<<std::endl;
    success = true;
  }
#ifdef HAVE_MPI
  MPI_Finalize();
#endif
  return ( success ? EXIT_SUCCESS : EXIT_FAILURE );
}