mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Ghost entries skipped for ilu apply and GL operator in AMG/CPR hierarchy.
This works since the ghost entries are the last entries
This commit is contained in:
parent
302503e172
commit
6c62753803
@ -273,9 +273,11 @@ using SeqOpW = Opm::WellModelMatrixAdapter<OBM<N>, BV<N>, BV<N>, false>;
|
|||||||
// Parallel communicator and operators.
|
// Parallel communicator and operators.
|
||||||
using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
|
using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
|
||||||
template <int N>
|
template <int N>
|
||||||
using ParOpM = Dune::OverlappingSchwarzOperator<OBM<N>, BV<N>, BV<N>, Comm>;
|
using ParOpM = Opm::GhostLastMatrixAdapter<OBM<N>, BV<N>, BV<N>, Comm>;
|
||||||
template <int N>
|
template <int N>
|
||||||
using ParOpW = Opm::WellModelGhostLastMatrixAdapter<OBM<N>, BV<N>, BV<N>, true>;
|
using ParOpW = Opm::WellModelGhostLastMatrixAdapter<OBM<N>, BV<N>, BV<N>, true>;
|
||||||
|
template <int N>
|
||||||
|
using ParOpD = Dune::OverlappingSchwarzOperator<OBM<N>, BV<N>, BV<N>, Comm>;
|
||||||
|
|
||||||
// Note: we must instantiate the constructor that is a template.
|
// Note: we must instantiate the constructor that is a template.
|
||||||
// This is only needed in the parallel case, since otherwise the Comm type is
|
// This is only needed in the parallel case, since otherwise the Comm type is
|
||||||
@ -292,7 +294,8 @@ template Dune::FlexibleSolver<Operator>::FlexibleSolver(Operator& op,
|
|||||||
INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpM<N>); \
|
INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpM<N>); \
|
||||||
INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpW<N>); \
|
INSTANTIATE_FLEXIBLESOLVER_OP(SeqOpW<N>); \
|
||||||
INSTANTIATE_FLEXIBLESOLVER_OP(ParOpM<N>); \
|
INSTANTIATE_FLEXIBLESOLVER_OP(ParOpM<N>); \
|
||||||
INSTANTIATE_FLEXIBLESOLVER_OP(ParOpW<N>);
|
INSTANTIATE_FLEXIBLESOLVER_OP(ParOpW<N>); \
|
||||||
|
INSTANTIATE_FLEXIBLESOLVER_OP(ParOpD<N>);
|
||||||
|
|
||||||
#else // HAVE_MPI
|
#else // HAVE_MPI
|
||||||
|
|
||||||
|
@ -111,7 +111,7 @@ void FlexibleSolverInfo<Matrix,Vector,Comm>::create(const Matrix& matrix,
|
|||||||
if (parallel) {
|
if (parallel) {
|
||||||
#if HAVE_MPI
|
#if HAVE_MPI
|
||||||
if (!wellOperator_) {
|
if (!wellOperator_) {
|
||||||
using ParOperatorType = Dune::OverlappingSchwarzOperator<Matrix, Vector, Vector, Comm>;
|
using ParOperatorType = Opm::GhostLastMatrixAdapter<Matrix, Vector, Vector, Comm>;
|
||||||
auto pop = std::make_unique<ParOperatorType>(matrix, comm);
|
auto pop = std::make_unique<ParOperatorType>(matrix, comm);
|
||||||
using FlexibleSolverType = Dune::FlexibleSolver<ParOperatorType>;
|
using FlexibleSolverType = Dune::FlexibleSolver<ParOperatorType>;
|
||||||
auto sol = std::make_unique<FlexibleSolverType>(*pop, comm, prm,
|
auto sol = std::make_unique<FlexibleSolverType>(*pop, comm, prm,
|
||||||
|
@ -185,6 +185,35 @@ void convertToCRS(const M& A, CRS& lower, CRS& upper, InvVector& inv)
|
|||||||
assert(colcount == numUpper);
|
assert(colcount == numUpper);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
template <class PI>
|
||||||
|
size_t set_interiorSize( [[maybe_unused]] size_t N, size_t interiorSize, [[maybe_unused]] const PI& comm)
|
||||||
|
{
|
||||||
|
return interiorSize;
|
||||||
|
}
|
||||||
|
|
||||||
|
#if HAVE_MPI
|
||||||
|
template<>
|
||||||
|
size_t set_interiorSize(size_t N, size_t interiorSize, const Dune::OwnerOverlapCopyCommunication<int,int>& comm)
|
||||||
|
{
|
||||||
|
if (interiorSize<N)
|
||||||
|
return interiorSize;
|
||||||
|
auto indexSet = comm.indexSet();
|
||||||
|
|
||||||
|
size_t new_is = 0;
|
||||||
|
for (auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx)
|
||||||
|
{
|
||||||
|
if (idx->local().attribute()==1)
|
||||||
|
{
|
||||||
|
auto loc = idx->local().local();
|
||||||
|
if (loc > new_is) {
|
||||||
|
new_is = loc;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
return new_is + 1;
|
||||||
|
}
|
||||||
|
#endif
|
||||||
|
|
||||||
} // end namespace detail
|
} // end namespace detail
|
||||||
|
|
||||||
|
|
||||||
@ -411,6 +440,11 @@ update()
|
|||||||
{
|
{
|
||||||
OPM_TIMEBLOCK(iluDecomposition);
|
OPM_TIMEBLOCK(iluDecomposition);
|
||||||
if (iluIteration_ == 0) {
|
if (iluIteration_ == 0) {
|
||||||
|
|
||||||
|
if (comm_) {
|
||||||
|
interiorSize_ = detail::set_interiorSize(A_->N(), interiorSize_, *comm_);
|
||||||
|
}
|
||||||
|
|
||||||
// create ILU-0 decomposition
|
// create ILU-0 decomposition
|
||||||
if (ordering_.empty())
|
if (ordering_.empty())
|
||||||
{
|
{
|
||||||
|
@ -203,7 +203,8 @@ struct StandardPreconditioners {
|
|||||||
// is the overlapping schwarz operator. This could be extended
|
// is the overlapping schwarz operator. This could be extended
|
||||||
// later, but at this point no other operators are compatible
|
// later, but at this point no other operators are compatible
|
||||||
// with the AMG hierarchy construction.
|
// with the AMG hierarchy construction.
|
||||||
if constexpr (std::is_same_v<O, Dune::OverlappingSchwarzOperator<M, V, V, C>>) {
|
if constexpr (std::is_same_v<O, Dune::OverlappingSchwarzOperator<M, V, V, C>> ||
|
||||||
|
std::is_same_v<O, Opm::GhostLastMatrixAdapter<M, V, V, C>>) {
|
||||||
F::addCreator("amg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
|
F::addCreator("amg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
|
||||||
using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<V, V>>;
|
using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<V, V>>;
|
||||||
const std::string smoother = prm.get<std::string>("smoother", "ParOverILU0");
|
const std::string smoother = prm.get<std::string>("smoother", "ParOverILU0");
|
||||||
@ -785,14 +786,28 @@ using OpBPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<MatrixBlock<dou
|
|||||||
Dune::BlockVector<Dune::FieldVector<double, Dim>>,
|
Dune::BlockVector<Dune::FieldVector<double, Dim>>,
|
||||||
Dune::BlockVector<Dune::FieldVector<double, Dim>>,
|
Dune::BlockVector<Dune::FieldVector<double, Dim>>,
|
||||||
CommPar>;
|
CommPar>;
|
||||||
|
template<int Dim>
|
||||||
|
using OpGLFPar = Opm::GhostLastMatrixAdapter<Dune::BCRSMatrix<Dune::FieldMatrix<double,Dim,Dim>>,
|
||||||
|
Dune::BlockVector<Dune::FieldVector<double,Dim>>,
|
||||||
|
Dune::BlockVector<Dune::FieldVector<double,Dim>>,
|
||||||
|
CommPar>;
|
||||||
|
|
||||||
|
template<int Dim>
|
||||||
|
using OpGLBPar = Opm::GhostLastMatrixAdapter<Dune::BCRSMatrix<MatrixBlock<double,Dim,Dim>>,
|
||||||
|
Dune::BlockVector<Dune::FieldVector<double,Dim>>,
|
||||||
|
Dune::BlockVector<Dune::FieldVector<double,Dim>>,
|
||||||
|
CommPar>;
|
||||||
|
|
||||||
#define INSTANCE_PF_PAR(Dim) \
|
#define INSTANCE_PF_PAR(Dim) \
|
||||||
template class PreconditionerFactory<OpBSeq<Dim>, CommPar>; \
|
template class PreconditionerFactory<OpBSeq<Dim>, CommPar>; \
|
||||||
template class PreconditionerFactory<OpFPar<Dim>, CommPar>; \
|
template class PreconditionerFactory<OpFPar<Dim>, CommPar>; \
|
||||||
template class PreconditionerFactory<OpBPar<Dim>, CommPar>; \
|
template class PreconditionerFactory<OpBPar<Dim>, CommPar>; \
|
||||||
|
template class PreconditionerFactory<OpGLFPar<Dim>,CommPar>; \
|
||||||
|
template class PreconditionerFactory<OpGLBPar<Dim>,CommPar>; \
|
||||||
template class PreconditionerFactory<OpW<Dim, false>, CommPar>; \
|
template class PreconditionerFactory<OpW<Dim, false>, CommPar>; \
|
||||||
template class PreconditionerFactory<OpWG<Dim, true>, CommPar>; \
|
template class PreconditionerFactory<OpWG<Dim, true>, CommPar>; \
|
||||||
template class PreconditionerFactory<OpBPar<Dim>, CommSeq>;
|
template class PreconditionerFactory<OpBPar<Dim>,CommSeq>; \
|
||||||
|
template class PreconditionerFactory<OpGLBPar<Dim>,CommSeq>;
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
#define INSTANCE_PF_SEQ(Dim) \
|
#define INSTANCE_PF_SEQ(Dim) \
|
||||||
|
@ -27,6 +27,7 @@
|
|||||||
#include <opm/simulators/linalg/twolevelmethodcpr.hh>
|
#include <opm/simulators/linalg/twolevelmethodcpr.hh>
|
||||||
|
|
||||||
#include <dune/istl/paamg/pinfo.hh>
|
#include <dune/istl/paamg/pinfo.hh>
|
||||||
|
#include <opm/simulators/linalg/WellOperators.hpp>
|
||||||
|
|
||||||
#include <cstddef>
|
#include <cstddef>
|
||||||
|
|
||||||
@ -83,10 +84,10 @@ namespace Opm
|
|||||||
PressureVectorType<Scalar>>;
|
PressureVectorType<Scalar>>;
|
||||||
template<class Scalar, class Comm>
|
template<class Scalar, class Comm>
|
||||||
using ParCoarseOperatorType
|
using ParCoarseOperatorType
|
||||||
= Dune::OverlappingSchwarzOperator<PressureMatrixType<Scalar>,
|
= Opm::GhostLastMatrixAdapter<PressureMatrixType<Scalar>,
|
||||||
PressureVectorType<Scalar>,
|
PressureVectorType<Scalar>,
|
||||||
PressureVectorType<Scalar>,
|
PressureVectorType<Scalar>,
|
||||||
Comm>;
|
Comm>;
|
||||||
template<class Scalar, class Comm>
|
template<class Scalar, class Comm>
|
||||||
using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
|
using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
|
||||||
SeqCoarseOperatorType<Scalar>,
|
SeqCoarseOperatorType<Scalar>,
|
||||||
|
@ -25,6 +25,7 @@
|
|||||||
#include <opm/simulators/linalg/twolevelmethodcpr.hh>
|
#include <opm/simulators/linalg/twolevelmethodcpr.hh>
|
||||||
#include <opm/simulators/linalg/PropertyTree.hpp>
|
#include <opm/simulators/linalg/PropertyTree.hpp>
|
||||||
#include <opm/simulators/linalg/matrixblock.hh>
|
#include <opm/simulators/linalg/matrixblock.hh>
|
||||||
|
#include <opm/simulators/linalg/WellOperators.hpp>
|
||||||
|
|
||||||
#include <cstddef>
|
#include <cstddef>
|
||||||
|
|
||||||
@ -36,10 +37,10 @@ namespace Opm { namespace Details {
|
|||||||
PressureVectorType<Scalar>>;
|
PressureVectorType<Scalar>>;
|
||||||
template<class Scalar, class Comm>
|
template<class Scalar, class Comm>
|
||||||
using ParCoarseOperatorType
|
using ParCoarseOperatorType
|
||||||
= Dune::OverlappingSchwarzOperator<PressureMatrixType<Scalar>,
|
= Opm::GhostLastMatrixAdapter<PressureMatrixType<Scalar>,
|
||||||
PressureVectorType<Scalar>,
|
PressureVectorType<Scalar>,
|
||||||
PressureVectorType<Scalar>,
|
PressureVectorType<Scalar>,
|
||||||
Comm>;
|
Comm>;
|
||||||
template<class Scalar, class Comm>
|
template<class Scalar, class Comm>
|
||||||
using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
|
using CoarseOperatorType = std::conditional_t<std::is_same<Comm, Dune::Amg::SequentialInformation>::value,
|
||||||
SeqCoarseOperatorType<Scalar>,
|
SeqCoarseOperatorType<Scalar>,
|
||||||
|
@ -29,6 +29,8 @@
|
|||||||
#include <opm/common/TimingMacros.hpp>
|
#include <opm/common/TimingMacros.hpp>
|
||||||
|
|
||||||
#include <opm/simulators/linalg/matrixblock.hh>
|
#include <opm/simulators/linalg/matrixblock.hh>
|
||||||
|
#include <dune/common/shared_ptr.hh>
|
||||||
|
#include <dune/istl/paamg/smoother.hh>
|
||||||
|
|
||||||
#include <cstddef>
|
#include <cstddef>
|
||||||
|
|
||||||
@ -325,7 +327,128 @@ protected:
|
|||||||
std::size_t interiorSize_;
|
std::size_t interiorSize_;
|
||||||
};
|
};
|
||||||
|
|
||||||
|
/*!
|
||||||
|
\brief Dune linear operator that assumes ghost rows are ordered after
|
||||||
|
interior rows. Avoids some computations because of this.
|
||||||
|
|
||||||
|
This is similar to WellModelGhostLastMatrixAdapter, with the difference that
|
||||||
|
here we do not have a well model, and also do calcilate the interiorSize
|
||||||
|
using the parallel index set. Created for use in AMG/CPR smoothers.
|
||||||
|
*/
|
||||||
|
template<class M, class X, class Y, class C>
|
||||||
|
class GhostLastMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
|
||||||
|
{
|
||||||
|
public:
|
||||||
|
typedef M matrix_type;
|
||||||
|
typedef X domain_type;
|
||||||
|
typedef Y range_type;
|
||||||
|
typedef typename X::field_type field_type;
|
||||||
|
|
||||||
|
|
||||||
|
typedef C communication_type;
|
||||||
|
|
||||||
|
Dune::SolverCategory::Category category() const override
|
||||||
|
{
|
||||||
|
return Dune::SolverCategory::overlapping;
|
||||||
|
}
|
||||||
|
|
||||||
|
//! constructor: just store a reference to a matrix
|
||||||
|
GhostLastMatrixAdapter (const M& A,
|
||||||
|
const communication_type& comm)
|
||||||
|
: A_( Dune::stackobject_to_shared_ptr(A) ), comm_(comm)
|
||||||
|
{
|
||||||
|
interiorSize_ = setInteriorSize(comm_);
|
||||||
|
}
|
||||||
|
|
||||||
|
GhostLastMatrixAdapter (const std::shared_ptr<M> A,
|
||||||
|
const communication_type& comm)
|
||||||
|
: A_( A ), comm_(comm)
|
||||||
|
{
|
||||||
|
interiorSize_ = setInteriorSize(comm_);
|
||||||
|
}
|
||||||
|
|
||||||
|
virtual void apply( const X& x, Y& y ) const override
|
||||||
|
{
|
||||||
|
for (auto row = A_->begin(); row.index() < interiorSize_; ++row)
|
||||||
|
{
|
||||||
|
y[row.index()]=0;
|
||||||
|
auto endc = (*row).end();
|
||||||
|
for (auto col = (*row).begin(); col != endc; ++col)
|
||||||
|
(*col).umv(x[col.index()], y[row.index()]);
|
||||||
|
}
|
||||||
|
|
||||||
|
ghostLastProject( y );
|
||||||
|
}
|
||||||
|
|
||||||
|
// y += \alpha * A * x
|
||||||
|
virtual void applyscaleadd (field_type alpha, const X& x, Y& y) const override
|
||||||
|
{
|
||||||
|
for (auto row = A_->begin(); row.index() < interiorSize_; ++row)
|
||||||
|
{
|
||||||
|
auto endc = (*row).end();
|
||||||
|
for (auto col = (*row).begin(); col != endc; ++col)
|
||||||
|
(*col).usmv(alpha, x[col.index()], y[row.index()]);
|
||||||
|
}
|
||||||
|
|
||||||
|
ghostLastProject( y );
|
||||||
|
}
|
||||||
|
|
||||||
|
virtual const matrix_type& getmat() const override { return *A_; }
|
||||||
|
|
||||||
|
size_t getInteriorSize() const { return interiorSize_;}
|
||||||
|
|
||||||
|
private:
|
||||||
|
void ghostLastProject(Y& y) const
|
||||||
|
{
|
||||||
|
size_t end = y.size();
|
||||||
|
for (size_t i = interiorSize_; i < end; ++i)
|
||||||
|
y[i] = 0; //project to interiorsize, i.e. ignore ghost
|
||||||
|
}
|
||||||
|
|
||||||
|
size_t setInteriorSize(const communication_type& comm) const
|
||||||
|
{
|
||||||
|
auto indexSet = comm.indexSet();
|
||||||
|
|
||||||
|
size_t is = 0;
|
||||||
|
// Loop over index set
|
||||||
|
for (auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx) {
|
||||||
|
//Only take "owner" indices
|
||||||
|
if (idx->local().attribute()==1) {
|
||||||
|
//get local index
|
||||||
|
auto loc = idx->local().local();
|
||||||
|
// if loc is higher than "old interior size", update it
|
||||||
|
if (loc > is) {
|
||||||
|
is = loc;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
return is + 1; //size is plus 1 since we start at 0
|
||||||
|
}
|
||||||
|
const std::shared_ptr<const matrix_type> A_ ;
|
||||||
|
const communication_type& comm_;
|
||||||
|
size_t interiorSize_;
|
||||||
|
};
|
||||||
|
|
||||||
} // namespace Opm
|
} // namespace Opm
|
||||||
|
|
||||||
|
namespace Dune {
|
||||||
|
namespace Amg {
|
||||||
|
|
||||||
|
template<class M, class X, class Y, class C>
|
||||||
|
class ConstructionTraits<Opm::GhostLastMatrixAdapter<M,X,Y,C> >
|
||||||
|
{
|
||||||
|
public:
|
||||||
|
typedef ParallelOperatorArgs<M,C> Arguments;
|
||||||
|
|
||||||
|
static inline std::shared_ptr<Opm::GhostLastMatrixAdapter<M,X,Y,C>> construct(const Arguments& args)
|
||||||
|
{
|
||||||
|
return std::make_shared<Opm::GhostLastMatrixAdapter<M,X,Y,C>>
|
||||||
|
(args.matrix_, args.comm_);
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
} // end namespace Amg
|
||||||
|
} // end namespace Dune
|
||||||
|
|
||||||
#endif // OPM_WELLOPERATORS_HEADER_INCLUDED
|
#endif // OPM_WELLOPERATORS_HEADER_INCLUDED
|
||||||
|
|
||||||
|
Loading…
Reference in New Issue
Block a user