#include <TpetraExt_MatrixMatrix.hpp>
#include <Tpetra_Core.hpp>
#include <Tpetra_CrsMatrix_fwd.hpp>
#include <Tpetra_Map_fwd.hpp>
#include <Tpetra_Vector_fwd.hpp>
#include <Teuchos_Comm.hpp>
#include <Teuchos_CommHelpers.hpp>
#include <Teuchos_CommandLineProcessor.hpp>
#include <Teuchos_DefaultComm.hpp>
#include <Teuchos_FancyOStream.hpp>
#include <Teuchos_ParameterList.hpp>
#include <Teuchos_RCP.hpp>
#include <Teuchos_StandardCatchMacros.hpp>
#include <Teuchos_Tuple.hpp>
#include <Teuchos_VerboseObject.hpp>
#include "BelosTpetraAdapter.hpp"
template <typename ScalarType>
int run(int argc, char *argv[]) {
using Teuchos::ParameterList;
using Teuchos::RCP;
using Teuchos::rcp;
using Teuchos::tuple;
using ST = typename Tpetra::MultiVector<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 V = typename Tpetra::Vector<ST, LO, GO, NT>;
using MV = typename Tpetra::MultiVector<ST, LO, GO, NT>;
using OP = typename Tpetra::Operator<ST, LO, GO, NT>;
using MAP = typename Tpetra::Map<LO, GO, NT>;
using MAT = typename Tpetra::CrsMatrix<ST, LO, GO, NT>;
using SCT = typename Teuchos::ScalarTraits<ST>;
using MT = typename SCT::magnitudeType;
Teuchos::GlobalMPISession session(&argc, &argv, NULL);
RCP<const Teuchos::Comm<int>> comm = Tpetra::getDefaultComm();
bool verbose = false;
bool success = true;
try {
bool procVerbose = false;
int frequency = -1;
int blocksize = 1;
int numRhs = 1;
int maxIters = 30;
int maxDeflate = 4;
int maxSave = 8;
std::string ortho("ICGS");
MT tol = sqrt(std::numeric_limits<ST>::epsilon());
Teuchos::CommandLineProcessor cmdp(false, true);
cmdp.setOption("verbose", "quiet", &verbose, "Print messages and results");
cmdp.setOption("frequency", &frequency, "Solvers frequency for printing residuals (#iters)");
cmdp.setOption("tol", &tol, "Relative residual tolerance used by PCPG solver");
cmdp.setOption("num-rhs", &numRhs, "Number of right-hand sides to be solved for");
cmdp.setOption("max-iters", &maxIters,
"Maximum number of iterations per linear system (-1 = "
"adapted to problem/block size)");
cmdp.setOption("num-deflate", &maxDeflate, "Number of vectors deflated from the linear system");
cmdp.setOption("num-save", &maxSave, "Number of vectors saved from old Krylov subspaces");
cmdp.setOption("ortho-type", &ortho, "Orthogonalization type, either DGKS, ICGS or IMGS");
if (cmdp.parse(argc, argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
return -1;
}
if (!verbose)
frequency = -1;
int numTimeStep = 4;
GO numElePerDirection = 14 * comm->getSize();
size_t numNodes = (numElePerDirection - 1) * (numElePerDirection - 1);
RCP<MAP> map = rcp(new MAP(numNodes, 0, comm));
RCP<MAT> stiff = rcp(new MAT(map, numNodes));
RCP<MAT> mass = rcp(new MAT(map, numNodes));
RCP<V> vecLHS = rcp(new V(map));
RCP<V> vecRHS = rcp(new V(map));
RCP<MV> LHS, RHS;
ST ko = 8.0 / 3.0, k1 = -1.0 / 3.0;
ST h = 1.0 / static_cast<ST>(numElePerDirection);
ST mo = h * h * 4.0 / 9.0, m1 = h * h / 9.0, m2 = h * h / 36.0;
ST pi = 4.0 * atan(1.0), valueLHS;
GO node, iX, iY;
for (LO lid = map->getMinLocalIndex(); lid <= map->getMaxLocalIndex(); lid++) {
node = map->getGlobalElement(lid);
iX = node % (numElePerDirection - 1);
iY = (node - iX) / (numElePerDirection - 1);
GO pos = node;
stiff->insertGlobalValues(node, tuple(pos), tuple(ko));
mass->insertGlobalValues(node, tuple(pos),
tuple(mo));
valueLHS = sin(pi * h * ((ST)iX + 1)) * cos(2.0 * pi * h * ((ST)iY + 1));
vecLHS->replaceGlobalValue(1, valueLHS);
if (iY > 0) {
pos = iX + (iY - 1) * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m1));
}
if (iY < numElePerDirection - 2) {
pos = iX + (iY + 1) * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m1));
}
if (iX > 0) {
pos = iX - 1 + iY * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m1));
if (iY > 0) {
pos = iX - 1 + (iY - 1) * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m2));
}
if (iY < numElePerDirection - 2) {
pos = iX - 1 + (iY + 1) * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m2));
}
}
if (iX < numElePerDirection - 2) {
pos = iX + 1 + iY * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m1));
if (iY > 0) {
pos = iX + 1 + (iY - 1) * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m2));
}
if (iY < numElePerDirection - 2) {
pos = iX + 1 + (iY + 1) * (numElePerDirection - 1);
stiff->insertGlobalValues(node, tuple(pos), tuple(k1));
mass->insertGlobalValues(node, tuple(pos), tuple(m2));
}
}
}
stiff->fillComplete();
mass->fillComplete();
const ST one = SCT::one();
ST hdt = .00005;
RCP<MAT> A = Tpetra::MatrixMatrix::add(one, false, *mass, hdt, false, *stiff);
hdt = -hdt;
RCP<MAT> B = Tpetra::MatrixMatrix::add(one, false, *mass, hdt, false, *stiff);
B->apply(*vecLHS, *vecRHS);
procVerbose = verbose && (comm->getRank() == 0);
LHS = Teuchos::rcp_implicit_cast<MV>(vecLHS);
RHS = Teuchos::rcp_implicit_cast<MV>(vecRHS);
const int numGlobalElements = RHS->getGlobalLength();
if (maxIters == -1)
maxIters = numGlobalElements / blocksize - 1;
RCP<ParameterList> belosList = rcp(new ParameterList());
belosList->set("Block Size",
blocksize);
belosList->set("Maximum Iterations",
maxIters);
belosList->set("Convergence Tolerance",
tol);
belosList->set("Num Deflated Blocks",
maxDeflate);
belosList->set("Num Saved Blocks",
maxSave);
belosList->set("Orthogonalization", ortho);
if (numRhs > 1) {
belosList->set("Show Maximum Residual Norm Only",
true);
}
if (verbose) {
if (frequency > 0)
belosList->set("Output Frequency", frequency);
} else
RCP<LinearProblem> problem = rcp(new LinearProblem(A, LHS, RHS));
bool set = problem->setProblem();
if (set == false) {
if (procVerbose)
std::cout << std::endl
<< "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
return -1;
}
RCP<PCPGSolMgr> solver = rcp(new PCPGSolMgr(problem, belosList));
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 << "Maximum number of iterations allowed: " << maxIters << std::endl;
std::cout << "Relative residual tolerance: " << tol << std::endl;
std::cout << std::endl;
}
bool badRes;
for (int timeStep = 0; timeStep < numTimeStep; timeStep++) {
if (timeStep) {
B->apply(*LHS, *RHS);
set = problem->setProblem(LHS, RHS);
if (set == false) {
if (procVerbose)
std::cout << std::endl
<< "ERROR: Belos::LinearProblem failed to set up correctly!" << std::endl;
return -1;
}
}
std::vector<ST> rhs_norm(numRhs);
MVT::MvNorm(*RHS, rhs_norm);
std::cout << " RHS norm is ... " << rhs_norm[0] << std::endl;
badRes = false;
std::vector<ST> actual_resids(numRhs);
MV resid(map, numRhs);
OPT::Apply(*A, *LHS, resid);
MVT::MvAddMv(-1.0, resid, 1.0, *RHS, resid);
MVT::MvNorm(resid, actual_resids);
MVT::MvNorm(*RHS, rhs_norm);
std::cout << " RHS norm is ... " << rhs_norm[0] << std::endl;
if (procVerbose) {
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;
break;
}
}
if (procVerbose) {
if (success)
std::cout << std::endl << "SUCCESS: Belos converged!" << std::endl;
else
std::cout << std::endl << "ERROR: Belos did not converge!" << 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);
}
Class which describes the linear problem to be solved by the iterative solver.
Declaration and definition of Belos::PCPGSolMgr (PCPG iterative linear solver).
Traits class which defines basic operations on multivectors.
Class which defines basic traits for the operator type.
PCPG iterative linear solver.
ReturnType
Whether the Belos solve converged for all linear systems.