www.mooseframework.org
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
PolynomialFit Class Reference

This class applies the Least Squares algorithm to a set of points to provide a smooth curve for sampling values. More...

#include <PolynomialFit.h>

Public Member Functions

 PolynomialFit (std::vector< Real > X, std::vector< Real > Y, unsigned int order, bool truncate_order=false)
 
virtual ~PolynomialFit ()=default
 
void generate ()
 This function generates the polynomial fit. More...
 
Real sample (Real x)
 This function will take an independent variable input and will return the dependent variable based on the generated fit. More...
 
void dumpSampleFile (std::string base_name, std::string x_label="X", std::string y_label="Y", float xmin=0, float xmax=0, float ymin=0, float ymax=0)
 This function will dump GNUPLOT input files that can be run to show the data points and function fits. More...
 
unsigned int getSampleSize ()
 This function returns the size of the array holding the points, i.e. More...
 
const std::vector< Real > & getCoefficients ()
 Get a const reference to the coefficients of the least squares fit. More...
 

Private Member Functions

void fillMatrix ()
 This is a helper function that creates the matrix necessary for the Least Squares algorithm. More...
 
void doLeastSquares ()
 This function is the wrapper for the LAPACK dgels function and is called by generate. More...
 

Private Attributes

std::vector< Real > _x
 
std::vector< Real > _y
 
std::vector< Real > _matrix
 
std::vector< Real > _coeffs
 
unsigned int _order
 
bool _truncate_order
 

Static Private Attributes

static int _file_number = 0
 

Detailed Description

This class applies the Least Squares algorithm to a set of points to provide a smooth curve for sampling values.

Requires: LAPACK

Definition at line 29 of file PolynomialFit.h.

Constructor & Destructor Documentation

PolynomialFit::PolynomialFit ( std::vector< Real >  X,
std::vector< Real >  Y,
unsigned int  order,
bool  truncate_order = false 
)

Definition at line 24 of file PolynomialFit.C.

28  : _x(x), _y(y), _order(order), _truncate_order(truncate_order)
29 {
30  if (_truncate_order) // && (_x.size() / 10) < _order)
31  {
32  if (_x.size() == 1)
33  _order = 0;
34  else
35  {
36  _order = (_x.size() / 10) + 1;
37 
38  if (_order > order)
39  _order = order;
40  }
41  }
42  else if (_x.size() < order)
43  throw std::domain_error(
44  "Polynomial Fit requires an order less than the size of the input vector");
45 }
static PetscErrorCode Vec x
std::vector< Real > _x
Definition: PolynomialFit.h:94
bool _truncate_order
Definition: PolynomialFit.h:99
std::vector< Real > _y
Definition: PolynomialFit.h:95
unsigned int _order
Definition: PolynomialFit.h:98
virtual PolynomialFit::~PolynomialFit ( )
virtualdefault

Member Function Documentation

void PolynomialFit::doLeastSquares ( )
private

This function is the wrapper for the LAPACK dgels function and is called by generate.

Definition at line 75 of file PolynomialFit.C.

Referenced by generate().

76 {
77  char mode = 'N';
78  int num_rows = _x.size();
79  int num_coeff = _order + 1;
80  int num_rhs = 1;
81  int buffer_size = -1;
82  Real opt_buffer_size;
83  Real * buffer;
84  int return_value = 0;
85 
86  // Must copy _y because the call to dgels destroys the original values
87  std::vector<Real> rhs = _y;
88 
89  FORTRAN_CALL(dgels)
90  (&mode,
91  &num_rows,
92  &num_coeff,
93  &num_rhs,
94  &_matrix[0],
95  &num_rows,
96  &rhs[0],
97  &num_rows,
98  &opt_buffer_size,
99  &buffer_size,
100  &return_value);
101  if (return_value)
102  throw std::runtime_error("Call to Fortran routine 'dgels' returned non-zero exit code");
103 
104  buffer_size = (int)opt_buffer_size;
105 
106  buffer = new Real[buffer_size];
107  FORTRAN_CALL(dgels)
108  (&mode,
109  &num_rows,
110  &num_coeff,
111  &num_rhs,
112  &_matrix[0],
113  &num_rows,
114  &rhs[0],
115  &num_rows,
116  buffer,
117  &buffer_size,
118  &return_value);
119  delete[] buffer;
120 
121  if (return_value)
122  throw std::runtime_error("Call to Fortran routine 'dgels' returned non-zero exit code");
123 
124  _coeffs.resize(num_coeff);
125  for (int i = 0; i < num_coeff; ++i)
126  _coeffs[i] = rhs[i];
127 }
std::vector< Real > _matrix
Definition: PolynomialFit.h:96
void FORTRAN_CALL() dgels(...)
std::vector< Real > _x
Definition: PolynomialFit.h:94
std::vector< Real > _y
Definition: PolynomialFit.h:95
std::vector< Real > _coeffs
Definition: PolynomialFit.h:97
unsigned int _order
Definition: PolynomialFit.h:98
void PolynomialFit::dumpSampleFile ( std::string  base_name,
std::string  x_label = "X",
std::string  y_label = "Y",
float  xmin = 0,
float  xmax = 0,
float  ymin = 0,
float  ymax = 0 
)

This function will dump GNUPLOT input files that can be run to show the data points and function fits.

Definition at line 145 of file PolynomialFit.C.

152 {
153  std::stringstream filename, filename_pts;
154  const unsigned char fill_character = '0';
155  const unsigned int field_width = 4;
156 
157  filename.fill(fill_character);
158  filename << base_name;
159  filename.width(field_width);
160  filename << _file_number << ".plt";
161 
162  filename_pts.fill(fill_character);
163  filename_pts << base_name << "_pts";
164  filename_pts.width(field_width);
165  filename_pts << _file_number << ".dat";
166 
167  /* First dump the GNUPLOT file with the Least Squares Equations */
168  std::ofstream out(filename.str().c_str());
169  out.precision(15);
170  out.fill(fill_character);
171 
172  out << "set terminal postscript color enhanced\n"
173  << "set output \"" << base_name;
174  out.width(field_width);
175  out << _file_number << ".eps\"\n"
176  << "set xlabel \"" << x_label << "\"\n"
177  << "set ylabel \"" << y_label << "\"\n";
178  if (xmin != 0 && xmax != 0)
179  out << "set xrange [" << xmin << ":" << xmax << "]\n";
180  if (ymin != 0 && ymax != 0)
181  out << "set yrange [" << ymin << ":" << ymax << "]\n";
182  out << "set key left top\n"
183  << "f(x)=";
184 
185  for (unsigned int i = 0; i < _coeffs.size(); ++i)
186  {
187  if (i)
188  out << "+";
189 
190  out << _coeffs[i];
191  for (unsigned int j = 0; j < i; ++j)
192  out << "*x";
193  }
194  out << "\nplot f(x) with lines, '" << filename_pts.str() << "' using 1:2 title \"Points\"\n";
195  out.close();
196 
197  libmesh_assert(_x.size() == _y.size());
198 
199  out.open(filename_pts.str().c_str());
200  /* Next dump the data points into a seperate file */
201  for (unsigned int i = 0; i < _x.size(); ++i)
202  out << _x[i] << " " << _y[i] << "\n";
203  out << std::endl;
204 
205  ++_file_number;
206  out.close();
207 }
static int _file_number
std::vector< Real > _x
Definition: PolynomialFit.h:94
std::vector< Real > _y
Definition: PolynomialFit.h:95
std::vector< Real > _coeffs
Definition: PolynomialFit.h:97
void PolynomialFit::fillMatrix ( )
private

This is a helper function that creates the matrix necessary for the Least Squares algorithm.

Definition at line 55 of file PolynomialFit.C.

Referenced by generate().

56 {
57  unsigned int num_rows = _x.size();
58  unsigned int num_cols = _order + 1;
59  _matrix.resize(num_rows * num_cols);
60 
61  for (unsigned int col = 0; col <= _order; ++col)
62  {
63  for (unsigned int row = 0; row < num_rows; ++row)
64  {
65  Real value = 1;
66  for (unsigned int i = 0; i < col; ++i)
67  value *= _x[row];
68 
69  _matrix[(col * num_rows) + row] = value;
70  }
71  }
72 }
std::vector< Real > _matrix
Definition: PolynomialFit.h:96
std::vector< Real > _x
Definition: PolynomialFit.h:94
unsigned int _order
Definition: PolynomialFit.h:98
void PolynomialFit::generate ( )

This function generates the polynomial fit.

This function must be called prior to using sample or dumpSample File. Note: If you pass a vectors that contain duplicate independent measures the call to LAPACK will fail

Definition at line 48 of file PolynomialFit.C.

Referenced by LeastSquaresFit::execute().

49 {
50  fillMatrix();
52 }
void doLeastSquares()
This function is the wrapper for the LAPACK dgels function and is called by generate.
Definition: PolynomialFit.C:75
void fillMatrix()
This is a helper function that creates the matrix necessary for the Least Squares algorithm...
Definition: PolynomialFit.C:55
const std::vector< Real > & PolynomialFit::getCoefficients ( )

Get a const reference to the coefficients of the least squares fit.

Definition at line 216 of file PolynomialFit.C.

Referenced by LeastSquaresFit::execute().

217 {
218  return _coeffs;
219 }
std::vector< Real > _coeffs
Definition: PolynomialFit.h:97
unsigned int PolynomialFit::getSampleSize ( )

This function returns the size of the array holding the points, i.e.

the number of sample points

Definition at line 210 of file PolynomialFit.C.

211 {
212  return _x.size();
213 }
std::vector< Real > _x
Definition: PolynomialFit.h:94
Real PolynomialFit::sample ( Real  x)

This function will take an independent variable input and will return the dependent variable based on the generated fit.

Definition at line 130 of file PolynomialFit.C.

Referenced by LeastSquaresFit::execute().

131 {
132  unsigned int size = _coeffs.size();
133  Real value = 0;
134 
135  Real curr_x = 1;
136  for (unsigned int i = 0; i < size; ++i)
137  {
138  value += _coeffs[i] * curr_x;
139  curr_x *= x;
140  }
141  return value;
142 }
static PetscErrorCode Vec x
std::vector< Real > _coeffs
Definition: PolynomialFit.h:97

Member Data Documentation

std::vector<Real> PolynomialFit::_coeffs
private

Definition at line 97 of file PolynomialFit.h.

Referenced by doLeastSquares(), dumpSampleFile(), getCoefficients(), and sample().

int PolynomialFit::_file_number = 0
staticprivate

Definition at line 101 of file PolynomialFit.h.

Referenced by dumpSampleFile().

std::vector<Real> PolynomialFit::_matrix
private

Definition at line 96 of file PolynomialFit.h.

Referenced by doLeastSquares(), and fillMatrix().

unsigned int PolynomialFit::_order
private

Definition at line 98 of file PolynomialFit.h.

Referenced by doLeastSquares(), fillMatrix(), and PolynomialFit().

bool PolynomialFit::_truncate_order
private

Definition at line 99 of file PolynomialFit.h.

Referenced by PolynomialFit().

std::vector<Real> PolynomialFit::_x
private
std::vector<Real> PolynomialFit::_y
private

Definition at line 95 of file PolynomialFit.h.

Referenced by doLeastSquares(), and dumpSampleFile().


The documentation for this class was generated from the following files: