Reference documentation for deal.II version 9.6.0
\(\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\}}\)
Loading...
Searching...
No Matches
function_cspline.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_function_cspline_h
16#define dealii_function_cspline_h
17
18#include <deal.II/base/config.h>
19
20#ifdef DEAL_II_WITH_GSL
22# include <deal.II/base/mutex.h>
23# include <deal.II/base/point.h>
24
25# include <gsl/gsl_spline.h>
26
27# include <memory>
28
30
31namespace Functions
32{
34 int,
35 << "Interpolation points vector size can not be <" << arg1
36 << ">.");
37
39 int,
40 int,
41 << "The size of interpolation points <" << arg1
42 << "> is different from the size of interpolation values <"
43 << arg2 << ">.");
44
45
47 int,
48 double,
49 double,
50 << "The input interpolation points are not strictly ordered : "
51 << std::endl
52 << "x[" << arg1 << "] = " << arg2 << " >= x[" << (arg1 + 1)
53 << "] = " << arg3 << '.');
54
57 double,
58 double,
59 double,
60 << "Spline function can not be evaluated outside of the interpolation range: "
61 << std::endl
62 << arg1 << " is not in [" << arg2 << ';' << arg3 << "].");
63
74 template <int dim>
75 class CSpline : public Function<dim>
76 {
77 public:
83 CSpline(const std::vector<double> &interpolation_points,
84 const std::vector<double> &interpolation_values);
85
86 virtual double
87 value(const Point<dim> &point,
88 const unsigned int component = 0) const override;
89
90 virtual Tensor<1, dim>
91 gradient(const Point<dim> &p,
92 const unsigned int component = 0) const override;
93
95 hessian(const Point<dim> &p,
96 const unsigned int component = 0) const override;
97
98 virtual double
99 laplacian(const Point<dim> &p,
100 const unsigned int component = 0) const override;
101
105 virtual std::size_t
106 memory_consumption() const override;
107
108 private:
112 const std::vector<double> interpolation_points;
113
117 const std::vector<double> interpolation_values;
118
122 std::unique_ptr<gsl_interp_accel, void (*)(gsl_interp_accel *)> acc;
123
127 std::unique_ptr<gsl_spline, void (*)(gsl_spline *)> cspline;
128
133 };
134} // namespace Functions
135
137
138#endif
139
140#endif
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const override
Threads::Mutex acc_mutex
virtual double value(const Point< dim > &point, const unsigned int component=0) const override
CSpline(const std::vector< double > &interpolation_points, const std::vector< double > &interpolation_values)
virtual std::size_t memory_consumption() const override
std::unique_ptr< gsl_interp_accel, void(*)(gsl_interp_accel *)> acc
const std::vector< double > interpolation_points
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const override
const std::vector< double > interpolation_values
std::unique_ptr< gsl_spline, void(*)(gsl_spline *)> cspline
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
static ::ExceptionBase & ExcCSplineSizeMismatch(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
static ::ExceptionBase & ExcCSplineRange(double arg1, double arg2, double arg3)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition exceptions.h:562
static ::ExceptionBase & ExcCSplineOrder(int arg1, double arg2, double arg3)
static ::ExceptionBase & ExcCSplineEmpty(int arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516