tmdlib is hosted by Hepforge, IPPP Durham
TMDlib  1.0.X
TMDlib_InterpolationKS.h
Go to the documentation of this file.
1 //------------------------------------------------------------------------------
2 //
3 // Copyright (c) 2012, Sebastian Sapeta
4 //
5 // This is part of KSgluon package.
6 //
7 // KSgluon is free software; you can redistribute it and/or modify
8 // it under the terms of the GNU General Public License as published by
9 // the Free Software Foundation; either version 2 of the License, or
10 // (at your option) any later version.
11 //
12 // KSgluon is distributed in the hope that it will be useful,
13 // but WITHOUT ANY WARRANTY; without even the implied warranty of
14 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 // GNU General Public License for more details.
16 //
17 // You should have received a copy of the GNU General Public License
18 // along with FastJet; if not, write to the Free Software
19 // Foundation, Inc.:
20 // 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
21 //
22 //------------------------------------------------------------------------------
23 #ifndef __TMDlib_INTERPOLATIONKS_H__
24 #define __TMDlib_INTERPOLATIONKS_H__
25 #include <cstdlib>
26 #include <cstdio>
27 #include <cmath>
28 #include <iostream>
29 #include <memory>
30 #include <vector>
31 #include <boost/multi_array/base.hpp>
32 #include <boost/shared_ptr.hpp>
33 #include <boost/multi_array.hpp>
34 #include <gsl/gsl_errno.h>
35 #include <gsl/gsl_interp.h>
36 #include <gsl/gsl_spline.h>
37 
38 
39 #include <gsl/gsl_interp2d.h>
40 #include <gsl/gsl_spline2d.h>
41 using namespace std;
42 
43 //------------------------------------------------------------------------------
46 //------------------------------------------------------------------------------
48  public:
52  virtual ~BaseInterpolation() {}
53 
55  virtual vector<pair <double, double> > grid_limits() const {
56  return vector<pair <double, double> >(); }
57 
59  virtual double interp (const vector<double> &x) { return -1.0; }
60 
71  unsigned int eval_status() {return _eval_status;}
72 
73  protected:
74  unsigned int _eval_status;
75 
76 };
77 
78 
79 //------------------------------------------------------------------------------
111 
112  public:
121  Interpolation2D(const vector<double> &x1vv, const vector<double> &x2vv,
122  boost::multi_array<double, 2> ym,
123  gsl_interp_type int_type = *gsl_interp_linear);
124 
126  ~Interpolation2D();
127 
133  vector<pair <double, double> > grid_limits() const;
134 
136  double interp (const vector<double> &x);
137 
138  private:
140  gsl_interp_type _int_type; //this must stay - we leave old constructor
141  gsl_interp2d_type _int2d_type; //possib: gsl_interp2d_bilinear or bicubic
142  unsigned int n1, n2;
143  boost::multi_array<double, 2> y;
144  gsl_interp_accel* x1acc;
145  gsl_interp_accel* x2acc;
146  gsl_spline2d* spline2d;
147  vector<double> x1v, x2v;
148 
149 };
150 //------------------------------------------------------------------------------
151 //
160 //------------------------------------------------------------------------------
162 
163  public:
166  CubicSpline3D(const vector<double> &x1vv, const vector<double> &x2vv,
167  const vector<double> &x3vv, boost::multi_array<double, 3> ym);
168 
169  ~CubicSpline3D();
170 
171  //double interp (double x3) {
172  // int i = 0, j = 0;
173  // return gsl_spline_eval (spline[i][j], x3, acc[i][j]);
174  //}
175 
179 
182  vector<pair <double, double> > grid_limits() const;
183  double interp (const vector<double> &x);
184 
185 
186  private:
187  unsigned int n1, n2, n3;
188  boost::multi_array<double, 3> y;
189  vector<gsl_interp_accel*> x1acc;
190  vector<gsl_interp_accel*> x2acc;
191  vector<gsl_spline2d*> spline2d;
192  vector<double> x1v, x2v, x3v;
193 
194 };
195 
196 
197 
198 #endif // __TMDlib_INTERPOLATIONKS_H__
Definition: TMDlib_InterpolationKS.h:161
Definition: TMDlib_InterpolationKS.h:110
BaseInterpolation()
default constructor
Definition: TMDlib_InterpolationKS.h:50
virtual ~BaseInterpolation()
default destructor
Definition: TMDlib_InterpolationKS.h:52
Definition: TMDlib_InterpolationKS.h:47
virtual double interp(const vector< double > &x)
return the interpolated value corresponding to the vector (x1,x2,x3,..)
Definition: TMDlib_InterpolationKS.h:59
virtual vector< pair< double, double > > grid_limits() const
grid limits in form of vector of pairs: [[x1min,x1max], [x2min,x2max],..]
Definition: TMDlib_InterpolationKS.h:55
unsigned int eval_status()
Definition: TMDlib_InterpolationKS.h:71