mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Added CUILU0 to the PreconditionerFactory.
This commit is contained in:
parent
dfc5331696
commit
dea49a5406
@ -40,6 +40,14 @@
|
||||
#include <dune/istl/paamg/kamg.hh>
|
||||
#include <dune/istl/paamg/fastamg.hh>
|
||||
|
||||
#include <config.h>
|
||||
#if HAVE_CUDA
|
||||
#include <opm/simulators/linalg/cuistl/CuSeqILU0.hpp>
|
||||
#include <opm/simulators/linalg/cuistl/PreconditionerAdapter.hpp>
|
||||
#include <opm/simulators/linalg/cuistl/PreconditionerConvertFieldTypeAdapter.hpp>
|
||||
#endif
|
||||
|
||||
|
||||
namespace Opm {
|
||||
|
||||
template<class Smoother>
|
||||
@ -398,6 +406,30 @@ struct StandardPreconditioners<Operator,Dune::Amg::SequentialInformation>
|
||||
using LevelTransferPolicy = Opm::PressureTransferPolicy<O, Dune::Amg::SequentialInformation, true>;
|
||||
return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(op, prm, weightsCalculator, pressureIndex);
|
||||
});
|
||||
|
||||
#if HAVE_CUDA
|
||||
F::addCreator("CUILU0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
|
||||
const double w = prm.get<double>("relaxation", 1.0);
|
||||
using field_type = typename V::field_type;
|
||||
using CuILU0 = typename Opm::cuistl::CuSeqILU0<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
|
||||
return std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CuILU0>>(std::make_shared<CuILU0>(op.getmat(), w));
|
||||
});
|
||||
|
||||
F::addCreator("CUILU0Float", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
|
||||
const double w = prm.get<double>("relaxation", 1.0);
|
||||
using block_type = typename V::block_type;
|
||||
using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
|
||||
using matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
|
||||
using CuILU0 = typename Opm::cuistl::CuSeqILU0<matrix_type_to, Opm::cuistl::CuVector<float>, Opm::cuistl::CuVector<float>>;
|
||||
using Adapter = typename Opm::cuistl::PreconditionerAdapter<VTo, VTo, CuILU0>;
|
||||
using Converter = typename Opm::cuistl::PreconditionerConvertFieldTypeAdapter<Adapter, M, V, V>;
|
||||
auto converted = std::make_shared<Converter>(op.getmat());
|
||||
auto adapted = std::make_shared<Adapter>(std::make_shared<CuILU0>(converted->getConvertedMatrix(), w));
|
||||
converted->setUnderlyingPreconditioner(adapted);
|
||||
return converted;
|
||||
|
||||
});
|
||||
#endif
|
||||
}
|
||||
};
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user