SplineFit.h
1//-*-C++-*-
2/***************************************************************************
3 *
4 * Copyright (C) 2010 by Paul Demorest
5 * Licensed under the Academic Free License version 2.1
6 *
7 ***************************************************************************/
8
9#ifndef __Pulsar_SplineFit_h
10#define __Pulsar_SplineFit_h
11
12#include "Estimate.h"
13#include "Reference.h"
14
15#include <gsl/gsl_bspline.h>
16#include <gsl/gsl_matrix.h>
17
20class SplineFit : public Reference::Able {
21
22public:
23
25 SplineFit ();
26
28 virtual ~SplineFit ();
29
31 void reset();
32
34 void set_order(int n) { order=n; calculated=false; }
35
37 int get_order() { return order; }
38
40 void set_uniform_breaks(int nint);
41
43 void add_data(double x, Estimate<double> y);
44
46 void compute();
47
49 double evaluate(double x);
50
52 double evaluate_deriv(double x);
53
55 double get_rchi2();
56
57protected:
58
60 std::vector<double> x;
61
63 std::vector< Estimate<double> > y;
64
66 std::vector<double> bp;
67
69 int order;
70
72 double chi2;
73
75 int ndof;
76
79
81 gsl_vector *coeffs;
82
84 gsl_matrix *cov;
85
87 bool check_range(double x);
88
90 void interval_check(bool fix=false);
91
93 void free_workspaces();
94
96 gsl_bspline_workspace *bwork;
97
98private:
99
100};
101
102#endif
Manages Reference::To references to the instance.
Definition ReferenceAble.h:35
int order
Spline order (0=const, 3=cubic, etc)
Definition SplineFit.h:69
double chi2
The fit chi2.
Definition SplineFit.h:72
void set_order(int n)
Set the degree of the fit.
Definition SplineFit.h:34
bool calculated
Has the fit been calculated?
Definition SplineFit.h:78
std::vector< Estimate< double > > y
The y values/errors for the fit.
Definition SplineFit.h:63
void compute()
Compute the fit using current data.
Definition SplineFit.C:198
virtual ~SplineFit()
Destructor.
Definition SplineFit.C:29
double evaluate_deriv(double x)
Evaluate the fit solution's derivative at the given x.
Definition SplineFit.C:117
bool check_range(double x)
Check if a requested x val is in the fit range.
Definition SplineFit.C:55
void reset()
Clear all current data, results.
Definition SplineFit.C:34
int get_order()
Get the current degree.
Definition SplineFit.h:37
void interval_check(bool fix=false)
Check if spline intervals and data make sense.
Definition SplineFit.C:141
std::vector< double > bp
The spline breakpoints.
Definition SplineFit.h:66
void set_uniform_breaks(int nint)
Set uniform breakpoints to span the data.
Definition SplineFit.C:75
gsl_matrix * cov
The fit cov matrix.
Definition SplineFit.h:84
gsl_bspline_workspace * bwork
bspline temp space
Definition SplineFit.h:96
int ndof
Fit NDOF.
Definition SplineFit.h:75
void free_workspaces()
Free spline workspaces.
Definition SplineFit.C:45
SplineFit()
Default constructor.
Definition SplineFit.C:18
double get_rchi2()
Get the reduced chi2 of the fit.
Definition SplineFit.C:133
void add_data(double x, Estimate< double > y)
Add a data point.
Definition SplineFit.C:68
double evaluate(double x)
Evaluate the fit solution at the given x.
Definition SplineFit.C:89
gsl_vector * coeffs
The fitted coeffs.
Definition SplineFit.h:81
std::vector< double > x
The x values for the fit.
Definition SplineFit.h:60

Generated using doxygen 1.14.0