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