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 INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_DUAL_SPLITTING_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_DUAL_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 TimeIntBDFDualSplitting : 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 TimeIntBDFDualSplitting(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 ~TimeIntBDFDualSplitting()
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 read_restart_vectors(BoostInputArchiveType & ia) final;
87
88 void
89 write_restart_vectors(BoostOutputArchiveType & oa) const final;
90
91 void
92 do_timestep_solve() final;
93
94 void
95 prepare_vectors_for_next_timestep() final;
96
97 void
98 convective_step();
99
100 void
101 evaluate_convective_term();
102
103 void
104 update_time_integrator_constants() final;
105
106 void
107 initialize_current_solution() final;
108
109 void
110 initialize_former_multistep_dof_vectors() final;
111
112 void
113 initialize_velocity_dbc();
114
115 void
116 pressure_step();
117
118 void
119 rhs_pressure(VectorType & rhs) const;
120
121 void
122 projection_step();
123
124 void
125 rhs_projection(VectorType & rhs) const;
126
127 void
128 penalty_step();
129
130 void
131 viscous_step();
132
133 void
134 rhs_viscous(VectorType & rhs,
135 VectorType const & velocity_mass_operator,
136 VectorType const & transport_velocity) const;
137
138 void
139 solve_steady_problem() final;
140
141 double
142 evaluate_residual();
143
144 VectorType const &
145 get_velocity(unsigned int i /* t_{n-i} */) const final;
146
147 VectorType const &
148 get_pressure(unsigned int i /* t_{n-i} */) const final;
149
150 void
151 set_velocity(VectorType const & velocity, unsigned int const i /* t_{n-i} */) final;
152
153 void
154 set_pressure(VectorType const & pressure, unsigned int const i /* t_{n-i} */) final;
155
156 std::shared_ptr<Operator> pde_operator;
157
158 std::vector<VectorType> velocity;
159
160 VectorType velocity_np;
161
162 std::vector<VectorType> pressure;
163
164 VectorType pressure_np;
165
166 std::vector<VectorType> velocity_dbc;
167 VectorType velocity_dbc_np;
168
169 // required for strongly-coupled partitioned FSI
170 VectorType pressure_last_iter;
171 VectorType velocity_projection_last_iter;
172 VectorType velocity_viscous_last_iter;
173
174 // iteration counts
175 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */>
176 iterations_pressure;
177 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */>
178 iterations_projection;
179 std::pair<
180 unsigned int /* calls */,
181 std::tuple<unsigned long long, unsigned long long> /* iteration counts {Newton, linear} */>
182 iterations_viscous;
183
184 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_penalty;
185 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_mass;
186
187 // time integrator constants: extrapolation scheme
188 ExtrapolationConstants extra_pressure_nbc;
189};
190
191} // namespace IncNS
192} // namespace ExaDG
193
194#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_DUAL_SPLITTING_H_ \
195 */
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:34
Definition driver.cpp:33