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

Skip to content

Commit cf1f762

Browse files
committed
BDF-2 integrator now fully supports variable time stepping. Now defaults are BDF-2 and A-SATS.
1 parent 8a69c0a commit cf1f762

File tree

6 files changed

+93
-68
lines changed

6 files changed

+93
-68
lines changed

examples/diffusion_2d/diffusion_2d.cpp

Lines changed: 4 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -80,10 +80,9 @@ int main()
8080
// parameters defined in solver_options.h
8181
dae::SolverOptions opt;
8282

83-
opt.dt_init = 0.0001; // Change initial time step
84-
opt.dt_increase_factor = 1.0; // Do not increase the time step
85-
opt.fact_every_iter = false; // Gain some speed. The matrices will be
86-
// factorized only once each time step.
83+
opt.dt_init = 5.0e-5; // Change initial time step
84+
opt.fact_every_iter = false; // Gain some speed. The matrices will be
85+
// factorized only once each time step.
8786

8887
// We can override Jacobian class from dae-cpp library and provide
8988
// analytical Jacobian. But we will use numerically estimated one.
@@ -212,7 +211,7 @@ int solution_check(dae::state_type &x, MKL_INT N, double t, double D)
212211
<< "% deviation from the analytical value)\n";
213212
std::cout << "Maximum relative error: " << err_max << "%\n";
214213

215-
if(err_max < 2.0 && err_conc < 1.0e-6)
214+
if(err_max < 1.0 && err_conc < 1.0e-6)
216215
return 0;
217216
else
218217
return 1;

examples/perovskite/perovskite.cpp

Lines changed: 3 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -99,10 +99,9 @@ int main()
9999
// parameters defined in solver_options.h
100100
dae::SolverOptions opt;
101101

102-
opt.dt_increase_factor = 1.4; // Reduce time step increase factor to get
103-
// better accuracy (default value is 2.0)
104-
opt.fact_every_iter = false; // Gain some speed. The matrices will be
105-
// factorized only once each time step.
102+
opt.time_stepping = 1; // Choose Stability-based time stepper
103+
opt.fact_every_iter = false; // Gain some speed. The matrices will be
104+
// factorized only once each time step.
106105

107106
// Create an instance of the solver with particular RHS, Mass matrix,
108107
// Jacobian and solver options

src/solver.cpp

Lines changed: 59 additions & 44 deletions
Original file line numberDiff line numberDiff line change
@@ -42,9 +42,13 @@ void Solver::operator()(state_type &x)
4242
// Initialise time integrator
4343
TimeIntegrator ti(m_rhs, m_jac, m_mass, m_opt, size);
4444

45-
// Initial time and time step
45+
// Initial time
4646
double t = 0.0;
47-
double dt = m_opt.dt_init;
47+
48+
// Time step
49+
double dt[2];
50+
dt[0] = m_opt.dt_init;
51+
dt[1] = 0.0;
4852

4953
// Initial output
5054
if(m_opt.verbosity > 0)
@@ -118,6 +122,9 @@ void Solver::operator()(state_type &x)
118122
double ddum; // Double dummy
119123
MKL_INT idum; // Integer dummy
120124

125+
// Memory control variables
126+
int peak_mem1 = 0, peak_mem2 = 0, peak_mem3 = 0;
127+
121128
// Intel MKL PARDISO iparm parameter
122129
MKL_INT iparm[64];
123130

@@ -131,16 +138,13 @@ void Solver::operator()(state_type &x)
131138

132139
// TODO: Start timer here
133140

134-
// Memory control variables
135-
int peak_mem1 = 0, peak_mem2 = 0, peak_mem3 = 0;
141+
t += dt[0];
136142

137143
bool final_time_step = false;
138144
int step_counter = 0;
139145

140-
while(t < (m_t1 + dt * 0.5))
146+
while(t < (m_t1 + dt[0] * 0.5))
141147
{
142-
t += dt; // Time step lapse
143-
144148
step_counter++;
145149

146150
if(m_opt.verbosity > 0)
@@ -267,24 +271,25 @@ void Solver::operator()(state_type &x)
267271

268272
// Decrease the time step, scrape the current time iteration and
269273
// carry out it again.
270-
t -= dt;
274+
t -= dt[0];
271275
step_counter--;
272276
final_time_step = false;
273-
dt /= m_opt.dt_decrease_factor;
274-
current_scheme = 1;
275-
if(dt < m_opt.dt_min)
277+
dt[0] /= m_opt.dt_decrease_factor;
278+
current_scheme = 1; // Fall back to BDF-1 for better stability
279+
if(dt[0] < m_opt.dt_min)
276280
{
277281
// This actually means solution error
278-
// TODO: stop the solver
279-
std::cout << "\nERROR: The time step was reduced to " << dt
282+
// TODO: Correctly stop the solver
283+
std::cout << "\nERROR: The time step was reduced to " << dt[0]
280284
<< " but the Newton method failed to converge\n";
281285
exit(4);
282286
}
283287
x = x_prev[0];
288+
t += dt[0];
284289
continue;
285290
}
286291

287-
// The solver reached the target time t1
292+
// The solver has reached the target time t1
288293
if(final_time_step)
289294
{
290295
break;
@@ -293,23 +298,25 @@ void Solver::operator()(state_type &x)
293298
// Simple yet efficient adaptive time stepping
294299
if(m_opt.time_stepping == 1) // S-SATS
295300
{
301+
dt[1] = dt[0];
302+
296303
if(iter < m_opt.dt_increase_threshold)
297304
{
298-
dt *= m_opt.dt_increase_factor;
305+
dt[0] *= m_opt.dt_increase_factor;
299306
if(m_opt.dt_increase_factor != 1.0)
300-
current_scheme = 1;
301-
if(dt > m_opt.dt_max)
302-
dt = m_opt.dt_max;
307+
current_scheme = 2;
308+
if(dt[0] > m_opt.dt_max)
309+
dt[0] = m_opt.dt_max;
303310
if(m_opt.verbosity > 0)
304311
std::cout << '>';
305312
}
306313
else if(iter >= m_opt.dt_decrease_threshold - 1)
307314
{
308-
dt /= m_opt.dt_decrease_factor;
315+
dt[0] /= m_opt.dt_decrease_factor;
309316
if(m_opt.dt_decrease_factor != 1.0)
310-
current_scheme = 1;
311-
if(dt < m_opt.dt_min)
312-
dt = m_opt.dt_min;
317+
current_scheme = 2;
318+
if(dt[0] < m_opt.dt_min)
319+
dt[0] = m_opt.dt_min;
313320
if(m_opt.verbosity > 0)
314321
std::cout << '<';
315322
}
@@ -331,49 +338,55 @@ void Solver::operator()(state_type &x)
331338
// Monitor function
332339
double eta = norm1 / (norm2 + m_opt.dt_eps_m);
333340

334-
// The time step can be increased
335-
if(eta < m_opt.dt_eta_min)
336-
{
337-
dt *= m_opt.dt_increase_factor;
338-
if(m_opt.dt_increase_factor != 1.0)
339-
current_scheme = 1;
340-
if(dt > m_opt.dt_max)
341-
dt = m_opt.dt_max;
342-
if(m_opt.verbosity > 0)
343-
std::cout << '>';
344-
}
345-
346341
// The time step should be reduced, scrape the current time
347342
// iteration
348343
if(eta > m_opt.dt_eta_max)
349344
{
350345
if(m_opt.verbosity > 0)
351346
std::cout << " <- redo: eta = " << eta;
352347

353-
t -= dt;
348+
t -= dt[0];
354349
step_counter--;
355350
final_time_step = false;
356-
dt /= m_opt.dt_decrease_factor;
357-
current_scheme = 1;
358-
if(dt < m_opt.dt_min)
351+
dt[0] /= m_opt.dt_decrease_factor;
352+
if(step_counter && m_opt.bdf_order == 2)
353+
current_scheme = 2;
354+
else
355+
current_scheme = 1;
356+
if(dt[0] < m_opt.dt_min)
359357
{
360358
// This actually means solution error
361-
// TODO: stop the solver
362-
std::cout << "\nERROR: The time step was reduced to " << dt
363-
<< " but the relative error is still above the "
364-
"threshold\n";
359+
// TODO: Correctly stop the solver
360+
std::cout << "\nERROR: The time step was reduced to "
361+
<< dt[0] << " but the relative error is still above the threshold\n";
365362
exit(5);
366363
}
367364
x = x_prev[0];
365+
t += dt[0];
368366
continue;
369367
}
368+
369+
dt[1] = dt[0];
370+
371+
// The time step can be increased
372+
if(eta < m_opt.dt_eta_min)
373+
{
374+
dt[0] *= m_opt.dt_increase_factor;
375+
if(m_opt.dt_increase_factor != 1.0)
376+
current_scheme = 2;
377+
if(dt[0] > m_opt.dt_max)
378+
dt[0] = m_opt.dt_max;
379+
if(m_opt.verbosity > 0)
380+
std::cout << '>';
381+
}
370382
}
371383

372-
if(t + dt > m_t1)
384+
// Looks like the solver has reached the target time t1
385+
if(t + dt[0] > m_t1)
373386
{
374387
final_time_step = true;
375388
// Adjust the last time step size
376-
dt = m_t1 - t;
389+
dt[0] = m_t1 - t;
377390
}
378391

379392
// Rewrite solution history
@@ -383,6 +396,8 @@ void Solver::operator()(state_type &x)
383396
}
384397
x_prev[0] = x;
385398

399+
t += dt[0]; // Time step lapse
400+
386401
} // while t
387402

388403
if(m_opt.verbosity > 0)

src/solver_options.h

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -30,13 +30,14 @@ class SolverOptions
3030

3131
// Order of BDF implicit numerical integration method:
3232
// 1 - first order BDF, 2 - BDF-2, ..., 6 - BDF-6
33-
int bdf_order = 6;
33+
// Default is BDF-2 since it fully supports variable time stepping
34+
int bdf_order = 2;
3435

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

4142
// Maximum number of Newton iterations. If the Newton method fails to
4243
// converge after max_Newton_iter iterations, the solver reduces time step

src/time_integrator.cpp

Lines changed: 23 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -71,30 +71,41 @@ TimeIntegrator::TimeIntegrator(RHS &rhs, Jacobian &jac, MassMatrix &mass,
7171
void TimeIntegrator::operator()(sparse_matrix_holder &Jt, state_type &b,
7272
sparse_matrix_holder &J, state_type &x,
7373
const state_type_matrix &x_prev, const double t,
74-
const double dt)
74+
const double dt[])
7575
{
7676
const MKL_INT size = (MKL_INT)(x.size());
7777

78-
const double invdt = 1.0 / dt;
79-
8078
double alpha = 0;
8179
int scheme = m_scheme - 1;
8280

8381
state_type dxdt(size);
8482

85-
for(MKL_INT i = 0; i < size; i++)
83+
if(m_scheme == 2 && m_opt.bdf_order == 2) // Variable time stepper for BDF-2
8684
{
87-
double val = 0;
88-
val += BDF_COEF[scheme][0] * x[i];
89-
for(int d = 1; d <= m_scheme; d++)
85+
for(MKL_INT i = 0; i < size; i++)
9086
{
91-
val += BDF_COEF[scheme][d] * x_prev[d - 1][i];
87+
dxdt[i] = (2.0*dt[0] + dt[1])/(dt[0]*(dt[0] + dt[1])) * x[i] - (dt[0] + dt[1])/(dt[0]*dt[1]) * x_prev[0][i] + dt[0]/(dt[1]*(dt[0] + dt[1])) * x_prev[1][i];
9288
}
93-
val *= invdt / BDF_COEF[scheme][7];
94-
dxdt[i] = val;
89+
90+
alpha = (2.0*dt[0] + dt[1])/(dt[0]*(dt[0] + dt[1]));
9591
}
92+
else // Constant time stepper
93+
{
94+
const double invdt = 1.0 / dt[0];
9695

97-
alpha = invdt * ALPHA_COEF[scheme];
96+
for(MKL_INT i = 0; i < size; i++)
97+
{
98+
double val = BDF_COEF[scheme][0] * x[i];
99+
for(int d = 1; d <= m_scheme; d++)
100+
{
101+
val += BDF_COEF[scheme][d] * x_prev[d - 1][i];
102+
}
103+
val *= invdt / BDF_COEF[scheme][7];
104+
dxdt[i] = val;
105+
}
106+
107+
alpha = invdt * ALPHA_COEF[scheme];
108+
}
98109

99110
// Calculate RHS
100111
m_rhs(x, b, t);
@@ -132,7 +143,7 @@ void TimeIntegrator::operator()(sparse_matrix_holder &Jt, state_type &b,
132143
Jt.ja.resize(nzmax);
133144

134145
// Replaces deprecated mkl_dcsradd()
135-
// Jt: = J - M*invdt
146+
// Jt: = J - M*alpha
136147
matrix_add(-alpha, m_M, J, Jt);
137148
}
138149

src/time_integrator.h

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -73,7 +73,7 @@ class TimeIntegrator
7373
void operator()(sparse_matrix_holder &Jt, state_type &b,
7474
sparse_matrix_holder &J, state_type &x,
7575
const state_type_matrix &x_prev, const double t,
76-
const double dt);
76+
const double dt[]);
7777
};
7878

7979
} // namespace daecpp_namespace_name

0 commit comments

Comments
 (0)