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\}}\)
function_cspline.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2016 - 2020 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_function_cspline_h
17#define dealii_function_cspline_h
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_GSL
23# include <deal.II/base/point.h>
25
26# include <gsl/gsl_spline.h>
27
29
30namespace Functions
31{
33 int,
34 << "Interpolation points vector size can not be <" << arg1
35 << ">.");
36
38 int,
39 int,
40 << "The size of interpolation points <" << arg1
41 << "> is different from the size of interpolation values <"
42 << arg2 << ">.");
43
44
46 int,
47 double,
48 double,
49 << "The input interpolation points are not strictly ordered : "
50 << std::endl
51 << "x[" << arg1 << "] = " << arg2 << " >= x[" << (arg1 + 1)
52 << "] = " << arg3 << ".");
53
56 double,
57 double,
58 double,
59 << "Spline function can not be evaluated outside of the interpolation range: "
60 << std::endl
61 << arg1 << " is not in [" << arg2 << ";" << arg3 << "].");
62
73 template <int dim>
74 class CSpline : public Function<dim>
75 {
76 public:
82 CSpline(const std::vector<double> &interpolation_points,
83 const std::vector<double> &interpolation_values);
84
88 virtual ~CSpline() override;
89
90 virtual double
91 value(const Point<dim> & point,
92 const unsigned int component = 0) const override;
93
94 virtual Tensor<1, dim>
95 gradient(const Point<dim> & p,
96 const unsigned int component = 0) const override;
97
99 hessian(const Point<dim> & p,
100 const unsigned int component = 0) const override;
101
102 virtual double
103 laplacian(const Point<dim> & p,
104 const unsigned int component = 0) const override;
105
109 virtual std::size_t
110 memory_consumption() const override;
111
112 private:
116 const std::vector<double> interpolation_points;
117
121 const std::vector<double> interpolation_values;
122
126 gsl_interp_accel *acc;
127
131 gsl_spline *cspline;
132
137 };
138} // namespace Functions
139
141
142#endif
143
144#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
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
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
gsl_interp_accel * acc
virtual ~CSpline() override
Definition: point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcCSplineSizeMismatch(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
static ::ExceptionBase & ExcCSplineRange(double arg1, double arg2, double arg3)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:561
static ::ExceptionBase & ExcCSplineOrder(int arg1, double arg2, double arg3)
static ::ExceptionBase & ExcCSplineEmpty(int arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188