Change order of constructor arguments.

This allows simplification, and having just two constructors.
This commit is contained in:
Atgeirr Flø Rasmussen 2020-06-17 09:31:19 +02:00
parent a6d186551b
commit 2dc2e053d1
6 changed files with 53 additions and 84 deletions

View File

@ -36,32 +36,23 @@ namespace Dune
/// Create a sequential solver.
template <class MatrixType, class VectorType>
FlexibleSolver<MatrixType, VectorType>::
FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix,
FlexibleSolver(const MatrixType& matrix,
const boost::property_tree::ptree& prm,
const std::function<VectorType()>& weightsCalculator)
{
init(prm, matrix, weightsCalculator, Dune::Amg::SequentialInformation());
init(matrix, Dune::Amg::SequentialInformation(), prm, weightsCalculator);
}
/// Create a parallel solver (if Comm is e.g. OwnerOverlapCommunication).
template <class MatrixType, class VectorType>
template <class Comm>
FlexibleSolver<MatrixType, VectorType>::
FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
// const Comm& comm)
const typename std::enable_if<IsComm<Comm>::value, Comm>::type& comm)
FlexibleSolver(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<VectorType()>& weightsCalculator)
{
init(prm, matrix, std::function<VectorType()>(), comm);
}
/// Create a parallel solver (if Comm is e.g. OwnerOverlapCommunication).
template <class MatrixType, class VectorType>
template <class Comm>
FlexibleSolver<MatrixType, VectorType>::
FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function<VectorType()>& weightsCalculator, const Comm& comm)
{
init(prm, matrix, weightsCalculator, comm);
init(matrix, comm, prm, weightsCalculator);
}
template <class MatrixType, class VectorType>
@ -183,8 +174,10 @@ namespace Dune
template <class Comm>
void
FlexibleSolver<MatrixType, VectorType>::
init(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function<VectorType()> weightsCalculator, const Comm& comm)
init(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<VectorType()> weightsCalculator)
{
initOpPrecSp(matrix, prm, weightsCalculator, comm);
initSolver(prm, comm.communicator().rank()==0);
@ -203,45 +196,36 @@ template <int N>
using OBM = Dune::BCRSMatrix<Opm::MatrixBlock<double, N, N>>;
// Variants using Dune::FieldMatrix blocks.
// template class Dune::FlexibleSolver<BM<1>, BV<1>>;
// template class Dune::FlexibleSolver<BM<2>, BV<2>>;
template class Dune::FlexibleSolver<BM<1>, BV<1>>;
template class Dune::FlexibleSolver<BM<2>, BV<2>>;
template class Dune::FlexibleSolver<BM<3>, BV<3>>;
// template class Dune::FlexibleSolver<BM<4>, BV<4>>;
/*
template class Dune::FlexibleSolver<BM<4>, BV<4>>;
// Variants using Opm::MatrixBlock blocks.
template class Dune::FlexibleSolver<OBM<1>, BV<1>>;
template class Dune::FlexibleSolver<OBM<2>, BV<2>>;
template class Dune::FlexibleSolver<OBM<3>, BV<3>>;
template class Dune::FlexibleSolver<OBM<4>, BV<4>>;
using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
template Dune::FlexibleSolver<OBM<1>, BV<1>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const OBM<1>& matrix,
const Comm&);
template Dune::FlexibleSolver<OBM<1>, BV<1>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
const std::function<BV<1>()>& weightsCalculator,
const Comm& comm);
template Dune::FlexibleSolver<OBM<2>, BV<2>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const OBM<2>& matrix,
const Comm&);
template Dune::FlexibleSolver<OBM<2>, BV<2>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
const std::function<BV<2>()>& weightsCalculator,
const Comm& comm);
template Dune::FlexibleSolver<OBM<3>, BV<3>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const OBM<3>& matrix,
const Comm&);
template Dune::FlexibleSolver<OBM<3>, BV<3>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
const std::function<BV<3>()>& weightsCalculator,
const Comm& comm);
template Dune::FlexibleSolver<OBM<4>, BV<4>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const OBM<4>& matrix,
const Comm&);
template Dune::FlexibleSolver<OBM<4>, BV<4>>::FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
const std::function<BV<4>()>& weightsCalculator,
const Comm& comm);
*/
template Dune::FlexibleSolver<OBM<1>, BV<1>>::FlexibleSolver(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<BV<1>()>& weightsCalculator);
template Dune::FlexibleSolver<OBM<2>, BV<2>>::FlexibleSolver(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<BV<2>()>& weightsCalculator);
template Dune::FlexibleSolver<OBM<3>, BV<3>>::FlexibleSolver(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<BV<3>()>& weightsCalculator);
template Dune::FlexibleSolver<OBM<4>, BV<4>>::FlexibleSolver(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<BV<4>()>& weightsCalculator);

View File

@ -35,19 +35,6 @@
namespace Dune
{
template<class C>
struct IsComm : std::false_type
{};
template<>
struct IsComm<Dune::Amg::SequentialInformation> : std::true_type
{};
#if HAVE_MPI
template<class Index>
struct IsComm<Dune::OwnerOverlapCopyCommunication<Index>> : std::true_type
{};
#endif
/// A solver class that encapsulates all needed objects for a linear solver
/// (operator, scalar product, iterative solver and preconditioner) and sets
/// them up based on runtime parameters, using the PreconditionerFactory for
@ -60,20 +47,16 @@ public:
using VectorType = VectorTypeT;
/// Create a sequential solver.
FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix,
FlexibleSolver(const MatrixType& matrix,
const boost::property_tree::ptree& prm,
const std::function<VectorTypeT()>& weightsCalculator = std::function<VectorTypeT()>());
/// Create a parallel solver (if Comm is e.g. OwnerOverlapCommunication).
template <class Comm>
FlexibleSolver(const boost::property_tree::ptree& prm,
const MatrixType& matrix,
// const Comm& comm);
const typename std::enable_if<IsComm<Comm>::value, Comm>::type& comm);
/// Create a parallel solver (if Comm is e.g. OwnerOverlapCommunication).
template <class Comm>
FlexibleSolver(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function<VectorTypeT()>& weightsCalculator, const Comm& comm);
FlexibleSolver(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<VectorTypeT()>& weightsCalculator = std::function<VectorTypeT()>());
virtual void apply(VectorType& x, VectorType& rhs, Dune::InverseOperatorResult& res) override;
@ -105,8 +88,10 @@ private:
// Main initialization routine.
// Call with Comm == Dune::Amg::SequentialInformation to get a serial solver.
template <class Comm>
void init(const boost::property_tree::ptree& prm, const MatrixType& matrix,
const std::function<VectorTypeT()> weightsCalculator, const Comm& comm);
void init(const MatrixType& matrix,
const Comm& comm,
const boost::property_tree::ptree& prm,
const std::function<VectorTypeT()> weightsCalculator);
std::shared_ptr<AbstractOperatorType> linearoperator_;
std::shared_ptr<AbstractPrecondType> preconditioner_;

View File

@ -905,10 +905,10 @@ protected:
if (isParallel()) {
#if HAVE_MPI
assert(noGhostMat_);
flexibleSolver_.reset(new FlexibleSolverType(prm_, *noGhostMat_, weightsCalculator, *comm_));
flexibleSolver_.reset(new FlexibleSolverType(*noGhostMat_, *comm_, prm_, weightsCalculator));
#endif
} else {
flexibleSolver_.reset(new FlexibleSolverType(prm_, *matrix_, weightsCalculator));
flexibleSolver_.reset(new FlexibleSolverType(*matrix_, prm_, weightsCalculator));
}
}
else

View File

@ -199,11 +199,11 @@ public:
if (isParallel()) {
#if HAVE_MPI
matrix_ = &mat.istlMatrix();
solver_.reset(new SolverType(prm_, mat.istlMatrix(), weightsCalculator, *comm_));
solver_.reset(new SolverType(mat.istlMatrix(), *comm_, prm_, weightsCalculator));
#endif
} else {
matrix_ = &mat.istlMatrix();
solver_.reset(new SolverType(prm_, mat.istlMatrix(), weightsCalculator));
solver_.reset(new SolverType(mat.istlMatrix(), prm_, weightsCalculator));
}
rhs_ = b;
} else {

View File

@ -61,9 +61,9 @@ namespace Amg
: linsolver_()
{
if (op.category() == Dune::SolverCategory::overlapping) {
linsolver_.reset(new Solver(prm, op.getmat(), std::function<X()>(), comm));
linsolver_.reset(new Solver(op.getmat(), comm, prm, std::function<X()>()));
} else {
linsolver_.reset(new Solver(prm, op.getmat(), std::function<X()>()));
linsolver_.reset(new Solver(op.getmat(), prm, std::function<X()>()));
}
}

View File

@ -70,7 +70,7 @@ testSolver(const boost::property_tree::ptree& prm, const std::string& matrix_fil
prm.get<int>("preconditioner.pressure_var_index"),
transpose);
};
Dune::FlexibleSolver<Matrix, Vector> solver(prm, matrix, wc);
Dune::FlexibleSolver<Matrix, Vector> solver(matrix, prm, wc);
Vector x(rhs.size());
Dune::InverseOperatorResult res;
solver.apply(x, rhs, res);