Logo ROOT   6.13/01
Reference Guide
GSLInterpolator.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 GSLInterpolator
26 //
27 // Created by: moneta at Sun Nov 28 08:54:48 2004
28 //
29 // Last update: Sun Nov 28 08:54:48 2004
30 //
31 
32 #include "GSLInterpolator.h"
33 
34 
35 #include <cassert>
36 
37 namespace ROOT {
38 namespace Math {
39 
40 
42  fResetNErrors(true),
43  fAccel(0),
44  fSpline(0)
45 {
46  // constructor given type and vectors of (x,y) points
47 
48  switch ( type )
49  {
51  fInterpType = gsl_interp_linear;
52  break ;
54  fInterpType = gsl_interp_polynomial;
55  break ;
56  // depened on GSL linear algebra
58  fInterpType = gsl_interp_cspline ;
59  break ;
61  fInterpType = gsl_interp_cspline_periodic ;
62  break ;
64  fInterpType = gsl_interp_akima;
65  break ;
67  fInterpType = gsl_interp_akima_periodic;
68  break ;
69  default :
70  // cspline
71  fInterpType = gsl_interp_cspline;
72  break ;
73  }
74  // allocate objects
75 
76  if (size >= fInterpType->min_size)
77  fSpline = gsl_spline_alloc( fInterpType, size);
78 
79 }
80 
81 bool GSLInterpolator::Init(unsigned int size, const double *x, const double * y) {
82  // initialize interpolation object with the given data
83  // if given size is different a new interpolator object is created
84  if (fSpline == 0)
85  fSpline = gsl_spline_alloc( fInterpType, size);
86 
87  else {
88  gsl_interp * interp = fSpline->interp;
89  if (size != interp->size) {
90  // free and reallocate a new object
91  gsl_spline_free(fSpline);
92  fSpline = gsl_spline_alloc( fInterpType, size);
93 
94  }
95  }
96  if (!fSpline) return false;
97 
98  int iret = gsl_spline_init( fSpline , x , y , size );
99  if (iret != 0) return false;
100 
101  if(fAccel==0)
102  fAccel = gsl_interp_accel_alloc() ;
103  else
104  gsl_interp_accel_reset(fAccel);
105 
106  // if (fSpline == 0 || fAccel == 0)
107  // throw std::exception();
108  assert (fSpline != 0);
109  assert (fAccel != 0);
110  // reset counter for error messages
111  fResetNErrors = true;
112  return true;
113 }
114 
116 {
117  // free gsl objects
118  if (fSpline != 0) gsl_spline_free(fSpline);
119  if (fAccel != 0) gsl_interp_accel_free( fAccel);
120 }
121 
123 {
124  // dummy copy ctr
125 }
126 
128 {
129  // dummy assignment operator
130  if (this == &rhs) return *this; // time saving self-test
131 
132  return *this;
133 }
134 
135 
136 } // namespace Math
137 } // namespace ROOT
Namespace for new ROOT classes and functions.
Definition: TFoamSampler.h:19
you should not use this method at all Int_t y
Definition: TRolke.cxx:630
GSLInterpolator & operator=(const GSLInterpolator &)
* x
Deprecated and error prone model selection interface.
Definition: TRolke.cxx:630
bool Init(unsigned int ndata, const double *x, const double *y)
Interpolation class based on GSL interpolation functions.
Type
Enumeration defining the types of interpolation methods availables.
Namespace for new Math classes and functions.
const gsl_interp_type * fInterpType
GSLInterpolator(unsigned int ndata, Interpolation::Type type)