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