22#ifndef INCLUDE_EXADG_TIME_INTEGRATION_BDF_CONSTANTS_H_
23#define INCLUDE_EXADG_TIME_INTEGRATION_BDF_CONSTANTS_H_
30#include <deal.II/base/conditional_ostream.h>
31#include <deal.II/base/exceptions.h>
34#include <exadg/time_integration/time_integration_constants_base.h>
41class BDFTimeIntegratorConstants :
public TimeIntegratorConstantsBase
44 BDFTimeIntegratorConstants(
unsigned int const order,
bool const start_with_low_order);
50 get_alpha(
unsigned int const i)
const;
53 print(dealii::ConditionalOStream & pcout)
const final;
57 set_constant_time_step(
unsigned int const current_order)
final;
60 set_adaptive_time_step(
unsigned int const current_order,
61 std::vector<double>
const & time_steps)
final;
70 std::vector<double> alpha;
41class BDFTimeIntegratorConstants :
public TimeIntegratorConstantsBase {
…};
79template<
typename VectorType>
81compute_bdf_time_derivative(
VectorType & derivative,
83 std::vector<VectorType const *>
const & previous_solutions,
85 double const & time_step_size)
87 derivative.equ(bdf.get_gamma0() / time_step_size, solution_np);
89 for(
unsigned int i = 0; i < previous_solutions.size(); ++i)
90 derivative.add(-bdf.get_alpha(i) / time_step_size, *previous_solutions[i]);
93template<
typename VectorType>
95compute_bdf_time_derivative(VectorType & derivative,
96 VectorType
const & solution_np,
97 std::vector<VectorType>
const & previous_solutions,
98 BDFTimeIntegratorConstants
const & bdf,
99 double const & time_step_size)
101 std::vector<VectorType const *> previous_solutions_ptrs(previous_solutions.size());
103 std::transform(previous_solutions.begin(),
104 previous_solutions.end(),
105 previous_solutions_ptrs.begin(),
106 [](VectorType
const & t) { return &t; });
108 compute_bdf_time_derivative(
109 derivative, solution_np, previous_solutions_ptrs, bdf, time_step_size);
112template<
typename VectorType>
114compute_bdf_time_derivative(
VectorType & derivative,
115 std::vector<VectorType const *>
const & previous_solutions_np,
116 std::vector<double>
const & times_np)
119 previous_solutions_np.size() == times_np.size() and times_np.size() > 0,
121 "times and previous_solutions_np handed to compute_bdf_time_derivative() have different sizes or size 0."));
123 std::vector<VectorType const *> previous_solutions(previous_solutions_np.begin() + 1,
124 previous_solutions_np.end());
127 std::vector<double> time_steps(times_np.size() - 1);
128 for(
unsigned int i = 0; i < time_steps.size(); ++i)
129 time_steps[i] = times_np[i] - times_np[i + 1];
133 unsigned int const order = previous_solutions.size();
135 bdf.update(order,
true, time_steps);
138 compute_bdf_time_derivative(
139 derivative, *previous_solutions_np[0], previous_solutions, bdf, time_steps[0]);
Definition bdf_constants.h:42