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