ExaDG
Loading...
Searching...
No Matches
laplace_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 LAPLACE_OPERATOR_H
23#define LAPLACE_OPERATOR_H
24
25#include <exadg/grid/grid_data.h>
26#include <exadg/operators/interior_penalty_parameter.h>
27#include <exadg/operators/operator_base.h>
28#include <exadg/operators/operator_type.h>
29#include <exadg/poisson/user_interface/boundary_descriptor.h>
30
31namespace ExaDG
32{
33namespace Poisson
34{
35namespace Operators
36{
38{
39 LaplaceKernelData() : IP_factor(1.0)
40 {
41 }
42
43 double IP_factor;
44};
45
46template<int dim, typename Number, int n_components = 1>
48{
49private:
50 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
51
52 typedef dealii::VectorizedArray<Number> scalar;
53
54 typedef FaceIntegrator<dim, n_components, Number> IntegratorFace;
55
56public:
57 LaplaceKernel() : degree(1), tau(dealii::make_vectorized_array<Number>(0.0))
58 {
59 }
60
61 void
62 reinit(dealii::MatrixFree<dim, Number> const & matrix_free,
63 LaplaceKernelData const & data_in,
64 unsigned int const dof_index)
65 {
66 data = data_in;
67
68 dealii::FiniteElement<dim> const & fe = matrix_free.get_dof_handler(dof_index).get_fe();
69 degree = fe.degree;
70
71 calculate_penalty_parameter(matrix_free, dof_index);
72 }
73
74 void
75 calculate_penalty_parameter(dealii::MatrixFree<dim, Number> const & matrix_free,
76 unsigned int const dof_index)
77 {
78 IP::calculate_penalty_parameter<dim, Number>(array_penalty_parameter, matrix_free, dof_index);
79 }
80
82 get_integrator_flags(bool const is_dg) const
83 {
84 IntegratorFlags flags;
85
86 flags.cell_evaluate = dealii::EvaluationFlags::gradients;
87 flags.cell_integrate = dealii::EvaluationFlags::gradients;
88
89 if(is_dg)
90 {
91 flags.face_evaluate = dealii::EvaluationFlags::values | dealii::EvaluationFlags::gradients;
92 flags.face_integrate = dealii::EvaluationFlags::values | dealii::EvaluationFlags::gradients;
93 }
94 else
95 {
96 // evaluation of Neumann BCs for continuous elements
97 flags.face_evaluate = dealii::EvaluationFlags::nothing;
98 flags.face_integrate = dealii::EvaluationFlags::values;
99 }
100
101 return flags;
102 }
103
104 static MappingFlags
105 get_mapping_flags(bool const compute_interior_face_integrals,
106 bool const compute_boundary_face_integrals)
107 {
108 MappingFlags flags;
109
110 flags.cells = dealii::update_gradients | dealii::update_JxW_values;
111
112 if(compute_interior_face_integrals)
113 {
114 flags.inner_faces =
115 dealii::update_gradients | dealii::update_JxW_values | dealii::update_normal_vectors;
116 }
117
118 if(compute_boundary_face_integrals)
119 {
120 flags.boundary_faces = dealii::update_gradients | dealii::update_JxW_values |
121 dealii::update_normal_vectors | dealii::update_quadrature_points;
122 }
123
124 return flags;
125 }
126
127 void
128 reinit_face(IntegratorFace & integrator_m,
129 IntegratorFace & integrator_p,
130 unsigned int const dof_index) const
131 {
132 tau = std::max(integrator_m.read_cell_data(array_penalty_parameter),
133 integrator_p.read_cell_data(array_penalty_parameter)) *
134 IP::get_penalty_factor<dim, Number>(
135 degree,
137 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
138 data.IP_factor);
139 }
140
141 void
142 reinit_boundary_face(IntegratorFace & integrator_m, unsigned int const dof_index) const
143 {
144 tau = integrator_m.read_cell_data(array_penalty_parameter) *
145 IP::get_penalty_factor<dim, Number>(
146 degree,
148 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
149 data.IP_factor);
150 }
151
152 void
153 reinit_face_cell_based(dealii::types::boundary_id const boundary_id,
154 IntegratorFace & integrator_m,
155 IntegratorFace & integrator_p,
156 unsigned int const dof_index) const
157 {
158 if(boundary_id == dealii::numbers::internal_face_boundary_id) // internal face
159 {
160 tau = std::max(integrator_m.read_cell_data(array_penalty_parameter),
161 integrator_p.read_cell_data(array_penalty_parameter)) *
162 IP::get_penalty_factor<dim, Number>(
163 degree,
165 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
166 data.IP_factor);
167 }
168 else // boundary face
169 {
170 tau = integrator_m.read_cell_data(array_penalty_parameter) *
171 IP::get_penalty_factor<dim, Number>(
172 degree,
174 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
175 data.IP_factor);
176 }
177 }
178
179 template<typename T>
180 inline DEAL_II_ALWAYS_INLINE //
181 T
182 calculate_gradient_flux(T const & value_m, T const & value_p) const
183 {
184 return -0.5 * (value_m - value_p);
185 }
186
187 template<typename T>
188 inline DEAL_II_ALWAYS_INLINE //
189 T
190 calculate_value_flux(T const & normal_gradient_m,
191 T const & normal_gradient_p,
192 T const & value_m,
193 T const & value_p) const
194 {
195 return 0.5 * (normal_gradient_m + normal_gradient_p) - tau * (value_m - value_p);
196 }
197
198private:
200
201 unsigned int degree;
202
203 dealii::AlignedVector<scalar> array_penalty_parameter;
204
205 mutable scalar tau;
206};
207
208} // namespace Operators
209
210template<int rank, int dim>
212{
213 LaplaceOperatorData() : OperatorBaseData(), quad_index_gauss_lobatto(0)
214 {
215 }
216
218
219 // continuous FE:
220 // for DirichletCached boundary conditions, another quadrature rule
221 // is needed to set the constrained DoFs.
222 unsigned int quad_index_gauss_lobatto;
223
224 std::shared_ptr<BoundaryDescriptor<rank, dim> const> bc;
225};
226
227template<int dim, typename Number, int n_components>
228class LaplaceOperator : public OperatorBase<dim, Number, n_components>
229{
230private:
231 static unsigned int const rank =
232 (n_components == 1) ? 0 : ((n_components == dim) ? 1 : dealii::numbers::invalid_unsigned_int);
233
236
237 typedef typename Base::IntegratorCell IntegratorCell;
238 typedef typename Base::IntegratorFace IntegratorFace;
239
240 typedef typename Base::Range Range;
241
242 typedef dealii::Tensor<rank, dim, dealii::VectorizedArray<Number>> value;
243
244 typedef typename Base::VectorType VectorType;
245
246public:
247 typedef Number value_type;
248
249 void
250 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
251 dealii::AffineConstraints<Number> const & affine_constraints,
252 LaplaceOperatorData<rank, dim> const & data);
253
255 get_data() const
256 {
257 return operator_data;
258 }
259
260 void
261 calculate_penalty_parameter(dealii::MatrixFree<dim, Number> const & matrix_free,
262 unsigned int const dof_index);
263
264 void
265 update_penalty_parameter();
266
267 // continuous FE: This function sets the inhomogeneous Dirichlet boundary values for Dirichlet
268 // degrees of freedom.
269 void
270 set_inhomogeneous_boundary_values(VectorType & solution) const final;
271
272 // only relevant for discontinuous Galerkin discretization (DG):
273 // Some more functionality on top of what is provided by the base class.
274 // This function evaluates the inhomogeneous boundary face integrals in DG where the
275 // Dirichlet boundary condition is extracted from a dof vector instead of a dealii::Function<dim>.
276 void
277 rhs_add_dirichlet_bc_from_dof_vector(VectorType & dst, VectorType const & src) const;
278
279private:
280 void
281 reinit_face_derived(IntegratorFace & integrator_m,
282 IntegratorFace & integrator_p,
283 unsigned int const face) const final;
284
285 void
286 reinit_boundary_face_derived(IntegratorFace & integrator_m, unsigned int const face) const final;
287
288 void
289 reinit_face_cell_based_derived(IntegratorFace & integrator_m,
290 IntegratorFace & integrator_p,
291 unsigned int const cell,
292 unsigned int const face,
293 dealii::types::boundary_id const boundary_id) const final;
294
295 void
296 do_cell_integral(IntegratorCell & integrator) const final;
297
298 void
299 do_face_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
300
301 void
302 do_face_int_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
303
304 void
305 do_face_ext_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
306
307 void
308 do_boundary_integral(IntegratorFace & integrator_m,
309 OperatorType const & operator_type,
310 dealii::types::boundary_id const & boundary_id) const final;
311
312 void
313 cell_loop_empty(dealii::MatrixFree<dim, Number> const & matrix_free,
314 VectorType & dst,
315 VectorType const & src,
316 Range const & range) const;
317
318 void
319 face_loop_empty(dealii::MatrixFree<dim, Number> const & matrix_free,
320 VectorType & dst,
321 VectorType const & src,
322 Range const & range) const;
323
324 // only relevant for discontinuous Galerkin discretization (DG)
325 void
326 boundary_face_loop_inhom_operator_dirichlet_bc_from_dof_vector(
327 dealii::MatrixFree<dim, Number> const & matrix_free,
328 VectorType & dst,
329 VectorType const & src,
330 Range const & range) const;
331
332 // only relevant for discontinuous Galerkin discretization (DG)
333 void
334 do_boundary_integral_dirichlet_bc_from_dof_vector(
335 IntegratorFace & integrator_m,
336 OperatorType const & operator_type,
337 dealii::types::boundary_id const & boundary_id) const;
338
339 // continuous FE: calculates Neumann boundary integral
340 void
341 do_boundary_integral_continuous(IntegratorFace & integrator_m,
342 dealii::types::boundary_id const & boundary_id) const final;
343
344 LaplaceOperatorData<rank, dim> operator_data;
345
347};
348
349} // namespace Poisson
350} // namespace ExaDG
351
352#endif
Definition operator_base.h:98
Definition laplace_operator.h:229
Definition laplace_operator.h:48
Definition driver.cpp:33
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70
Definition integrator_flags.h:30
Definition mapping_flags.h:31
Definition operator_base.h:59
Definition laplace_operator.h:212
Definition laplace_operator.h:38