Skip to content
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

Core library cleanup #1457

Merged
merged 4 commits into from
Mar 14, 2025
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
26 changes: 12 additions & 14 deletions include/core/bdf.h
Original file line number Diff line number Diff line change
Expand Up @@ -27,7 +27,7 @@ using namespace dealii;
* should be in decreasing order.
*
* For example, if the method is a BDF2 (\f$n=2\f$), it uses three values for
* the time: \f$t\f$, \f$t-\Delta t_1\f$ and \f$t-\Delta t_2\f$. Thus the @p
* the time: \f$t\f$, \f$t-\Delta t_1\f$ and \f$t-\Delta t_2\f$. Thus, the @p
* time_steps vector should contain \f$\Delta t_1\f$ and \f$\Delta t_2\f$.
*
* @return Vector containing BDF integration coefficient values.
Expand All @@ -36,8 +36,8 @@ using namespace dealii;
*/
Vector<double>
calculate_bdf_coefficients(
const Parameters::SimulationControl::TimeSteppingMethod method,
const std::vector<double> &time_steps);
Parameters::SimulationControl::TimeSteppingMethod method,
const std::vector<double> &time_steps);


/**
Expand Down Expand Up @@ -65,13 +65,11 @@ calculate_bdf_coefficients(
* calculation of BDF coefficients.
*/
Vector<double>
delta(const unsigned int order,
const unsigned int n,
const unsigned int j,
delta(unsigned int order,
unsigned int n,
unsigned int j,
const Vector<double> &times);



/**
* @brief Get the maximum number of previous time steps supposed by the BDF
* schemes implemented in Lethe.
Expand All @@ -97,17 +95,17 @@ maximum_number_of_previous_solutions()
*/
inline unsigned int
number_of_previous_solutions(
Parameters::SimulationControl::TimeSteppingMethod method)
const Parameters::SimulationControl::TimeSteppingMethod method)
Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Shouldn't the const be absent here, given the link in the description?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Very good question friend. Because this is an inline function the function is defined afterwards. So here the header does not contain the function prototype but the function itself. Consequently, the const has to be there because the function is defined straight up afterwards :)

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Oh that makes sense! Thank you for the explanation

{
if (method == Parameters::SimulationControl::TimeSteppingMethod::bdf1 ||
method == Parameters::SimulationControl::TimeSteppingMethod::steady_bdf)
return 1;
else if (method == Parameters::SimulationControl::TimeSteppingMethod::bdf2)
if (method == Parameters::SimulationControl::TimeSteppingMethod::bdf2)
return 2;
else if (method == Parameters::SimulationControl::TimeSteppingMethod::bdf3)
if (method == Parameters::SimulationControl::TimeSteppingMethod::bdf3)
return 3;
else
return 0;

return 0;
}

/**
Expand All @@ -131,7 +129,7 @@ number_of_previous_solutions(
*/

template <typename DataType>
inline void
void
bdf_extrapolate(const std::vector<double> &time_vector,
const std::vector<std::vector<DataType>> &solution_vector,
const unsigned int number_of_previous_solutions,
Expand Down
Loading