16 #include <deal.II/base/polynomials_piecewise.h> 19 DEAL_II_NAMESPACE_OPEN
25 template <
typename number>
28 const unsigned int n_intervals,
29 const unsigned int interval,
30 const bool spans_next_interval)
31 : polynomial(coefficients_on_interval)
32 , n_intervals(n_intervals)
34 , spans_two_intervals(spans_next_interval)
42 template <
typename number>
45 std::vector<number> &values)
const 49 value(x, values.size() - 1, values.data());
54 template <
typename number>
57 const unsigned int n_derivatives,
58 number * values)
const 62 double derivative_change_sign = 1.;
65 const number step = 1. / n_intervals;
67 if (spans_two_intervals)
69 const double offset = step * interval;
70 if (x < offset || x > offset + step + step)
72 for (
unsigned int k = 0; k <= n_derivatives; ++k)
76 else if (x < offset + step)
80 derivative_change_sign = -1.;
81 y = offset + step + step - x;
86 const double offset = step * interval;
87 if (x < offset || x > offset + step)
89 for (
unsigned int k = 0; k <= n_derivatives; ++k)
99 if ((std::abs(y) < 1e-14 &&
100 (interval > 0 || derivative_change_sign == -1.)) ||
101 (std::abs(y - step) < 1e-14 &&
102 (interval < n_intervals - 1 || derivative_change_sign == -1.)))
104 values[0] = value(x);
105 for (
unsigned int d = 1; d <= n_derivatives; ++d)
111 polynomial.value(y, n_derivatives, values);
114 for (
unsigned int j = 1; j <= n_derivatives; j += 2)
115 values[j] *= derivative_change_sign;
120 std::vector<PiecewisePolynomial<double>>
122 const unsigned int n_subdivisions,
123 const unsigned int base_degree)
125 std::vector<Polynomial<double>> p_base =
126 LagrangeEquidistant::generate_complete_basis(base_degree);
127 for (
auto &polynomial : p_base)
128 polynomial.scale(n_subdivisions);
130 std::vector<PiecewisePolynomial<double>> p;
131 p.reserve(n_subdivisions * base_degree + 1);
133 p.emplace_back(p_base[0], n_subdivisions, 0,
false);
134 for (
unsigned int s = 0; s < n_subdivisions; ++s)
135 for (
unsigned int i = 0; i < base_degree; ++i)
136 p.emplace_back(p_base[i + 1],
139 i == (base_degree - 1) && s < n_subdivisions - 1);
149 template class PiecewisePolynomial<float>;
150 template class PiecewisePolynomial<double>;
151 template class PiecewisePolynomial<long double>;
154 DEAL_II_NAMESPACE_CLOSE
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
std::vector< PiecewisePolynomial< double > > generate_complete_Lagrange_basis_on_subdivisions(const unsigned int n_subdivisions, const unsigned int base_degree)
PiecewisePolynomial(const Polynomial< number > &coefficients_on_interval, const unsigned int n_intervals, const unsigned int interval, const bool spans_next_interval)
static ::ExceptionBase & ExcZero()