Isis 3 Programmer Reference
ZeroBufferFit.h
1 #ifndef ZeroBufferFit_h
2 #define ZeroBufferFit_h
3 
10 /* SPDX-License-Identifier: CC0-1.0 */
11 
12 #include <string>
13 #include <vector>
14 
15 #include "HiCalTypes.h"
16 #include "HiCalUtil.h"
17 #include "HiCalConf.h"
18 #include "NonLinearLSQ.h"
19 #include "Module.h"
20 #include "IException.h"
21 
22 namespace Isis {
23 
42  class ZeroBufferFit : public NonLinearLSQ, public Module {
43 
44  public:
45  ZeroBufferFit(const HiCalConf &conf);
47  virtual ~ZeroBufferFit() { }
48 
54  void setBin(int bin) { _timet.setBin(bin); }
55 
61  void setLineTime(double ltime) { _timet.setLineTime(ltime); }
62 
71  inline int size() const { return (_data.dim()); }
72 
84  int nSize() const { return (_b2.dim()); }
85 
94  int nParms() const { return (4); }
95 
97  void setabsErr(double absError) { _absErr = absError; }
99  void setrelErr(double relError) { _relErr = relError; }
101  double absErr() const { return (_absErr); }
103  double relErr() const { return (_relErr); }
104 
105  HiVector Solve(const HiVector &d);
106  NLVector guess();
107  int checkIteration(const int Iter, const NLVector &fitcoefs,
108  const NLVector &uncerts, double cplxconj,
109  int Istatus);
110 
111  NLVector f_x(const NLVector &a);
112  NLMatrix df_x(const NLVector &a);
113 
115  double Chisq() const { return (_chisq); }
117  int DoF() const { return (nSize() - nParms()); }
118 
119  HiVector Yfit() const;
120  HiVector Normalize(const HiVector &v);
121 
122  private:
123  HiLineTimeEqn _timet; // This is the X data set
124  HiVector _data; // Typically will be the HiRISE buffer data
125  HiVector _b2; // Data buffer used in fitting
126  double _absErr; // Absolute error convergence test
127  double _relErr; // Relative error convergence test
128  double _maxLog; // Maximum log value to constrain
129  int _badLines; // Exclude lines at end of buffer
130  int _sWidth; // Width of guestimate filter
131  int _sIters; // Filter interations
132  bool _skipFit; // Skip fitting and pass input through
133  bool _useLinFit; // Use linear fit on failure of LM, else Zf
134  int _minLines; // Minimum number of lines to fit (Default: 100)
135  HiVector _cc; // Parameter of 2-D fit
136  HiVector _guess; // Initial guestimate of solutions
137  HiVector _coefs; // Coefficients of solution
138  HiVector _uncert; // Uncertanties
139  double _chisq; // ChiSq of NonLinear equation
140 
141  HiVector poly_fit(const HiVector &d, const double line0 = 0.0) const;
142  virtual void printOn(std::ostream &o) const;
144  int goodLines(const HiVector &d) const { return (d.dim() - _badLines); }
146  bool gotGoodLines(const HiVector &d) const {
147  return (goodLines(d) >= _minLines);
148  }
149  };
150 
151 } // namespace Isis
152 #endif
Isis::ZeroBufferFit::setrelErr
void setrelErr(double relError)
Sets the relative error parameter
Definition: ZeroBufferFit.h:99
Isis::ZeroBufferFit::Solve
HiVector Solve(const HiVector &d)
Compute non-linear fit to (typically) ZeroBufferSmooth module.
Definition: ZeroBufferFit.cpp:91
Isis::ZeroBufferFit::size
int size() const
Returns the size of the data buffer.
Definition: ZeroBufferFit.h:71
Isis::Module
Module manages HiRISE calibration vectors from various sources.
Definition: Module.h:39
Isis::ZeroBufferFit::goodLines
int goodLines(const HiVector &d) const
Returns the number of good lines in the image.
Definition: ZeroBufferFit.h:144
Isis::ZeroBufferFit::ZeroBufferFit
ZeroBufferFit(const HiCalConf &conf)
Compute second level drift correction (Zf module)
Definition: ZeroBufferFit.cpp:38
Isis::ZeroBufferFit::setabsErr
void setabsErr(double absError)
Sets the absolute error parameter.
Definition: ZeroBufferFit.h:97
Isis::ZeroBufferFit
Computes non-linear lsq fit of HiRISE Drift (Zd module)
Definition: ZeroBufferFit.h:42
Isis::ZeroBufferFit::absErr
double absErr() const
Returns the current value of the absolute error.
Definition: ZeroBufferFit.h:101
Isis::ZeroBufferFit::gotGoodLines
bool gotGoodLines(const HiVector &d) const
Determines if the vector contains any valid lines.
Definition: ZeroBufferFit.h:146
Isis::ZeroBufferFit::setLineTime
void setLineTime(double ltime)
Set scan line time.
Definition: ZeroBufferFit.h:61
Isis::ZeroBufferFit::DoF
int DoF() const
Returns the Degrees of Freedom.
Definition: ZeroBufferFit.h:117
Isis::ZeroBufferFit::Normalize
HiVector Normalize(const HiVector &v)
Compute normalized solution vector from result.
Definition: ZeroBufferFit.cpp:282
Isis::NonLinearLSQ
NonLinearLSQ Computes a fit using a Levenberg-Marquardt algorithm.
Definition: NonLinearLSQ.h:42
Isis::HiLineTimeEqn
Compute HiRISE line times.
Definition: HiCalUtil.h:361
Isis::ZeroBufferFit::printOn
virtual void printOn(std::ostream &o) const
Provides virtualized dump of data from this module.
Definition: ZeroBufferFit.cpp:319
Isis::ZeroBufferFit::df_x
NLMatrix df_x(const NLVector &a)
Computes the first derivative of the function at the current iteration.
Definition: ZeroBufferFit.cpp:245
Isis::HiVector
TNT::Array1D< double > HiVector
1-D Buffer
Definition: HiCalTypes.h:27
Isis::ZeroBufferFit::setBin
void setBin(int bin)
Set binning/summing mode.
Definition: ZeroBufferFit.h:54
Isis::ZeroBufferFit::nParms
int nParms() const
Number of parameter to be fitted.
Definition: ZeroBufferFit.h:94
Isis::ZeroBufferFit::Yfit
HiVector Yfit() const
Computes the solution vector using current coefficents.
Definition: ZeroBufferFit.cpp:266
Isis::ZeroBufferFit::relErr
double relErr() const
Returns the current value of the relative error.
Definition: ZeroBufferFit.h:103
Isis::ZeroBufferFit::poly_fit
HiVector poly_fit(const HiVector &d, const double line0=0.0) const
Compute a polyonomial fit using multivariate statistics.
Definition: ZeroBufferFit.cpp:305
Isis::ZeroBufferFit::checkIteration
int checkIteration(const int Iter, const NLVector &fitcoefs, const NLVector &uncerts, double cplxconj, int Istatus)
Computes the interation check for convergence.
Definition: ZeroBufferFit.cpp:219
Isis::ZeroBufferFit::nSize
int nSize() const
Returns the size of the fitted buffer.
Definition: ZeroBufferFit.h:84
Isis::ZeroBufferFit::f_x
NLVector f_x(const NLVector &a)
Computes the function value at the current iteration.
Definition: ZeroBufferFit.cpp:227
Isis::ZeroBufferFit::guess
NLVector guess()
Compute the initial guess of the fit.
Definition: ZeroBufferFit.cpp:162
Isis::ZeroBufferFit::~ZeroBufferFit
virtual ~ZeroBufferFit()
Destructor.
Definition: ZeroBufferFit.h:47
Isis::ZeroBufferFit::Chisq
double Chisq() const
Returns the Chi-Square value of the fit solution.
Definition: ZeroBufferFit.h:115
Isis
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16