ryujin 2.1.1 revision a15074459a388761bd8df6bd4ef7e6abe9d8077b
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 "flux_library.h"
9
10#include <convenience_macros.h>
11#include <deal.II/base/config.h>
12#include <discretization.h>
14#include <openmp.h>
15#include <patterns_conversion.h>
16#include <simd.h>
17
18#include <deal.II/base/parameter_acceptor.h>
19#include <deal.II/base/tensor.h>
20
21#include <array>
22
23namespace ryujin
24{
25 namespace ScalarConservation
26 {
27 template <int dim, typename Number>
28 class HyperbolicSystemView;
29
36 class HyperbolicSystem final : public dealii::ParameterAcceptor
37 {
38 public:
42 static inline std::string problem_name = "Scalar conservation equation";
43
47 HyperbolicSystem(const std::string &subsection = "/HyperbolicSystem");
48
55 template <int dim, typename Number>
56 auto view() const
57 {
59 }
60
61 private:
66 std::string flux_;
67
69 using Flux = FluxLibrary::Flux;
70 std::shared_ptr<Flux> selected_flux_;
71
72 template <int dim, typename Number>
75 }; /* HyperbolicSystem */
76
77
83 template <int dim, typename Number>
85 {
86 public:
91 HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
92 : hyperbolic_system_(hyperbolic_system)
93 {
94 }
95
99 template <int dim2, typename Number2>
100 auto view() const
101 {
102 return HyperbolicSystemView<dim2, Number2>{hyperbolic_system_};
103 }
104
109
114
115 DEAL_II_ALWAYS_INLINE inline const std::string &flux() const
116 {
117 return hyperbolic_system_.flux_;
118 }
119
120 DEAL_II_ALWAYS_INLINE inline ScalarNumber
122 {
123 const auto &flux = hyperbolic_system_.selected_flux_;
124 return ScalarNumber(flux->derivative_approximation_delta());
125 }
126
128
132
136 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
137 flux_function(const Number &u) const;
138
142 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
143 flux_gradient_function(const Number &u) const;
144
146
150
151 private:
152 const HyperbolicSystem &hyperbolic_system_;
153
154 public:
156
160
164 static constexpr unsigned int problem_dimension = 1;
165
170 using state_type = dealii::Tensor<1, problem_dimension, Number>;
171
176
181 static inline const auto component_names =
182 std::array<std::string, problem_dimension>{"u"};
183
188 static inline const auto primitive_component_names =
189 std::array<std::string, problem_dimension>{"u"};
190
194 using flux_type =
195 dealii::Tensor<1, problem_dimension, dealii::Tensor<1, dim, Number>>;
196
201
203
207
211 static constexpr unsigned int n_precomputed_initial_values = 0;
212
217 std::array<Number, n_precomputed_initial_values>;
218
225
229 static inline const auto precomputed_initial_names =
230 std::array<std::string, n_precomputed_initial_values>{};
231
235 static constexpr unsigned int n_precomputed_values = 2. * dim;
236
240 using precomputed_state_type = std::array<Number, n_precomputed_values>;
241
247
251 static inline const auto precomputed_names =
252 []() -> std::array<std::string, n_precomputed_values> {
253 if constexpr (dim == 1)
254 return {"f", "df"};
255 else if constexpr (dim == 2)
256 return {"f_1", "f_2", "df_1", "df_2"};
257 else if constexpr (dim == 3)
258 return {"f_1", "f_2", "f_3", "df_1", "df_2", "df_3"};
259 __builtin_trap();
260 }();
261
265 static constexpr unsigned int n_precomputation_cycles = 1;
266
271 template <typename DISPATCH, typename SPARSITY>
272 void precomputation_loop(unsigned int /*cycle*/,
273 const DISPATCH &dispatch_check,
274 precomputed_vector_type & /*precomputed_values*/,
275 const SPARSITY & /*sparsity_simd*/,
276 const vector_type & /*U*/,
277 unsigned int /*left*/,
278 unsigned int /*right*/) const;
279
281
285
292 static Number state(const state_type &U);
293
300 Number square_entropy(const Number &u) const;
301
309 Number square_entropy_derivative(const Number &u) const;
310
317 Number kruzkov_entropy(const Number &k, const Number &u) const;
318
326 Number kruzkov_entropy_derivative(const Number &k, const Number &u) const;
327
333 bool is_admissible(const state_type & /*U*/) const
334 {
335 return true;
336 }
337
339
343
347 template <typename Lambda>
349 const dealii::types::boundary_id /*id*/,
350 const state_type &U,
351 const dealii::Tensor<1, dim, Number> & /*normal*/,
352 const Lambda & /*get_dirichlet_data*/) const;
353
355
359
364 dealii::Tensor<1, dim, Number> construct_flux_tensor(
365 const precomputed_state_type &precomputed_state) const;
366
372 dealii::Tensor<1, dim, Number> construct_flux_gradient_tensor(
373 const precomputed_state_type &precomputed_state) const;
374
396 const precomputed_initial_vector_type & /*piv*/,
397 const unsigned int i,
398 const state_type & /*U_i*/) const;
399
402 const precomputed_initial_vector_type & /*piv*/,
403 const unsigned int *js,
404 const state_type & /*U_j*/) const;
405
412 const flux_contribution_type &flux_j,
413 const dealii::Tensor<1, dim, Number> &c_ij) const;
414
416 static constexpr bool have_high_order_flux = false;
417
421 const dealii::Tensor<1, dim, Number> &c_ij) const = delete;
422
424
428
430 static constexpr bool have_source_terms = false;
431
433 const unsigned int i,
434 const state_type &U_i,
435 const ScalarNumber tau) const = delete;
436
438 const unsigned int *js,
439 const state_type &U_j,
440 const ScalarNumber tau) const = delete;
441
443
447
458 template <typename ST>
460 {
461 return state;
462 }
463
468 state_type from_primitive_state(const state_type &primitive_state) const
469 {
470 return primitive_state;
471 }
472
478 {
479 return state;
480 }
481
487 template <typename Lambda>
489 const Lambda & /*lambda*/) const
490 {
491 return state;
492 }
493
494 }; /* HyperbolicSystemView */
495
496
497 /*
498 * -------------------------------------------------------------------------
499 * Inline definitions
500 * -------------------------------------------------------------------------
501 */
502
503
504 inline HyperbolicSystem::HyperbolicSystem(const std::string &subsection)
505 : ParameterAcceptor(subsection)
506 {
507 flux_ = "burgers";
508 add_parameter("flux",
509 flux_,
510 "The scalar flux. Valid names are given by any of the "
511 "subsections defined below");
512
513 /*
514 * And finally populate the flux list with all flux configurations
515 * defined in the FluxLibrary namespace:
516 */
517 FluxLibrary::populate_flux_list(flux_list_, subsection);
518
519 const auto populate_functions = [this]() {
520 bool initialized = false;
521 for (auto &it : flux_list_)
522
523 /* Populate flux functions: */
524 if (it->name() == flux_) {
525 selected_flux_ = it;
526 it->parse_parameters_call_back();
527 problem_name = "Scalar conservation equation (" + it->name() +
528 ": " + it->flux_formula() + ")";
529 initialized = true;
530 break;
531 }
532
533 AssertThrow(initialized,
534 dealii::ExcMessage(
535 "Could not find a flux description with name \"" +
536 flux_ + "\""));
537 };
538
539 ParameterAcceptor::parse_parameters_call_back.connect(populate_functions);
540 populate_functions();
541 }
542
543
544 template <int dim, typename Number>
545 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
547 {
548 const auto &flux = hyperbolic_system_.selected_flux_;
549 dealii::Tensor<1, dim, Number> result;
550
551 /* This access by calling into value() repeatedly is terrible: */
552 for (unsigned int k = 0; k < dim; ++k) {
553 if constexpr (std::is_same_v<ScalarNumber, Number>) {
554 result[k] = flux->value(u, k);
555 } else {
556 for (unsigned int s = 0; s < Number::size(); ++s) {
557 result[k][s] = flux->value(u[s], k);
558 }
559 }
560 }
561
562 return result;
563 }
564
565
566 template <int dim, typename Number>
567 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
569 const Number &u) const
570 {
571 const auto &flux = hyperbolic_system_.selected_flux_;
572 dealii::Tensor<1, dim, Number> result;
573
574 /* This access by calling into value() repeatedly is terrible: */
575 for (unsigned int k = 0; k < dim; ++k) {
576 if constexpr (std::is_same_v<ScalarNumber, Number>) {
577 result[k] = flux->gradient(u, k);
578 } else {
579 for (unsigned int s = 0; s < Number::size(); ++s) {
580 result[k][s] = flux->gradient(u[s], k);
581 }
582 }
583 }
584
585 return result;
586 }
587
588
589 template <int dim, typename Number>
590 template <typename DISPATCH, typename SPARSITY>
591 DEAL_II_ALWAYS_INLINE inline void
593 unsigned int cycle [[maybe_unused]],
594 const DISPATCH &dispatch_check,
595 precomputed_vector_type &precomputed_values,
596 const SPARSITY &sparsity_simd,
597 const vector_type &U,
598 unsigned int left,
599 unsigned int right) const
600 {
601 Assert(cycle == 0, dealii::ExcInternalError());
602
603 /* We are inside a thread parallel context */
604
605 unsigned int stride_size = get_stride_size<Number>;
606
608 for (unsigned int i = left; i < right; i += stride_size) {
609
610 /* Skip constrained degrees of freedom: */
611 const unsigned int row_length = sparsity_simd.row_length(i);
612 if (row_length == 1)
613 continue;
614
615 dispatch_check(i);
616
617 const auto U_i = U.template get_tensor<Number>(i);
618 const auto u_i = state(U_i);
619
620 const auto f_i = flux_function(u_i);
621 const auto df_i = flux_gradient_function(u_i);
622
624
625 for (unsigned int k = 0; k < n_precomputed_values / 2; ++k) {
626 prec_i[k] = f_i[k];
627 prec_i[dim + k] = df_i[k];
628 }
629
630 precomputed_values.template write_tensor<Number>(prec_i, i);
631 }
632 }
633
634
635 template <int dim, typename Number>
636 DEAL_II_ALWAYS_INLINE inline Number
638 {
639 return U[0];
640 }
641
642
643 template <int dim, typename Number>
644 DEAL_II_ALWAYS_INLINE inline Number
646 {
647 return ScalarNumber(0.5) * u * u;
648 }
649
650
651 template <int dim, typename Number>
652 DEAL_II_ALWAYS_INLINE inline Number
654 const Number &u) const
655 {
656 return u;
657 }
658
659
660 template <int dim, typename Number>
661 DEAL_II_ALWAYS_INLINE inline Number
663 const Number &u) const
664 {
665 return std::abs(k - u);
666 }
667
668
669 template <int dim, typename Number>
670 DEAL_II_ALWAYS_INLINE inline Number
672 const Number &k, const Number &u) const
673 {
674 constexpr auto gte = dealii::SIMDComparison::greater_than_or_equal;
675 // return sgn(u-k):
676 return dealii::compare_and_apply_mask<gte>(u, k, Number(1.), Number(-1.));
677 }
678
679
680 template <int dim, typename Number>
681 template <typename Lambda>
682 DEAL_II_ALWAYS_INLINE inline auto
684 dealii::types::boundary_id id,
685 const state_type &U,
686 const dealii::Tensor<1, dim, Number> & /*normal*/,
687 const Lambda &get_dirichlet_data) const -> state_type
688 {
689 state_type result = U;
690
691 if (id == Boundary::dirichlet) {
692 result = get_dirichlet_data();
693
694 } else if (id == Boundary::slip) {
695 AssertThrow(
696 false,
697 dealii::ExcMessage("Invalid boundary ID »Boundary::slip«, slip "
698 "boundary conditions are unavailable for scalar "
699 "conservation equations."));
700 __builtin_trap();
701
702 } else if (id == Boundary::no_slip) {
703 AssertThrow(
704 false,
705 dealii::ExcMessage("Invalid boundary ID »Boundary::no_slip«, "
706 "no-slip boundary conditions are unavailable "
707 "for scalar conservation equations."));
708 __builtin_trap();
709
710 } else if (id == Boundary::dynamic) {
711 AssertThrow(
712 false,
713 dealii::ExcMessage("Invalid boundary ID »Boundary::dynamic«, "
714 "dynamic boundary conditions are unavailable "
715 "for scalar conservation equations."));
716 __builtin_trap();
717 }
718
719 return result;
720 }
721
722
723 template <int dim, typename Number>
724 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
726 const precomputed_state_type &precomputed_state) const
727 {
728 dealii::Tensor<1, dim, Number> result;
729
730 if constexpr (dim == 1) {
731 const auto &[f, df] = precomputed_state;
732 result[0] = f;
733
734 } else if constexpr (dim == 2) {
735 const auto &[f_1, f_2, df_1, df_2] = precomputed_state;
736 result[0] = f_1;
737 result[1] = f_2;
738
739 } else if constexpr (dim == 3) {
740 const auto &[f_1, f_2, f_3, df_1, df_2, df_3] = precomputed_state;
741 result[0] = f_1;
742 result[1] = f_2;
743 result[2] = f_3;
744 }
745
746 return result;
747 }
748
749
750 template <int dim, typename Number>
751 DEAL_II_ALWAYS_INLINE inline dealii::Tensor<1, dim, Number>
753 const precomputed_state_type &precomputed_state) const
754 {
755 dealii::Tensor<1, dim, Number> result;
756
757 if constexpr (dim == 1) {
758 const auto &[f, df] = precomputed_state;
759 result[0] = df;
760
761 } else if constexpr (dim == 2) {
762 const auto &[f_1, f_2, df_1, df_2] = precomputed_state;
763 result[0] = df_1;
764 result[1] = df_2;
765
766 } else if constexpr (dim == 3) {
767 const auto &[f_1, f_2, f_3, df_1, df_2, df_3] = precomputed_state;
768 result[0] = df_1;
769 result[1] = df_2;
770 result[2] = df_3;
771 }
772
773 return result;
774 }
775
776
777 template <int dim, typename Number>
778 DEAL_II_ALWAYS_INLINE inline auto
780 const precomputed_vector_type &pv,
781 const precomputed_initial_vector_type & /*piv*/,
782 const unsigned int i,
783 const state_type & /*U_i*/) const -> flux_contribution_type
784 {
785 /* The flux contribution is a rank 2 tensor, thus a little dance: */
787 result[0] = construct_flux_tensor(
788 pv.template get_tensor<Number, precomputed_state_type>(i));
789 return result;
790 }
791
792
793 template <int dim, typename Number>
794 DEAL_II_ALWAYS_INLINE inline auto
796 const precomputed_vector_type &pv,
797 const precomputed_initial_vector_type & /*piv*/,
798 const unsigned int *js,
799 const state_type & /*U_j*/) const -> flux_contribution_type
800 {
801 /* The flux contribution is a rank 2 tensor, thus a little dance: */
803 result[0] = construct_flux_tensor(
804 pv.template get_tensor<Number, precomputed_state_type>(js));
805 return result;
806 }
807
808
809 template <int dim, typename Number>
810 DEAL_II_ALWAYS_INLINE inline auto
812 const flux_contribution_type &flux_i,
813 const flux_contribution_type &flux_j,
814 const dealii::Tensor<1, dim, Number> &c_ij) const -> state_type
815 {
816 return -contract(add(flux_i, flux_j), c_ij);
817 }
818
819 } // namespace ScalarConservation
820} // namespace ryujin
std::array< Number, n_precomputed_values > precomputed_state_type
dealii::Tensor< 1, problem_dimension, Number > state_type
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
dealii::Tensor< 1, dim, Number > construct_flux_gradient_tensor(const precomputed_state_type &precomputed_state) const
DEAL_II_ALWAYS_INLINE const std::string & flux() const
static constexpr unsigned int n_precomputed_initial_values
HyperbolicSystemView(const HyperbolicSystem &hyperbolic_system)
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_function(const Number &u) const
flux_contribution_type flux_contribution(const precomputed_vector_type &pv, const precomputed_initial_vector_type &, const unsigned int i, const state_type &) const
dealii::Tensor< 1, problem_dimension, dealii::Tensor< 1, dim, Number > > flux_type
state_type nodal_source(const precomputed_vector_type &pv, const unsigned int *js, const state_type &U_j, const ScalarNumber tau) const =delete
state_type from_primitive_state(const state_type &primitive_state) const
state_type apply_boundary_conditions(const dealii::types::boundary_id, const state_type &U, const dealii::Tensor< 1, dim, Number > &, const Lambda &) const
Number square_entropy_derivative(const Number &u) const
state_type nodal_source(const precomputed_vector_type &pv, const unsigned int i, const state_type &U_i, const ScalarNumber tau) const =delete
state_type to_primitive_state(const state_type &state) const
DEAL_II_ALWAYS_INLINE ScalarNumber derivative_approximation_delta() const
dealii::Tensor< 1, dim, Number > construct_flux_tensor(const precomputed_state_type &precomputed_state) 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
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim, Number > flux_gradient_function(const Number &u) const
static constexpr unsigned int n_precomputed_values
dealii::Tensor< 1, problem_dimension, Number > state_type
std::array< Number, n_precomputed_values > precomputed_state_type
std::array< Number, n_precomputed_initial_values > precomputed_initial_state_type
void precomputation_loop(unsigned int, const DISPATCH &dispatch_check, precomputed_vector_type &, const SPARSITY &, const vector_type &, unsigned int, unsigned int) const
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
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")
#define RYUJIN_OMP_FOR
Definition: openmp.h:70
std::set< std::shared_ptr< Flux > > flux_list_type
Definition: flux_library.h:14
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)