Reference documentation for deal.II version 9.0.0
function_cspline.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2017 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 at
12 // the top level of the deal.II distribution.
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
22 #include <deal.II/base/function.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/thread_management.h>
25 #include <gsl/gsl_spline.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 namespace Functions
30 {
32  int,
33  << "Interpolation points vector size can not be <"<<arg1<<">."
34  );
35 
37  int,
38  int,
39  << "The size of interpolation points <"<<arg1<<"> is different from the size of interpolation values <" << arg2 <<">."
40  );
41 
42 
44  int,
45  double,
46  double,
47  << "The input interpolation points are not strictly ordered : " << std::endl << "x[" << arg1 << "] = "<< arg2 <<" >= x["<<(arg1+1)<<"] = "<<arg3 <<"."
48  );
49 
51  double,
52  double,
53  double,
54  << "Spline function can not be evaluated outside of the interpolation range: "<< std::endl << arg1 << " is not in ["<< arg2<<";"<<arg3<<"]."
55  );
56 
68  template <int dim>
69  class CSpline : public Function<dim>
70  {
71  public:
77  CSpline(const std::vector<double> &interpolation_points,
78  const std::vector<double> &interpolation_values);
79 
83  virtual ~CSpline();
84 
85  virtual double value (const Point<dim> &point,
86  const unsigned int component = 0) const;
87 
88  virtual Tensor<1,dim> gradient (const Point<dim> &p,
89  const unsigned int component = 0) const;
90 
91  virtual SymmetricTensor<2,dim> hessian (const Point<dim> &p,
92  const unsigned int component = 0) const;
93 
94  virtual double laplacian(const Point< dim > &p,
95  const unsigned int component = 0) const;
96 
97  std::size_t memory_consumption () const;
98 
99  private:
103  const std::vector<double> interpolation_points;
104 
108  const std::vector<double> interpolation_values;
109 
113  gsl_interp_accel *acc;
114 
118  gsl_spline *cspline;
119 
124  };
125 }
126 
127 DEAL_II_NAMESPACE_CLOSE
128 
129 #endif
130 
131 #endif
132 
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const
virtual double laplacian(const Point< dim > &p, const unsigned int component=0) const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:346
static ::ExceptionBase & ExcCSplineEmpty(int arg1)
virtual double value(const Point< dim > &point, const unsigned int component=0) const
const std::vector< double > interpolation_points
static ::ExceptionBase & ExcCSplineOrder(int arg1, double arg2, double arg3)
Threads::Mutex acc_mutex
gsl_interp_accel * acc
const std::vector< double > interpolation_values
static ::ExceptionBase & ExcCSplineRange(double arg1, double arg2, double arg3)
static ::ExceptionBase & ExcCSplineSizeMismatch(int arg1, int arg2)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:370
CSpline(const std::vector< double > &interpolation_points, const std::vector< double > &interpolation_values)