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