File failed to load: https://isis.astrogeology.usgs.gov/6.0.0/Object/assets/jax/output/NativeMML/config.js
 |
Isis Developer Reference
|
Go to the documentation of this file.
11 #include "tnt/tnt_array2d.h"
103 int sparseRows = 0,
int sparseCols = 0,
bool jigsaw =
false);
105 void AddKnown(
const std::vector<double> &input,
double expected,
106 double weight = 1.0);
108 std::vector<double>
GetInput(
int row)
const;
118 double Evaluate(
const std::vector<double> &input);
121 void Weight(
int index,
double weight);
130 return p_expected.size();
138 std::vector<double>
GetEpsilons ()
const {
return p_epsilonsSparse; }
145 void SolveCholesky () {}
148 void FillSparseA(
const std::vector<double> &data);
149 bool ApplyParameterWeights();
152 std::vector<double> p_epsilonsSparse;
153 std::vector<double> p_parameterWeights;
155 arma::SpMat<double> p_sparseA;
156 arma::SpMat<double> p_normals;
158 arma::mat p_SLU_Factor;
164 int p_currentFillRow;
167 int p_constrainedParameters;
168 int p_degreesOfFreedom;
172 std::vector<std::vector<double> > p_input;
174 std::vector<double> p_expected;
176 std::vector<double> p_sqrtWeight;
179 std::vector<double> p_residuals;
LeastSquares(Isis::BasisFunction &basis, bool sparse=false, int sparseRows=0, int sparseCols=0, bool jigsaw=false)
Creates a LeastSquares Object.
Definition: LeastSquares.cpp:23
int Solve(Isis::LeastSquares::SolveMethod method=SVD)
After all the data has been registered through AddKnown, invoke this method to solve the system of eq...
Definition: LeastSquares.cpp:205
int GetDegreesOfFreedom()
Definition: LeastSquares.h:134
std::vector< double > GetInput(int row) const
This method returns the data at the given row.
Definition: LeastSquares.cpp:158
void SetCoefficients(const std::vector< double > &coefs)
Set the coefficients for the equation.
Definition: BasisFunction.cpp:42
int Knowns() const
The number of knowns (or times AddKnown was invoked) linear combination of the variables.
Definition: LeastSquares.h:129
Generic linear equation class.
Definition: BasisFunction.h:48
@ Unknown
A type of error that cannot be classified as any of the other error types.
Definition: IException.h:118
@ SPARSE
Sparse.
Definition: LeastSquares.h:114
double Evaluate(const std::vector< double > &input)
Invokes the BasisFunction Evaluate method.
Definition: LeastSquares.cpp:553
double Term(int c) const
Returns the cth term.
Definition: BasisFunction.h:97
QString Name() const
Returns the name of the equation.
Definition: BasisFunction.h:80
SolveMethod
Definition: LeastSquares.h:112
int Rows() const
This methods returns the number of rows in the matrix.
Definition: LeastSquares.cpp:187
double Residual(int i) const
Returns the ith residual.
Definition: LeastSquares.cpp:590
QString toString(bool boolToConvert)
Global function to convert a boolean to a string.
Definition: IString.cpp:211
void Reset()
Definition: LeastSquares.cpp:523
~LeastSquares()
Destroys the LeastSquares object.
Definition: LeastSquares.cpp:62
Generic least square fitting class.
Definition: LeastSquares.h:99
#define _FILEINFO_
Macro for the filename and line number.
Definition: IException.h:24
virtual void Expand(const std::vector< double > &vars)
This is the function you should replace depending on your needs.
Definition: BasisFunction.cpp:123
void SetParameterWeights(const std::vector< double > weights)
Definition: LeastSquares.h:139
std::vector< double > GetEpsilons() const
Definition: LeastSquares.h:138
double GetExpected(int row) const
This method returns the expected value at the given row.
Definition: LeastSquares.cpp:173
Isis exception class.
Definition: IException.h:91
@ QRD
QR Decomposition.
Definition: LeastSquares.h:113
int Variables() const
Returns the number of variables in the equation.
Definition: BasisFunction.h:72
@ Programmer
This error is for when a programmer made an API call that was illegal.
Definition: IException.h:146
void SetNumberOfConstrainedParameters(int n)
Definition: LeastSquares.h:140
void ResetSparse()
Definition: LeastSquares.h:137
void AddKnown(const std::vector< double > &input, double expected, double weight=1.0)
Invoke this method for each set of knowns.
Definition: LeastSquares.cpp:96
@ SVD
Singular Value Decomposition.
Definition: LeastSquares.h:112
double Evaluate(const std::vector< double > &vars)
Compute the equation using the input variables.
Definition: BasisFunction.cpp:64
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16
int Coefficients() const
Returns the number of coefficients for the equation.
Definition: BasisFunction.h:64
void Weight(int index, double weight)
Reset the weight for the ith known.
Definition: LeastSquares.cpp:614
double GetSigma0()
Definition: LeastSquares.h:133
std::vector< double > Residuals() const
Returns a vector of residuals (errors).
Definition: LeastSquares.cpp:570