Skip to content
Snippets Groups Projects
Commit dcdcdc5e authored by Santiago Ospina De Los Ríos's avatar Santiago Ospina De Los Ríos
Browse files

Add newton parameters to the time stepper

* This change requires the parameter tree to contain new parameters
* Updates the tests to always include a default set of newton parameters
parent 3c18ba75
No related branches found
No related tags found
1 merge request!40Resolve "strange simulation results"
Pipeline #30580 passed
......@@ -183,10 +183,13 @@ public:
* @brief Construct a new object
*
* @param rk_method PDELab time stepping Runge-Kutta paramameters
* @param solver_parameters parameter options for the newton solver
*/
RKStepper(
std::unique_ptr<PDELab::TimeSteppingParameterInterface<Time>> rk_method)
std::unique_ptr<PDELab::TimeSteppingParameterInterface<Time>> rk_method,
const ParameterTree& solver_parameters)
: _rk_method{ std::move(rk_method) }
, _solver_parameters{ solver_parameters }
, _logger{ Logging::Logging::componentLogger({}, "stepper") }
{
_logger.detail("Setting up time stepper"_fmt);
......@@ -197,9 +200,10 @@ public:
* @brief Construct a new object
*
* @param rk_method Keyword with the RK method type
* @param solver_parameters parameter options for the newton solver
*/
RKStepper( std::string rk_method)
: RKStepper(make_rk_method(rk_method))
RKStepper(std::string rk_method, const ParameterTree& solver_parameters)
: RKStepper(make_rk_method(rk_method), solver_parameters)
{}
/**
......@@ -294,7 +298,7 @@ public:
const Logging::Logger& logger() const { return _logger; }
private:
//! Setup and cache solver for a system
//! Setup and cache solver for a given system
template<class System>
auto& get_solver(const System& system) const
{
......@@ -329,6 +333,51 @@ private:
auto non_linear_operator =
std::make_unique<NonLinearOperator>(grid_operator, *linear_solver);
// Add settings to the newton solver
auto reduction = _solver_parameters.template get<double>("reduction");
non_linear_operator->setReduction(reduction);
auto min_reduction =
_solver_parameters.template get<double>("min_linear_reduction");
non_linear_operator->setMinLinearReduction(min_reduction);
auto fixed_reduction =
_solver_parameters.template get<bool>("fixed_linear_reduction");
non_linear_operator->setFixedLinearReduction(fixed_reduction);
auto max_it = _solver_parameters.template get<uint>("max_iterations");
non_linear_operator->setMaxIterations(max_it);
auto abs_limit = _solver_parameters.template get<double>("absolute_limit");
non_linear_operator->setAbsoluteLimit(abs_limit);
auto reassemble =
_solver_parameters.template get<double>("reassemble_threshold");
non_linear_operator->setReassembleThreshold(reassemble);
auto keep_matrix = _solver_parameters.template get<bool>("keep_matrix");
non_linear_operator->setKeepMatrix(keep_matrix);
auto force_iteration =
_solver_parameters.template get<bool>("force_iteration");
non_linear_operator->setForceIteration(force_iteration);
const auto& ls = _solver_parameters.sub("linear_search", true);
auto ls_strategy = ls.template get<std::string>("strategy");
try {
non_linear_operator->setLineSearchStrategy(ls_strategy);
} catch (const Dune::Exception& e) {
_logger.error("Not valid linear search strategy: {}"_fmt, ls_strategy);
DUNE_THROW(IOError, "Not valid linear search strategy: " << ls_strategy);
}
auto ls_max_it = ls.template get<uint>("max_iterations");
non_linear_operator->setLineSearchMaxIterations(ls_max_it);
auto ls_damping = ls.template get<double>("damping_factor");
non_linear_operator->setLineSearchDampingFactor(ls_damping);
_logger.trace("Get one step operator"_fmt);
auto one_step_operator = std::make_unique<OneStepOperator>(
*_rk_method, grid_operator, *non_linear_operator);
......@@ -346,6 +395,7 @@ private:
private:
std::unique_ptr<const PDELab::TimeSteppingParameterInterface<Time>>
_rk_method;
const ParameterTree _solver_parameters;
const Logging::Logger _logger;
mutable std::any _internal_state;
};
......@@ -459,8 +509,10 @@ make_default_stepper(const ParameterTree& config)
log.trace("Decrease factor: {}"_fmt,decrease_factor);
log.trace("Runge-Kutta method: {}"_fmt,rk_type);
const auto& newton_parameters = config.sub("newton", true);
using RKStepper = Dune::Copasi::RKStepper<double>;
auto rk_stepper = RKStepper{ rk_type };
auto rk_stepper = RKStepper{ rk_type, newton_parameters };
using Stepper = Dune::Copasi::SimpleAdaptiveStepper<RKStepper>;
return Stepper{
std::move(rk_stepper), min_step, max_step, decrease_factor, increase_factor
......
File moved
######################## Time stepping settings #############################
# Common options for the time stepping section including its newton solver
[model.time_stepping]
rk_method = alexander_2
decrease_factor = 0.5
increase_factor = 1.5
[model.time_stepping.newton]
reduction = 1e-8
min_linear_reduction = 1e-3
fixed_linear_reduction = false
max_iterations = 40
absolute_limit = 1e-12
reassemble_threshold = 0.0
keep_matrix = true
force_iteration = false
[model.time_stepping.newton.linear_search]
strategy = hackbuschReusken
max_iterations = 10
damping_factor = 0.5
# Test: Check that complicated problems have an approximately same analytic and
# numeric jacobian
import logging.ini
import default_logging.ini
import default_time_stepping.ini
__name = test_NFkappaB
......@@ -14,14 +15,11 @@ initial_level = 0
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = 0.
end = 365
initial_step = 1e-0
min_step = 1e-3
max_step = 1e-0
decrease_factor = 0.9
increase_factor = 1.1
[model.compartments]
cytoplasm = 0
......
# Test: Compare a simple multidomain diffusion-reaction problem with respect
# to a reference solution produced by stable versions of the program
import logging.ini
import default_logging.ini
import default_time_stepping.ini
__name = test_cell
......@@ -14,7 +15,6 @@ initial_level = 0
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = 0.
end = 0.01
step = 2e-3
......
import logging.ini
# Test: Compare a simple exponential ODE with its analytic result
import default_logging.ini
import default_time_stepping.ini
_grow_rate = 0.2
_initial = 1
......@@ -14,14 +17,11 @@ initial_level = 0
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = 0
end = 10
initial_step = 1e-2
min_step = 1e-3
max_step = 1e-1
decrease_factor = 0.5
increase_factor = 1.5
[model.writer]
file_path = {__name}
......
# Test: Check that the first steps of a pure diffusive problem evolve as a
# gaussian process
import logging.ini
import default_logging.ini
import default_time_stepping.ini
_diffusion = 0.005
_t0 = 1.
......@@ -19,14 +20,11 @@ initial_level = 4
order = 0, 1 | expand geometry_type
[model.time_stepping]
rk_method = alexander_2
begin = {_t0}
end = 2
initial_step = 1e-2
min_step = 1e-3
max_step = 1e-1
decrease_factor = 0.5
increase_factor = 1.5
[model.writer]
file_path = {__name}
......
# Test: Check that the first steps of a pure diffusive problem evolve as a
# gaussian process (mixed geometry version)
import default_logging.ini
import default_time_stepping.ini
__name = test_gauss_mixed_geometry
_diffusion = 0.005
......@@ -16,14 +19,11 @@ initial_level = 4
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = {_t0}
end = 2
initial_step = 1e-2
min_step = 1e-3
max_step = 1e-1
decrease_factor = 0.5
increase_factor = 1.5
[model.writer]
file_path = {__name}
......
# Test that checks fluxes on membranes with a reference vtk file produced with
# previous versions of dune_copasi_md
import logging.ini
import default_logging.ini
import default_time_stepping.ini
__name = test_in_out
......@@ -14,15 +15,11 @@ initial_level = 1
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = 0
end = 1
initial_step = 0.3
min_step = 1e-3
max_step = 0.35
decrease_factor = 0.5
increase_factor = 1.5
[model.writer]
file_path = {__name}
......
# Test: Check the interpolation output for math expression initialization
import logging.ini
import default_logging.ini
import default_time_stepping.ini
_diffusion = 0.01
_t0 = 1.
......@@ -18,15 +19,11 @@ initial_level = 6
order = 0, 1 | expand geometry_type
[model.time_stepping]
rk_method = alexander_2
begin = {_t0}
end = {_t0}
step = 0.1
initial_step = 0.1
min_step = 1e-3
max_step = 1e-0
decrease_factor = 0.9
increase_factor = 1.1
[model.compartments]
domain = 0
......
# Test: An inifile without `writer` section must be valid
import logging.ini
import default_logging.ini
import default_time_stepping.ini
__name = test_no_writer
......@@ -13,14 +14,11 @@ initial_level = 0
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = 0
end = 1
initial_step = 0.1
min_step = 1e-3
max_step = 1e-0
decrease_factor = 0.9
increase_factor = 1.1
[model.compartments]
domain = 0
......
# Test: Checks interpolation of TIFF files with a reference vtk file
# produced with previous versions of dune_copasi_sd/md
import logging.ini
import default_logging.ini
import default_time_stepping.ini
__name = test_tiff
......@@ -14,15 +15,11 @@ initial_level = 3
order = 1
[model.time_stepping]
rk_method = alexander_2
begin = 0.
end = 0.
step = 0.1
initial_step = 0.1
min_step = 1e-3
max_step = 1e-0
decrease_factor = 0.9
increase_factor = 1.1
[model.writer]
file_path = {__name}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment