ExaDG
Loading...
Searching...
No Matches
convective_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 CONV_DIFF_CONVECTION_OPERATOR
23#define CONV_DIFF_CONVECTION_OPERATOR
24
25#include <exadg/convection_diffusion/user_interface/boundary_descriptor.h>
26#include <exadg/convection_diffusion/user_interface/parameters.h>
27#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
28#include <exadg/matrix_free/integrators.h>
29#include <exadg/operators/operator_base.h>
30
31#include <memory>
32
33namespace ExaDG
34{
35namespace ConvDiff
36{
37namespace Operators
38{
39template<int dim>
41{
43 : formulation(FormulationConvectiveTerm::DivergenceFormulation),
44 velocity_type(TypeVelocityField::Function),
45 dof_index_velocity(1),
46 numerical_flux_formulation(NumericalFluxConvectiveOperator::Undefined)
47 {
48 }
49
50 // formulation
51 FormulationConvectiveTerm formulation;
52
53 // analytical vs. numerical velocity field
54 TypeVelocityField velocity_type;
55
56 // TypeVelocityField::DoFVector
57 unsigned int dof_index_velocity;
58
59 // TypeVelocityField::Function
60 std::shared_ptr<dealii::Function<dim>> velocity;
61
62 // numerical flux (e.g., central flux vs. Lax-Friedrichs flux)
63 NumericalFluxConvectiveOperator numerical_flux_formulation;
64};
65
66template<int dim, typename Number>
68{
69private:
70 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
71
72 typedef CellIntegrator<dim, dim, Number> CellIntegratorVelocity;
73 typedef FaceIntegrator<dim, dim, Number> FaceIntegratorVelocity;
74
75 typedef dealii::VectorizedArray<Number> scalar;
76 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
77
78 typedef CellIntegrator<dim, 1, Number> IntegratorCell;
79 typedef FaceIntegrator<dim, 1, Number> IntegratorFace;
80
81public:
82 void
83 reinit(dealii::MatrixFree<dim, Number> const & matrix_free,
84 ConvectiveKernelData<dim> const & data_in,
85 unsigned int const quad_index,
86 bool const use_own_velocity_storage)
87 {
88 data = data_in;
89
90 if(data.velocity_type == TypeVelocityField::DoFVector)
91 {
92 integrator_velocity =
93 std::make_shared<CellIntegratorVelocity>(matrix_free, data.dof_index_velocity, quad_index);
94
95 integrator_velocity_m = std::make_shared<FaceIntegratorVelocity>(matrix_free,
96 true,
97 data.dof_index_velocity,
98 quad_index);
99
100 integrator_velocity_p = std::make_shared<FaceIntegratorVelocity>(matrix_free,
101 false,
102 data.dof_index_velocity,
103 quad_index);
104
105 if(use_own_velocity_storage)
106 {
107 velocity.reset();
108 matrix_free.initialize_dof_vector(velocity.own(), data.dof_index_velocity);
109 }
110 }
111 }
112
114 get_integrator_flags() const
115 {
116 IntegratorFlags flags;
117
118 if(data.formulation == FormulationConvectiveTerm::DivergenceFormulation)
119 {
120 flags.cell_evaluate = dealii::EvaluationFlags::values;
121 flags.cell_integrate = dealii::EvaluationFlags::gradients;
122 }
123 else if(data.formulation == FormulationConvectiveTerm::ConvectiveFormulation)
124 {
125 flags.cell_evaluate = dealii::EvaluationFlags::gradients;
126 flags.cell_integrate = dealii::EvaluationFlags::values;
127 }
128 else
129 {
130 AssertThrow(false, dealii::ExcMessage("Not implemented."));
131 }
132
133 flags.face_evaluate = dealii::EvaluationFlags::values;
134 flags.face_integrate = dealii::EvaluationFlags::values;
135
136 return flags;
137 }
138
139 static MappingFlags
140 get_mapping_flags()
141 {
142 MappingFlags flags;
143
144 flags.cells = dealii::update_gradients | dealii::update_JxW_values |
145 dealii::update_quadrature_points; // q-points due to analytical velocity field
146 flags.inner_faces = dealii::update_JxW_values | dealii::update_quadrature_points |
147 dealii::update_normal_vectors; // q-points due to analytical velocity field
148 flags.boundary_faces =
149 dealii::update_JxW_values | dealii::update_quadrature_points | dealii::update_normal_vectors;
150
151 return flags;
152 }
153
154 dealii::LinearAlgebra::distributed::Vector<Number> const &
155 get_velocity() const
156 {
157 return *velocity;
158 }
159
160 void
161 set_velocity_copy(VectorType const & velocity_in) const
162 {
163 AssertThrow(data.velocity_type == TypeVelocityField::DoFVector,
164 dealii::ExcMessage("Invalid parameter velocity_type."));
165
166 velocity.own() = velocity_in;
167
168 velocity->update_ghost_values();
169 }
170
171 void
172 set_velocity_ptr(VectorType const & velocity_in) const
173 {
174 AssertThrow(data.velocity_type == TypeVelocityField::DoFVector,
175 dealii::ExcMessage("Invalid parameter velocity_type."));
176
177 velocity.reset(velocity_in);
178
179 velocity->update_ghost_values();
180 }
181
182 void
183 reinit_cell(unsigned int const cell) const
184 {
185 if(data.velocity_type == TypeVelocityField::DoFVector)
186 {
187 integrator_velocity->reinit(cell);
188 integrator_velocity->gather_evaluate(*velocity, dealii::EvaluationFlags::values);
189 }
190 }
191
192 void
193 reinit_face(unsigned int const face) const
194 {
195 if(data.velocity_type == TypeVelocityField::DoFVector)
196 {
197 integrator_velocity_m->reinit(face);
198 integrator_velocity_m->gather_evaluate(*velocity, dealii::EvaluationFlags::values);
199
200 integrator_velocity_p->reinit(face);
201 integrator_velocity_p->gather_evaluate(*velocity, dealii::EvaluationFlags::values);
202 }
203 }
204
205 void
206 reinit_boundary_face(unsigned int const face) const
207 {
208 if(data.velocity_type == TypeVelocityField::DoFVector)
209 {
210 integrator_velocity_m->reinit(face);
211 integrator_velocity_m->gather_evaluate(*velocity, dealii::EvaluationFlags::values);
212 }
213 }
214
215 void
216 reinit_face_cell_based(unsigned int const cell,
217 unsigned int const face,
218 dealii::types::boundary_id const boundary_id) const
219 {
220 if(data.velocity_type == TypeVelocityField::DoFVector)
221 {
222 integrator_velocity_m->reinit(cell, face);
223 integrator_velocity_m->gather_evaluate(*velocity, dealii::EvaluationFlags::values);
224
225 if(boundary_id == dealii::numbers::internal_face_boundary_id) // internal face
226 {
227 // TODO: Matrix-free implementation in deal.II does currently not allow to access data of
228 // the neighboring element in case of cell-based face loops.
229 // integrator_velocity_p->reinit(cell, face);
230 // integrator_velocity_p->gather_evaluate(*velocity, dealii::EvaluationFlags::values);
231 }
232 }
233 }
234
235 /*
236 * This function calculates the numerical flux using the central flux.
237 */
238 inline DEAL_II_ALWAYS_INLINE //
239 scalar
240 calculate_central_flux(scalar const & value_m,
241 scalar const & value_p,
242 scalar const & normal_velocity) const
243 {
244 scalar average_value = 0.5 * (value_m + value_p);
245
246 return normal_velocity * average_value;
247 }
248
249 /*
250 * The same as above, but with discontinuous velocity field.
251 */
252 inline DEAL_II_ALWAYS_INLINE //
253 scalar
254 calculate_central_flux(scalar const & value_m,
255 scalar const & value_p,
256 scalar const & normal_velocity_m,
257 scalar const & normal_velocity_p) const
258 {
259 return 0.5 * (normal_velocity_m * value_m + normal_velocity_p * value_p);
260 }
261
262 /*
263 * This function calculates the numerical flux using the Lax-Friedrichs flux.
264 */
265 inline DEAL_II_ALWAYS_INLINE //
266 scalar
267 calculate_lax_friedrichs_flux(scalar const & value_m,
268 scalar const & value_p,
269 scalar const & normal_velocity) const
270 {
271 scalar average_value = 0.5 * (value_m + value_p);
272 scalar jump_value = value_m - value_p;
273 scalar lambda = std::abs(normal_velocity);
274
275 return normal_velocity * average_value + 0.5 * lambda * jump_value;
276 }
277
278 /*
279 * The same as above, but with discontinuous velocity field.
280 */
281 inline DEAL_II_ALWAYS_INLINE //
282 scalar
283 calculate_lax_friedrichs_flux(scalar const & value_m,
284 scalar const & value_p,
285 scalar const & normal_velocity_m,
286 scalar const & normal_velocity_p) const
287 {
288 scalar jump_value = value_m - value_p;
289 scalar lambda = std::max(std::abs(normal_velocity_m), std::abs(normal_velocity_p));
290
291 return 0.5 * (normal_velocity_m * value_m + normal_velocity_p * value_p) +
292 0.5 * lambda * jump_value;
293 }
294
295 /*
296 * This function calculates the numerical flux where the type of the numerical flux depends on the
297 * specified parameter. This function handles both analytical and numerical velocity fields.
298 */
299 inline DEAL_II_ALWAYS_INLINE //
300 scalar
301 calculate_flux(unsigned int const q,
302 IntegratorFace & integrator,
303 scalar const & value_m,
304 scalar const & value_p,
305 vector const & normal_m,
306 Number const & time,
307 bool const exterior_velocity_available) const
308 {
309 scalar flux = dealii::make_vectorized_array<Number>(0.0);
310
311 if(data.velocity_type == TypeVelocityField::Function)
312 {
313 dealii::Point<dim, scalar> q_points = integrator.quadrature_point(q);
314
315 vector velocity = FunctionEvaluator<1, dim, Number>::value(*(data.velocity), q_points, time);
316
317 scalar normal_velocity = velocity * normal_m;
318
319 // flux functions are the same for DivergenceFormulation and ConvectiveFormulation
320
321 if(data.numerical_flux_formulation == NumericalFluxConvectiveOperator::CentralFlux)
322 {
323 flux = calculate_central_flux(value_m, value_p, normal_velocity);
324 }
325 else if(data.numerical_flux_formulation == NumericalFluxConvectiveOperator::LaxFriedrichsFlux)
326 {
327 flux = calculate_lax_friedrichs_flux(value_m, value_p, normal_velocity);
328 }
329 }
330 else if(data.velocity_type == TypeVelocityField::DoFVector)
331 {
332 vector velocity_m = integrator_velocity_m->get_value(q);
333 vector velocity_p =
334 exterior_velocity_available ? integrator_velocity_p->get_value(q) : velocity_m;
335
336 scalar normal_velocity_m = velocity_m * normal_m;
337 scalar normal_velocity_p = velocity_p * normal_m;
338
339 if(data.formulation == FormulationConvectiveTerm::DivergenceFormulation)
340 {
341 if(data.numerical_flux_formulation == NumericalFluxConvectiveOperator::CentralFlux)
342 {
343 flux = calculate_central_flux(value_m, value_p, normal_velocity_m, normal_velocity_p);
344 }
345 else if(data.numerical_flux_formulation ==
346 NumericalFluxConvectiveOperator::LaxFriedrichsFlux)
347 {
348 flux =
349 calculate_lax_friedrichs_flux(value_m, value_p, normal_velocity_m, normal_velocity_p);
350 }
351 }
352 else if(data.formulation == FormulationConvectiveTerm::ConvectiveFormulation)
353 {
354 scalar normal_velocity = 0.5 * (normal_velocity_m + normal_velocity_p);
355
356 if(data.numerical_flux_formulation == NumericalFluxConvectiveOperator::CentralFlux)
357 {
358 flux = calculate_central_flux(value_m, value_p, normal_velocity);
359 }
360 else if(data.numerical_flux_formulation ==
361 NumericalFluxConvectiveOperator::LaxFriedrichsFlux)
362 {
363 flux = calculate_lax_friedrichs_flux(value_m, value_p, normal_velocity);
364 }
365 }
366 else
367 {
368 AssertThrow(false, dealii::ExcMessage("Not implemented."));
369 }
370 }
371 else
372 {
373 AssertThrow(false, dealii::ExcMessage("Not implemented."));
374 }
375
376 return flux;
377 }
378
379 inline DEAL_II_ALWAYS_INLINE //
380 vector
381 calculate_average_velocity(unsigned int const q,
382 IntegratorFace & integrator,
383 Number const & time,
384 bool const exterior_velocity_available) const
385 {
386 vector velocity;
387
388 if(data.velocity_type == TypeVelocityField::Function)
389 {
390 dealii::Point<dim, scalar> q_points = integrator.quadrature_point(q);
391
392 velocity = FunctionEvaluator<1, dim, Number>::value(*(data.velocity), q_points, time);
393 }
394 else if(data.velocity_type == TypeVelocityField::DoFVector)
395 {
396 vector velocity_m = integrator_velocity_m->get_value(q);
397 vector velocity_p =
398 exterior_velocity_available ? integrator_velocity_p->get_value(q) : velocity_m;
399
400 velocity = 0.5 * (velocity_m + velocity_p);
401 }
402 else
403 {
404 AssertThrow(false, dealii::ExcMessage("Not implemented."));
405 }
406
407 return velocity;
408 }
409
410 inline DEAL_II_ALWAYS_INLINE //
411 std::tuple<scalar, scalar>
412 calculate_flux_interior_and_neighbor(unsigned int const q,
413 IntegratorFace & integrator,
414 scalar const & value_m,
415 scalar const & value_p,
416 vector const & normal_m,
417 Number const & time,
418 bool const exterior_velocity_available) const
419 {
420 scalar fluxM =
421 calculate_flux(q, integrator, value_m, value_p, normal_m, time, exterior_velocity_available);
422 scalar fluxP = -fluxM;
423
424 if(data.formulation == FormulationConvectiveTerm::ConvectiveFormulation)
425 {
426 vector velocity =
427 calculate_average_velocity(q, integrator, time, exterior_velocity_available);
428 scalar normal_velocity = velocity * normal_m;
429
430 // second term appears since the strong formulation is implemented (integration by parts
431 // is performed twice)
432 fluxM = fluxM - normal_velocity * value_m;
433 // opposite signs since n⁺ = - n⁻
434 fluxP = fluxP + normal_velocity * value_p;
435 }
436
437 return std::make_tuple(fluxM, fluxP);
438 }
439
440 inline DEAL_II_ALWAYS_INLINE //
441 scalar
442 calculate_flux_interior(unsigned int const q,
443 IntegratorFace & integrator,
444 scalar const & value_m,
445 scalar const & value_p,
446 vector const & normal_m,
447 Number const & time,
448 bool const exterior_velocity_available) const
449 {
450 scalar flux =
451 calculate_flux(q, integrator, value_m, value_p, normal_m, time, exterior_velocity_available);
452
453 if(data.formulation == FormulationConvectiveTerm::ConvectiveFormulation)
454 {
455 vector velocity =
456 calculate_average_velocity(q, integrator, time, exterior_velocity_available);
457 scalar normal_velocity = velocity * normal_m;
458
459 // second term appears since the strong formulation is implemented (integration by parts
460 // is performed twice)
461 flux = flux - normal_velocity * value_m;
462 }
463
464 return flux;
465 }
466
467
468 /*
469 * Volume flux, i.e., the term occurring in the volume integral
470 */
471 inline DEAL_II_ALWAYS_INLINE //
472 vector
473 get_volume_flux_divergence_form(scalar const & value,
474 IntegratorCell & integrator,
475 unsigned int const q,
476 Number const & time) const
477 {
478 vector velocity;
479
480 if(data.velocity_type == TypeVelocityField::Function)
481 {
482 velocity = FunctionEvaluator<1, dim, Number>::value(*(data.velocity),
483 integrator.quadrature_point(q),
484 time);
485 }
486 else if(data.velocity_type == TypeVelocityField::DoFVector)
487 {
488 velocity = integrator_velocity->get_value(q);
489 }
490 else
491 {
492 AssertThrow(false, dealii::ExcMessage("Not implemented."));
493 }
494
495 return (-value * velocity);
496 }
497
498 inline DEAL_II_ALWAYS_INLINE //
499 scalar
500 get_volume_flux_convective_form(vector const & gradient,
501 IntegratorCell & integrator,
502 unsigned int const q,
503 Number const & time) const
504 {
505 vector velocity;
506
507 if(data.velocity_type == TypeVelocityField::Function)
508 {
509 velocity = FunctionEvaluator<1, dim, Number>::value(*(data.velocity),
510 integrator.quadrature_point(q),
511 time);
512 }
513 else if(data.velocity_type == TypeVelocityField::DoFVector)
514 {
515 velocity = integrator_velocity->get_value(q);
516 }
517 else
518 {
519 AssertThrow(false, dealii::ExcMessage("Not implemented."));
520 }
521
522 return (velocity * gradient);
523 }
524
525private:
527
528 mutable lazy_ptr<VectorType> velocity;
529
530 std::shared_ptr<CellIntegratorVelocity> integrator_velocity;
531 std::shared_ptr<FaceIntegratorVelocity> integrator_velocity_m;
532 std::shared_ptr<FaceIntegratorVelocity> integrator_velocity_p;
533};
534
535} // namespace Operators
536
537
538template<int dim>
540{
542 {
543 }
544
546
547 std::shared_ptr<BoundaryDescriptor<dim> const> bc;
548};
549
550template<int dim, typename Number>
551class ConvectiveOperator : public OperatorBase<dim, Number, 1>
552{
553private:
555
556 typedef typename Base::IntegratorCell IntegratorCell;
557 typedef typename Base::IntegratorFace IntegratorFace;
558
559 typedef typename Base::VectorType VectorType;
560
561 typedef dealii::VectorizedArray<Number> scalar;
562 typedef dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> vector;
563
564public:
565 void
566 initialize(dealii::MatrixFree<dim, Number> const & matrix_free,
567 dealii::AffineConstraints<Number> const & affine_constraints,
568 ConvectiveOperatorData<dim> const & data,
569 std::shared_ptr<Operators::ConvectiveKernel<dim, Number>> kernel);
570
571 dealii::LinearAlgebra::distributed::Vector<Number> const &
572 get_velocity() const;
573
574 void
575 set_velocity_copy(VectorType const & velocity) const;
576
577 void
578 set_velocity_ptr(VectorType const & velocity) const;
579
580private:
581 void
582 reinit_cell_derived(IntegratorCell & integrator, unsigned int const cell) const final;
583
584 void
585 reinit_face_derived(IntegratorFace & integrator_m,
586 IntegratorFace & integrator_p,
587 unsigned int const face) const final;
588
589 void
590 reinit_boundary_face_derived(IntegratorFace & integrator_m, unsigned int const face) const final;
591
592 void
593 reinit_face_cell_based_derived(IntegratorFace & integrator_m,
594 IntegratorFace & integrator_p,
595 unsigned int const cell,
596 unsigned int const face,
597 dealii::types::boundary_id const boundary_id) const final;
598
599 void
600 do_cell_integral(IntegratorCell & integrator) const final;
601
602 void
603 do_face_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
604
605 void
606 do_face_int_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
607
608 void
609 do_face_ext_integral(IntegratorFace & integrator_m, IntegratorFace & integrator_p) const final;
610
611 void
612 do_boundary_integral(IntegratorFace & integrator_m,
613 OperatorType const & operator_type,
614 dealii::types::boundary_id const & boundary_id) const final;
615
616 // TODO can be removed later once matrix-free evaluation allows accessing neighboring data for
617 // cell-based face loops
618 void
619 do_face_int_integral_cell_based(IntegratorFace & integrator_m,
620 IntegratorFace & integrator_p) const final;
621
622 ConvectiveOperatorData<dim> operator_data;
623
624 std::shared_ptr<Operators::ConvectiveKernel<dim, Number>> kernel;
625};
626} // namespace ConvDiff
627} // namespace ExaDG
628
629#endif
Definition convective_operator.h:552
Definition convective_operator.h:68
Definition operator_base.h:98
Definition lazy_ptr.h:29
Definition driver.cpp:33
Definition convective_operator.h:540
Definition convective_operator.h:41
Definition evaluate_functions.h:40
Definition integrator_flags.h:30
Definition mapping_flags.h:31
Definition operator_base.h:59