Isis 3 Programmer Reference
|
Computes non-linear lsq fit of HiRISE Drift (Zd module) More...
#include <ZeroBufferFit.h>
Public Types | |
typedef TNT::Array1D< double > | NLVector |
typedef TNT::Array2D< double > | NLMatrix |
Public Member Functions | |
ZeroBufferFit (const HiCalConf &conf) | |
Compute second level drift correction (Zf module) More... | |
virtual | ~ZeroBufferFit () |
Destructor. More... | |
void | setBin (int bin) |
Set binning/summing mode. More... | |
void | setLineTime (double ltime) |
Set scan line time. More... | |
int | size () const |
Returns the size of the data buffer. More... | |
int | nSize () const |
Returns the size of the fitted buffer. More... | |
int | nParms () const |
Number of parameter to be fitted. More... | |
void | setabsErr (double absError) |
Sets the absolute error parameter. More... | |
void | setrelErr (double relError) |
Sets the relative error parameter. More... | |
double | absErr () const |
Returns the current value of the absolute error. More... | |
double | relErr () const |
Returns the current value of the relative error. More... | |
HiVector | Solve (const HiVector &d) |
Compute non-linear fit to (typically) ZeroBufferSmooth module. More... | |
NLVector | guess () |
Compute the initial guess of the fit. More... | |
int | checkIteration (const int Iter, const NLVector &fitcoefs, const NLVector &uncerts, double cplxconj, int Istatus) |
Computes the interation check for convergence. More... | |
NLVector | f_x (const NLVector &a) |
Computes the function value at the current iteration. More... | |
NLMatrix | df_x (const NLVector &a) |
Computes the first derivative of the function at the current iteration. More... | |
double | Chisq () const |
Returns the Chi-Square value of the fit solution. More... | |
int | DoF () const |
Returns the Degrees of Freedom. More... | |
HiVector | Yfit () const |
Computes the solution vector using current coefficents. More... | |
HiVector | Normalize (const HiVector &v) |
Compute normalized solution vector from result. More... | |
void | setMaxIters (int m) |
Sets the maximum number of iterations. More... | |
int | maxIters () const |
Maximum number iterations for valid solution. More... | |
int | curvefit () |
int | status () const |
Return status of last fit processing. More... | |
bool | success () const |
Determine success from last fit processing. More... | |
bool | success (int status) const |
Check status for success of the given condition. More... | |
std::string | statusstr () const |
Return error message pertaining to last fit procesing. More... | |
std::string | statusstr (int status) const |
Return error message given status condition. More... | |
NLVector | coefs () const |
Return coefficients from last fit processing. More... | |
NLVector | uncert () const |
Return uncertainties from last fit processing. More... | |
int | nIterations () const |
Return number of iterations from last fit processing. More... | |
QString | name () const |
Returns name of component. More... | |
QString | getcsvFile () const |
Returns expanded name of last CSV file loaded by loadCsv. More... | |
virtual void | Process (const Module &c) |
Invokes the process method on the Module vector. More... | |
virtual void | Process (const HiVector &v) |
Default processing behavior makes a reference copy of data array. More... | |
HiVector | loadCsv (const QString &csvBase, const HiCalConf &conf, const DbProfile &prof, const int &elements=0) |
Provide generic loading of CSV file for all modules. More... | |
const HiVector & | ref () const |
Return data via a const reference. More... | |
double | operator() (int index) const |
Return individual element of the data vector. More... | |
const HiHistory & | History () const |
Return recorded history of events. More... | |
virtual void | record (PvlGroup &pvl, const QString keyname="ModuleHistory") const |
Record history in Pvl group object. More... | |
void | Dump (const QString &fname) const |
Dumps the component to a specified file. More... | |
Protected Types | |
enum | { DefaultWidth = 10, DefaultPrecision = 6 } |
Protected Member Functions | |
void | Terminate (const std::string &message="") |
void | Abort (const std::string &reason="") |
bool | doContinue () const |
QString | formatDbl (const double &value) const |
Properly format values that could be special pixels. More... | |
Protected Attributes | |
QString | _name |
Name of component. More... | |
QString | _csvFile |
Fully expanded name of CSV file if present. More... | |
HiHistory | _history |
Hierarchial component history. More... | |
int | _fmtWidth |
Default field with of double. More... | |
int | _fmtPrecision |
Default field with of double. More... | |
Private Member Functions | |
HiVector | poly_fit (const HiVector &d, const double line0=0.0) const |
Compute a polyonomial fit using multivariate statistics. More... | |
virtual void | printOn (std::ostream &o) const |
Provides virtualized dump of data from this module. More... | |
int | goodLines (const HiVector &d) const |
Returns the number of good lines in the image. More... | |
bool | gotGoodLines (const HiVector &d) const |
Determines if the vector contains any valid lines. More... | |
Private Attributes | |
HiLineTimeEqn | _timet |
HiVector | _data |
HiVector | _b2 |
double | _absErr |
double | _relErr |
double | _maxLog |
int | _badLines |
int | _sWidth |
int | _sIters |
bool | _skipFit |
bool | _useLinFit |
int | _minLines |
HiVector | _cc |
HiVector | _guess |
HiVector | _coefs |
HiVector | _uncert |
double | _chisq |
Computes non-linear lsq fit of HiRISE Drift (Zd module)
This class is best used with individual HiRISE images as the number of lines is critical to proper use. It is best applied by getting the buffer as a reference and applying it during systematic processing.
This class models the drift correction using the Levenberg-Marquardt algorithm.
2008-05-23 Kris Becker Added the ZdOnFailUseLinear option.
2010-10-28 Kris Becker Updated parameter names removing "Zd" and replacing with "ZeroBufferFit".
Definition at line 56 of file ZeroBufferFit.h.
Isis::ZeroBufferFit::ZeroBufferFit | ( | const HiCalConf & | conf | ) |
Compute second level drift correction (Zf module)
This class provides the second level drift correction that is
conf |
Definition at line 53 of file ZeroBufferFit.cpp.
References Isis::Module::_history, Isis::ConfKey(), Isis::DbProfile::exists(), Isis::IsEqual(), Isis::IsTrueValue(), Isis::NonLinearLSQ::maxIters(), Isis::DbProfile::Name(), Isis::NonLinearLSQ::setMaxIters(), Isis::toDouble(), Isis::ToDouble(), Isis::toInt(), Isis::ToInteger(), and Isis::ToString().
|
inlinevirtual |
Destructor.
Definition at line 61 of file ZeroBufferFit.h.
|
inlinevirtual |
Returns the current value of the absolute error.
Reimplemented from Isis::NonLinearLSQ.
Definition at line 115 of file ZeroBufferFit.h.
|
virtual |
Computes the interation check for convergence.
Iter | Current iteration |
fitcoefs | Vector of current fit coefficients |
uncerts | Uncertainties |
cplxconj | Complex conjugate of the current iteration |
Istatus | State of current iteration |
Reimplemented from Isis::NonLinearLSQ.
Definition at line 234 of file ZeroBufferFit.cpp.
|
inline |
Returns the Chi-Square value of the fit solution.
Definition at line 129 of file ZeroBufferFit.h.
Referenced by Solve().
|
inlineinherited |
Return coefficients from last fit processing.
Definition at line 118 of file NonLinearLSQ.h.
Referenced by Solve().
|
virtual |
Computes the first derivative of the function at the current iteration.
Implements Isis::NonLinearLSQ.
Definition at line 260 of file ZeroBufferFit.cpp.
References Isis::MIN().
|
inline |
Returns the Degrees of Freedom.
Definition at line 131 of file ZeroBufferFit.h.
References nParms(), and nSize().
Referenced by Solve().
|
inlineinherited |
Dumps the component to a specified file.
fname | Name of file to dump contents to |
Definition at line 146 of file Module.h.
References _FILEINFO_, Isis::FileName::expanded(), and Isis::IException::User.
|
virtual |
Computes the function value at the current iteration.
Implements Isis::NonLinearLSQ.
Definition at line 242 of file ZeroBufferFit.cpp.
References Isis::MIN().
|
inlineprotectedinherited |
Properly format values that could be special pixels.
This method applies ISIS special pixel value conventions to properly print pixel values.
[in] | (double) | value Input value to test for specialness and print as requested by caller |
[in] | (int) | width Width of field in which to print the value |
[in] | (int) | prec Precision used to format the value |
Definition at line 184 of file Module.h.
References Isis::Module::_fmtPrecision, Isis::Module::_fmtWidth, Isis::IsSpecial(), and Isis::PixelToString().
Referenced by Isis::ZeroBufferSmooth::printOn(), printOn(), Isis::ZeroDark::printOn(), Isis::ZeroReverse::printOn(), and Isis::Module::printOn().
|
inlineinherited |
Returns expanded name of last CSV file loaded by loadCsv.
Definition at line 86 of file Module.h.
References Isis::Module::_csvFile.
|
inlineprivate |
Returns the number of good lines in the image.
Definition at line 158 of file ZeroBufferFit.h.
Referenced by gotGoodLines(), and Solve().
|
inlineprivate |
Determines if the vector contains any valid lines.
Definition at line 160 of file ZeroBufferFit.h.
References goodLines().
|
virtual |
Compute the initial guess of the fit.
This method provides the non-linear fit with an initial guess of the solution. It involves a linear fit to the latter half of the data to provide the first two coefficents, the difference of the averages of the residuals at both ends of the data set and 5 times the last line time as the final (fourth) element...a bit involved really.
Implements Isis::NonLinearLSQ.
Definition at line 177 of file ZeroBufferFit.cpp.
References Isis::Module::_history, Isis::Statistics::AddData(), Isis::Statistics::Average(), Isis::MIN(), poly_fit(), Isis::Module::ref(), Isis::Statistics::Reset(), and Isis::ToString().
|
inlineinherited |
Return recorded history of events.
Definition at line 131 of file Module.h.
References Isis::Module::_history.
Referenced by Isis::ZeroReverse::init(), and Isis::ZeroBufferSmooth::init().
|
inlineinherited |
Provide generic loading of CSV file for all modules.
This method provides generalized access to CSV files through the standardized format.
cvsBase | Name of base keyword for CSV file |
conf | Configuration parameters |
prof | Module profile parameters |
samples | Number of expect elements to be read from CSV file |
Definition at line 116 of file Module.h.
References Isis::Module::_csvFile, and Isis::Module::_history.
|
inlineinherited |
Maximum number iterations for valid solution.
Definition at line 85 of file NonLinearLSQ.h.
Referenced by ZeroBufferFit().
|
inlineinherited |
|
inlineinherited |
Return number of iterations from last fit processing.
Definition at line 122 of file NonLinearLSQ.h.
Referenced by Solve().
Compute normalized solution vector from result.
Definition at line 297 of file ZeroBufferFit.cpp.
References Isis::Module::_history, and Isis::ToString().
|
inlinevirtual |
Number of parameter to be fitted.
This is the number of parameters that ZeroBufferFit needs to fit. This method is a requirement of the NonLinearLSQ class.
Implements Isis::NonLinearLSQ.
Definition at line 108 of file ZeroBufferFit.h.
Referenced by DoF().
|
inlinevirtual |
Returns the size of the fitted buffer.
Important: This returns the size of the buffer being fitted and not the size of original data buffer. This is a requirement of the NonLinearLSQ class. USE WITH CAUTION!
Implements Isis::NonLinearLSQ.
Definition at line 98 of file ZeroBufferFit.h.
Referenced by DoF().
|
inlineinherited |
Return individual element of the data vector.
Definition at line 128 of file Module.h.
References Isis::Module::_data.
|
private |
Compute a polyonomial fit using multivariate statistics.
Used as a fallback solution, this method computes a linear statistical solution from the linear regression analysis of the multivariate statistics of the data.
d | Data vector to fit |
line0 | Current line number in the image |
Definition at line 320 of file ZeroBufferFit.cpp.
References Isis::MultivariateStatistics::AddData(), and Isis::MultivariateStatistics::LinearRegression().
|
privatevirtual |
Provides virtualized dump of data from this module.
Reimplemented from Isis::Module.
Definition at line 334 of file ZeroBufferFit.cpp.
References Isis::Module::_fmtWidth, Isis::Module::_history, Isis::Module::formatDbl(), and Yfit().
|
inlinevirtualinherited |
Invokes the process method on the Module vector.
Definition at line 91 of file Module.h.
References Isis::Module::ref().
|
inlinevirtualinherited |
Default processing behavior makes a reference copy of data array.
Reimplemented in Isis::SplineFill, and Isis::LowPassFilter.
Definition at line 96 of file Module.h.
References Isis::Module::_data.
|
inlinevirtualinherited |
Record history in Pvl group object.
Definition at line 134 of file Module.h.
References Isis::Module::_history.
|
inlineinherited |
Return data via a const reference.
Definition at line 126 of file Module.h.
References Isis::Module::_data.
Referenced by guess(), Isis::ZeroReverse::init(), Isis::ZeroBufferSmooth::init(), and Isis::Module::Process().
|
inlinevirtual |
Returns the current value of the relative error.
Reimplemented from Isis::NonLinearLSQ.
Definition at line 117 of file ZeroBufferFit.h.
|
inline |
Sets the absolute error parameter.
Definition at line 111 of file ZeroBufferFit.h.
|
inline |
Set binning/summing mode.
bin | Summing mode of observatio |
Definition at line 68 of file ZeroBufferFit.h.
|
inline |
|
inlineinherited |
Sets the maximum number of iterations.
m | User provides the maximum number iterations |
Definition at line 77 of file NonLinearLSQ.h.
Referenced by ZeroBufferFit().
|
inline |
Sets the relative error parameter.
Definition at line 113 of file ZeroBufferFit.h.
|
inline |
Returns the size of the data buffer.
This is the size of the original data buffer.
Definition at line 85 of file ZeroBufferFit.h.
Compute non-linear fit to (typically) ZeroBufferSmooth module.
This method computes a non-linear fit to the result of the ZeroBufferSmooth module. There are several things that can go wrong and some config conditions that dictate behavior of this process.
Should the image be a short exposure (i.e., not many lines) the fit will not succeed so it simply skips this entire module providing the input result (d) as the solution. This wil also occur when the user has selected the skip option for the module.
A fit is attempted on the ZeroBufferSmooth data. The non-linear solution must converge within the specifed number of iterations (MaximumIterations) or a polynomial fit will be used in leui of a valid solution.
d | ZeroBufferSmooth data solution as input to this method |
Definition at line 106 of file ZeroBufferFit.cpp.
References Isis::Module::_history, Chisq(), Isis::NonLinearLSQ::coefs(), DoF(), goodLines(), gotGoodLines(), Isis::NonLinearLSQ::nIterations(), poly_fit(), Isis::NonLinearLSQ::statusstr(), Isis::NonLinearLSQ::success(), Isis::ToString(), Isis::NonLinearLSQ::uncert(), and Yfit().
|
inlineinherited |
Return status of last fit processing.
Definition at line 96 of file NonLinearLSQ.h.
Referenced by Isis::NonLinearLSQ::statusstr(), and Isis::NonLinearLSQ::success().
|
inlineinherited |
Return error message pertaining to last fit procesing.
Definition at line 102 of file NonLinearLSQ.h.
Referenced by Solve().
|
inlineinherited |
Return error message given status condition.
Definition at line 106 of file NonLinearLSQ.h.
References Isis::NonLinearLSQ::status().
|
inlineinherited |
Determine success from last fit processing.
Definition at line 98 of file NonLinearLSQ.h.
Referenced by Solve().
|
inlineinherited |
Check status for success of the given condition.
Definition at line 100 of file NonLinearLSQ.h.
References Isis::NonLinearLSQ::status().
|
inlineinherited |
Return uncertainties from last fit processing.
Definition at line 120 of file NonLinearLSQ.h.
Referenced by Solve().
HiVector Isis::ZeroBufferFit::Yfit | ( | ) | const |
Computes the solution vector using current coefficents.
Definition at line 281 of file ZeroBufferFit.cpp.
References gotGoodLines().
|
protectedinherited |
Fully expanded name of CSV file if present.
Definition at line 165 of file Module.h.
Referenced by Isis::Module::getcsvFile(), and Isis::Module::loadCsv().
|
protectedinherited |
Default field with of double.
Definition at line 169 of file Module.h.
Referenced by Isis::Module::formatDbl().
|
protectedinherited |
Default field with of double.
Definition at line 168 of file Module.h.
Referenced by Isis::Module::formatDbl(), Isis::ZeroBufferSmooth::printOn(), printOn(), Isis::ZeroDark::printOn(), and Isis::ZeroReverse::printOn().
|
protectedinherited |
Hierarchial component history.
Definition at line 167 of file Module.h.
Referenced by Isis::GainUnitConversion::getTempDepGain(), guess(), Isis::Module::History(), Isis::ZeroReverse::init(), Isis::ZeroBufferSmooth::init(), Isis::Module::loadCsv(), Normalize(), Isis::ZeroBufferSmooth::printOn(), printOn(), Isis::ZeroDark::printOn(), Isis::ZeroReverse::printOn(), Isis::Module::printOn(), Isis::LowPassFilter::Process(), Isis::SplineFill::Process(), Isis::Module::record(), Solve(), and ZeroBufferFit().
|
protectedinherited |