casacore
|
#include <LSQFit.h>
Classes | |
struct | AsReal |
struct | Complex |
struct | Conjugate |
struct | Real |
Simple classes to overload templated memberfunctions. More... | |
struct | Separable |
Public Types | |
enum | ReadyCode { NONREADY , SOLINCREMENT , DERIVLEVEL , MAXITER , NOREDUCTION , SINGULAR , N_ReadyCode } |
State of the non-linear solution. More... | |
enum | ErrorField { NC , SUMWEIGHT , SUMLL , CHI2 , N_ErrorField } |
Offset of fields in error_p data area. More... | |
Public Member Functions | |
LSQFit (uInt nUnknowns, uInt nConstraints=0) | |
Construct an object with the number of unknowns and constraints, using the default collinearity factor and the default Levenberg-Marquardt adjustment factor. More... | |
LSQFit (uInt nUnknowns, const LSQReal &, uInt nConstraints=0) | |
Allow explicit Real specification. More... | |
LSQFit (uInt nUnknowns, const LSQComplex &, uInt nConstraints=0) | |
Allow explicit Complex specification. More... | |
LSQFit () | |
Default constructor (empty, only usable after a set(nUnknowns) ) More... | |
LSQFit (const LSQFit &other) | |
Copy constructor (deep copy) More... | |
LSQFit & | operator= (const LSQFit &other) |
Assignment (deep copy) More... | |
~LSQFit () | |
Bool | invert (uInt &nRank, Bool doSVD=False) |
Triangularize the normal equations and determine the rank nRank of the normal equations and, in the case of an SVD solution, the constraint equations. More... | |
template<class U > | |
void | copy (const Double *beg, const Double *end, U &sol, LSQReal) |
Copy date from beg to end; converting if necessary to complex data. More... | |
template<class U > | |
void | copy (const Double *beg, const Double *end, U &sol, LSQComplex) |
template<class U > | |
void | copy (const Double *beg, const Double *end, U *sol, LSQReal) |
template<class U > | |
void | copy (const Double *beg, const Double *end, U *sol, LSQComplex) |
template<class U > | |
void | uncopy (Double *beg, const Double *end, U &sol, LSQReal) |
template<class U > | |
void | uncopy (Double *beg, const Double *end, U &sol, LSQComplex) |
template<class U > | |
void | uncopy (Double *beg, const Double *end, U *sol, LSQReal) |
template<class U > | |
void | uncopy (Double *beg, const Double *end, U *sol, LSQComplex) |
template<class U > | |
void | copyDiagonal (U &errors, LSQReal) |
template<class U > | |
void | copyDiagonal (U &errors, LSQComplex) |
template<class U > | |
void | solve (U *sol) |
Solve normal equations. More... | |
template<class U > | |
void | solve (std::complex< U > *sol) |
template<class U > | |
void | solve (U &sol) |
template<class U > | |
Bool | solveLoop (uInt &nRank, U *sol, Bool doSVD=False) |
Solve a loop in a non-linear set. More... | |
template<class U > | |
Bool | solveLoop (uInt &nRank, std::complex< U > *sol, Bool doSVD=False) |
template<class U > | |
Bool | solveLoop (uInt &nRank, U &sol, Bool doSVD=False) |
template<class U > | |
Bool | solveLoop (Double &fit, uInt &nRank, U *sol, Bool doSVD=False) |
template<class U > | |
Bool | solveLoop (Double &fit, uInt &nRank, std::complex< U > *sol, Bool doSVD=False) |
template<class U > | |
Bool | solveLoop (Double &fit, uInt &nRank, U &sol, Bool doSVD=False) |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const U &obs, Bool doNorm=True, Bool doKnown=True) |
Make normal equations using the cEq condition equation (cArray) (with nUnknowns elements) and a weight weight , given the known observed value obs . More... | |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const U &obs, LSQFit::Real, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const std::complex< U > &obs, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Complex, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Separable, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::AsReal, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Conjugate, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const U &obs, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const V &cEq2, const U &weight, const U &obs, const U &obs2, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const U &obs, LSQFit::Real, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const std::complex< U > &obs, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Complex, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Separable, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::AsReal, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNorm (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Conjugate, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const U &obs, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const U &obs, LSQFit::Real, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const std::complex< U > &obs, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Complex, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Separable, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const std::complex< U > &obs, LSQFit::AsReal, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V > | |
void | makeNorm (const std::vector< std::pair< uInt, V > > &cEq, const U &weight, const std::complex< U > &obs, LSQFit::Conjugate, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNormSorted (uInt nIndex, const W &cEqIndex, const V &cEq, const U &weight, const U &obs, Bool doNorm=True, Bool doKnown=True) |
template<class U , class V , class W > | |
void | makeNormSorted (uInt nIndex, const W &cEqIndex, const V &cEq, const V &cEq2, const U &weight, const U &obs, const U &obs2, Bool doNorm=True, Bool doKnown=True) |
template<class U > | |
Bool | getConstraint (uInt n, U *cEq) const |
Get the n-th (from 0 to the rank deficiency, or missing rank, see e.g. More... | |
template<class U > | |
Bool | getConstraint (uInt n, std::complex< U > *cEq) const |
template<class U > | |
Bool | getConstraint (uInt n, U &cEq) const |
template<class U , class V > | |
Bool | setConstraint (uInt n, const V &cEq, const U &obs) |
Add a new constraint equation (updating nConstraints); or set a numbered constraint equation (0..nConstraints-1). More... | |
template<class U , class V > | |
Bool | setConstraint (uInt n, const V &cEq, const std::complex< U > &obs) |
template<class U , class V , class W > | |
Bool | setConstraint (uInt n, uInt nIndex, const W &cEqIndex, const V &cEq, const U &obs) |
template<class U , class V , class W > | |
Bool | setConstraint (uInt n, uInt nIndex, const W &cEqIndex, const V &cEq, const std::complex< U > &obs) |
template<class U , class V > | |
Bool | addConstraint (const V &cEq, const U &obs) |
template<class U , class V > | |
Bool | addConstraint (const V &cEq, const std::complex< U > &obs) |
template<class U , class V , class W > | |
Bool | addConstraint (uInt nIndex, const W &cEqIndex, const V &cEq, const U &obs) |
template<class U , class V , class W > | |
Bool | addConstraint (uInt nIndex, const W &cEqIndex, const V &cEq, const std::complex< U > &obs) |
Bool | merge (const LSQFit &other) |
Merge other LSQFit object (i.e. More... | |
Bool | merge (const LSQFit &other, uInt nIndex, const uInt *nEqIndex) |
Bool | merge (const LSQFit &other, uInt nIndex, const std::vector< uInt > &nEqIndex) |
template<class W > | |
Bool | merge (const LSQFit &other, uInt nIndex, const W &nEqIndex) |
void | reset () |
Reset status to empty. More... | |
void | set (uInt nUnknowns, uInt nConstraints=0) |
Set new sizes (default is for Real) More... | |
void | set (Int nUnknowns, Int nConstraints=0) |
void | set (uInt nUnknowns, const LSQReal &, uInt nConstraints=0) |
void | set (Int nUnknowns, const LSQReal &, Int nConstraints=0) |
void | set (uInt nUnknowns, const LSQComplex &, uInt nConstraints=0) |
void | set (Int nUnknowns, const LSQComplex &, Int nConstraints=0) |
void | set (Double factor=1e-6, Double LMFactor=1e-3) |
Set new factors (collinearity factor , and Levenberg-Marquardt LMFactor ) More... | |
void | setEpsValue (Double epsval=1e-8) |
Set new value solution test. More... | |
void | setEpsDerivative (Double epsder=1e-8) |
Set new derivative test. More... | |
void | setMaxIter (uInt maxiter=0) |
Set maximum number of iterations. More... | |
uInt | nIterations () const |
Get number of iterations done. More... | |
void | setBalanced (Bool balanced=False) |
Set the expected form of the normal equations. More... | |
LSQFit::ReadyCode | isReady () const |
Ask the state of the non-linear solutions. More... | |
const std::string & | readyText () const |
template<class U > | |
Bool | getCovariance (U *covar) |
Get the covariance matrix (of size nUnknowns * nUnknowns ) More... | |
template<class U > | |
Bool | getCovariance (std::complex< U > *covar) |
template<class U > | |
Bool | getErrors (U *errors) |
Get main diagonal of covariance function (of size nUnknowns ) More... | |
template<class U > | |
Bool | getErrors (std::complex< U > *errors) |
template<class U > | |
Bool | getErrors (U &errors) |
uInt | nUnknowns () const |
Get the number of unknowns. More... | |
uInt | nConstraints () const |
Get the number of constraints. More... | |
uInt | getDeficiency () const |
Get the rank deficiency Warning: Note that the number is returned assuming real values; For complex values it has to be halved More... | |
Double | getChi () const |
Get chi^2 (both are identical); the standard deviation (per observation) and the standard deviation per weight unit. More... | |
Double | getChi2 () const |
Double | getSD () const |
Double | getWeightedSD () const |
void | debugIt (uInt &nun, uInt &np, uInt &ncon, uInt &ner, uInt &rank, Double *&nEq, Double *&known, Double *&constr, Double *&er, uInt *&piv, Double *&sEq, Double *&sol, Double &prec, Double &nonlin) const |
Debug: More... | |
Bool | fromRecord (String &error, const RecordInterface &in) |
Create an LSQFit object from a record. More... | |
Bool | toRecord (String &error, RecordInterface &out) const |
Create a record from an LSQFit object. More... | |
const String & | ident () const |
Get identification of record. More... | |
void | toAipsIO (AipsIO &) const |
Save or restore using AipsIO. More... | |
void | fromAipsIO (AipsIO &) |
Static Public Attributes | |
static Real | REAL |
And values to use. More... | |
static Complex | COMPLEX |
static Separable | SEPARABLE |
static AsReal | ASREAL |
static Conjugate | CONJUGATE |
Protected Types | |
enum | StateBit { INVERTED , TRIANGLE , NONLIN , N_StateBit } |
Bits that can be set/referenced. More... | |
Protected Member Functions | |
Double * | rowrt (uInt i) const |
Get pointer in rectangular array. More... | |
Double * | rowru (uInt i) const |
void | init () |
Initialise areas. More... | |
void | clear () |
Clear areas. More... | |
void | deinit () |
De-initialise area. More... | |
void | solveIt () |
Solve normal equations. More... | |
Bool | solveItLoop (Double &fit, uInt &nRank, Bool doSVD=False) |
One non-linear LM loop. More... | |
void | solveMR (uInt nin) |
Solve missing rank part. More... | |
Bool | invertRect () |
Invert rectangular matrix (i.e. More... | |
Double | normSolution (const Double *sol) const |
Get the norm of the current solution vector. More... | |
Double | normInfKnown (const Double *known) const |
Get the infinite norm of the known vector. More... | |
Bool | mergeIt (const LSQFit &other, uInt nIndex, const uInt *nEqIndex) |
Merge sparse normal equations. More... | |
void | save (Bool all=True) |
Save current status (or part) More... | |
void | restore (Bool all=True) |
Restore current status. More... | |
void | copy (const LSQFit &other, Bool all=True) |
Copy data. More... | |
void | extendConstraints (uInt n) |
Extend the constraint equation area to the specify number of equations. More... | |
void | createNCEQ () |
Create the solution equation area nceq_p and fill it. More... | |
void | getWorkSOL () |
Get work areas for solutions, covariance. More... | |
void | getWorkCOV () |
Static Protected Member Functions | |
static Double | realMC (const std::complex< Double > &x, const std::complex< Double > &y) |
Calculate the real or imag part of x*conj(y) More... | |
static Double | imagMC (const std::complex< Double > &x, const std::complex< Double > &y) |
static Float | realMC (const std::complex< Float > &x, const std::complex< Float > &y) |
static Float | imagMC (const std::complex< Float > &x, const std::complex< Float > &y) |
Protected Attributes | |
uInt | state_p |
Bits set to indicate state. More... | |
uInt | nun_p |
Number of unknowns. More... | |
uInt | ncon_p |
Number of constraints. More... | |
uInt | n_p |
Matrix size (will be n_p = nun_p + ncon_p) More... | |
uInt | r_p |
Rank of normal equations (normally n_p) More... | |
Double | prec_p |
Collinearity precision. More... | |
Double | startnon_p |
Levenberg start factor. More... | |
Double | nonlin_p |
Levenberg current factor. More... | |
Double | stepfactor_p |
Levenberg step factor. More... | |
Double | epsval_p |
Test value for [incremental] solution in non-linear loop. More... | |
Double | epsder_p |
Test value for known vector in non-linear loop. More... | |
Bool | balanced_p |
Indicator for a well balanced normal equation. More... | |
uInt | maxiter_p |
Maximum number of iterations for non-linear solution. More... | |
uInt | niter_p |
Iteration count for non-linear solution. More... | |
ReadyCode | ready_p |
Indicate the non-linear state. More... | |
uInt * | piv_p |
Pivot table (n_p) More... | |
LSQMatrix * | norm_p |
Normal equations (triangular nun_p * nun_p) More... | |
uInt | nnc_p |
Current length nceq_p. More... | |
LSQMatrix * | nceq_p |
Normal combined with constraint equations for solutions (triangular nnc_p*nnc_p) More... | |
Double * | known_p |
Known part equations (n_p) More... | |
Double * | error_p |
Counts for errors (N_ErrorField) More... | |
Double * | constr_p |
Constraint equation area (nun_p*ncon_p)) More... | |
Double * | sol_p |
Solution area (n_p) More... | |
LSQFit * | nar_p |
Save area for non-linear case (size determined internally) More... | |
Double * | lar_p |
Save area for non-symmetric (i.e. More... | |
Double * | wsol_p |
Work areas for interim solutions and covariance. More... | |
Double * | wcov_p |
Static Protected Attributes | |
static const String | recid |
Record field names. More... | |
static const String | state |
static const String | nun |
static const String | ncon |
static const String | prec |
static const String | startnon |
static const String | nonlin |
static const String | rank |
static const String | nnc |
static const String | piv |
static const String | constr |
static const String | known |
static const String | errors |
static const String | sol |
static const String | lar |
static const String | wsol |
static const String | wcov |
static const String | nceq |
static const String | nar |
Basic class for the least squares fitting
From Least SQuares and Fitting
The LSQFit class contains the basic functions to do all the fitting described in the Note about fitting. It handles real, and complex equations;
linear and non-linear (Levenberg-Marquardt) solutions;
regular (with optional constraints) or Singular Value Decomposition (SVD
).
In essence they are a set of routines to generate normal equations (makeNorm()
) in triangular form from a set of condition equations;
to do a Cholesky-type decomposition of the normal equations (either regular or SVD
) and test its rank (invert()
);
to do a quasi inversion of the decomposed equations (solve()
) to obtain the solution and/or the errors.
All calculations are done in place. Methods to obtain additional information about the fitting process are available.
This class can be used as a stand-alone class outside of the Casacore environment. In that case the aips.h include file can be replaced if necessary by appropriate typedefs for Double, Float and uInt.
The interface to the methods have standard data or standard STL iterator arguments only. They can be used with any container having an STL random-access iterator interface. Especially they can be used with carrays
(necessary templates provided), Casacore Vectors (necessary templates provided in LSQaips
), standard random access STL containers (like std::vector
).
The normal operation of the class consists of the following steps:
Create an LSQFit object. The information that can be provided in the constructor of the object, either directly, or indirectly using the set()
commands, is (see Note 224):
Separately settable are:
Create the normal equations used in solving the set of condition equations of the user, by using the makeNorm()
methods. Separate makenorm()
methods are provided for sparse condition equations (e.g. if data for 3 antennas are provided, rather than for all 64)
If there are user provided constraints, either limiting constraints like the sum of the angles in a triangle is 180 degrees, or constraints to add missing information if e.g. only differences between parameters have been measured, they can be added to the normal equations with the setConstraint()
or the addConstraint()
methods. Lagrange multipliers will be used to solve the extended normal equations.
The normal equations are triangu;arised (using the collinearity factor as a check for solvability) with the invert()
method. If the normal equations are non-solvable an error is returned, or a switch to an SVD solution is made if indicated in the invert
call.
The solutions and adjustment errors are obtained with the solve()
method. A non-linear loop in a Levenberg-Marquardt adjustment can be obtained (together with convergence information), with the solveLoop()
method (see below) replacing the combination of invert
and solve
.
Non-linear loops are done by looping through the data using makeNorm()
calls, and upgrade the solution with the solveLoop()
method. The normal equations are upgraded by changing LM factors. Upgrade depends on the 'balanced' factor. The LM factor is either added in some way to all diagonal elements (if balanced) or all diagonal elements are multiplied by (1+factor)
After each loop convergence can be tested by the isReady()
call; which will return False
or a non-zero code indicating the ready reason. Reasons for stopping can be:
setEpsValue()
, default 1e-8) value. setEpsDerivative()
, default 1e-8, value. solveLoop
call, which is cost-free getCovariance(), getErrors()
, getChi()
(or getChi2
), getSD
and getWeightedSD
methods after a solve()
or after the final loop in a non-linear solution (of course, when necessary only). An LSQFit object can be re-used by issuing the reset()
command, or set()
of new values. If an unknown has not been used in the condition equations at all, the doDiagonal()
will make sure a proper solution is obtained, with missing unknowns zeroed.
Most of the calculations are done in place; however, enough data is saved that it is possible to continue with the same (partial) normal equations after e.g. an interim solution.
If the normal equations are produced in separate partial sets (e.g. in a multi-processor environment) a merge()
method can combine them.
Tip: It is suggested to add any possible constraint equations after the merge;
A debugIt()
method provides read access to all internal information.
The member definitions are split over three files. The second one contains the templated member function definitions, to bypass the problem of duplicate definitions of non-templated members when pre-compiling them. The third contains methods for saving objects as Records or through AipsIO.
Warning: No boundary checks on input and output containers is done for faster execution; In general these tests should be done at the higher level routines, like the LinearFit and NonLinearFit classes which should be checked for usage of LSQFit;
The contents can be saved in a record (toRecord
), and an object can be created from a record (fromRecord
). The record identifier is 'lfit'.
The object can also be saved or restored using AipsIO.
See the tLSQFit.cc and tLSQaips.cc program for extensive examples.
The following example will first create 2 condition equations for 3 unknowns (the third is degenerate). It will first create normal equations for a 2 unknown solution and solve; then it will create normal equations for a 3 unknown solution, and solve (note that the degenerate will be set to 0. The last one will use SVD and one condition equation.r
Which will produce the output:
The class was written to be able to do complex, real standard and SVD solutions in a simple and fast way.
|
protected |
Construct an object with the number of unknowns and constraints, using the default collinearity factor and the default Levenberg-Marquardt adjustment factor.
Assume real
Allow explicit Real specification.
casacore::LSQFit::LSQFit | ( | uInt | nUnknowns, |
const LSQComplex & | , | ||
uInt | nConstraints = 0 |
||
) |
Allow explicit Complex specification.
casacore::LSQFit::LSQFit | ( | ) |
Default constructor (empty, only usable after a set(nUnknowns)
)
casacore::LSQFit::LSQFit | ( | const LSQFit & | other | ) |
Copy constructor (deep copy)
casacore::LSQFit::~LSQFit | ( | ) |
Bool casacore::LSQFit::addConstraint | ( | const V & | cEq, |
const std::complex< U > & | obs | ||
) |
Bool casacore::LSQFit::addConstraint | ( | const V & | cEq, |
const U & | obs | ||
) |
Bool casacore::LSQFit::addConstraint | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const std::complex< U > & | obs | ||
) |
Bool casacore::LSQFit::addConstraint | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | obs | ||
) |
|
protected |
Clear areas.
void casacore::LSQFit::copy | ( | const Double * | beg, |
const Double * | end, | ||
U & | sol, | ||
LSQComplex | |||
) |
void casacore::LSQFit::copy | ( | const Double * | beg, |
const Double * | end, | ||
U & | sol, | ||
LSQReal | |||
) |
Copy date from beg to end; converting if necessary to complex data.
void casacore::LSQFit::copy | ( | const Double * | beg, |
const Double * | end, | ||
U * | sol, | ||
LSQComplex | |||
) |
void casacore::LSQFit::copy | ( | const Double * | beg, |
const Double * | end, | ||
U * | sol, | ||
LSQReal | |||
) |
Copy data.
If all False, only the relevant data for non-linear solution are copied (normal equations, knows and errors).
void casacore::LSQFit::copyDiagonal | ( | U & | errors, |
LSQComplex | |||
) |
void casacore::LSQFit::copyDiagonal | ( | U & | errors, |
LSQReal | |||
) |
|
protected |
Create the solution equation area nceq_p and fill it.
void casacore::LSQFit::debugIt | ( | uInt & | nun, |
uInt & | np, | ||
uInt & | ncon, | ||
uInt & | ner, | ||
uInt & | rank, | ||
Double *& | nEq, | ||
Double *& | known, | ||
Double *& | constr, | ||
Double *& | er, | ||
uInt *& | piv, | ||
Double *& | sEq, | ||
Double *& | sol, | ||
Double & | prec, | ||
Double & | nonlin | ||
) | const |
Debug:
nun =
number of unknowns np =
total number of solved unknowns (nun+ncon) ncon =
number of constraint equations ner =
number of elements in chi2 vector rank =
rank) nEq =
normal equation (nun*nun as triangular matrix) known =
known vector (np) constr =
constraint matrix (ncon*nun) er =
error info vector (ner) piv =
pivot vector (np) sEq =
normal solution equation (np*np triangular) sol =
internal solution vector (np) prec =
collinearity precision nonlin =
current Levenberg factor-1 Note that all pointers may be 0.
|
protected |
De-initialise area.
|
protected |
Extend the constraint equation area to the specify number of equations.
void casacore::LSQFit::fromAipsIO | ( | AipsIO & | ) |
Bool casacore::LSQFit::fromRecord | ( | String & | error, |
const RecordInterface & | in | ||
) |
Create an LSQFit object from a record.
An error message is generated, and False returned if an invalid record is given. A valid record will return True. Error messages are postfixed to error.
Double casacore::LSQFit::getChi | ( | ) | const |
Get chi^2 (both are identical); the standard deviation (per observation) and the standard deviation per weight unit.
Referenced by casacore::GenericL2Fit< T >::chiSquare(), and getChi2().
|
inline |
Get the n-th
(from 0 to the rank deficiency, or missing rank, see e.g.
getDeficiency()
) constraint equation as determined by invert()
in SVD-mode in cEq[nUnknown]
. False returned for illegal n. Note that nMissing will be equal to the number of unknowns (nUnknowns
, or double that for the complex case) minus the rank as returned from the invert()
method.
Bool casacore::LSQFit::getCovariance | ( | std::complex< U > * | covar | ) |
Bool casacore::LSQFit::getCovariance | ( | U * | covar | ) |
Get the covariance matrix (of size nUnknowns * nUnknowns
)
Referenced by casacore::LSQaips::getCovariance().
|
inline |
Get the rank deficiency
Warning: Note that the number is returned assuming real values; For complex values it has to be halved
Definition at line 775 of file LSQFit.h.
Referenced by casacore::GenericL2Fit< T >::getRank().
Bool casacore::LSQFit::getErrors | ( | std::complex< U > * | errors | ) |
Bool casacore::LSQFit::getErrors | ( | U & | errors | ) |
Bool casacore::LSQFit::getErrors | ( | U * | errors | ) |
Get main diagonal of covariance function (of size nUnknowns
)
Referenced by casacore::LSQaips::getErrors().
Double casacore::LSQFit::getSD | ( | ) | const |
Double casacore::LSQFit::getWeightedSD | ( | ) | const |
|
protected |
|
protected |
Get work areas for solutions, covariance.
const String& casacore::LSQFit::ident | ( | ) | const |
Get identification of record.
|
protected |
Initialise areas.
Triangularize the normal equations and determine the rank nRank
of the normal equations and, in the case of an SVD
solution, the constraint equations.
The collinearity factor is used to determine if the system can be solved (in essence it is the square of the sine of the angle between a column in the normal equations and the plane suspended by the other columns: if too parallel, the equations are degenerate). If doSVD
is given as False, False is returned if rank not maximal, else an SVD
solution is done.
|
protected |
Invert rectangular matrix (i.e.
when constraints present)
|
inline |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::AsReal | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Complex | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Conjugate | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Separable | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const U & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const std::vector< std::pair< uInt, V > > & | cEq, |
const U & | weight, | ||
const U & | obs, | ||
LSQFit::Real | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::AsReal | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Complex | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Conjugate | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Separable | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const U & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
Make normal equations using the cEq
condition equation (cArray) (with nUnknowns
elements) and a weight weight
, given the known observed value obs
.
doNorm
and doKnown
can be used to e.g. re-use existing normal equations, i.e. the condition equations, but make a new known side (i.e. new observations).
The versions with cEqIndex[]
indicate which of the nUnknowns
are actually present in the condition equation (starting indexing at 0); the other terms are supposed to be zero. E.g. if a 12-telescope array has an equation only using telescopes 2 and 4, the lengths of cEqIndex
and cEq
will be both 2, and the index will contain 1 and 3 (when telescope numbering starts at 1) or 2 and 4 (when telescope numbering starts at 0. The index is given as an iterator (and hence can be a raw pointer)
The complex versions can have different interpretation of the inputs, where the complex number can be seen either as a complex number; as two real numbers, or as coefficients of equations with complex conjugates. See Note 224) for the details.
Versions with pair assume that the pairs are created by the SparseDiff automatic differentiation class. The pair is an index and a value. The indices are assumed to be sorted.
Special (makeNormSorted) indexed versions exist which assume that the given indices are sorted (which is the case for the LOFAR BBS environment).
Some versions exist with two sets of equations (cEq2, obs2). If two simultaneous equations are created they will be faster.
Note that the use of const U &
is due to a Float->Double conversion problem on Solaris. Linux was ok.
void casacore::LSQFit::makeNorm | ( | const V & | cEq, |
const U & | weight, | ||
const U & | obs, | ||
LSQFit::Real | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const std::complex< U > & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::AsReal | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Complex | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Conjugate | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const std::complex< U > & | obs, | ||
LSQFit::Separable | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const U & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const U & | obs, | ||
LSQFit::Real | , | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNorm | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const V & | cEq2, | ||
const U & | weight, | ||
const U & | obs, | ||
const U & | obs2, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNormSorted | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | weight, | ||
const U & | obs, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
void casacore::LSQFit::makeNormSorted | ( | uInt | nIndex, |
const W & | cEqIndex, | ||
const V & | cEq, | ||
const V & | cEq2, | ||
const U & | weight, | ||
const U & | obs, | ||
const U & | obs2, | ||
Bool | doNorm = True , |
||
Bool | doKnown = True |
||
) |
Merge other LSQFit
object (i.e.
the normal equation and related information) into this
. Both objects must have the same number of unknowns, and be pure normal equations (i.e. no invert(), solve(), solveLoop()
or statistics calls should have been made). If merging cannot be done, False
is returned. The index case (the index is an iterator) assumes that the normal equations to be merged are a sparse subset of the complete matrix. The index 'vector' specifies which unknowns are present. An index outside the scope of the final equations will be skipped.
Tip: For highest numerical precision in the case of a larger number of partial normal equations to be merged, it is best to merge them in pairs (and repeat);
|
protected |
Merge sparse normal equations.
Referenced by merge().
|
inline |
|
inline |
Get the infinite norm of the known vector.
Get the norm of the current solution vector.
|
inline |
Get the number of unknowns.
Definition at line 769 of file LSQFit.h.
References nun_p.
Referenced by casacore::LSQaips::getErrors(), casacore::GenericL2Fit< T >::getRank(), set(), and casacore::LSQaips::solve().
Assignment (deep copy)
Referenced by casacore::LSQaips::operator=().
const std::string& casacore::LSQFit::readyText | ( | ) | const |
void casacore::LSQFit::reset | ( | ) |
Reset status to empty.
Set new factors (collinearity factor
, and Levenberg-Marquardt LMFactor
)
|
inline |
Definition at line 729 of file LSQFit.h.
References nConstraints(), nUnknowns(), and set().
Definition at line 725 of file LSQFit.h.
References nConstraints(), nUnknowns(), and set().
Definition at line 719 of file LSQFit.h.
References nConstraints(), nUnknowns(), and set().
void casacore::LSQFit::set | ( | uInt | nUnknowns, |
const LSQComplex & | , | ||
uInt | nConstraints = 0 |
||
) |
Definition at line 722 of file LSQFit.h.
References nConstraints(), nUnknowns(), and set().
Set the expected form of the normal equations.
Definition at line 746 of file LSQFit.h.
References balanced_p.
Bool casacore::LSQFit::setConstraint | ( | uInt | n, |
const V & | cEq, | ||
const std::complex< U > & | obs | ||
) |
Bool casacore::LSQFit::setConstraint | ( | uInt | n, |
const V & | cEq, | ||
const U & | obs | ||
) |
Add a new constraint equation (updating nConstraints); or set a numbered constraint equation (0..nConstraints-1).
False if illegal number n. The constraints are equations with nUnknowns
terms, and a constant value. E.g. measuring three angles of a triangle could lead to equation [1,1,1]
with obs as 3.1415
. Note that each complex constraint will be converted into two real constraints (see Note 224).
Bool casacore::LSQFit::setConstraint | ( | uInt | n, |
uInt | nIndex, | ||
const W & | cEqIndex, | ||
const V & | cEq, | ||
const std::complex< U > & | obs | ||
) |
Bool casacore::LSQFit::setConstraint | ( | uInt | n, |
uInt | nIndex, | ||
const W & | cEqIndex, | ||
const V & | cEq, | ||
const U & | obs | ||
) |
|
inline |
|
inline |
|
inline |
void casacore::LSQFit::solve | ( | std::complex< U > * | sol | ) |
void casacore::LSQFit::solve | ( | U & | sol | ) |
void casacore::LSQFit::solve | ( | U * | sol | ) |
Solve normal equations.
The solution will be given in sol
.
Referenced by casacore::LSQaips::solve().
|
protected |
Solve normal equations.
One non-linear LM loop.
Bool casacore::LSQFit::solveLoop | ( | Double & | fit, |
uInt & | nRank, | ||
std::complex< U > * | sol, | ||
Bool | doSVD = False |
||
) |
Bool casacore::LSQFit::solveLoop | ( | Double & | fit, |
uInt & | nRank, | ||
U & | sol, | ||
Bool | doSVD = False |
||
) |
Bool casacore::LSQFit::solveLoop | ( | Double & | fit, |
uInt & | nRank, | ||
U * | sol, | ||
Bool | doSVD = False |
||
) |
Bool casacore::LSQFit::solveLoop | ( | uInt & | nRank, |
std::complex< U > * | sol, | ||
Bool | doSVD = False |
||
) |
Solve a loop in a non-linear set.
The methods with the fit
argument are deprecated. Use the combination without the 'fit' parameter, and the isReady()
call. The 'fit' parameter returns for each loop a goodness of fit indicator. If it is >0; more loops are necessary. If it is negative, and has an absolute value of say less than.001, it is probably ok, and the iterations can be stopped. Other arguments are as for solve()
and invert()
. The sol
is used for both input (parameter guess) and output.
Referenced by casacore::LSQaips::solveLoop().
|
protected |
Solve missing rank part.
Bool casacore::LSQFit::toRecord | ( | String & | error, |
RecordInterface & | out | ||
) | const |
Create a record from an LSQFit object.
The return will be False and an error message generated only if the object does not contain a valid object. Error messages are postfixed to error.
void casacore::LSQFit::uncopy | ( | Double * | beg, |
const Double * | end, | ||
U & | sol, | ||
LSQComplex | |||
) |
void casacore::LSQFit::uncopy | ( | Double * | beg, |
const Double * | end, | ||
U & | sol, | ||
LSQReal | |||
) |
void casacore::LSQFit::uncopy | ( | Double * | beg, |
const Double * | end, | ||
U * | sol, | ||
LSQComplex | |||
) |
void casacore::LSQFit::uncopy | ( | Double * | beg, |
const Double * | end, | ||
U * | sol, | ||
LSQReal | |||
) |
|
protected |
Indicator for a well balanced normal equation.
A balanced equation is one with similar values in the main diagonal.
Definition at line 892 of file LSQFit.h.
Referenced by setBalanced().
|
protected |
|
protected |
Test value for known vector in non-linear loop.
||known||inf is tested
Definition at line 889 of file LSQFit.h.
Referenced by setEpsDerivative().
|
protected |
Test value for [incremental] solution in non-linear loop.
The ||sol increment||/||sol||
is tested
Definition at line 886 of file LSQFit.h.
Referenced by setEpsValue().
|
protected |
|
staticprotected |
Definition at line 856 of file LSQFit.h.
Referenced by casacore::LSQaips::getErrors().
|
protected |
|
protected |
|
protected |
Maximum number of iterations for non-linear solution.
If a non-zero maximum number of iterations is set, the value is tested in non-linear loops
Definition at line 896 of file LSQFit.h.
Referenced by nIterations(), and setMaxIter().
|
protected |
Matrix size (will be n_p = nun_p + ncon_p)
Definition at line 873 of file LSQFit.h.
Referenced by getDeficiency(), and rowrt().
|
protected |
|
protected |
|
protected |
|
protected |
Iteration count for non-linear solution.
Definition at line 898 of file LSQFit.h.
Referenced by nIterations().
|
protected |
|
protected |
|
protected |
|
protected |
Number of unknowns.
Definition at line 869 of file LSQFit.h.
Referenced by nUnknowns(), and rowru().
|
protected |
|
protected |
Rank of normal equations (normally n_p)
Definition at line 875 of file LSQFit.h.
Referenced by getDeficiency().
|
protected |
|
staticprotected |
|
staticprotected |
Definition at line 857 of file LSQFit.h.
Referenced by casacore::LSQaips::solve(), and casacore::LSQaips::solveLoop().
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |