#include <Tpetra_Core.hpp>
#include <Tpetra_CrsMatrix.hpp>
#include <Tpetra_MultiVector.hpp>
#include <Galeri_XpetraMaps.hpp>
#include <Galeri_XpetraMatrixTypes.hpp>
#include <Galeri_XpetraProblemFactory.hpp>
#include <Teuchos_RCP.hpp>
#include <Teuchos_Comm.hpp>
#include <Teuchos_CommHelpers.hpp>
#include <Teuchos_DefaultComm.hpp>
#include "Teuchos_ParameterList.hpp"
#include <Teuchos_TypeNameTraits.hpp>
#include <Teuchos_oblackholestream.hpp>
#include "Teuchos_StandardCatchMacros.hpp"
#include "Teuchos_CommandLineProcessor.hpp"
#include "BelosTpetraAdapter.hpp"
template <typename ScalarType>
int run(int argc, char *argv[]) {
using ST = typename Tpetra::MultiVector<ScalarType>::scalar_type;
using LO = typename Tpetra::MultiVector<>::local_ordinal_type;
using GO = typename Tpetra::MultiVector<>::global_ordinal_type;
using NT = typename Tpetra::MultiVector<>::node_type;
using OP = typename Tpetra::Operator<ST,LO,GO,NT>;
using MV = typename Tpetra::MultiVector<ST,LO,GO,NT>;
using MT = typename Teuchos::ScalarTraits<ST>::magnitudeType;
using tmap_t = Tpetra::Map<LO,GO,NT>;
using tvector_t = Tpetra::Vector<ST,LO,GO,NT>;
using trowmatrix_t = Tpetra::RowMatrix<ST,LO,GO,NT>;
using tcrsmatrix_t = Tpetra::CrsMatrix<ST,LO,GO,NT>;
using Teuchos::RCP;
using Teuchos::rcp;
using Teuchos::ParameterList;
Teuchos::GlobalMPISession mpiSession (&argc, &argv, &std::cout);
const auto comm = Tpetra::getDefaultComm();
const int myPID = comm->getRank();
bool verbose = false;
bool success = true;
try {
bool procVerbose = 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 GMRES solver.");
cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
cmdp.setOption("block-size",&blockSize,"Block size used by GMRES.");
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 GMRES solver.");
cmdp.setOption("nx",&nx,"Number of discretization points in each direction of 3D Laplacian.");
if (cmdp.parse(argc,argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
return -1;
}
if (!verbose)
frequency = -1;
procVerbose = ( verbose && (myPID==0) );
if (procVerbose) {
}
Teuchos::ParameterList GaleriList;
GaleriList.set ("n", nx * nx );
GaleriList.set ("nx", nx);
GaleriList.set ("ny", nx);
auto Map = RCP{Galeri::Xpetra::CreateMap<ST,GO,tmap_t>("Cartesian2D", comm, GaleriList)};
auto GaleriProblem = Galeri::Xpetra::BuildProblem<ST,LO,GO,tmap_t,tcrsmatrix_t,MV>("Laplace2D", Map, GaleriList);
auto A = GaleriProblem->BuildMatrix();
RCP<MV> B = rcp (new MV (Map, numrhs));
RCP<MV> X = rcp (new MV (Map, numrhs));
RCP<MV> Xexact = rcp (new MV (Map, numrhs));
MVT::MvRandom(*Xexact);
OPT::Apply(*A, *Xexact, *B );
const int numGlobalElements = B->getGlobalLength();
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 (procVerbose)
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 (procVerbose) {
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 (procVerbose)
std::cout << "Number of iterations performed for this solve: " << numIters << std::endl;
bool badRes = false;
std::vector<ST> actualResids( numrhs );
std::vector<ST> rhsNorm( numrhs );
MV resid(Map, numrhs);
OPT::Apply( *A, *X, resid );
MVT::MvAddMv( -1.0, resid, 1.0, *B, resid );
MVT::MvNorm( resid, actualResids );
MVT::MvNorm( *B, rhsNorm );
if (procVerbose) {
std::cout<< "---------- Actual Residuals (normalized) ----------"<<std::endl<<std::endl;
for ( int i=0; i<numrhs; i++) {
ST actRes = actualResids[i]/rhsNorm[i];
std::cout<<"Problem "<<i<<" : \t"<< actRes <<std::endl;
if (actRes > tol) badRes = true;
}
}
success = false;
if (procVerbose)
std::cout << "End Result: TEST FAILED" << std::endl;
} else {
if (procVerbose)
std::cout << "End Result: TEST PASSED" << std::endl;
}
}
TEUCHOS_STANDARD_CATCH_STATEMENTS(verbose, std::cerr, success);
return success ? EXIT_SUCCESS : EXIT_FAILURE;
}
int main(int argc, char *argv[]) {
return run<double>(argc,argv);
}
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
std::string Belos_Version()