#include <Tpetra_Map.hpp>
#include <Tpetra_Core.hpp>
#include <Tpetra_Vector.hpp>
#include <Tpetra_CrsMatrix.hpp>
#include <Teuchos_Comm.hpp>
#include <Teuchos_CommHelpers.hpp>
#include <Teuchos_DefaultComm.hpp>
#include "Teuchos_StandardCatchMacros.hpp"
#include "BelosTpetraTestFramework.hpp"
template <typename ScalarType>
int run(int argc, char *argv[]) {
using ST = typename Tpetra::Vector<ScalarType>::scalar_type;
using LO = typename Tpetra::Vector<>::local_ordinal_type;
using GO = typename Tpetra::Vector<>::global_ordinal_type;
using NT = typename Tpetra::Vector<>::node_type;
using OP = typename Tpetra::Operator<ST,LO,GO,NT>;
using MV = typename Tpetra::MultiVector<ST,LO,GO,NT>;
using tmap_t = Tpetra::Map<LO,GO,NT>;
using tvector_t = Tpetra::Vector<ST,LO,GO,NT>;
using tcrsmatrix_t = Tpetra::CrsMatrix<ST,LO,GO,NT>;
const auto comm = Tpetra::getDefaultComm();
const int myPID = comm->getRank();
bool verbose = false;
bool success = true;
try {
bool procVerbose = false;
int frequency = -1;
int blockSize = 1;
int numrhs = 1;
int maxIters = -1;
int maxSubspace = 50;
int maxRestarts = 15;
std::string filename("orsirr1.hb");
MT tol = 1.0e-5;
cmdp.
setOption(
"verbose",
"quiet",&verbose,
"Print messages and results.");
cmdp.
setOption(
"frequency",&frequency,
"Solvers frequency for printing residuals (#iters).");
cmdp.
setOption(
"filename",&filename,
"Filename for test matrix. Acceptable file extensions: *.hb,*.mtx,*.triU,*.triS");
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.");
return -1;
}
if (!verbose)
frequency = -1;
procVerbose = ( verbose && (myPID==0) );
if (procVerbose) {
}
Belos::Tpetra::HarwellBoeingReader<tcrsmatrix_t> reader( comm );
RCP<tcrsmatrix_t> A = reader.readFromFile( filename );
RCP<const tmap_t> Map = A->getDomainMap();
RCP<MV> B, X;
X =
rcp(
new MV(Map,numrhs) );
MVT::MvRandom( *X );
B =
rcp(
new MV(Map,numrhs) );
OPT::Apply( *A, *X, *B );
MVT::MvInit( *X, 0.0 );
const int numGlobalElements = B->getGlobalLength();
if (maxIters == -1)
maxIters = numGlobalElements - 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 );
}
else
if (set == false) {
if (procVerbose)
std::cout << std::endl << "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
return -1;
}
RCP< Belos::SolverManager<ST,MV,OP> > newSolver
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 restart cycle: " << maxIters << std::endl;
std::cout << "Relative residual tolerance: " << tol << std::endl;
std::cout << 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 << std::endl << "ERROR: Belos did not converge!" << std::endl;
} else {
success = true;
if (procVerbose)
std::cout << std::endl << "SUCCESS: Belos converged!" << std::endl;
}
}
return success ? EXIT_SUCCESS : EXIT_FAILURE;
}
int main(int argc, char *argv[]) {
return run<double>(argc,argv);
}