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>
33class OperatorDualSplitting;
34
35template<int dim, typename Number>
36class TimeIntBDFDualSplitting : public TimeIntBDF<dim, Number>
37{
38private:
40
41 typedef typename Base::VectorType VectorType;
42
44
45public:
46 TimeIntBDFDualSplitting(std::shared_ptr<Operator> operator_in,
47 std::shared_ptr<HelpersALE<dim, Number> const> helpers_ale_in,
48 std::shared_ptr<PostProcessorInterface<Number>> postprocessor_in,
49 Parameters const & param_in,
50 MPI_Comm const & mpi_comm_in,
51 bool const is_test_in);
52
54 {
55 }
56
57 void
58 postprocessing_stability_analysis();
59
60 void
61 print_iterations() const final;
62
63 VectorType const &
64 get_velocity() const final;
65
66 VectorType const &
67 get_velocity_np() const final;
68
69 VectorType const &
70 get_pressure() const final;
71
72 VectorType const &
73 get_pressure_np() const final;
74
75private:
76 void
77 allocate_vectors() final;
78
79 void
80 setup_derived() final;
81
82 void
83 read_restart_vectors(boost::archive::binary_iarchive & ia) final;
84
85 void
86 write_restart_vectors(boost::archive::binary_oarchive & oa) const final;
87
88 void
89 do_timestep_solve() final;
90
91 void
92 prepare_vectors_for_next_timestep() final;
93
94 void
95 convective_step();
96
97 void
98 evaluate_convective_term();
99
100 void
101 update_time_integrator_constants() final;
102
103 void
104 initialize_current_solution() final;
105
106 void
107 initialize_former_multistep_dof_vectors() final;
108
109 void
110 initialize_velocity_dbc();
111
112 void
113 pressure_step();
114
115 void
116 rhs_pressure(VectorType & rhs) const;
117
118 void
119 projection_step();
120
121 void
122 rhs_projection(VectorType & rhs) const;
123
124 void
125 penalty_step();
126
127 void
128 viscous_step();
129
130 void
131 rhs_viscous(VectorType & rhs) const;
132
133 void
134 solve_steady_problem() final;
135
136 double
137 evaluate_residual();
138
139 VectorType const &
140 get_velocity(unsigned int i /* t_{n-i} */) const final;
141
142 VectorType const &
143 get_pressure(unsigned int i /* t_{n-i} */) const final;
144
145 void
146 set_velocity(VectorType const & velocity, unsigned int const i /* t_{n-i} */) final;
147
148 void
149 set_pressure(VectorType const & pressure, unsigned int const i /* t_{n-i} */) final;
150
151 std::shared_ptr<Operator> pde_operator;
152
153 std::vector<VectorType> velocity;
154
155 VectorType velocity_np;
156
157 std::vector<VectorType> pressure;
158
159 VectorType pressure_np;
160
161 std::vector<VectorType> velocity_dbc;
162 VectorType velocity_dbc_np;
163
164 // required for strongly-coupled partitioned FSI
165 VectorType pressure_last_iter;
166 VectorType velocity_projection_last_iter;
167 VectorType velocity_viscous_last_iter;
168
169 // iteration counts
170 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */>
171 iterations_pressure;
172 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */>
173 iterations_projection;
174 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_viscous;
175 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_penalty;
176 std::pair<unsigned int /* calls */, unsigned long long /* iteration counts */> iterations_mass;
177
178 // time integrator constants: extrapolation scheme
179 ExtrapolationConstants extra_pressure_nbc;
180};
181
182} // namespace IncNS
183} // namespace ExaDG
184
185#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_TIME_INT_BDF_DUAL_SPLITTING_H_ \
186 */
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 time_int_bdf_dual_splitting.h:37
Definition time_int_bdf.h:46
Definition driver.cpp:33