ryujin 2.1.1 revision ef0fcd4010d109b860652ece4a7b8963fb7d46b1
hyperbolic_system.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 - 2024 by the ryujin authors
4//
5
6#pragma once
7
8#include <compile_time_options.h>
9
10#include "flux_library.h"
11
12#include <convenience_macros.h>
13#include <deal.II/base/config.h>
14#include <discretization.h>
16#include <openmp.h>
17#include <patterns_conversion.h>
18#include <simd.h>
19#include <state_vector.h>
20
21#include <deal.II/base/parameter_acceptor.h>
22#include <deal.II/base/tensor.h>
23
24#include <array>
25
26namespace ryujin
27{
28 namespace ScalarConservation
29 {
30 template <int dim, typename Number>
31 class HyperbolicSystemView;
32
39 class HyperbolicSystem final : public dealii::ParameterAcceptor
40 {
41 public:
45 static inline std::string problem_name = "Scalar conservation equation";
46
50 HyperbolicSystem(const std::string &subsection = "/HyperbolicSystem");
51
58 template <int dim, typename Number>
59 auto view() const
60 {
62 }
63
64 private:
69 std::string flux_;
70
72 using Flux = FluxLibrary::Flux;
73 std::shared_ptr<Flux> selected_flux_;
74
75 template <int dim, typename Number>
78 }; /* HyperbolicSystem */
79
80
86 template <int dim, typename Number>
88 {
89 public:
94 HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
95 : hyperbolic_system_(hyperbolic_system)
96 {
97 }
98
102 template <int dim2, typename Number2>
103 auto view() const
104 {
105 return HyperbolicSystemView<dim2, Number2>{hyperbolic_system_};
106 }
107
112
117
118 DEAL_II_ALWAYS_INLINE inline const std::string &flux() const
119 {
120 return hyperbolic_system_.flux_;
121 }
122
123 DEAL_II_ALWAYS_INLINE inline ScalarNumber
125 {
126 const auto &flux = hyperbolic_system_.selected_flux_;
127 return ScalarNumber(flux->derivative_approximation_delta());
128 }
129
131
135
139 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
140 flux_function(const Number &u) const;
141
145 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
146 flux_gradient_function(const Number &u) const;
147
149
153
154 private:
155 const HyperbolicSystem &hyperbolic_system_;
156
157 public:
159
163
167 static constexpr unsigned int problem_dimension = 1;
168
172 using state_type = dealii::Tensor<1, problem_dimension, Number>;
173
177 using flux_type =
178 dealii::Tensor<1, problem_dimension, dealii::Tensor<1, dim, Number>>;
179
184
189 static inline const auto component_names =
190 std::array<std::string, problem_dimension>{"u"};
191
196 static inline const auto primitive_component_names =
197 std::array<std::string, problem_dimension>{"u"};
198
202 static constexpr unsigned int n_precomputed_values = 2 * dim;
203
207 using precomputed_type = std::array<Number, n_precomputed_values>;
208
212 static inline const auto precomputed_names =
213 []() -> std::array<std::string, n_precomputed_values> {
214 if constexpr (dim == 1)
215 return {"f", "df"};
216 else if constexpr (dim == 2)
217 return {"f_1", "f_2", "df_1", "df_2"};
218 else if constexpr (dim == 3)
219 return {"f_1", "f_2", "f_3", "df_1", "df_2", "df_3"};
220 __builtin_trap();
221 }();
222
226 static constexpr unsigned int n_initial_precomputed_values = 0;
227
232 std::array<Number, n_initial_precomputed_values>;
233
237 static inline const auto initial_precomputed_names =
238 std::array<std::string, n_initial_precomputed_values>{};
239
243 using StateVector = Vectors::
244 StateVector<ScalarNumber, problem_dimension, n_precomputed_values>;
245
251
257
265
267
271
272
276 static constexpr unsigned int n_precomputation_cycles = 1;
277
282 template <typename DISPATCH, typename SPARSITY>
283 void precomputation_loop(unsigned int cycle,
284 const DISPATCH &dispatch_check,
285 const SPARSITY &sparsity_simd,
286 StateVector &state_vector,
287 unsigned int left,
288 unsigned int right,
289 const bool skip_constrained_dofs = true) const;
290
292
296
303 static Number state(const state_type &U);
304
311 Number square_entropy(const Number &u) const;
312
320 Number square_entropy_derivative(const Number &u) const;
321
328 Number kruzkov_entropy(const Number &k, const Number &u) const;
329
337 Number kruzkov_entropy_derivative(const Number &k, const Number &u) const;
338
344 bool is_admissible(const state_type & /*U*/) const
345 {
346 return true;
347 }
348
350
354
358 template <typename Lambda>
360 apply_boundary_conditions(const dealii::types::boundary_id id,
361 const state_type &U,
362 const dealii::Tensor<1, dim, Number> &normal,
363 const Lambda &get_dirichlet_data) const;
364
366
370
375 dealii::Tensor<1, dim, Number>
376 construct_flux_tensor(const precomputed_type &precomputed_state) const;
377
383 dealii::Tensor<1, dim, Number> construct_flux_gradient_tensor(
384 const precomputed_type &precomputed_state) const;
385
407 const InitialPrecomputedVector & /*piv*/,
408 const unsigned int i,
409 const state_type & /*U_i*/) const;
410
413 const InitialPrecomputedVector & /*piv*/,
414 const unsigned int *js,
415 const state_type & /*U_j*/) const;
416
423 const flux_contribution_type &flux_j,
424 const dealii::Tensor<1, dim, Number> &c_ij) const;
425
427 static constexpr bool have_high_order_flux = false;
428
432 const dealii::Tensor<1, dim, Number> &c_ij) const = delete;
433
435
439
441 static constexpr bool have_source_terms = false;
442
444 const unsigned int i,
445 const state_type &U_i,
446 const ScalarNumber tau) const = delete;
447
449 const unsigned int *js,
450 const state_type &U_j,
451 const ScalarNumber tau) const = delete;
452
454
458
469 template <typename ST>
471 {
472 return state;
473 }
474
479 state_type from_primitive_state(const state_type &primitive_state) const
480 {
481 return primitive_state;
482 }
483
489 {
490 return state;
491 }
492
498 template <typename Lambda>
500 const Lambda & /*lambda*/) const
501 {
502 return state;
503 }
504
505 }; /* HyperbolicSystemView */
506
507
508 /*
509 * -------------------------------------------------------------------------
510 * Inline definitions
511 * -------------------------------------------------------------------------
512 */
513
514
515 inline HyperbolicSystem::HyperbolicSystem(const std::string &subsection)
516 : ParameterAcceptor(subsection)
517 {
518 flux_ = "burgers";
519 add_parameter("flux",
520 flux_,
521 "The scalar flux. Valid names are given by any of the "
522 "subsections defined below");
523
524 /*
525 * And finally populate the flux list with all flux configurations
526 * defined in the FluxLibrary namespace:
527 */
528 FluxLibrary::populate_flux_list(flux_list_, subsection);
529
530 const auto populate_functions = [this]() {
531 bool initialized = false;
532 for (auto &it : flux_list_)
533
534 /* Populate flux functions: */
535 if (it->name() == flux_) {
536 selected_flux_ = it;
537 it->parse_parameters_call_back();
538 problem_name = "Scalar conservation equation (" + it->name() +
539 ": " + it->flux_formula() + ")";
540 initialized = true;
541 break;
542 }
543
544 AssertThrow(initialized,
545 dealii::ExcMessage(
546 "Could not find a flux description with name \"" +
547 flux_ + "\""));
548 };
549
550 ParameterAcceptor::parse_parameters_call_back.connect(populate_functions);
551 populate_functions();
552 }
553
554
555 template <int dim, typename Number>
556 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
558 {
559 const auto &flux = hyperbolic_system_.selected_flux_;
560 dealii::Tensor<1, dim, Number> result;
561
562 /* This access by calling into value() repeatedly is terrible: */
563 for (unsigned int k = 0; k < dim; ++k) {
564 if constexpr (std::is_same_v<ScalarNumber, Number>) {
565 result[k] = flux->value(u, k);
566 } else {
567 for (unsigned int s = 0; s < Number::size(); ++s) {
568 result[k][s] = flux->value(u[s], k);
569 }
570 }
571 }
572
573 return result;
574 }
575
576
577 template <int dim, typename Number>
578 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
580 const Number &u) const
581 {
582 const auto &flux = hyperbolic_system_.selected_flux_;
583 dealii::Tensor<1, dim, Number> result;
584
585 /* This access by calling into value() repeatedly is terrible: */
586 for (unsigned int k = 0; k < dim; ++k) {
587 if constexpr (std::is_same_v<ScalarNumber, Number>) {
588 result[k] = flux->gradient(u, k);
589 } else {
590 for (unsigned int s = 0; s < Number::size(); ++s) {
591 result[k][s] = flux->gradient(u[s], k);
592 }
593 }
594 }
595
596 return result;
597 }
598
599
600 template <int dim, typename Number>
601 template <typename DISPATCH, typename SPARSITY>
602 DEAL_II_ALWAYS_INLINE inline void
604 unsigned int cycle [[maybe_unused]],
605 const DISPATCH &dispatch_check,
606 const SPARSITY &sparsity_simd,
607 StateVector &state_vector,
608 unsigned int left,
609 unsigned int right,
610 const bool skip_constrained_dofs /*= true*/) const
611 {
612 Assert(cycle == 0, dealii::ExcInternalError());
613
614 const auto &U = std::get<0>(state_vector);
615 auto &precomputed = std::get<1>(state_vector);
616
617 /* We are inside a thread parallel context */
618
619 unsigned int stride_size = get_stride_size<Number>;
620
622 for (unsigned int i = left; i < right; i += stride_size) {
623
624 /* Skip constrained degrees of freedom: */
625 const unsigned int row_length = sparsity_simd.row_length(i);
626 if (skip_constrained_dofs && row_length == 1)
627 continue;
628
629 dispatch_check(i);
630
631 const auto U_i = U.template get_tensor<Number>(i);
632 const auto u_i = state(U_i);
633
634 const auto f_i = flux_function(u_i);
635 const auto df_i = flux_gradient_function(u_i);
636
637 precomputed_type prec_i;
638
639 for (unsigned int k = 0; k < n_precomputed_values / 2; ++k) {
640 prec_i[k] = f_i[k];
641 prec_i[dim + k] = df_i[k];
642 }
643
644 precomputed.template write_tensor<Number>(prec_i, i);
645 }
646 }
647
648
649 template <int dim, typename Number>
650 DEAL_II_ALWAYS_INLINE inline Number
652 {
653 return U[0];
654 }
655
656
657 template <int dim, typename Number>
658 DEAL_II_ALWAYS_INLINE inline Number
660 {
661 return ScalarNumber(0.5) * u * u;
662 }
663
664
665 template <int dim, typename Number>
666 DEAL_II_ALWAYS_INLINE inline Number
668 const Number &u) const
669 {
670 return u;
671 }
672
673
674 template <int dim, typename Number>
675 DEAL_II_ALWAYS_INLINE inline Number
677 const Number &u) const
678 {
679 return std::abs(k - u);
680 }
681
682
683 template <int dim, typename Number>
684 DEAL_II_ALWAYS_INLINE inline Number
686 const Number &k, const Number &u) const
687 {
688 constexpr auto gte = dealii::SIMDComparison::greater_than_or_equal;
689 // return sgn(u-k):
690 return dealii::compare_and_apply_mask<gte>(u, k, Number(1.), Number(-1.));
691 }
692
693
694 template <int dim, typename Number>
695 template <typename Lambda>
696 DEAL_II_ALWAYS_INLINE inline auto
698 dealii::types::boundary_id id,
699 const state_type &U,
700 const dealii::Tensor<1, dim, Number> & /*normal*/,
701 const Lambda &get_dirichlet_data) const -> state_type
702 {
703 state_type result = U;
704
705 if (id == Boundary::dirichlet) {
706 result = get_dirichlet_data();
707
708 } else if (id == Boundary::dirichlet_momentum) {
709 AssertThrow(false,
710 dealii::ExcMessage(
711 "Invalid boundary ID »Boundary::dirichlet_momentum«, "
712 "enforcing Dirichlet boundary conditions on a momentum "
713 "is not possible for scalar conservation equations."));
714
715 } else if (id == Boundary::slip) {
716 AssertThrow(
717 false,
718 dealii::ExcMessage("Invalid boundary ID »Boundary::slip«, slip "
719 "boundary conditions are unavailable for scalar "
720 "conservation equations."));
721 __builtin_trap();
722
723 } else if (id == Boundary::no_slip) {
724 AssertThrow(
725 false,
726 dealii::ExcMessage("Invalid boundary ID »Boundary::no_slip«, "
727 "no-slip boundary conditions are unavailable "
728 "for scalar conservation equations."));
729 __builtin_trap();
730
731 } else if (id == Boundary::dynamic) {
732 AssertThrow(
733 false,
734 dealii::ExcMessage("Invalid boundary ID »Boundary::dynamic«, "
735 "dynamic boundary conditions are unavailable "
736 "for scalar conservation equations."));
737 __builtin_trap();
738
739 } else {
740 AssertThrow(false, dealii::ExcNotImplemented());
741 }
742
743 return result;
744 }
745
746
747 template <int dim, typename Number>
748 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
750 const precomputed_type &precomputed) const
751 {
752 dealii::Tensor<1, dim, Number> result;
753
754 if constexpr (dim == 1) {
755 const auto &[f, df] = precomputed;
756 result[0] = f;
757
758 } else if constexpr (dim == 2) {
759 const auto &[f_1, f_2, df_1, df_2] = precomputed;
760 result[0] = f_1;
761 result[1] = f_2;
762
763 } else if constexpr (dim == 3) {
764 const auto &[f_1, f_2, f_3, df_1, df_2, df_3] = precomputed;
765 result[0] = f_1;
766 result[1] = f_2;
767 result[2] = f_3;
768 }
769
770 return result;
771 }
772
773
774 template <int dim, typename Number>
775 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
777 const precomputed_type &precomputed) const
778 {
779 dealii::Tensor<1, dim, Number> result;
780
781 if constexpr (dim == 1) {
782 const auto &[f, df] = precomputed;
783 result[0] = df;
784
785 } else if constexpr (dim == 2) {
786 const auto &[f_1, f_2, df_1, df_2] = precomputed;
787 result[0] = df_1;
788 result[1] = df_2;
789
790 } else if constexpr (dim == 3) {
791 const auto &[f_1, f_2, f_3, df_1, df_2, df_3] = precomputed;
792 result[0] = df_1;
793 result[1] = df_2;
794 result[2] = df_3;
795 }
796
797 return result;
798 }
799
800
801 template <int dim, typename Number>
802 DEAL_II_ALWAYS_INLINE inline auto
804 const PrecomputedVector &pv,
805 const InitialPrecomputedVector & /*piv*/,
806 const unsigned int i,
807 const state_type & /*U_i*/) const -> flux_contribution_type
808 {
809 /* The flux contribution is a rank 2 tensor, thus a little dance: */
811 result[0] = construct_flux_tensor(
812 pv.template get_tensor<Number, precomputed_type>(i));
813 return result;
814 }
815
816
817 template <int dim, typename Number>
818 DEAL_II_ALWAYS_INLINE inline auto
820 const PrecomputedVector &pv,
821 const InitialPrecomputedVector & /*piv*/,
822 const unsigned int *js,
823 const state_type & /*U_j*/) const -> flux_contribution_type
824 {
825 /* The flux contribution is a rank 2 tensor, thus a little dance: */
827 result[0] = construct_flux_tensor(
828 pv.template get_tensor<Number, precomputed_type>(js));
829 return result;
830 }
831
832
833 template <int dim, typename Number>
834 DEAL_II_ALWAYS_INLINE inline auto
836 const flux_contribution_type &flux_i,
837 const flux_contribution_type &flux_j,
838 const dealii::Tensor<1, dim, Number> &c_ij) const -> state_type
839 {
840 return -contract(add(flux_i, flux_j), c_ij);
841 }
842
843 } // namespace ScalarConservation
844} // namespace ryujin
dealii::Tensor< 1, problem_dimension, Number > state_type
std::array< Number, n_precomputed_values > precomputed_type
Vectors::StateVector< ScalarNumber, problem_dimension, n_precomputed_values > StateVector
flux_contribution_type flux_contribution(const PrecomputedVector &pv, const InitialPrecomputedVector &, const unsigned int i, const state_type &) const
state_type apply_galilei_transform(const state_type &state, const Lambda &) const
static constexpr unsigned int n_precomputation_cycles
Number kruzkov_entropy_derivative(const Number &k, const Number &u) const
DEAL_II_ALWAYS_INLINE const std::string & flux() const
HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_function(const Number &u) const
dealii::Tensor< 1, problem_dimension, dealii::Tensor< 1, dim, Number > > flux_type
std::array< Number, n_initial_precomputed_values > initial_precomputed_type
state_type from_primitive_state(const state_type &primitive_state) const
std::array< Number, n_precomputed_values > precomputed_type
Number square_entropy_derivative(const Number &u) const
dealii::Tensor< 1, dim, Number > construct_flux_gradient_tensor(const precomputed_type &precomputed_state) const
state_type to_primitive_state(const state_type &state) const
DEAL_II_ALWAYS_INLINE ScalarNumber derivative_approximation_delta() const
state_type nodal_source(const PrecomputedVector &pv, const unsigned int *js, const state_type &U_j, const ScalarNumber tau) const =delete
void precomputation_loop(unsigned int cycle, const DISPATCH &dispatch_check, const SPARSITY &sparsity_simd, StateVector &state_vector, unsigned int left, unsigned int right, const bool skip_constrained_dofs=true) const
typename get_value_type< Number >::type ScalarNumber
static Number state(const state_type &U)
Number kruzkov_entropy(const Number &k, const Number &u) const
state_type nodal_source(const PrecomputedVector &pv, const unsigned int i, const state_type &U_i, const ScalarNumber tau) const =delete
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_gradient_function(const Number &u) const
static constexpr unsigned int n_initial_precomputed_values
static constexpr unsigned int n_precomputed_values
dealii::Tensor< 1, problem_dimension, Number > state_type
state_type expand_state(const ST &state) const
state_type high_order_flux_divergence(const flux_contribution_type &, const flux_contribution_type &, const dealii::Tensor< 1, dim, Number > &c_ij) const =delete
Vectors::StateVector< ScalarNumber, problem_dimension, n_precomputed_values > StateVector
dealii::Tensor< 1, dim, Number > construct_flux_tensor(const precomputed_type &precomputed_state) const
state_type apply_boundary_conditions(const dealii::types::boundary_id id, const state_type &U, const dealii::Tensor< 1, dim, Number > &normal, const Lambda &get_dirichlet_data) const
state_type flux_divergence(const flux_contribution_type &flux_i, const flux_contribution_type &flux_j, const dealii::Tensor< 1, dim, Number > &c_ij) const
HyperbolicSystem(const std::string &subsection="/HyperbolicSystem")
@ dirichlet_momentum
#define RYUJIN_OMP_FOR
Definition: openmp.h:70
std::set< std::shared_ptr< Flux > > flux_list_type
Definition: flux_library.h:16
void populate_flux_list(flux_list_type &flux_list, const std::string &subsection)
DEAL_II_ALWAYS_INLINE FT add(const FT &flux_left_ij, const FT &flux_right_ij)
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, problem_dim, T > contract(const FT &flux_ij, const TT &c_ij)