mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
move NonLinearSolver parameters to TypeTag-free parameter system
This commit is contained in:
parent
dd1bc8d75b
commit
94c30a74b8
@ -76,7 +76,7 @@ namespace Opm::Properties {
|
|||||||
namespace TTag {
|
namespace TTag {
|
||||||
struct FlowProblem {
|
struct FlowProblem {
|
||||||
using InheritsFrom = std::tuple<FlowTimeSteppingParameters,
|
using InheritsFrom = std::tuple<FlowTimeSteppingParameters,
|
||||||
FlowNonLinearSolver, FlowBaseProblem, BlackOilModel>;
|
FlowBaseProblem, BlackOilModel>;
|
||||||
};
|
};
|
||||||
}
|
}
|
||||||
// default in flow is to formulate the equations in surface volumes
|
// default in flow is to formulate the equations in surface volumes
|
||||||
|
@ -39,40 +39,13 @@
|
|||||||
|
|
||||||
#include <memory>
|
#include <memory>
|
||||||
|
|
||||||
namespace Opm::Properties::TTag {
|
|
||||||
|
|
||||||
struct FlowNonLinearSolver {};
|
|
||||||
|
|
||||||
}
|
|
||||||
|
|
||||||
namespace Opm::Parameters {
|
namespace Opm::Parameters {
|
||||||
|
|
||||||
template<class TypeTag, class MyTypeTag>
|
template<class Scalar>
|
||||||
struct NewtonMaxRelax { using type = Properties::UndefinedProperty; };
|
struct NewtonMaxRelax { static constexpr Scalar value = 0.5; };
|
||||||
|
|
||||||
// we are reusing NewtonMaxIterations from opm-models
|
struct NewtonMinIterations { static constexpr int value = 2; };
|
||||||
// See opm/models/nonlinear/newtonmethodproperties.hh
|
struct NewtonRelaxationType { static constexpr auto value = "dampen"; };
|
||||||
|
|
||||||
template<class TypeTag, class MyTypeTag>
|
|
||||||
struct NewtonMinIterations { using type = Properties::UndefinedProperty; };
|
|
||||||
|
|
||||||
template<class TypeTag, class MyTypeTag>
|
|
||||||
struct NewtonRelaxationType { using type = Properties::UndefinedProperty; };
|
|
||||||
|
|
||||||
template<class TypeTag>
|
|
||||||
struct NewtonMaxRelax<TypeTag, Properties::TTag::FlowNonLinearSolver>
|
|
||||||
{
|
|
||||||
using type = GetPropType<TypeTag, Properties::Scalar>;
|
|
||||||
static constexpr type value = 0.5;
|
|
||||||
};
|
|
||||||
|
|
||||||
template<class TypeTag>
|
|
||||||
struct NewtonMinIterations<TypeTag, Properties::TTag::FlowNonLinearSolver>
|
|
||||||
{ static constexpr int value = 2; };
|
|
||||||
|
|
||||||
template<class TypeTag>
|
|
||||||
struct NewtonRelaxationType<TypeTag, Properties::TTag::FlowNonLinearSolver>
|
|
||||||
{ static constexpr auto value = "dampen"; };
|
|
||||||
|
|
||||||
} // namespace Opm::Parameters
|
} // namespace Opm::Parameters
|
||||||
|
|
||||||
@ -125,11 +98,11 @@ void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
|
|||||||
reset();
|
reset();
|
||||||
|
|
||||||
// overload with given parameters
|
// overload with given parameters
|
||||||
relaxMax_ = Parameters::get<TypeTag, Parameters::NewtonMaxRelax>();
|
relaxMax_ = Parameters::Get<Parameters::NewtonMaxRelax<Scalar>>();
|
||||||
maxIter_ = Parameters::Get<Parameters::NewtonMaxIterations>();
|
maxIter_ = Parameters::Get<Parameters::NewtonMaxIterations>();
|
||||||
minIter_ = Parameters::get<TypeTag, Parameters::NewtonMinIterations>();
|
minIter_ = Parameters::Get<Parameters::NewtonMinIterations>();
|
||||||
|
|
||||||
const auto& relaxationTypeString = Parameters::get<TypeTag, Parameters::NewtonRelaxationType>();
|
const auto& relaxationTypeString = Parameters::Get<Parameters::NewtonRelaxationType>();
|
||||||
if (relaxationTypeString == "dampen") {
|
if (relaxationTypeString == "dampen") {
|
||||||
relaxType_ = NonlinearRelaxType::Dampen;
|
relaxType_ = NonlinearRelaxType::Dampen;
|
||||||
} else if (relaxationTypeString == "sor") {
|
} else if (relaxationTypeString == "sor") {
|
||||||
@ -142,13 +115,13 @@ void stabilizeNonlinearUpdate(BVector& dx, BVector& dxOld,
|
|||||||
|
|
||||||
static void registerParameters()
|
static void registerParameters()
|
||||||
{
|
{
|
||||||
Parameters::registerParam<TypeTag, Parameters::NewtonMaxRelax>
|
Parameters::Register<Parameters::NewtonMaxRelax<Scalar>>
|
||||||
("The maximum relaxation factor of a Newton iteration");
|
("The maximum relaxation factor of a Newton iteration");
|
||||||
Parameters::Register<Parameters::NewtonMaxIterations>
|
Parameters::Register<Parameters::NewtonMaxIterations>
|
||||||
("The maximum number of Newton iterations per time step");
|
("The maximum number of Newton iterations per time step");
|
||||||
Parameters::registerParam<TypeTag, Parameters::NewtonMinIterations>
|
Parameters::Register<Parameters::NewtonMinIterations>
|
||||||
("The minimum number of Newton iterations per time step");
|
("The minimum number of Newton iterations per time step");
|
||||||
Parameters::registerParam<TypeTag, Parameters::NewtonRelaxationType>
|
Parameters::Register<Parameters::NewtonRelaxationType>
|
||||||
("The type of relaxation used by Newton method");
|
("The type of relaxation used by Newton method");
|
||||||
|
|
||||||
Parameters::SetDefault<Parameters::NewtonMaxIterations>(20);
|
Parameters::SetDefault<Parameters::NewtonMaxIterations>(20);
|
||||||
|
Loading…
Reference in New Issue
Block a user