ExaDG
Loading...
Searching...
No Matches
divergence_operator.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_SPATIAL_DISCRETIZATION_OPERATORS_DIVERGENCE_OPERATOR_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_DIVERGENCE_OPERATOR_H_
24
25
26#include <exadg/incompressible_navier_stokes/spatial_discretization/operators/weak_boundary_conditions.h>
27#include <exadg/incompressible_navier_stokes/user_interface/parameters.h>
28#include <exadg/matrix_free/integrators.h>
29#include <exadg/operators/mapping_flags.h>
30
31namespace ExaDG
32{
33namespace IncNS
34{
35namespace Operators
36{
37template<int dim, typename Number>
39{
40private:
41 typedef CellIntegrator<dim, dim, Number> CellIntegratorU;
42
43 typedef dealii::VectorizedArray<Number> scalar;
44 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
45
46public:
47 static MappingFlags
48 get_mapping_flags()
49 {
50 MappingFlags flags;
51
52 flags.cells = dealii::update_JxW_values | dealii::update_gradients;
53 flags.inner_faces = dealii::update_JxW_values | dealii::update_normal_vectors;
54 flags.boundary_faces =
55 dealii::update_JxW_values | dealii::update_quadrature_points | dealii::update_normal_vectors;
56
57 return flags;
58 }
59
60 /*
61 * This function implements the central flux as numerical flux function.
62 */
63 inline DEAL_II_ALWAYS_INLINE //
64 vector
65 calculate_flux(vector const & value_m, vector const & value_p) const
66 {
67 return 0.5 * (value_m + value_p);
68 }
69
70 /*
71 * Volume flux, i.e., the term occurring in the volume integral for
72 * weak formulation (performing integration-by-parts)
73 */
74 inline DEAL_II_ALWAYS_INLINE //
75 vector
76 get_volume_flux_weak(CellIntegratorU & velocity, unsigned int const q) const
77 {
78 // minus sign due to integration by parts
79 return -velocity.get_value(q);
80 }
81
82 /*
83 * Volume flux, i.e., the term occurring in the volume integral for
84 * strong formulation (no integration-by-parts, or integration-by-parts performed twice)
85 */
86 inline DEAL_II_ALWAYS_INLINE //
87 scalar
88 get_volume_flux_strong(CellIntegratorU & velocity, unsigned int const q) const
89 {
90 return velocity.get_divergence(q);
91 }
92};
93} // namespace Operators
94
95template<int dim>
96struct DivergenceOperatorData
97{
98 DivergenceOperatorData()
99 : dof_index_velocity(0),
100 dof_index_pressure(1),
101 quad_index(0),
102 integration_by_parts(true),
103 use_boundary_data(true),
104 formulation(FormulationVelocityDivergenceTerm::Weak)
105 {
106 }
107
108 unsigned int dof_index_velocity;
109 unsigned int dof_index_pressure;
110
111 unsigned int quad_index;
112
113 bool integration_by_parts;
114 bool use_boundary_data;
115
116 FormulationVelocityDivergenceTerm formulation;
117
118 std::shared_ptr<BoundaryDescriptorU<dim> const> bc;
119};
120
121template<int dim, typename Number>
122class DivergenceOperator
123{
124public:
125 typedef DivergenceOperator<dim, Number> This;
126
127 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
128
129 typedef dealii::VectorizedArray<Number> scalar;
130 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
131
132 typedef std::pair<unsigned int, unsigned int> Range;
133
134 typedef CellIntegrator<dim, dim, Number> CellIntegratorU;
135 typedef CellIntegrator<dim, 1, Number> CellIntegratorP;
136
137 typedef FaceIntegrator<dim, dim, Number> FaceIntegratorU;
138 typedef FaceIntegrator<dim, 1, Number> FaceIntegratorP;
139
140 DivergenceOperator();
141
142 void
143 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
144 DivergenceOperatorData<dim> const & data);
145
147 get_operator_data() const;
148
149 // homogeneous operator
150 void
151 apply(VectorType & dst, VectorType const & src) const;
152
153 void
154 apply_add(VectorType & dst, VectorType const & src) const;
155
156 // inhomogeneous operator
157 void
158 rhs(VectorType & dst, Number const evaluation_time) const;
159
160 void
161 rhs_bc_from_dof_vector(VectorType & dst, VectorType const & src) const;
162
163 void
164 rhs_add(VectorType & dst, Number const evaluation_time) const;
165
166 // full operator, i.e., homogeneous and inhomogeneous contributions
167 void
168 evaluate(VectorType & dst, VectorType const & src, Number const evaluation_time) const;
169
170 void
171 evaluate_add(VectorType & dst, VectorType const & src, Number const evaluation_time) const;
172
173private:
174 void
175 do_cell_integral_weak(CellIntegratorP & pressure, CellIntegratorU & velocity) const;
176
177 void
178 do_cell_integral_strong(CellIntegratorP & pressure, CellIntegratorU & velocity) const;
179
180 void
181 do_face_integral(FaceIntegratorU & velocity_m,
182 FaceIntegratorU & velocity_p,
183 FaceIntegratorP & pressure_m,
184 FaceIntegratorP & pressure_p) const;
185
186 void
187 do_boundary_integral(FaceIntegratorU & velocity,
188 FaceIntegratorP & pressure,
189 OperatorType const & operator_type,
190 dealii::types::boundary_id const & boundary_id) const;
191
192 void
193 do_boundary_integral_from_dof_vector(FaceIntegratorU & velocity,
194 FaceIntegratorU & velocity_exterior,
195 FaceIntegratorP & pressure,
196 OperatorType const & operator_type,
197 dealii::types::boundary_id const & boundary_id) const;
198
199 void
200 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
201 VectorType & dst,
202 VectorType const & src,
203 Range const & cell_range) const;
204
205 void
206 face_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
207 VectorType & dst,
208 VectorType const & src,
209 Range const & face_range) const;
210
211 void
212 boundary_face_loop_hom_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
213 VectorType & dst,
214 VectorType const & src,
215 Range const & face_range) const;
216
217 void
218 boundary_face_loop_full_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
219 VectorType & dst,
220 VectorType const & src,
221 Range const & face_range) const;
222
223 void
224 cell_loop_inhom_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
225 VectorType & dst,
226 VectorType const & src,
227 Range const & cell_range) const;
228
229 void
230 face_loop_inhom_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
231 VectorType & dst,
232 VectorType const & src,
233 Range const & face_range) const;
234
235 void
236 boundary_face_loop_inhom_operator(dealii::MatrixFree<dim, Number> const & matrix_free,
237 VectorType & dst,
238 VectorType const & src,
239 std::pair<unsigned int, unsigned int> const & face_range) const;
240
241 void
242 boundary_face_loop_inhom_operator_bc_from_dof_vector(
243 dealii::MatrixFree<dim, Number> const & matrix_free,
244 VectorType & dst,
245 VectorType const & src,
246 Range const & face_range) const;
247
248 dealii::MatrixFree<dim, Number> const * matrix_free;
249
251
252 mutable double time;
253
255
256 // needed if Dirichlet boundary condition is evaluated from dof vector
257 mutable VectorType const * velocity_bc;
258};
259
260} // namespace IncNS
261} // namespace ExaDG
262
263
264#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_DIVERGENCE_OPERATOR_H_ \
265 */
Definition divergence_operator.h:39
Definition driver.cpp:33
Definition divergence_operator.h:97
Definition mapping_flags.h:31