ExaDG
Loading...
Searching...
No Matches
bdf_constants.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_TIME_INTEGRATION_BDF_CONSTANTS_H_
23#define INCLUDE_EXADG_TIME_INTEGRATION_BDF_CONSTANTS_H_
24
25// C/C++
26#include <algorithm>
27#include <vector>
28
29// deal.II
30#include <deal.II/base/conditional_ostream.h>
31#include <deal.II/base/exceptions.h>
32
33// ExaDG
34#include <exadg/time_integration/time_integration_constants_base.h>
35
36namespace ExaDG
37{
42{
43public:
44 BDFTimeIntegratorConstants(unsigned int const order, bool const start_with_low_order);
45
46 double
47 get_gamma0() const;
48
49 double
50 get_alpha(unsigned int const i) const;
51
52 void
53 print(dealii::ConditionalOStream & pcout) const final;
54
55private:
56 void
57 set_constant_time_step(unsigned int const current_order) final;
58
59 void
60 set_adaptive_time_step(unsigned int const current_order,
61 std::vector<double> const & time_steps) final;
62
63 /*
64 * BDF time integrator constants:
65 *
66 * du/dt = (gamma_0 u^{n+1} - alpha_0 u^{n} - alpha_1 u^{n-1} ... - alpha_{J-1} u^{n-J+1})/dt
67 */
68 double gamma0;
69
70 std::vector<double> alpha;
71};
72
73/*
74 * Calculates the time derivative
75 *
76 * derivative = du/dt = (gamma_0 u^{n+1} - alpha_0 u^{n} - alpha_1 u^{n-1} ... - alpha_{J-1}
77 * u^{n-J+1})/dt
78 */
79template<typename VectorType>
80void
81compute_bdf_time_derivative(VectorType & derivative,
82 VectorType const & solution_np,
83 std::vector<VectorType const *> const & previous_solutions,
85 double const & time_step_size)
86{
87 derivative.equ(bdf.get_gamma0() / time_step_size, solution_np);
88
89 for(unsigned int i = 0; i < previous_solutions.size(); ++i)
90 derivative.add(-bdf.get_alpha(i) / time_step_size, *previous_solutions[i]);
91}
92
93template<typename VectorType>
94void
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)
100{
101 std::vector<VectorType const *> previous_solutions_ptrs(previous_solutions.size());
102
103 std::transform(previous_solutions.begin(),
104 previous_solutions.end(),
105 previous_solutions_ptrs.begin(),
106 [](VectorType const & t) { return &t; });
107
108 compute_bdf_time_derivative(
109 derivative, solution_np, previous_solutions_ptrs, bdf, time_step_size);
110}
111
112template<typename VectorType>
113void
114compute_bdf_time_derivative(VectorType & derivative,
115 std::vector<VectorType const *> const & previous_solutions_np,
116 std::vector<double> const & times_np)
117{
118 AssertThrow(
119 previous_solutions_np.size() == times_np.size() and times_np.size() > 0,
120 dealii::ExcMessage(
121 "times and previous_solutions_np handed to compute_bdf_time_derivative() have different sizes or size 0."));
122
123 std::vector<VectorType const *> previous_solutions(previous_solutions_np.begin() + 1,
124 previous_solutions_np.end());
125
126 // compute time step sizes from times
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];
130
131
132 // construct BDF constants
133 unsigned int const order = previous_solutions.size();
134 BDFTimeIntegratorConstants bdf(order, true);
135 bdf.update(order, true, time_steps);
136
137 // compute temporal derivative
138 compute_bdf_time_derivative(
139 derivative, *previous_solutions_np[0], previous_solutions, bdf, time_steps[0]);
140}
141
142
143} // namespace ExaDG
144
145#endif /* INCLUDE_EXADG_TIME_INTEGRATION_BDF_CONSTANTS_H_ */
Definition bdf_constants.h:42
Definition time_integration_constants_base.h:30
Definition driver.cpp:33