MUQ  0.4.3
Regression.h
Go to the documentation of this file.
1 #ifndef REGRESSION_H_
2 #define REGRESSION_H_
3 
4 #include <boost/property_tree/ptree.hpp>
5 
6 #include <Eigen/QR>
7 
9 
11 #include "MUQ/Modeling/WorkPiece.h"
12 
14 
16 
17 namespace muq {
18  namespace Approximation {
19  class Regression : public muq::Modeling::WorkPiece, public std::enable_shared_from_this<Regression> {
20  public:
21 
29  Regression(boost::property_tree::ptree const& pt);
30 
32 
37  void Fit(std::vector<Eigen::VectorXd> xs, std::vector<Eigen::VectorXd> const& ys, Eigen::VectorXd const& center);
38 
40 
44  void Fit(std::vector<Eigen::VectorXd> const& xs, std::vector<Eigen::VectorXd> const& ys);
45 
46  int NumInterpolationPoints() const;
47 
54  std::pair<Eigen::VectorXd, double> PoisednessConstant(std::vector<Eigen::VectorXd> xs, Eigen::VectorXd const& center, int kn = -1) const;
55 
57  const unsigned int order;
58 
59  private:
60 
61  virtual void EvaluateImpl(muq::Modeling::ref_vector<boost::any> const& inputs) override;
62 
63  void ComputeBasisDerivatives(Eigen::VectorXd const& point, std::vector<Eigen::VectorXd>& gradient) const;
64 
66 
72  Eigen::MatrixXd ComputeCoefficients(std::vector<Eigen::VectorXd> const& xs, std::vector<Eigen::VectorXd> const& ys) const;
73 
75 
79  Eigen::MatrixXd ComputeCoefficientsRHS(Eigen::MatrixXd const& vand, std::vector<Eigen::VectorXd> const& ys_data) const;
80 
82 
86  Eigen::MatrixXd VandermondeMatrix(std::vector<Eigen::VectorXd> const& xs) const;
87 
89 
93  Eigen::ArrayXd CenterPoints(std::vector<Eigen::VectorXd>& xs);
94 
96 
100  Eigen::ArrayXd CenterPoints(std::vector<Eigen::VectorXd>& xs, Eigen::VectorXd const& center) const;
101 
103 
109  Eigen::ArrayXd CenterPoints(std::vector<Eigen::VectorXd>& xs, Eigen::VectorXd const& center, unsigned int const kn) const;
110 
112  public:
113 
114  PoisednessCost(std::shared_ptr<Regression const> parent, std::vector<Eigen::RowVectorXd> const& lagrangeCoeff, unsigned int const inDim);
115 
116  virtual ~PoisednessCost() = default;
117 
118  virtual double Cost() override;
119 
120  private:
121 
122  virtual void GradientImpl(unsigned int outWrt, unsigned int inWrt, muq::Modeling::ref_vector<Eigen::VectorXd> const& input, Eigen::VectorXd const& sensitivity) override;
123 
124  std::shared_ptr<Regression const> parent;
125 
126  const std::vector<Eigen::RowVectorXd>& lagrangeCoeff;
127  };
128 
130  public:
131 
132  PoisednessConstraint(unsigned int const inDim, double const alpha);
133 
134  virtual ~PoisednessConstraint() = default;
135 
136  private:
137 
138  virtual void EvaluateImpl(muq::Modeling::ref_vector<Eigen::VectorXd> const& input) override;
139 
140  virtual void JacobianImpl(unsigned int const outputDimWrt,
141  unsigned int const inputDimWrt,
142  muq::Modeling::ref_vector<Eigen::VectorXd> const& input) override;
143 
144  const double alpha;
145  };
146 
148  const unsigned int inputDim;
149 
151 
154  const double alpha;
155 
157  std::shared_ptr<muq::Utilities::MultiIndexSet> multi;
158 
160  std::shared_ptr<IndexedScalarBasis> poly;
161 
163  Eigen::VectorXd currentCenter;
164 
166  Eigen::ArrayXd currentRadius;
167 
169  Eigen::MatrixXd coeff;
170 
172  boost::property_tree::ptree optPt;
173  };
174  } // namespace Approximation
175 } // namespace muq
176 
177 #endif
virtual void JacobianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, muq::Modeling::ref_vector< Eigen::VectorXd > const &input) override
Definition: Regression.cpp:347
PoisednessConstraint(unsigned int const inDim, double const alpha)
Definition: Regression.cpp:338
virtual void EvaluateImpl(muq::Modeling::ref_vector< Eigen::VectorXd > const &input) override
Definition: Regression.cpp:341
std::shared_ptr< Regression const > parent
Definition: Regression.h:124
const std::vector< Eigen::RowVectorXd > & lagrangeCoeff
Definition: Regression.h:126
virtual void GradientImpl(unsigned int outWrt, unsigned int inWrt, muq::Modeling::ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &sensitivity) override
Compute the gradient of the cost function.
Definition: Regression.cpp:305
PoisednessCost(std::shared_ptr< Regression const > parent, std::vector< Eigen::RowVectorXd > const &lagrangeCoeff, unsigned int const inDim)
Definition: Regression.cpp:293
Eigen::ArrayXd currentRadius
Current radius of inputs.
Definition: Regression.h:166
boost::property_tree::ptree optPt
Parameters for the poisedness cosntant optimization.
Definition: Regression.h:172
const double alpha
The radius of the poisedness constraint.
Definition: Regression.h:154
void ComputeBasisDerivatives(Eigen::VectorXd const &point, std::vector< Eigen::VectorXd > &gradient) const
Definition: Regression.cpp:263
virtual void EvaluateImpl(muq::Modeling::ref_vector< boost::any > const &inputs) override
User-implemented function that determines the behavior of this muq::Modeling::WorkPiece.
Definition: Regression.cpp:61
const unsigned int order
The order of the regression.
Definition: Regression.h:57
std::shared_ptr< IndexedScalarBasis > poly
The polynomial basis (in one variable) used to compute the Vandermonde matrix.
Definition: Regression.h:160
Eigen::MatrixXd ComputeCoefficients(std::vector< Eigen::VectorXd > const &xs, std::vector< Eigen::VectorXd > const &ys) const
Compute the coefficients for the basis functions.
Definition: Regression.cpp:117
Eigen::ArrayXd CenterPoints(std::vector< Eigen::VectorXd > &xs)
Center the input points.
Definition: Regression.cpp:216
const unsigned int inputDim
The input dimension.
Definition: Regression.h:148
Eigen::MatrixXd ComputeCoefficientsRHS(Eigen::MatrixXd const &vand, std::vector< Eigen::VectorXd > const &ys_data) const
Compute the right hand side given data to compute the polynomial coefficients.
Definition: Regression.cpp:140
void Fit(std::vector< Eigen::VectorXd > xs, std::vector< Eigen::VectorXd > const &ys, Eigen::VectorXd const &center)
Compute the coeffiecents of the polynomial given data.
Definition: Regression.cpp:95
Eigen::MatrixXd coeff
Coeffients for the polynomial basis.
Definition: Regression.h:169
Regression(boost::property_tree::ptree const &pt)
Definition: Regression.cpp:18
std::pair< Eigen::VectorXd, double > PoisednessConstant(std::vector< Eigen::VectorXd > xs, Eigen::VectorXd const &center, int kn=-1) const
Definition: Regression.cpp:222
Eigen::MatrixXd VandermondeMatrix(std::vector< Eigen::VectorXd > const &xs) const
Create the Vandermonde matrix.
Definition: Regression.cpp:157
Eigen::VectorXd currentCenter
Current center of the inputs.
Definition: Regression.h:163
std::shared_ptr< muq::Utilities::MultiIndexSet > multi
The multi-index to so we know the order of each term.
Definition: Regression.h:157
Provides an abstract interface for defining vector-valued model components.
Definition: ModPiece.h:148
Base class for MUQ's modelling envronment.
Definition: WorkPiece.h:40
The cost function for an optimization routine.
Definition: CostFunction.h:74
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37