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