#include "ml_config.h"
#include "ml_common.h"
#if defined(HAVE_ML_MLAPI) && defined(HAVE_ML_GALERI)
#include "ml_include.h"
using namespace Teuchos;
using namespace MLAPI;
public:
TwoLevelDDAdditive(
const Operator & FineMatrix,
FineMatrix_(FineMatrix),
R_(R),
P_(P),
FineSolver_(FineSolver),
CoarseSolver_(CoarseSolver)
{}
TwoLevelDDAdditive(const TwoLevelDDAdditive& rhs) :
FineMatrix_(rhs.GetFineMatrix()),
R_(rhs.GetR()),
P_(rhs.GetP()),
FineSolver_(rhs.GetFineSolver()),
CoarseSolver_(rhs.GetCoarseSolver())
{}
TwoLevelDDAdditive& operator=(const TwoLevelDDAdditive& rhs)
{
if (this != &rhs) {
FineMatrix_ = rhs.GetFineMatrix();
R_ = rhs.GetR();
P_ = rhs.GetP();
FineSolver_ = rhs.GetFineSolver();
CoarseSolver_ = rhs.GetCoarseSolver();
}
return(*this);
}
{
x_f = FineSolver_ * r_f;
r_c = R_ * r_f;
r_c = CoarseSolver_ * r_c;
x_f = x_f + P_ * r_c;
return(0);
}
const Space GetOperatorDomainSpace()
const {
return(FineMatrix_.GetDomainSpace());
}
const Space GetOperatorRangeSpace()
const {
return(FineMatrix_.GetRangeSpace());
}
const Space GetDomainSpace()
const {
return(FineMatrix_.GetDomainSpace());
}
const Space GetRangeSpace()
const {
return(FineMatrix_.GetRangeSpace());
}
{
return(FineMatrix_);
}
{
return(R_);
}
{
return(P_);
}
{
return(FineSolver_);
}
{
return(CoarseSolver_);
}
std::ostream& Print(std::ostream& os, const bool verbose = true) const
{
os << "*** MLAPI::TwoLevelDDAdditive ***" << std::endl;
}
return(os);
}
private:
};
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#endif
try {
int NumGlobalElements = 10000;
Space FineSpace(NumGlobalElements);
LHS = 0.0;
RHS.Random();
#ifdef HAVE_ML_METIS
MLList.
set(
"aggregation: type",
"METIS");
MLList.
set(
"aggregation: nodes per aggregate",64);
#else
MLList.
set(
"aggregation: type",
"Uncoupled");
#endif
CoarseMatrix =
GetRAP(R,FineMatrix,P);
FineSolver.
Reshape(FineMatrix,
"symmetric Gauss-Seidel", MLList);
CoarseSolver.
Reshape(CoarseMatrix,
"Amesos-KLU", MLList);
TwoLevelDDAdditive MLAPIPrec(FineMatrix,FineSolver,CoarseSolver,R,P);
MLList.
set(
"krylov: max iterations", 1550);
MLList.
set(
"krylov: tolerance", 1e-9);
MLList.
set(
"krylov: type",
"gmres");
MLList.
set(
"krylov: output", 16);
Krylov(FineMatrix, LHS, RHS, MLAPIPrec, MLList);
}
catch (const int e) {
std::cerr << "Caught integer exception, code = " << e << std::endl;
}
catch (...) {
std::cerr << "Caught exception..." << std::endl;
}
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(EXIT_SUCCESS);
}
#else
#include "ml_include.h"
int main(int argc, char *argv[])
{
#ifdef HAVE_MPI
MPI_Init(&argc,&argv);
#endif
puts("This MLAPI example requires the following configuration options:");
puts("\t--enable-epetra");
puts("\t--enable-teuchos");
puts("\t--enable-ifpack");
puts("\t--enable-amesos");
puts("\t--enable-galeri");
puts("Please check your configure line.");
#ifdef HAVE_MPI
MPI_Finalize();
#endif
return(0);
}
#endif // #if defined(HAVE_ML_MLAPI)