ExaDG
Loading...
Searching...
No Matches
divergence_penalty_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_PENALTY_OPERATOR_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_DIVERGENCE_PENALTY_OPERATOR_H_
24
25#include <exadg/grid/calculate_characteristic_element_length.h>
26#include <exadg/incompressible_navier_stokes/user_interface/parameters.h>
27#include <exadg/matrix_free/integrators.h>
28#include <exadg/operators/integrator_flags.h>
29#include <exadg/operators/mapping_flags.h>
30
31namespace ExaDG
32{
33namespace IncNS
34{
35/*
36 * Divergence penalty operator:
37 *
38 * ( div(v_h) , tau_div * div(u_h) )_Omega^e
39 *
40 * where
41 *
42 * v_h : test function
43 * u_h : solution
44 * tau_div: divergence penalty factor
45 *
46 * use convective term: tau_div_conv = K * ||U||_mean * h_eff
47 *
48 * where h_eff = h / (k_u+1) with a characteristic
49 * element length h derived from the element volume V_e
50 *
51 * use viscous term: tau_div_viscous = K * nu
52 *
53 * use both terms: tau_div = tau_div_conv + tau_div_viscous
54 *
55 */
56
57namespace Operators
58{
59struct DivergencePenaltyKernelData
60{
61 DivergencePenaltyKernelData()
62 : type_penalty_parameter(TypePenaltyParameter::ConvectiveTerm),
63 viscosity(0.0),
64 degree(1),
65 penalty_factor(1.0)
66 {
67 }
68
69 // type of penalty parameter (viscous and/or convective terms)
70 TypePenaltyParameter type_penalty_parameter;
71
72 // viscosity, needed for computation of penalty factor
73 double viscosity;
74
75 // degree of finite element shape functions
76 unsigned int degree;
77
78 // the penalty term can be scaled by 'penalty_factor'
79 double penalty_factor;
80};
81
82template<int dim, typename Number>
83class DivergencePenaltyKernel
84{
85private:
86 typedef CellIntegrator<dim, dim, Number> IntegratorCell;
87
88 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
89
90 typedef dealii::VectorizedArray<Number> scalar;
91
92public:
93 DivergencePenaltyKernel()
94 : matrix_free(nullptr), dof_index(0), quad_index(0), array_penalty_parameter(0)
95 {
96 }
97
98 void
99 reinit(dealii::MatrixFree<dim, Number> const & matrix_free,
100 unsigned int const dof_index,
101 unsigned int const quad_index,
102 DivergencePenaltyKernelData const & data)
103 {
104 this->matrix_free = &matrix_free;
105
106 this->dof_index = dof_index;
107 this->quad_index = quad_index;
108
109 this->data = data;
110
111 unsigned int n_cells = matrix_free.n_cell_batches() + matrix_free.n_ghost_cell_batches();
112 array_penalty_parameter.resize(n_cells);
113 }
114
116 get_data()
117 {
118 return this->data;
119 }
120
122 get_integrator_flags() const
123 {
124 IntegratorFlags flags;
125
126 flags.cell_evaluate = dealii::EvaluationFlags::gradients;
127 flags.cell_integrate = dealii::EvaluationFlags::gradients;
128
129 // no face integrals
130
131 return flags;
132 }
133
134 static MappingFlags
135 get_mapping_flags()
136 {
137 MappingFlags flags;
138
139 flags.cells = dealii::update_JxW_values | dealii::update_gradients;
140
141 // no face integrals
142
143 return flags;
144 }
145
146 void
147 calculate_penalty_parameter(VectorType const & velocity)
148 {
149 velocity.update_ghost_values();
150
151 IntegratorCell integrator(*matrix_free, dof_index, quad_index);
152
153 dealii::AlignedVector<scalar> JxW_values(integrator.n_q_points);
154
155 ElementType const element_type =
156 get_element_type(matrix_free->get_dof_handler(dof_index).get_triangulation());
157
158 unsigned int n_cells = matrix_free->n_cell_batches() + matrix_free->n_ghost_cell_batches();
159 for(unsigned int cell = 0; cell < n_cells; ++cell)
160 {
161 scalar tau_convective = dealii::make_vectorized_array<Number>(0.0);
162 scalar tau_viscous = dealii::make_vectorized_array<Number>(data.viscosity);
163
164 if(data.type_penalty_parameter == TypePenaltyParameter::ConvectiveTerm or
165 data.type_penalty_parameter == TypePenaltyParameter::ViscousAndConvectiveTerms)
166 {
167 integrator.reinit(cell);
168 integrator.read_dof_values(velocity);
169 integrator.evaluate(dealii::EvaluationFlags::values);
170
171 scalar volume = dealii::make_vectorized_array<Number>(0.0);
172 scalar norm_U_mean = dealii::make_vectorized_array<Number>(0.0);
173 for(unsigned int q = 0; q < integrator.n_q_points; ++q)
174 {
175 volume += integrator.JxW(q);
176 norm_U_mean += integrator.JxW(q) * integrator.get_value(q).norm();
177 }
178 norm_U_mean /= volume;
179
180 scalar h = calculate_characteristic_element_length(volume, dim, element_type);
181 scalar h_eff = calculate_high_order_element_length(h, data.degree, true);
182
183 tau_convective = norm_U_mean * h_eff;
184 }
185
186 if(data.type_penalty_parameter == TypePenaltyParameter::ConvectiveTerm)
187 {
188 array_penalty_parameter[cell] = data.penalty_factor * tau_convective;
189 }
190 else if(data.type_penalty_parameter == TypePenaltyParameter::ViscousTerm)
191 {
192 array_penalty_parameter[cell] = data.penalty_factor * tau_viscous;
193 }
194 else if(data.type_penalty_parameter == TypePenaltyParameter::ViscousAndConvectiveTerms)
195 {
196 array_penalty_parameter[cell] = data.penalty_factor * (tau_convective + tau_viscous);
197 }
198 }
199 }
200
201 void
202 reinit_cell(IntegratorCell & integrator) const
203 {
204 tau = integrator.read_cell_data(array_penalty_parameter);
205 }
206
207 /*
208 * Volume flux, i.e., the term occurring in the volume integral
209 */
210 inline DEAL_II_ALWAYS_INLINE //
211 scalar
212 get_volume_flux(IntegratorCell const & integrator, unsigned int const q) const
213 {
214 return tau * integrator.get_divergence(q);
215 }
216
217private:
218 dealii::MatrixFree<dim, Number> const * matrix_free;
219
220 unsigned int dof_index;
221 unsigned int quad_index;
222
224
225 dealii::AlignedVector<scalar> array_penalty_parameter;
226
227 mutable scalar tau;
228};
229
230} // namespace Operators
231
232struct DivergencePenaltyData
233{
234 DivergencePenaltyData() : dof_index(0), quad_index(0)
235 {
236 }
237
238 unsigned int dof_index;
239 unsigned int quad_index;
240};
241
242template<int dim, typename Number>
243class DivergencePenaltyOperator
244{
245private:
246 typedef DivergencePenaltyOperator<dim, Number> This;
247
248 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
249
250 typedef std::pair<unsigned int, unsigned int> Range;
251
252 typedef CellIntegrator<dim, dim, Number> IntegratorCell;
253
255
256public:
257 DivergencePenaltyOperator();
258
259 void
260 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
261 DivergencePenaltyData const & data,
262 std::shared_ptr<Kernel> const kernel);
263
264 void
265 update(VectorType const & velocity);
266
267 void
268 apply(VectorType & dst, VectorType const & src) const;
269
270 void
271 apply_add(VectorType & dst, VectorType const & src) const;
272
273private:
274 void
275 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
276 VectorType & dst,
277 VectorType const & src,
278 Range const & range) const;
279
280 void
281 do_cell_integral(IntegratorCell & integrator) const;
282
283 dealii::MatrixFree<dim, Number> const * matrix_free;
284
286
287 std::shared_ptr<Operators::DivergencePenaltyKernel<dim, Number>> kernel;
288};
289
290} // namespace IncNS
291} // namespace ExaDG
292
293#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_DIVERGENCE_PENALTY_OPERATOR_H_ \
294 */
Definition divergence_penalty_operator.h:84
Definition driver.cpp:33
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70
Number calculate_high_order_element_length(Number const element_length, unsigned int const fe_degree, bool const is_dg)
Definition calculate_characteristic_element_length.h:108
Number calculate_characteristic_element_length(Number const element_volume, unsigned int const dim, ElementType const &element_type)
Definition calculate_characteristic_element_length.h:68
Definition divergence_penalty_operator.h:233
Definition divergence_penalty_operator.h:60
Definition integrator_flags.h:30
Definition mapping_flags.h:31