ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
flux_function.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include "flux.h"
9
10#include <deal.II/base/function_parser.h>
11
12#include <boost/algorithm/string/classification.hpp>
13#include <boost/algorithm/string/split.hpp>
14
15#include <numeric>
16
17namespace ryujin
18{
19 namespace FluxLibrary
20 {
26 class Function : public Flux
27 {
28 public:
29 Function(const std::string &subsection)
30 : Flux("function", subsection)
31 {
32 expression_ = "0.5*u*u";
33 add_parameter("expression",
34 expression_,
35 "A function expression of the flux as a function of "
36 "state used to create a muparser object to evaluate the "
37 "flux. For two, or three dimensional fluxes, components "
38 "are separated with a semicolon (;).");
39
40 this->derivative_approximation_delta_ = 1.0e-10;
41 add_parameter("derivative approximation delta",
43 "Step size of the central difference quotient to compute "
44 "an approximation of the flux derivative");
45
46 /*
47 * Set up the muparser object with the final flux description from
48 * the parameter file:
49 */
50 const auto set_up_muparser = [this] {
51 std::vector<std::string> split_expressions;
52 boost::split(split_expressions, expression_, boost::is_any_of(";"));
53
54 const auto size = split_expressions.size();
55
56 Assert(0 < size && size <= 3,
57 dealii::ExcMessage(
58 "user specified flux description must be either one, two, "
59 "or three strings separated by a comma"));
60 flux_function_ = std::make_unique<dealii::FunctionParser<1>>(
61 size, 0.0, this->derivative_approximation_delta_);
62 flux_function_->initialize({"u"}, split_expressions, {});
63
64 flux_formula_ = "f(u)={" + expression_ + "}";
65 };
66
67 set_up_muparser();
68 ParameterAcceptor::parse_parameters_call_back.connect(set_up_muparser);
69 }
70
71
72 double value(const double state,
73 const unsigned int direction) const override
74 {
75 return flux_function_->value(dealii::Point<1>(state), direction);
76 }
77
78
79 double gradient(const double state,
80 const unsigned int direction) const override
81 {
82 return flux_function_->gradient(dealii::Point<1>(state), direction)[0];
83 }
84
85
86 private:
87 std::string expression_;
88
89 std::unique_ptr<dealii::FunctionParser<1>> flux_function_;
90 };
91 } // namespace FluxLibrary
92} // namespace ryujin
std::string flux_formula_
Definition: flux.h:87
double derivative_approximation_delta_
Definition: flux.h:86
Function(const std::string &subsection)
Definition: flux_function.h:29
double gradient(const double state, const unsigned int direction) const override
Definition: flux_function.h:79
double value(const double state, const unsigned int direction) const override
Definition: flux_function.h:72