ExaDG
Loading...
Searching...
No Matches
viscous_operator.h
1/*
2 * viscous_operator.h
3 *
4 * Created on: Nov 5, 2018
5 * Author: fehn
6 */
7
8#ifndef INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_VISCOUS_OPERATOR_H_
9#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_VISCOUS_OPERATOR_H_
10
11#include <exadg/grid/grid_data.h>
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/interior_penalty_parameter.h>
16#include <exadg/operators/operator_base.h>
17#include <exadg/operators/variable_coefficients.h>
18
19namespace ExaDG
20{
21namespace IncNS
22{
23namespace Operators
24{
26{
28 : IP_factor(1.0),
29 viscosity(1.0),
30 formulation_viscous_term(FormulationViscousTerm::DivergenceFormulation),
31 penalty_term_div_formulation(PenaltyTermDivergenceFormulation::Symmetrized),
32 IP_formulation(InteriorPenaltyFormulation::SIPG),
33 viscosity_is_variable(false),
34 variable_normal_vector(false)
35 {
36 }
37
38 double IP_factor;
39 double viscosity;
40 FormulationViscousTerm formulation_viscous_term;
41 PenaltyTermDivergenceFormulation penalty_term_div_formulation;
42 InteriorPenaltyFormulation IP_formulation;
43 bool viscosity_is_variable;
44 bool variable_normal_vector;
45};
46
47template<int dim, typename Number>
49{
50private:
51 typedef dealii::VectorizedArray<Number> scalar;
52 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
53 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
54
55 typedef CellIntegrator<dim, dim, Number> IntegratorCell;
56 typedef FaceIntegrator<dim, dim, Number> IntegratorFace;
57
58public:
59 ViscousKernel() : quad_index(0), degree(1), tau(dealii::make_vectorized_array<Number>(0.0))
60 {
61 }
62
63 void
64 reinit(dealii::MatrixFree<dim, Number> const & matrix_free,
65 ViscousKernelData const & data,
66 unsigned int const dof_index,
67 unsigned int const quad_index)
68 {
69 this->data = data;
70
71 this->quad_index = quad_index;
72
73 dealii::FiniteElement<dim> const & fe = matrix_free.get_dof_handler(dof_index).get_fe();
74 this->degree = fe.degree;
75
76 calculate_penalty_parameter(matrix_free, dof_index);
77
78 AssertThrow(data.viscosity >= 0.0, dealii::ExcMessage("Viscosity is not set!"));
79
80 if(data.viscosity_is_variable)
81 {
82 // allocate vectors for variable coefficients and initialize with constant viscosity
83 viscosity_coefficients.initialize(matrix_free, quad_index, true, false);
84 viscosity_coefficients.set_coefficients(data.viscosity);
85 }
86 }
87
88 void
89 calculate_penalty_parameter(dealii::MatrixFree<dim, Number> const & matrix_free,
90 unsigned int const dof_index)
91 {
92 IP::calculate_penalty_parameter<dim, Number>(array_penalty_parameter, matrix_free, dof_index);
93 }
94
95 ViscousKernelData const &
96 get_data() const
97 {
98 return this->data;
99 }
100
101 unsigned int
102 get_quad_index() const
103 {
104 return this->quad_index;
105 }
106
107 unsigned int
108 get_degree() const
109 {
110 return this->degree;
111 }
112
113 void
114 set_constant_coefficient(Number const & constant_coefficient)
115 {
116 viscosity_coefficients.set_coefficients(constant_coefficient);
117 }
118
119 void
120 set_coefficient_cell(unsigned int const cell, unsigned int const q, scalar const & value)
121 {
122 viscosity_coefficients.set_coefficient_cell(cell, q, value);
123 }
124
125 scalar
126 get_coefficient_face(unsigned int const face, unsigned int const q)
127 {
128 return viscosity_coefficients.get_coefficient_face(face, q);
129 }
130
131 void
132 set_coefficient_face(unsigned int const face, unsigned int const q, scalar const & value)
133 {
134 viscosity_coefficients.set_coefficient_face(face, q, value);
135 }
136
137 scalar
138 get_coefficient_face_neighbor(unsigned int const face, unsigned int const q)
139 {
140 return viscosity_coefficients.get_coefficient_face_neighbor(face, q);
141 }
142
143 void
144 set_coefficient_face_neighbor(unsigned int const face, unsigned int const q, scalar const & value)
145 {
146 viscosity_coefficients.set_coefficient_face_neighbor(face, q, value);
147 }
148
150 get_integrator_flags() const
151 {
152 IntegratorFlags flags;
153
154 flags.cell_evaluate = dealii::EvaluationFlags::gradients;
155 flags.cell_integrate = dealii::EvaluationFlags::gradients;
156
157 flags.face_evaluate = dealii::EvaluationFlags::values | dealii::EvaluationFlags::gradients;
158 flags.face_integrate = dealii::EvaluationFlags::values | dealii::EvaluationFlags::gradients;
159
160 return flags;
161 }
162
163 static MappingFlags
164 get_mapping_flags(bool const compute_interior_face_integrals,
165 bool const compute_boundary_face_integrals)
166 {
167 MappingFlags flags;
168
169 flags.cells = dealii::update_JxW_values | dealii::update_gradients;
170 if(compute_interior_face_integrals)
171 flags.inner_faces =
172 dealii::update_JxW_values | dealii::update_gradients | dealii::update_normal_vectors;
173 if(compute_boundary_face_integrals)
174 flags.boundary_faces = dealii::update_JxW_values | dealii::update_gradients |
175 dealii::update_normal_vectors | dealii::update_quadrature_points;
176
177 return flags;
178 }
179
180 void
181 reinit_face(IntegratorFace & integrator_m,
182 IntegratorFace & integrator_p,
183 unsigned int const dof_index) const
184 {
185 tau = std::max(integrator_m.read_cell_data(array_penalty_parameter),
186 integrator_p.read_cell_data(array_penalty_parameter)) *
187 IP::get_penalty_factor<dim, Number>(
188 degree,
190 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
191 data.IP_factor);
192 }
193
194 void
195 reinit_boundary_face(IntegratorFace & integrator_m, unsigned int const dof_index) const
196 {
197 tau = integrator_m.read_cell_data(array_penalty_parameter) *
198 IP::get_penalty_factor<dim, Number>(
199 degree,
201 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
202 data.IP_factor);
203 }
204
205 void
206 reinit_face_cell_based(dealii::types::boundary_id const boundary_id,
207 IntegratorFace & integrator_m,
208 IntegratorFace & integrator_p,
209 unsigned int const dof_index) const
210 {
211 if(boundary_id == dealii::numbers::internal_face_boundary_id) // internal face
212 {
213 tau = std::max(integrator_m.read_cell_data(array_penalty_parameter),
214 integrator_p.read_cell_data(array_penalty_parameter)) *
215 IP::get_penalty_factor<dim, Number>(
216 degree,
218 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
219 data.IP_factor);
220 }
221 else // boundary face
222 {
223 tau = integrator_m.read_cell_data(array_penalty_parameter) *
224 IP::get_penalty_factor<dim, Number>(
225 degree,
227 integrator_m.get_matrix_free().get_dof_handler(dof_index).get_triangulation()),
228 data.IP_factor);
229 }
230 }
231
232 /*
233 * This functions return the viscosity for a given cell in a given quadrature point
234 */
235 inline DEAL_II_ALWAYS_INLINE //
236 scalar
237 get_viscosity_cell(unsigned int const cell, unsigned int const q) const
238 {
239 scalar viscosity = dealii::make_vectorized_array<Number>(data.viscosity);
240
241 if(data.viscosity_is_variable)
242 {
243 viscosity = viscosity_coefficients.get_coefficient_cell(cell, q);
244 }
245
246 return viscosity;
247 }
248
249 /*
250 * This function calculates the average viscosity for interior faces.
251 */
252 inline DEAL_II_ALWAYS_INLINE //
253 scalar
254 calculate_average_viscosity(unsigned int const face, unsigned int const q) const
255 {
256 scalar average_viscosity = dealii::make_vectorized_array<Number>(0.0);
257
258 scalar coefficient_face = viscosity_coefficients.get_coefficient_face(face, q);
259 scalar coefficient_face_neighbor =
260 viscosity_coefficients.get_coefficient_face_neighbor(face, q);
261
262 // harmonic mean (harmonic weighting according to Schott and Rasthofer et al. (2015))
263 average_viscosity = 2.0 * coefficient_face * coefficient_face_neighbor /
264 (coefficient_face + coefficient_face_neighbor);
265
266 // arithmetic mean
267 // average_viscosity = 0.5 * (coefficient_face + coefficient_face_neighbor);
268
269 // maximum value
270 // average_viscosity = std::max(coefficient_face, coefficient_face_neighbor);
271
272 return average_viscosity;
273 }
274
275 /*
276 * This function returns the viscosity for interior faces.
277 */
278 inline DEAL_II_ALWAYS_INLINE //
279 scalar
280 get_viscosity_interior_face(unsigned int const face, unsigned int const q) const
281 {
282 scalar viscosity = dealii::make_vectorized_array<Number>(data.viscosity);
283
284 if(data.viscosity_is_variable)
285 {
286 viscosity = calculate_average_viscosity(face, q);
287 }
288
289 return viscosity;
290 }
291
292 /*
293 * This function returns the viscosity for boundary faces.
294 */
295 inline DEAL_II_ALWAYS_INLINE //
296 scalar
297 get_viscosity_boundary_face(unsigned int const face, unsigned int const q) const
298 {
299 scalar viscosity = dealii::make_vectorized_array<Number>(data.viscosity);
300
301 if(data.viscosity_is_variable)
302 {
303 viscosity = viscosity_coefficients.get_coefficient_face(face, q);
304 }
305
306 return viscosity;
307 }
308
309 /*
310 * Volume flux, i.e., the term occurring in the volume integral
311 */
312 inline DEAL_II_ALWAYS_INLINE //
313 tensor
314 get_volume_flux(tensor const & gradient, scalar const & viscosity) const
315 {
316 if(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation)
317 {
318 return viscosity * (gradient + transpose(gradient));
319 }
320 else if(data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation)
321 {
322 return viscosity * gradient;
323 }
324 else
325 {
326 AssertThrow(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation or
327 data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation,
328 dealii::ExcMessage("Specified formulation of viscous term is not implemented."));
329
330 return tensor();
331 }
332 }
333
338 inline DEAL_II_ALWAYS_INLINE //
339 tensor
340 calculate_gradient_flux(vector const & value_m,
341 vector const & value_p,
342 vector const & normal,
343 scalar const & viscosity) const
344 {
345 tensor flux;
346
347 vector jump_value = value_m - value_p;
348 tensor jump_tensor = outer_product(jump_value, normal);
349
350 if(data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation)
351 {
352 if(data.IP_formulation == InteriorPenaltyFormulation::NIPG)
353 {
354 flux = 0.5 * viscosity * jump_tensor;
355 }
356 else if(data.IP_formulation == InteriorPenaltyFormulation::SIPG)
357 {
358 flux = -0.5 * viscosity * jump_tensor;
359 }
360 else
361 {
362 AssertThrow(
363 false, dealii::ExcMessage("Specified interior penalty formulation is not implemented."));
364 }
365 }
366 else if(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation)
367 {
368 if(data.IP_formulation == InteriorPenaltyFormulation::NIPG)
369 {
370 flux = 0.5 * viscosity * (jump_tensor + transpose(jump_tensor));
371 }
372 else if(data.IP_formulation == InteriorPenaltyFormulation::SIPG)
373 {
374 flux = -0.5 * viscosity * (jump_tensor + transpose(jump_tensor));
375 }
376 else
377 {
378 AssertThrow(
379 false, dealii::ExcMessage("Specified interior penalty formulation is not implemented."));
380 }
381 }
382 else
383 {
384 AssertThrow(false,
385 dealii::ExcMessage("Specified formulation of viscous term is not implemented."));
386 }
387
388 return flux;
389 }
390
391 /*
392 * This function calculates the gradient in normal direction on element e
393 * depending on the formulation of the viscous term.
394 */
395 template<typename Integrator>
396 inline DEAL_II_ALWAYS_INLINE //
397 vector
398 calculate_normal_gradient(unsigned int const q, Integrator & integrator) const
399 {
400 tensor gradient;
401
402 if(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation)
403 {
404 /*
405 * F = 2 * nu * symmetric_gradient
406 * = 2.0 * nu * 1/2 (grad(u) + grad(u)^T)
407 */
408 gradient = dealii::make_vectorized_array<Number>(2.0) * integrator.get_symmetric_gradient(q);
409 }
410 else if(data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation)
411 {
412 /*
413 * F = nu * grad(u)
414 */
415 gradient = integrator.get_gradient(q);
416 }
417 else
418 {
419 AssertThrow(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation or
420 data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation,
421 dealii::ExcMessage("Specified formulation of viscous term is not implemented."));
422 }
423
424 vector normal_gradient = gradient * integrator.get_normal_vector(q);
425
426 return normal_gradient;
427 }
428
429
436 inline DEAL_II_ALWAYS_INLINE //
437 vector
438 calculate_value_flux(vector const & normal_gradient_m,
439 vector const & normal_gradient_p,
440 vector const & value_m,
441 vector const & value_p,
442 vector const & normal,
443 scalar const & viscosity) const
444 {
445 vector flux;
446
447 vector jump_value = value_m - value_p;
448 vector average_normal_gradient = 0.5 * (normal_gradient_m + normal_gradient_p);
449
450 if(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation)
451 {
452 if(data.penalty_term_div_formulation == PenaltyTermDivergenceFormulation::Symmetrized)
453 {
454 flux = viscosity * average_normal_gradient -
455 viscosity * tau * (jump_value + (jump_value * normal) * normal);
456 }
457 else if(data.penalty_term_div_formulation == PenaltyTermDivergenceFormulation::NotSymmetrized)
458 {
459 flux = viscosity * average_normal_gradient - viscosity * tau * jump_value;
460 }
461 else
462 {
463 AssertThrow(
464 data.penalty_term_div_formulation == PenaltyTermDivergenceFormulation::Symmetrized or
465 data.penalty_term_div_formulation == PenaltyTermDivergenceFormulation::NotSymmetrized,
466 dealii::ExcMessage("Specified formulation of viscous term is not implemented."));
467 }
468 }
469 else if(data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation)
470 {
471 flux = viscosity * average_normal_gradient - viscosity * tau * jump_value;
472 }
473 else
474 {
475 AssertThrow(data.formulation_viscous_term == FormulationViscousTerm::DivergenceFormulation or
476 data.formulation_viscous_term == FormulationViscousTerm::LaplaceFormulation,
477 dealii::ExcMessage("Specified formulation of viscous term is not implemented."));
478 }
479
480 return flux;
481 }
482
483 // clang-format off
484 /*
485 * This function calculates the interior velocity gradient in normal
486 * direction depending on the operator type.
487 *
488 * The variable normal_gradient_m has the meaning of F(u⁻)*n with
489 *
490 * Divergence formulation: F(u) = F_nu(u) / nu = nu * ( grad(u) + grad(u)^T )
491 * Laplace formulation: F(u) = F_nu(u) / nu = nu * grad(u)
492 *
493 */
494 // clang-format on
495 template<typename Integrator>
496 inline DEAL_II_ALWAYS_INLINE //
497 vector
498 calculate_interior_normal_gradient(unsigned int const q,
499 Integrator const & integrator,
500 OperatorType const & operator_type) const
501 {
502 vector normal_gradient_m;
503
504 if(operator_type == OperatorType::full or operator_type == OperatorType::homogeneous)
505 {
506 normal_gradient_m = calculate_normal_gradient(q, integrator);
507 }
508 else if(operator_type == OperatorType::inhomogeneous)
509 {
510 // do nothing, normal_gradient_m is already initialized with 0
511 }
512 else
513 {
514 AssertThrow(false, dealii::ExcMessage("Specified OperatorType is not implemented!"));
515 }
516
517 return normal_gradient_m;
518 }
519
520private:
521 ViscousKernelData data;
522
523 unsigned int quad_index;
524
525 unsigned int degree;
526
527 dealii::AlignedVector<scalar> array_penalty_parameter;
528
529 mutable scalar tau;
530
531 VariableCoefficients<dealii::VectorizedArray<Number>> viscosity_coefficients;
532};
533
534} // namespace Operators
535
536template<int dim>
538{
540 {
541 }
542
544
545 std::shared_ptr<BoundaryDescriptorU<dim> const> bc;
546};
547
548template<int dim, typename Number>
549class ViscousOperator : public OperatorBase<dim, Number, dim>
550{
551public:
552 typedef dealii::VectorizedArray<Number> scalar;
553 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
554 typedef dealii::Tensor<2, dim, dealii::VectorizedArray<Number>> tensor;
555
557
558 typedef typename Base::VectorType VectorType;
559 typedef typename Base::Range Range;
560 typedef typename Base::IntegratorCell IntegratorCell;
561 typedef typename Base::IntegratorFace IntegratorFace;
562
563 void
564 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
565 dealii::AffineConstraints<Number> const & affine_constraints,
566 ViscousOperatorData<dim> const & data,
567 std::shared_ptr<Operators::ViscousKernel<dim, Number>> viscous_kernel);
568
569 void
570 update();
571
572private:
573 void
574 reinit_face_derived(IntegratorFace & integrator_m,
575 IntegratorFace & integrator_p,
576 unsigned int const face) const final;
577
578 void
579 reinit_boundary_face_derived(IntegratorFace & integrator_m, unsigned int const face) const final;
580
581 void
582 reinit_face_cell_based_derived(IntegratorFace & integrator_m,
583 IntegratorFace & integrator_p,
584 unsigned int const cell,
585 unsigned int const face,
586 dealii::types::boundary_id const boundary_id) const final;
587
588 void
589 do_cell_integral(IntegratorCell & integrator) const final;
590
591 void
592 do_face_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
593
594 void
595 do_face_int_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
596
597 void
598 do_face_ext_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
599
600 void
601 do_boundary_integral(IntegratorFace & integrator,
602 OperatorType const & operator_type,
603 dealii::types::boundary_id const & boundary_id) const final;
604
605 ViscousOperatorData<dim> operator_data;
606
607 std::shared_ptr<Operators::ViscousKernel<dim, Number>> kernel;
608};
609
610} // namespace IncNS
611} // namespace ExaDG
612
613#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SPATIAL_DISCRETIZATION_OPERATORS_VISCOUS_OPERATOR_H_ \
614 */
Definition viscous_operator.h:49
DEAL_II_ALWAYS_INLINE vector calculate_value_flux(vector const &normal_gradient_m, vector const &normal_gradient_p, vector const &value_m, vector const &value_p, vector const &normal, scalar const &viscosity) const
Definition viscous_operator.h:438
DEAL_II_ALWAYS_INLINE tensor calculate_gradient_flux(vector const &value_m, vector const &value_p, vector const &normal, scalar const &viscosity) const
Definition viscous_operator.h:340
Definition viscous_operator.h:550
Definition operator_base.h:98
coefficient_type get_coefficient_cell(unsigned int const cell, unsigned int const q) const
Definition variable_coefficients.h:86
coefficient_type get_coefficient_face_neighbor(unsigned int const face, unsigned int const q) const
Definition variable_coefficients.h:127
void initialize(dealii::MatrixFree< dim, Number > const &matrix_free, unsigned int const quad_index, bool const store_face_data_in, bool const store_cell_based_face_data_in)
Definition variable_coefficients.h:57
void set_coefficient_cell(unsigned int const cell, unsigned int const q, coefficient_type const &coefficient)
Definition variable_coefficients.h:96
void set_coefficient_face(unsigned int const face, unsigned int const q, coefficient_type const &coefficient)
Definition variable_coefficients.h:116
void set_coefficients(coefficient_type const &constant_coefficient)
Definition variable_coefficients.h:77
void set_coefficient_face_neighbor(unsigned int const face, unsigned int const q, coefficient_type const &coefficient)
Definition variable_coefficients.h:137
coefficient_type get_coefficient_face(unsigned int const face, unsigned int const q) const
Definition variable_coefficients.h:107
Definition driver.cpp:33
ElementType get_element_type(dealii::Triangulation< dim > const &tria)
Definition grid_data.h:70
Definition viscous_operator.h:26
Definition viscous_operator.h:538
Definition integrator_flags.h:30
Definition mapping_flags.h:31
Definition operator_base.h:59