www.mooseframework.org
MonotoneCubicInterpolation.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include <vector>
13 
14 #include "libmesh/libmesh_common.h"
15 
16 using libMesh::Real;
17 
26 {
27 public:
32 
39  MonotoneCubicInterpolation(const std::vector<Real> & x, const std::vector<Real> & y);
40 
41  virtual ~MonotoneCubicInterpolation() = default;
42 
50  virtual void setData(const std::vector<Real> & x, const std::vector<Real> & y);
51 
56  virtual Real sample(const Real & x) const;
57 
63  virtual Real sampleDerivative(const Real & x) const;
64 
70  virtual Real sample2ndDerivative(const Real & x) const;
71 
77  virtual void dumpCSV(std::string filename, const std::vector<Real> & xnew);
78 
82  virtual unsigned int getSampleSize();
83 
84 protected:
85  // Error check routines run during initialization
86  virtual void errorCheck();
87  Real sign(const Real & x) const;
88 
89  // Building blocks of Hermite polynomials
90  Real phi(const Real & t) const;
91  Real psi(const Real & t) const;
92  Real phiPrime(const Real & t) const;
93  Real psiPrime(const Real & t) const;
94  Real phiDoublePrime(const Real & t) const;
95  Real psiDoublePrime(const Real & t) const;
96 
97  // Cubic Hermite polynomials
98  Real h1(const Real & xhi, const Real & xlo, const Real & x) const;
99  Real h2(const Real & xhi, const Real & xlo, const Real & x) const;
100  Real h3(const Real & xhi, const Real & xlo, const Real & x) const;
101  Real h4(const Real & xhi, const Real & xlo, const Real & x) const;
102  Real h1Prime(const Real & xhi, const Real & xlo, const Real & x) const;
103  Real h2Prime(const Real & xhi, const Real & xlo, const Real & x) const;
104  Real h3Prime(const Real & xhi, const Real & xlo, const Real & x) const;
105  Real h4Prime(const Real & xhi, const Real & xlo, const Real & x) const;
106  Real h1DoublePrime(const Real & xhi, const Real & xlo, const Real & x) const;
107  Real h2DoublePrime(const Real & xhi, const Real & xlo, const Real & x) const;
108  Real h3DoublePrime(const Real & xhi, const Real & xlo, const Real & x) const;
109  Real h4DoublePrime(const Real & xhi, const Real & xlo, const Real & x) const;
110 
111  // Interpolating cubic polynomial and derivatives
112  virtual Real p(const Real & xhi,
113  const Real & xlo,
114  const Real & fhi,
115  const Real & flo,
116  const Real & dhi,
117  const Real & dlo,
118  const Real & x) const;
119  virtual Real pPrime(const Real & xhi,
120  const Real & xlo,
121  const Real & fhi,
122  const Real & flo,
123  const Real & dhi,
124  const Real & dlo,
125  const Real & x) const;
126  virtual Real pDoublePrime(const Real & xhi,
127  const Real & xlo,
128  const Real & fhi,
129  const Real & flo,
130  const Real & dhi,
131  const Real & dlo,
132  const Real & x) const;
133 
134  // Algorithm routines
135  virtual void initialize_derivs();
136  virtual void modify_derivs(
137  const Real & alpha, const Real & beta, const Real & delta, Real & yp_lo, Real & yp_hi);
138  virtual void solve();
139  virtual void findInterval(const Real & x, unsigned int & klo, unsigned int & khi) const;
140 
141  std::vector<Real> _x;
142  std::vector<Real> _y;
143  std::vector<Real> _h;
144  std::vector<Real> _yp;
145  std::vector<Real> _delta;
146  std::vector<Real> _alpha;
147  std::vector<Real> _beta;
148 
149  unsigned int _n_knots;
150  unsigned int _n_intervals;
151  unsigned int _internal_knots;
152 };
Real h1Prime(const Real &xhi, const Real &xlo, const Real &x) const
virtual Real sample2ndDerivative(const Real &x) const
This function will take an independent variable input and will return the second derivative of the de...
Real h2Prime(const Real &xhi, const Real &xlo, const Real &x) const
Real h3Prime(const Real &xhi, const Real &xlo, const Real &x) const
virtual unsigned int getSampleSize()
This method returns the length of the independent variable vector.
Real h3DoublePrime(const Real &xhi, const Real &xlo, const Real &x) const
Real psiPrime(const Real &t) const
Real h4DoublePrime(const Real &xhi, const Real &xlo, const Real &x) const
This class interpolates values given a set of data pairs and an abscissa.
virtual void modify_derivs(const Real &alpha, const Real &beta, const Real &delta, Real &yp_lo, Real &yp_hi)
virtual Real pPrime(const Real &xhi, const Real &xlo, const Real &fhi, const Real &flo, const Real &dhi, const Real &dlo, const Real &x) const
virtual void setData(const std::vector< Real > &x, const std::vector< Real > &y)
Method generally used when MonotoneCubicInterpolation object was created using the empty constructor...
Real psiDoublePrime(const Real &t) const
Real phiPrime(const Real &t) const
virtual Real sampleDerivative(const Real &x) const
This function will take an independent variable input and will return the derivative of the dependent...
MonotoneCubicInterpolation()
Empty constructor.
Real h4(const Real &xhi, const Real &xlo, const Real &x) const
Real sign(const Real &x) const
virtual ~MonotoneCubicInterpolation()=default
virtual Real sample(const Real &x) const
This function will take an independent variable input and will return the dependent variable based on...
Real h1DoublePrime(const Real &xhi, const Real &xlo, const Real &x) const
Real h4Prime(const Real &xhi, const Real &xlo, const Real &x) const
Real h2(const Real &xhi, const Real &xlo, const Real &x) const
virtual Real pDoublePrime(const Real &xhi, const Real &xlo, const Real &fhi, const Real &flo, const Real &dhi, const Real &dlo, const Real &x) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void dumpCSV(std::string filename, const std::vector< Real > &xnew)
This function takes an array of independent variable values and writes a CSV file with values corresp...
Real h3(const Real &xhi, const Real &xlo, const Real &x) const
virtual void findInterval(const Real &x, unsigned int &klo, unsigned int &khi) const
Real h1(const Real &xhi, const Real &xlo, const Real &x) const
virtual Real p(const Real &xhi, const Real &xlo, const Real &fhi, const Real &flo, const Real &dhi, const Real &dlo, const Real &x) const
Real phiDoublePrime(const Real &t) const
Real h2DoublePrime(const Real &xhi, const Real &xlo, const Real &x) const