Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
Interpolator.cxx
Go to the documentation of this file.
1// @(#)root/mathmore:$Id$
2// Authors: L. Moneta, A. Zsenei 08/2005
3
4 /**********************************************************************
5 * *
6 * Copyright (c) 2004 ROOT Foundation, CERN/PH-SFT *
7 * *
8 * This library is free software; you can redistribute it and/or *
9 * modify it under the terms of the GNU General Public License *
10 * as published by the Free Software Foundation; either version 2 *
11 * of the License, or (at your option) any later version. *
12 * *
13 * This library is distributed in the hope that it will be useful, *
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
16 * General Public License for more details. *
17 * *
18 * You should have received a copy of the GNU General Public License *
19 * along with this library (see file COPYING); if not, write *
20 * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
21 * 330, Boston, MA 02111-1307 USA, or contact the author. *
22 * *
23 **********************************************************************/
24
25// Implementation file for class Interpolator using GSL
26//
27// Created by: moneta at Fri Nov 26 15:00:25 2004
28//
29// Last update: Fri Nov 26 15:00:25 2004
30//
31
32#include "Math/Interpolator.h"
33#include "GSLInterpolator.h"
34#include <algorithm>
35
36
37namespace ROOT {
38namespace Math {
39
41 // allocate GSL interpolation object
42 fInterp = new GSLInterpolator(ndata, type);
43}
44
45Interpolator::Interpolator(const std::vector<double> & x, const std::vector<double> & y, Interpolation::Type type)
46{
47 // allocate and initialize GSL interpolation object with data
48
49 size_t size = std::min( x.size(), y.size() );
50
52
53 fInterp->Init(size, &x.front(), &y.front() );
54
55}
56
57
59{
60 // destructor (delete underlined obj)
61 if (fInterp) delete fInterp;
62}
63
64bool Interpolator::SetData(unsigned int ndata, const double * x, const double *y) {
65 // set the interpolation data
66 return fInterp->Init(ndata, x, y);
67}
68bool Interpolator::SetData(const std::vector<double> & x, const std::vector<double> &y) {
69 // set the interpolation data
70 size_t size = std::min( x.size(), y.size() );
71 return fInterp->Init(size, &x.front(), &y.front());
72}
73
74
75double Interpolator::Eval( double x ) const
76{
77 // forward evaluation
78 return fInterp->Eval(x);
79}
80
81double Interpolator::Deriv( double x ) const
82{
83 // forward deriv evaluation
84 return fInterp->Deriv(x);
85}
86
87double Interpolator::Deriv2( double x ) const {
88 // forward deriv evaluation
89 return fInterp->Deriv2(x);
90}
91
92double Interpolator::Integ( double a, double b) const {
93 // forward integ evaluation
94 return fInterp->Integ(a,b);
95}
96
97std::string Interpolator::TypeGet() const {
98 // forward name request
99 return fInterp->Name();
100}
101std::string Interpolator::Type() const {
102 // forward name request
103 return fInterp->Name();
104}
105
106
107
108} // namespace Math
109} // namespace ROOT
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Interpolation class based on GSL interpolation functions.
double Deriv2(double x) const
double Eval(double x) const
bool Init(unsigned int ndata, const double *x, const double *y)
double Integ(double a, double b) const
double Deriv(double x) const
double Deriv2(double x) const
Return the second derivative of the interpolated function at point x.
std::string Type() const
Return the type of interpolation method.
bool SetData(const std::vector< double > &x, const std::vector< double > &y)
Set the data vector ( x[] and y[] ) To be efficient, the size of the data must be the same of the val...
double Integ(double a, double b) const
Return the Integral of the interpolated function over the range [a,b].
std::string TypeGet() const
double Deriv(double x) const
Return the derivative of the interpolated function at point x.
double Eval(double x) const
Return the interpolated value at point x.
GSLInterpolator * fInterp
Interpolator(unsigned int ndata=0, Interpolation::Type type=Interpolation::kCSPLINE)
Constructs an interpolator class from number of data points and with Interpolation::Type type.
Type
Enumeration defining the types of interpolation methods availables.
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Namespace for new Math classes and functions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...