Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
polynomials_piecewise.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2000 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_polynomials_piecewise_h
17#define dealii_polynomials_piecewise_h
18
19
20
21#include <deal.II/base/config.h>
22
24#include <deal.II/base/point.h>
27
28#include <vector>
29
31
41namespace Polynomials
42{
63 template <typename number>
65 {
66 public:
78 PiecewisePolynomial(const Polynomial<number> &coefficients_on_interval,
79 const unsigned int n_intervals,
80 const unsigned int interval,
81 const bool spans_next_interval);
82
89 number
90 value(const number x) const;
91
106 void
107 value(const number x, std::vector<number> &values) const;
108
124 void
125 value(const number x,
126 const unsigned int n_derivatives,
127 number * values) const;
128
133 unsigned int
134 degree() const;
135
141 template <class Archive>
142 void
143 serialize(Archive &ar, const unsigned int version);
144
148 virtual std::size_t
149 memory_consumption() const;
150
151 protected:
157
162 unsigned int n_intervals;
163
168 unsigned int interval;
169
175 };
176
177
178
184 std::vector<PiecewisePolynomial<double>>
186 const unsigned int n_subdivisions,
187 const unsigned int base_degree);
188
189} // namespace Polynomials
190
191
194/* -------------------------- inline functions --------------------- */
195
196namespace Polynomials
197{
198 template <typename number>
199 inline unsigned int
201 {
202 return polynomial.degree();
203 }
204
205
206
207 template <typename number>
208 inline number
210 {
211 AssertIndexRange(interval, n_intervals);
212 number y = x;
213 // shift polynomial if necessary
214 if (n_intervals > 1)
215 {
216 const number step = 1. / n_intervals;
217
218 // polynomial spans over two intervals
219 if (spans_two_intervals == true)
220 {
221 const number offset = step * interval;
222 if (x < offset)
223 return 0;
224 else if (x > offset + step + step)
225 return 0;
226 else if (x < offset + step)
227 y = x - offset;
228 else
229 y = offset + step + step - x;
230 }
231 else
232 {
233 const number offset = step * interval;
234 if (x < offset || x > offset + step)
235 return 0;
236 else
237 y = x - offset;
238 }
239
240 return polynomial.value(y);
241 }
242 else
243 return polynomial.value(x);
244 }
245
246
247
248 template <typename number>
249 template <class Archive>
250 inline void
251 PiecewisePolynomial<number>::serialize(Archive &ar, const unsigned int)
252 {
253 // forward to serialization function in the base class.
254 ar &static_cast<Subscriptor &>(*this);
255 ar &polynomial;
256 ar &n_intervals;
257 ar &interval;
258 ar &spans_two_intervals;
259 }
260
261} // namespace Polynomials
262
264
265#endif
PiecewisePolynomial(const Polynomial< number > &coefficients_on_interval, const unsigned int n_intervals, const unsigned int interval, const bool spans_next_interval)
void serialize(Archive &ar, const unsigned int version)
number value(const number x) const
virtual std::size_t memory_consumption() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
std::vector< PiecewisePolynomial< double > > generate_complete_Lagrange_basis_on_subdivisions(const unsigned int n_subdivisions, const unsigned int base_degree)