ExaDG
Loading...
Searching...
No Matches
time_int_bdf_consistent_splitting.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 EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_CONSISTENT_SPLITTING_H_
23#define EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_CONSISTENT_SPLITTING_H_
24
25#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf.h>
26
27namespace ExaDG
28{
29namespace IncNS
30{
31// forward declarations
32template<int dim, typename Number>
34
35template<int dim, typename Number>
36class TimeIntBDFConsistentSplitting : public TimeIntBDF<dim, Number>
37{
38private:
39 using BoostInputArchiveType = TimeIntBase::BoostInputArchiveType;
40 using BoostOutputArchiveType = TimeIntBase::BoostOutputArchiveType;
41
42 typedef TimeIntBDF<dim, Number> Base;
43
44 typedef typename Base::VectorType VectorType;
45
47
48public:
49 TimeIntBDFConsistentSplitting(std::shared_ptr<Operator> operator_in,
50 std::shared_ptr<HelpersALE<dim, Number> const> helpers_ale_in,
51 std::shared_ptr<PostProcessorInterface<Number>> postprocessor_in,
52 Parameters const & param_in,
53 MPI_Comm const & mpi_comm_in,
54 bool const is_test_in);
55
56 virtual ~TimeIntBDFConsistentSplitting()
57 {
58 }
59
60 void
61 postprocessing_stability_analysis();
62
63 void
64 print_iterations() const final;
65
66 VectorType const &
67 get_velocity() const final;
68
69 VectorType const &
70 get_velocity_np() const final;
71
72 VectorType const &
73 get_pressure() const final;
74
75 VectorType const &
76 get_pressure_np() const final;
77
78private:
79 void
80 allocate_vectors() final;
81
82 void
83 setup_derived() final;
84
85 void
86 do_timestep_solve() final;
87
88 void
89 prepare_vectors_for_next_timestep() final;
90
91 void
92 evaluate_convective_term();
93
94 void
95 update_time_integrator_constants() final;
96
97 void
98 initialize_current_solution() final;
99
100 void
101 initialize_former_multistep_dof_vectors() final;
102
103 void
104 initialize_velocity_dbc();
105
106 void
107 pressure_step();
108
109 void
110 rhs_pressure(VectorType & rhs) const;
111
112 void
113 momentum_step();
114
115 void
116 rhs_momentum(VectorType & rhs, VectorType const & transport_velocity) const;
117
118 void
119 penalty_step();
120
121 void
122 solve_steady_problem() final;
123
124 double
125 evaluate_residual();
126
127 VectorType const &
128 get_velocity(unsigned int i /* t_{n-i} */) const final;
129
130 VectorType const &
131 get_pressure(unsigned int i /* t_{n-i} */) const final;
132
133 void
134 set_velocity(VectorType const & velocity, unsigned int const i /* t_{n-i} */) final;
135
136 void
137 set_pressure(VectorType const & pressure, unsigned int const i /* t_{n-i} */) final;
138
139 std::shared_ptr<Operator> pde_operator;
140
141 std::vector<VectorType> velocity;
142
143 VectorType velocity_np;
144
145 std::vector<VectorType> pressure;
146
147 VectorType pressure_np;
148
149 std::vector<VectorType> velocity_divergence;
150
151 std::vector<VectorType> vec_convective_term_div;
152
153 std::vector<VectorType> velocity_dbc;
154 VectorType velocity_dbc_np;
155
156 // iteration counts
157 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */>
158 iterations_pressure;
159 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */>
160 iterations_projection;
161 std::pair<
162 unsigned int /* calls */,
163 std::tuple<unsigned long long, unsigned long long> /* iteration counts {Newton, linear} */>
164 iterations_viscous;
165
166 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_penalty;
167 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_mass;
168
169 // time integrator constants: extrapolation scheme
170 ExtrapolationConstants extra_pressure_nbc;
171 ExtrapolationConstants extra_pressure_rhs;
172};
173
174} // namespace IncNS
175} // namespace ExaDG
176
177#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_Consistent_SPLITTING_H_ \
178 */
Definition extrapolation_constants.h:37
Definition lambda_functions_ale.h:40
Definition operator_consistent_splitting.h:35
Definition parameters.h:46
Definition postprocessor_interface.h:37
Definition driver.cpp:33