#include "BelosEpetraAdapter.hpp"
#include "Epetra_CrsMatrix.h"
#include "Epetra_MultiVector.h"
#include "Galeri_Maps.h"
#include "Galeri_CrsMatrices.h"
#include <Teuchos_CommandLineProcessor.hpp>
#include <Teuchos_GlobalMPISession.hpp>
#include <Teuchos_oblackholestream.hpp>
#include <Teuchos_TypeNameTraits.hpp>
#include <Teuchos_StandardCatchMacros.hpp>
#ifdef EPETRA_MPI
# include "Epetra_MpiComm.h"
#else
# include "Epetra_SerialComm.h"
#endif
int
main (int argc, char *argv[])
{
int MyPID = 0;
typedef double ST;
typedef Teuchos::ScalarTraits<ST> SCT;
typedef SCT::magnitudeType MT;
typedef Epetra_MultiVector MV;
typedef Epetra_Operator OP;
using Teuchos::ParameterList;
using Teuchos::RCP;
using Teuchos::rcp;
#ifdef EPETRA_MPI
MPI_Init (&argc, &argv);
Epetra_MpiComm Comm (MPI_COMM_WORLD);
#else
Epetra_SerialComm Comm;
#endif
bool verbose = false;
bool success = true;
try {
bool proc_verbose = false;
bool debug = false;
int frequency = -1;
int blocksize = 1;
int numrhs = 1;
int maxiters = -1;
int maxsubspace = 50;
int maxrestarts = 15;
int nx = 10;
MT tol = 1.0e-5;
Teuchos::CommandLineProcessor cmdp(false,true);
cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
cmdp.setOption("debug","nondebug",&debug,"Print debugging information from solver.");
cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
cmdp.setOption("tol",&tol,"Relative residual tolerance used by solver.");
cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
cmdp.setOption("block-size",&blocksize,"Block size used by solver.");
cmdp.setOption("max-iters",&maxiters,"Maximum number of iterations per linear system (-1 = adapted to problem/block size).");
cmdp.setOption("max-subspace",&maxsubspace,"Maximum number of blocks the solver can use for the subspace.");
cmdp.setOption("max-restarts",&maxrestarts,"Maximum number of restarts allowed for solver.");
cmdp.setOption("nx",&nx,"Number of discretization points in each direction of 2D Laplacian.");
if (cmdp.parse(argc,argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
return -1;
}
if (!verbose)
frequency = -1;
Teuchos::ParameterList GaleriList;
GaleriList.set ("n", nx * nx);
GaleriList.set ("nx", nx);
GaleriList.set ("ny", nx);
RCP<Epetra_Map> Map = rcp (Galeri::CreateMap ("Linear", Comm, GaleriList));
RCP<Epetra_RowMatrix> A =
rcp (Galeri::CreateCrsMatrix ("Laplace2D", &*Map, GaleriList));
proc_verbose = verbose && (MyPID==0);
RCP<MV> B = rcp (new MV (*Map, numrhs));
RCP<MV> X = rcp (new MV (*Map, numrhs));
RCP<MV> Xexact = rcp (new MV (*Map, numrhs));
Xexact->Random ();
A->Apply( *Xexact, *B );
const int NumGlobalElements = B->GlobalLength();
if (maxiters == -1)
maxiters = NumGlobalElements/blocksize - 1;
ParameterList belosList;
belosList.set( "Num Blocks", maxsubspace);
belosList.set( "Block Size", blocksize );
belosList.set( "Maximum Iterations", maxiters );
belosList.set( "Maximum Restarts", maxrestarts );
belosList.set( "Convergence Tolerance", tol );
if (verbose) {
if (frequency > 0)
belosList.set( "Output Frequency", frequency );
}
if (debug) {
}
belosList.set( "Verbosity", verbosity );
if (set == false) {
if (proc_verbose)
std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
return -1;
}
std::string solverName = "Block GMRES";
RCP< Belos::SolverManager<double,MV,OP> > newSolver = factory.create (solverName, rcp(&belosList,false));
newSolver->setProblem( rcp(&problem,false) );
if (proc_verbose) {
std::cout << std::endl << std::endl;
std::cout << "Dimension of matrix: " << NumGlobalElements << std::endl;
std::cout << "Number of right-hand sides: " << numrhs << std::endl;
std::cout << "Block size used by solver: " << blocksize << std::endl;
std::cout << "Max number of restarts allowed: " << maxrestarts << std::endl;
std::cout << "Max number of Gmres iterations per linear system: " << maxiters << std::endl;
std::cout << "Relative residual tolerance: " << tol << std::endl;
std::cout << std::endl;
}
int numIters = newSolver->getNumIters();
if (proc_verbose)
std::cout << "Number of iterations performed for this solve: " << numIters << std::endl;
bool badRes = false;
std::vector<double> actual_resids( numrhs );
std::vector<double> rhs_norm( numrhs );
Epetra_MultiVector resid(*Map, numrhs);
OPT::Apply( *A, *X, resid );
MVT::MvAddMv( -1.0, resid, 1.0, *B, resid );
MVT::MvNorm( resid, actual_resids );
MVT::MvNorm( *B, rhs_norm );
if (proc_verbose) {
std::cout<< "---------- Actual Residuals (normalized) ----------"<<std::endl<<std::endl;
for ( int i=0; i<numrhs; i++) {
double actRes = actual_resids[i]/rhs_norm[i];
std::cout<<"Problem "<<i<<" : \t"<< actRes <<std::endl;
if (actRes > tol) badRes = true;
}
}
success = false;
if (proc_verbose)
std::cout << "End Result: TEST FAILED" << std::endl;
} else {
if (proc_verbose)
std::cout << "End Result: TEST PASSED" << std::endl;
}
}
TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
#ifdef EPETRA_MPI
MPI_Finalize();
#endif
return success ? EXIT_SUCCESS : EXIT_FAILURE;
}
virtual bool setProblem(const Teuchos::RCP< MV > &newX=Teuchos::null, const Teuchos::RCP< const MV > &newB=Teuchos::null)
Set up the linear problem manager.
Traits class which defines basic operations on multivectors.
Class which defines basic traits for the operator type.
ReturnType
Whether the Belos solve converged for all linear systems.
typename ::Belos::Impl::SolverFactorySelector< SC, MV, OP >::type SolverFactory