Thanks to visit codestin.com
Credit goes to github.com

Skip to content

Closed issue #3 - time stepping #4

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 3 commits into from
May 13, 2019
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion examples/diffusion_2d/diffusion_2d.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -212,7 +212,7 @@ int solution_check(dae::state_type &x, MKL_INT N, double t, double D)
<< "% deviation from the analytical value)\n";
std::cout << "Maximum relative error: " << err_max << "%\n";

if(err_max < 1.0 && err_conc < 1.0e-6)
if(err_max < 2.0 && err_conc < 1.0e-6)
return 0;
else
return 1;
Expand Down
11 changes: 4 additions & 7 deletions examples/perovskite/perovskite.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -99,13 +99,10 @@ int main()
// parameters defined in solver_options.h
dae::SolverOptions opt;

#ifdef DAE_SINGLE
opt.atol = 1.0e-3; // Absolute tolerance for single precision
#else
opt.atol = 1.0e-6; // Absolute tolerance for double precision
#endif
opt.fact_every_iter = false; // Gain some speed. The matrices will be
// factorized only once each time step.
opt.dt_increase_factor = 1.4; // Reduce time step increase factor to get
// better accuracy (default value is 2.0)
opt.fact_every_iter = false; // Gain some speed. The matrices will be
// factorized only once each time step.

// Create an instance of the solver with particular RHS, Mass matrix,
// Jacobian and solver options
Expand Down
96 changes: 87 additions & 9 deletions src/solver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -31,7 +31,8 @@ void Solver::operator()(state_type &x)
else if(m_opt.dt_init > m_t1)
{
double new_dt = m_t1 / 10.0;
std::cout << "WARNING: Initial time step is bigger than the integration time t1: "
std::cout << "WARNING: Initial time step is bigger than the "
"integration time t1: "
<< m_opt.dt_init
<< "\n The solver will use dt = t1/10 = " << new_dt
<< std::endl;
Expand Down Expand Up @@ -130,6 +131,7 @@ void Solver::operator()(state_type &x)

// TODO: Start timer here

// Memory control variables
int peak_mem1 = 0, peak_mem2 = 0, peak_mem3 = 0;

bool final_time_step = false;
Expand Down Expand Up @@ -254,6 +256,7 @@ void Solver::operator()(state_type &x)
{
break;
}

} // for iter

// Newton iterator failed to converge within max_Newton_iter iterations
Expand All @@ -268,27 +271,102 @@ void Solver::operator()(state_type &x)
step_counter--;
final_time_step = false;
dt /= m_opt.dt_decrease_factor;
current_scheme = 1;
if(dt < m_opt.dt_min)
{
// This actually means solution error
// TODO: stop the solver
std::cout << "\nERROR: The time step was reduced to " << dt
<< " but the Newton method failed to converge\n";
exit(4);
}
x = x_prev[0];
continue;
}

// The solver reached the target time t1
if(final_time_step)
{
break;
}

// Simple yet efficient adaptive time stepping
if(iter < m_opt.dt_increase_threshold)
if(m_opt.time_stepping == 1) // S-SATS
{
dt *= m_opt.dt_increase_factor;
if(m_opt.verbosity > 0)
std::cout << '>';
if(iter < m_opt.dt_increase_threshold)
{
dt *= m_opt.dt_increase_factor;
if(m_opt.dt_increase_factor != 1.0)
current_scheme = 1;
if(dt > m_opt.dt_max)
dt = m_opt.dt_max;
if(m_opt.verbosity > 0)
std::cout << '>';
}
else if(iter >= m_opt.dt_decrease_threshold - 1)
{
dt /= m_opt.dt_decrease_factor;
if(m_opt.dt_decrease_factor != 1.0)
current_scheme = 1;
if(dt < m_opt.dt_min)
dt = m_opt.dt_min;
if(m_opt.verbosity > 0)
std::cout << '<';
}
}
else if(iter >= m_opt.dt_decrease_threshold - 1)
else if(m_opt.time_stepping == 2) // A-SATS
{
dt /= m_opt.dt_decrease_factor;
if(m_opt.verbosity > 0)
std::cout << '<';
double norm1 = 0.0;
double norm2 = 0.0;

// Estimate NORM(C(n+1) - C(n)) and NORM(C(n))
for(MKL_INT i = 0; i < size; i++)
{
norm1 += (x[i] - x_prev[0][i]) * (x[i] - x_prev[0][i]);
norm2 += x_prev[0][i] * x_prev[0][i];
}
norm1 = sqrt(norm1);
norm2 = sqrt(norm2);

// Monitor function
double eta = norm1 / (norm2 + m_opt.dt_eps_m);

// The time step can be increased
if(eta < m_opt.dt_eta_min)
{
dt *= m_opt.dt_increase_factor;
if(m_opt.dt_increase_factor != 1.0)
current_scheme = 1;
if(dt > m_opt.dt_max)
dt = m_opt.dt_max;
if(m_opt.verbosity > 0)
std::cout << '>';
}

// The time step should be reduced, scrape the current time
// iteration
if(eta > m_opt.dt_eta_max)
{
if(m_opt.verbosity > 0)
std::cout << " <- redo: eta = " << eta;

t -= dt;
step_counter--;
final_time_step = false;
dt /= m_opt.dt_decrease_factor;
current_scheme = 1;
if(dt < m_opt.dt_min)
{
// This actually means solution error
// TODO: stop the solver
std::cout << "\nERROR: The time step was reduced to " << dt
<< " but the relative error is still above the "
"threshold\n";
exit(5);
}
x = x_prev[0];
continue;
}
}

if(t + dt > m_t1)
Expand Down
7 changes: 7 additions & 0 deletions src/solver_options.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -95,6 +95,13 @@ void SolverOptions::check_options()
bdf_order = 1;
}

if(time_stepping < 1 || time_stepping > 2)
{
// TODO: print warning
// fall back to S-SATS
time_stepping = 1;
}

// TODO: add more checks
}

Expand Down
34 changes: 28 additions & 6 deletions src/solver_options.h
Original file line number Diff line number Diff line change
Expand Up @@ -32,30 +32,52 @@ class SolverOptions
// 1 - first order BDF, 2 - BDF-2, ..., 6 - BDF-6
int bdf_order = 6;

// Time stepping algorithm:
// 1 - Stability-based Simple Adaptive Time Stepping (S-SATS),
// 2 - Accuracy-based Simple Adaptive Time Stepping (A-SATS) from
// https://www.sciencedirect.com/science/article/pii/S0377042705005534
int time_stepping = 1;

// Maximum number of Newton iterations. If the Newton method fails to
// converge after max_Newton_iter iterations, the solver reduces time step
// and tries to make the current step again.
int max_Newton_iter = 15;

// Absolute tolerance
// Absolute tolerance for the Newton algorithm
#ifdef DAE_SINGLE
double atol = 1.0e-3; // Absolute tolerance for single precision
#else
double atol = 1.0e-6; // Absolute tolerance for double precision
#endif

#ifdef DAE_SINGLE
double dt_eps_m =
1.0e-5; // The order of the rounding unit for single precision
#else
double dt_eps_m =
1.0e-10; // The order of the rounding unit for double precision
#endif

// Initial time step
double dt_init = 0.1;

// Minimum and maximum time steps
double dt_min = dt_eps_m;
double dt_max = 100.0;

// Verbosity level of the solver:
// 0 - be silent, 1 - prints some basic information, 2 - chatterbox
int verbosity = 1;

// Adaptive time stepping options
int dt_increase_threshold = 3;
int dt_decrease_threshold = 7;
double dt_increase_factor = 1.4;
double dt_decrease_factor = 1.4;
// Simple Adaptive Time Stepping options
int dt_increase_threshold = 3; // Time step amplification threshold
// (S-SATS only)
int dt_decrease_threshold = 7; // Time stepreduction threshold
// (S-SATS only)
double dt_increase_factor = 2.0; // Time step amplification factor
double dt_decrease_factor = 2.0; // Time step reduction factor
double dt_eta_min = 0.05; // Monitor function lower threshold (A-SATS only)
double dt_eta_max = 0.5; // Monitor function higher threshold (A-SATS only)

// Intel MKL PARDISO parameters (iparam). More about iparam:
// https://software.intel.com/en-us/mkl-developer-reference-c-pardiso-iparm-parameter
Expand Down