Documentation ¶
Index ¶
- func AsSymDense(m *mat64.Dense) (*mat64.SymDense, error)
- func DenseIdentity(n int) *mat64.Dense
- func HouseholderTransf(A *mat64.Dense, n, m int)
- func Identity(n int) *mat64.SymDense
- func IsNil(m mat64.Matrix) bool
- func NewChiSquare(kf LDKF, runs MonteCarloRuns, controls []*mat64.Vector, withNEES, withNIS bool) ([]float64, []float64, error)
- func NewHybridKF(x0 *mat64.Vector, P0 mat64.Symmetric, noise Noise, measSize int) (*HybridKF, *HybridKFEstimate, error)
- func NewInformation(i0 *mat64.Vector, I0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Information, *InformationEstimate, error)
- func NewInformationFromState(x0 *mat64.Vector, P0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Information, *InformationEstimate, error)
- func NewPurePredictorVanilla(x0 *mat64.Vector, Covar0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Vanilla, *VanillaEstimate, error)
- func NewSRIF(x0 *mat64.Vector, P0 mat64.Symmetric, measSize int, nonTriR bool, n Noise) (*SRIF, *SRIFEstimate, error)
- func NewSquareRoot(x0 *mat64.Vector, P0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*SquareRoot, *SquareRootEstimate, error)
- func NewVanilla(x0 *mat64.Vector, Covar0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Vanilla, *VanillaEstimate, error)
- func ScaledDenseIdentity(n int, s float64) *mat64.Dense
- func ScaledIdentity(n int, s float64) *mat64.SymDense
- func Sign(v float64) float64
- func VanLoan(A, Γ, W *mat64.Dense, Δt float64) (*mat64.Dense, *mat64.SymDense, error)
- type AWGN
- type BatchGroundTruth
- type BatchKF
- type BatchNoise
- type CSVExporter
- type DimensionAgreement
- type ErrorEstimate
- type Estimate
- type Exporter
- type FilterType
- type HybridKF
- func (kf *HybridKF) DisableEKF()
- func (kf *HybridKF) EKFEnabled() bool
- func (kf *HybridKF) EnableEKF()
- func (kf *HybridKF) GetNoise() Noise
- func (kf *HybridKF) Predict() (est Estimate, err error)
- func (kf *HybridKF) Prepare(Φ, Htilde *mat64.Dense)
- func (kf *HybridKF) PreparePNT(Γ *mat64.Dense)
- func (kf *HybridKF) SetNoise(n Noise)
- func (kf *HybridKF) SmoothAll(estimates []*HybridKFEstimate) (err error)
- func (kf *HybridKF) String() string
- func (kf *HybridKF) Update(realObservation, computedObservation *mat64.Vector) (est Estimate, err error)
- type HybridKFEstimate
- func (e HybridKFEstimate) Covariance() mat64.Symmetric
- func (e HybridKFEstimate) Gain() mat64.Matrix
- func (e HybridKFEstimate) Innovation() *mat64.Vector
- func (e HybridKFEstimate) IsWithin2σ() bool
- func (e HybridKFEstimate) IsWithinNσ(N float64) bool
- func (e HybridKFEstimate) Measurement() *mat64.Vector
- func (e HybridKFEstimate) ObservationDev() *mat64.Vector
- func (e HybridKFEstimate) PredCovariance() mat64.Symmetric
- func (e HybridKFEstimate) State() *mat64.Vector
- func (e HybridKFEstimate) String() string
- type Information
- func (kf *Information) GetInputControl() mat64.Matrix
- func (kf *Information) GetMeasurementMatrix() mat64.Matrix
- func (kf *Information) GetNoise() Noise
- func (kf *Information) GetStateTransition() mat64.Matrix
- func (kf *Information) Reset()
- func (kf *Information) SetInputControl(G mat64.Matrix)
- func (kf *Information) SetMeasurementMatrix(H mat64.Matrix)
- func (kf *Information) SetNoise(n Noise)
- func (kf *Information) SetStateTransition(F mat64.Matrix)
- func (kf *Information) String() string
- func (kf *Information) Update(measurement, control *mat64.Vector) (est Estimate, err error)
- type InformationEstimate
- func (e InformationEstimate) Covariance() mat64.Symmetric
- func (e InformationEstimate) Innovation() *mat64.Vector
- func (e InformationEstimate) IsWithin2σ() bool
- func (e InformationEstimate) IsWithinNσ(N float64) bool
- func (e InformationEstimate) Measurement() *mat64.Vector
- func (e InformationEstimate) PredCovariance() mat64.Symmetric
- func (e InformationEstimate) State() *mat64.Vector
- func (e InformationEstimate) String() string
- type LDKF
- type MonteCarloRun
- type MonteCarloRuns
- type NLDKF
- type Noise
- type Noiseless
- type SRIF
- func (kf *SRIF) DisableEKF()
- func (kf *SRIF) EKFEnabled() bool
- func (kf *SRIF) EnableEKF()
- func (kf *SRIF) Predict() (est Estimate, err error)
- func (kf *SRIF) Prepare(Φ, Htilde *mat64.Dense)
- func (kf *SRIF) PreparePNT(Γ *mat64.Dense)
- func (kf *SRIF) SetNoise(n Noise)
- func (kf *SRIF) SmoothAll(estimates []*SRIFEstimate) (err error)
- func (kf *SRIF) Update(realObservation, computedObservation *mat64.Vector) (est Estimate, err error)
- type SRIFEstimate
- func (e SRIFEstimate) Covariance() mat64.Symmetric
- func (e SRIFEstimate) Innovation() *mat64.Vector
- func (e SRIFEstimate) IsWithin2σ() bool
- func (e SRIFEstimate) IsWithinNσ(N float64) bool
- func (e SRIFEstimate) Measurement() *mat64.Vector
- func (e SRIFEstimate) ObservationDev() *mat64.Vector
- func (e SRIFEstimate) PredCovariance() mat64.Symmetric
- func (e SRIFEstimate) State() *mat64.Vector
- func (e SRIFEstimate) String() string
- type SquareRoot
- func (kf *SquareRoot) GetInputControl() mat64.Matrix
- func (kf *SquareRoot) GetMeasurementMatrix() mat64.Matrix
- func (kf *SquareRoot) GetNoise() Noise
- func (kf *SquareRoot) GetStateTransition() mat64.Matrix
- func (kf *SquareRoot) Reset()
- func (kf *SquareRoot) SetInputControl(G mat64.Matrix)
- func (kf *SquareRoot) SetMeasurementMatrix(H mat64.Matrix)
- func (kf *SquareRoot) SetNoise(n Noise)
- func (kf *SquareRoot) SetStateTransition(F mat64.Matrix)
- func (kf *SquareRoot) String() string
- func (kf *SquareRoot) Update(measurement, control *mat64.Vector) (est Estimate, err error)
- type SquareRootEstimate
- func (e SquareRootEstimate) Covariance() mat64.Symmetric
- func (e SquareRootEstimate) Gain() mat64.Matrix
- func (e SquareRootEstimate) Innovation() *mat64.Vector
- func (e SquareRootEstimate) IsWithin2σ() bool
- func (e SquareRootEstimate) IsWithinNσ(N float64) bool
- func (e SquareRootEstimate) Measurement() *mat64.Vector
- func (e SquareRootEstimate) PredCovariance() mat64.Symmetric
- func (e SquareRootEstimate) State() *mat64.Vector
- func (e SquareRootEstimate) String() string
- type Vanilla
- func (kf *Vanilla) GetInputControl() mat64.Matrix
- func (kf *Vanilla) GetMeasurementMatrix() mat64.Matrix
- func (kf *Vanilla) GetNoise() Noise
- func (kf *Vanilla) GetStateTransition() mat64.Matrix
- func (kf *Vanilla) Reset()
- func (kf *Vanilla) SetInputControl(G mat64.Matrix)
- func (kf *Vanilla) SetMeasurementMatrix(H mat64.Matrix)
- func (kf *Vanilla) SetNoise(n Noise)
- func (kf *Vanilla) SetStateTransition(F mat64.Matrix)
- func (kf *Vanilla) String() string
- func (kf *Vanilla) Update(measurement, control *mat64.Vector) (est Estimate, err error)
- type VanillaEstimate
- func (e VanillaEstimate) Covariance() mat64.Symmetric
- func (e VanillaEstimate) Gain() mat64.Matrix
- func (e VanillaEstimate) Innovation() *mat64.Vector
- func (e VanillaEstimate) IsWithin2σ() bool
- func (e VanillaEstimate) IsWithinNσ(N float64) bool
- func (e VanillaEstimate) Measurement() *mat64.Vector
- func (e VanillaEstimate) PredCovariance() mat64.Symmetric
- func (e VanillaEstimate) State() *mat64.Vector
- func (e VanillaEstimate) String() string
Constants ¶
This section is empty.
Variables ¶
This section is empty.
Functions ¶
func AsSymDense ¶
AsSymDense attempts return a SymDense from the provided Dense.
func DenseIdentity ¶
DenseIdentity returns an identity matrix of type Dense and of the provided size.
func HouseholderTransf ¶
HouseholderTransf performs the Householder transformation of the given A matrix. Changes are done directly in the provided matrix.
func NewChiSquare ¶
func NewChiSquare(kf LDKF, runs MonteCarloRuns, controls []*mat64.Vector, withNEES, withNIS bool) ([]float64, []float64, error)
NewChiSquare runs the Chi square tests from the MonteCarlo runs. These runs and the KF used are the ones tested via Chi square. The KF provided must be a pure predictor Vanilla KF and will be used to compute the intermediate steps of both the NEES and NIS tests. Returns NEESmeans, NISmeans and an error if applicable TODO: Change order of parameters.
func NewHybridKF ¶
func NewHybridKF(x0 *mat64.Vector, P0 mat64.Symmetric, noise Noise, measSize int) (*HybridKF, *HybridKFEstimate, error)
NewHybridKF returns a new hybrid Kalman Filter which can be used both as a CKF and EKF. Warning: there is a failsafe preventing any update prior to updating the matrices. Usage: ``` / kf.Prepare(Φ, Htilde)
estimate, err := kf.Update(realObs, computedObs)
/ ``` Parameters: - x0: initial state estimate - P0: initial covariance symmetric matrix - noise: Noise - measSize: number of rows of the measurement vector (not actually important)
func NewInformation ¶
func NewInformation(i0 *mat64.Vector, I0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Information, *InformationEstimate, error)
NewInformation returns a new Information KF. To get the next estimate, call Update() with the next measurement and the control vector. This will return a new InformationEstimate which contains everything of this step and an error if any. Parameters: - i0: initial information state (usually a zero vector) - I0: initial information matrix (usually a zero matrix) - F: state update matrix - G: control matrix (if all zeros, then control vector will not be used) - H: measurement update matrix - noise: Noise
func NewInformationFromState ¶
func NewInformationFromState(x0 *mat64.Vector, P0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Information, *InformationEstimate, error)
NewInformationFromState returns a new Information KF. To get the next estimate, call Update() with the next measurement and the control vector. This will return a new InformationEstimate which contains everything of this step and an error if any. Parameters: - i0: initial information state (usually a zero vector) - I0: initial information matrix (usually a zero matrix) - F: state update matrix - G: control matrix (if all zeros, then control vector will not be used) - H: measurement update matrix - noise: Noise
func NewPurePredictorVanilla ¶
func NewPurePredictorVanilla(x0 *mat64.Vector, Covar0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Vanilla, *VanillaEstimate, error)
NewPurePredictorVanilla returns a new Vanilla KF which only does prediction.
func NewSRIF ¶
func NewSRIF(x0 *mat64.Vector, P0 mat64.Symmetric, measSize int, nonTriR bool, n Noise) (*SRIF, *SRIFEstimate, error)
NewSRIF returns a new Square Root Information Filter. It uses the algorithms from "Statistical Orbit determination" by Tapley, Schutz & Born. Set nonTriR to `true` to NOT use the Householder transformation on \bar{R_k}.
func NewSquareRoot ¶
func NewSquareRoot(x0 *mat64.Vector, P0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*SquareRoot, *SquareRootEstimate, error)
NewSquareRoot returns a new Square Root KF. To get the next estimate, push to the MeasChan the next measurement and read from StateEst and MeasEst to get the next state estimate (\hat{x}{k+1}^{+}) and next measurement estimate (\hat{y}{k+1}). The Covar channel stores the next covariance of the system (P_{k+1}^{+}). Parameters: - x0: initial state - Covar0: initial covariance matrix - F: state update matrix - G: control matrix (if all zeros, then control vector will not be used) - H: measurement update matrix - noise: Noise
func NewVanilla ¶
func NewVanilla(x0 *mat64.Vector, Covar0 mat64.Symmetric, F, G, H mat64.Matrix, noise Noise) (*Vanilla, *VanillaEstimate, error)
NewVanilla returns a new Vanilla KF. To get the next estimate, simply push to the MeasChan the next measurement and read from StateEst and MeasEst to get the next state estimate (\hat{x}_{k+1}^{+}) and next measurement estimate (\hat{y}_{k+1}). The Covar channel stores the next covariance of the system (P_{k+1}^{+}). Parameters: - x0: initial state - Covar0: initial covariance matrix - F: state update matrix - G: control matrix (if all zeros, then control vector will not be used) - H: measurement update matrix - n: Noise
func ScaledDenseIdentity ¶
ScaledDenseIdentity returns an identity matrix time of type Dense a scaling factor of the provided size.
func ScaledIdentity ¶
ScaledIdentity returns an identity matrix time a scaling factor of the provided size.
Types ¶
type AWGN ¶
AWGN implements the Noise interface and generates an Additive white Gaussian noise.
func (*AWGN) Measurement ¶
Measurement implements the Noise interface.
func (*AWGN) MeasurementMatrix ¶
MeasurementMatrix implements the Noise interface.
func (*AWGN) ProcessMatrix ¶
ProcessMatrix implements the Noise interface.
type BatchGroundTruth ¶
type BatchGroundTruth struct {
// contains filtered or unexported fields
}
BatchGroundTruth computes the error of a given estimate from a known batch of states and measurements.
func NewBatchGroundTruth ¶
func NewBatchGroundTruth(states, measurements []*mat64.Vector) *BatchGroundTruth
NewBatchGroundTruth initializes a new batch ground truth.
func (*BatchGroundTruth) Error ¶
func (t *BatchGroundTruth) Error(k int, est Estimate) Estimate
Error returns an ErrorEstimate after comparing the provided state and measurements with the ground truths.
func (*BatchGroundTruth) ErrorWithOffset ¶
ErrorWithOffset returns an ErrorEstimate after comparing the provided state, adding the offset and measurements with the ground truths.
type BatchKF ¶
type BatchKF struct { Λ *mat64.Dense N *mat64.Vector Measurements []measurementInfo // contains filtered or unexported fields }
BatchKF defines a vanilla kalman filter. Use NewVanilla to initialize.
func NewBatchKF ¶
NewBatchKF returns a new hybrid Kalman Filter which can be used both as a CKF and EKF. Warning: there is a failsafe preventing any update prior to updating the matrices. Usage: ``` / kf.Prepare(Φ, Htilde)
estimate, err := kf.Update(realObs, computedObs)
/ ``` Parameters: - x0: initial state estimate - P0: initial covariance symmetric matrix - noise: Noise - measSize: number of rows of the measurement vector (not actually important)
func (*BatchKF) SetNextMeasurement ¶
SetNextMeasurement sets the next sequential measurement to the list of measurements to be taken into account for the filter.
type BatchNoise ¶
type BatchNoise struct {
// contains filtered or unexported fields
}
BatchNoise implements the Noise interface.
func (BatchNoise) Measurement ¶
func (n BatchNoise) Measurement(k int) *mat64.Vector
Measurement implements the Noise interface.
func (BatchNoise) MeasurementMatrix ¶
func (n BatchNoise) MeasurementMatrix() mat64.Symmetric
MeasurementMatrix implements the Noise interface.
func (BatchNoise) Process ¶
func (n BatchNoise) Process(k int) *mat64.Vector
Process implements the Noise interface.
func (BatchNoise) ProcessMatrix ¶
func (n BatchNoise) ProcessMatrix() mat64.Symmetric
ProcessMatrix implements the Noise interface.
func (BatchNoise) String ¶
func (n BatchNoise) String() string
String implements the Stringer interface.
type CSVExporter ¶
type CSVExporter struct {
// contains filtered or unexported fields
}
CSVExporter returns a new CSV exporter.
func NewCSVExporter ¶
func NewCSVExporter(headers []string, filepath, filename string) (e *CSVExporter, err error)
NewCSVExporter initializes a new CSV export.
func NewCustomCSVExporter ¶
func NewCustomCSVExporter(headers []string, filepath, filename string, covarBound float64) (e *CSVExporter, err error)
NewCustomCSVExporter initializes a new CSV export. NOTE: Prefix header with `_` to state that it does not relate to the estimation (and therefore won't have the added covariance bounds).
func (CSVExporter) Write ¶
func (e CSVExporter) Write(est Estimate) error
Write writes the estimate to the CSV file.
func (CSVExporter) WriteRaw ¶
func (e CSVExporter) WriteRaw(s string) error
WriteRaw writes a raw line to the CSV file.
func (CSVExporter) WriteRawLn ¶
func (e CSVExporter) WriteRawLn(s string) error
WriteRawLn writes a raw line to the CSV file.
type DimensionAgreement ¶
type DimensionAgreement uint8
DimensionAgreement defines how two matrices' dimensions should agree.
type ErrorEstimate ¶
type ErrorEstimate struct {
VanillaEstimate // This is effectively the same as a VanillaEstimate, so no change.
}
ErrorEstimate implements the Estimate interface and is used to show the error of an estimate.
type Estimate ¶
type Estimate interface { IsWithinNσ(N float64) bool // IsWithinNσ returns whether the estimation is within the N*σ bounds. State() *mat64.Vector // Returns \hat{x}_{k+1}^{+} Measurement() *mat64.Vector // Returns \hat{y}_{k}^{+} Innovation() *mat64.Vector // Returns y_{k} - H*\hat{x}_{k+1}^{-} Covariance() mat64.Symmetric // Return P_{k+1}^{+} PredCovariance() mat64.Symmetric // Return P_{k+1}^{-} String() string // Must implement the stringer interface. }
Estimate is returned from Update() in any KF. This allows to avoid some computations in other filters, e.g. in the Information filter.
type FilterType ¶
type FilterType uint8
FilterType allows for quick comparison of filters.
const ( // CKFType definition would be a tautology CKFType FilterType = iota + 1 // EKFType definition would be a tautology EKFType // UKFType definition would be a tautology UKFType // SRIFType definition would be a tautology SRIFType )
func (FilterType) String ¶
func (f FilterType) String() string
type HybridKF ¶
type HybridKF struct {
Φ, Htilde, Γ *mat64.Dense
Noise Noise
// contains filtered or unexported fields
}
HybridKF defines a hybrid kalman filter for non-linear dynamical systems. Use NewHybridKF to initialize.
func (*HybridKF) DisableEKF ¶
func (kf *HybridKF) DisableEKF()
DisableEKF switches this back to a CKF mode.
func (*HybridKF) EKFEnabled ¶
EKFEnabled returns whether the KF is in EKF mode.
func (*HybridKF) EnableEKF ¶
func (kf *HybridKF) EnableEKF()
EnableEKF switches this to an EKF mode.
func (*HybridKF) Predict ¶
Predict computes only the time update (or prediction). Will return an error if the KF is locked (call Prepare to unlock).
func (*HybridKF) PreparePNT ¶
PreparePNT prepares the process noise transition matrix and enabled the SNC for the next update. WARNING: If not called, the SNC *will not* be included.
func (*HybridKF) SmoothAll ¶
func (kf *HybridKF) SmoothAll(estimates []*HybridKFEstimate) (err error)
SmoothAll will smooth all the previous estimates using the provided data. Returns the smoothed estimates. Will return an error if there are more estimates than there should be. WARNING: overwrites the provided array of estimates.
type HybridKFEstimate ¶
type HybridKFEstimate struct {
Φ, Γ *mat64.Dense // Used for smoothing
Δobs *mat64.Vector
// contains filtered or unexported fields
}
HybridKFEstimate is the output of each update state of the Vanilla KF. It implements the Estimate interface.
func (HybridKFEstimate) Covariance ¶
func (e HybridKFEstimate) Covariance() mat64.Symmetric
Covariance implements the Estimate interface.
func (HybridKFEstimate) Gain ¶
func (e HybridKFEstimate) Gain() mat64.Matrix
Gain the Estimate interface.
func (HybridKFEstimate) Innovation ¶
func (e HybridKFEstimate) Innovation() *mat64.Vector
Innovation implements the Estimate interface.
func (HybridKFEstimate) IsWithin2σ ¶
func (e HybridKFEstimate) IsWithin2σ() bool
IsWithin2σ returns whether the estimation is within the 2σ bounds.
func (HybridKFEstimate) IsWithinNσ ¶
func (e HybridKFEstimate) IsWithinNσ(N float64) bool
IsWithinNσ returns whether the estimation is within the 2σ bounds.
func (HybridKFEstimate) Measurement ¶
func (e HybridKFEstimate) Measurement() *mat64.Vector
Measurement implements the Estimate interface.
func (HybridKFEstimate) ObservationDev ¶
func (e HybridKFEstimate) ObservationDev() *mat64.Vector
ObservationDev returns the observation deviation.
func (HybridKFEstimate) PredCovariance ¶
func (e HybridKFEstimate) PredCovariance() mat64.Symmetric
PredCovariance implements the Estimate interface.
func (HybridKFEstimate) State ¶
func (e HybridKFEstimate) State() *mat64.Vector
State implements the Estimate interface.
func (HybridKFEstimate) String ¶
func (e HybridKFEstimate) String() string
type Information ¶
type Information struct { Finv mat64.Matrix G mat64.Matrix H mat64.Matrix Qinv mat64.Matrix Rinv mat64.Matrix Noise Noise // contains filtered or unexported fields }
Information defines a vanilla kalman filter. Use NewVanilla to initialize.
func (*Information) GetInputControl ¶
func (kf *Information) GetInputControl() mat64.Matrix
GetInputControl returns the G matrix.
func (*Information) GetMeasurementMatrix ¶
func (kf *Information) GetMeasurementMatrix() mat64.Matrix
GetMeasurementMatrix returns the H matrix.
func (*Information) GetNoise ¶
func (kf *Information) GetNoise() Noise
GetNoise updates the F matrix.
func (*Information) GetStateTransition ¶
func (kf *Information) GetStateTransition() mat64.Matrix
GetStateTransition returns the F matrix. *WARNING:* Returns the *INVERSE* of F for the information filter.
func (*Information) Reset ¶
func (kf *Information) Reset()
Reset reinitializes the KF with its initial estimate.
func (*Information) SetInputControl ¶
func (kf *Information) SetInputControl(G mat64.Matrix)
SetInputControl updates the G matrix.
func (*Information) SetMeasurementMatrix ¶
func (kf *Information) SetMeasurementMatrix(H mat64.Matrix)
SetMeasurementMatrix updates the H matrix.
func (*Information) SetStateTransition ¶
func (kf *Information) SetStateTransition(F mat64.Matrix)
SetStateTransition updates the F matrix.
func (*Information) String ¶
func (kf *Information) String() string
type InformationEstimate ¶
type InformationEstimate struct {
// contains filtered or unexported fields
}
InformationEstimate is the output of each update state of the Information KF. It implements the Estimate interface.
func NewInformationEstimate ¶
func NewInformationEstimate(infoState, meas *mat64.Vector, infoMat, predInfoMat mat64.Symmetric) InformationEstimate
NewInformationEstimate initializes a new InformationEstimate.
func (InformationEstimate) Covariance ¶
func (e InformationEstimate) Covariance() mat64.Symmetric
Covariance implements the Estimate interface. *NOTE:* With the IF, one cannot view the covariance matrix until there is enough information.
func (InformationEstimate) Innovation ¶
func (e InformationEstimate) Innovation() *mat64.Vector
Innovation implements the Estimate interface.
func (InformationEstimate) IsWithin2σ ¶
func (e InformationEstimate) IsWithin2σ() bool
IsWithin2σ returns whether the estimation is within the 2σ bounds.
func (InformationEstimate) IsWithinNσ ¶
func (e InformationEstimate) IsWithinNσ(N float64) bool
IsWithinNσ returns whether the estimation is within the 2σ bounds.
func (InformationEstimate) Measurement ¶
func (e InformationEstimate) Measurement() *mat64.Vector
Measurement implements the Estimate interface.
func (InformationEstimate) PredCovariance ¶
func (e InformationEstimate) PredCovariance() mat64.Symmetric
PredCovariance implements the Estimate interface. *NOTE:* With the IF, one cannot view the prediction covariance matrix until there is enough information.
func (InformationEstimate) State ¶
func (e InformationEstimate) State() *mat64.Vector
State implements the Estimate interface.
func (InformationEstimate) String ¶
func (e InformationEstimate) String() string
type LDKF ¶
type LDKF interface { Update(measurement, control *mat64.Vector) (Estimate, error) GetNoise() Noise GetStateTransition() mat64.Matrix GetInputControl() mat64.Matrix GetMeasurementMatrix() mat64.Matrix SetStateTransition(mat64.Matrix) SetInputControl(mat64.Matrix) SetMeasurementMatrix(mat64.Matrix) SetNoise(Noise) Reset() String() string }
LDKF defines a linear dynamics Kalman Filter.
type MonteCarloRun ¶
type MonteCarloRun struct {
Estimates []Estimate
}
MonteCarloRun stores the results of an MC run.
type MonteCarloRuns ¶
type MonteCarloRuns struct { Runs []MonteCarloRun // contains filtered or unexported fields }
MonteCarloRuns stores MC runs.
func NewMonteCarloRuns ¶
func NewMonteCarloRuns(samples, steps, rowsH int, controls []*mat64.Vector, kf *Vanilla) MonteCarloRuns
NewMonteCarloRuns run monte carlos on the provided filter.
func (MonteCarloRuns) AsCSV ¶
func (mc MonteCarloRuns) AsCSV(headers []string) []string
AsCSV is used as a CSV serializer. Does not include the header.
func (MonteCarloRuns) Mean ¶
func (mc MonteCarloRuns) Mean(step int) (mean []float64)
Mean returns the mean of all the samples for the given time step.
func (MonteCarloRuns) StdDev ¶
func (mc MonteCarloRuns) StdDev(step int) (mean []float64)
StdDev returns the standard deviation of all the samples for the given time step.
type NLDKF ¶
type NLDKF interface { Prepare(Φ, Htilde *mat64.Dense) Predict() (est Estimate, err error) Update(realObservation, computedObservation *mat64.Vector) (est Estimate, err error) EKFEnabled() bool EnableEKF() DisableEKF() PreparePNT(Γ *mat64.Dense) SetNoise(n Noise) }
NLDKF defines a non-linear dynamics Kalman Filter. Operates and is architectured slightly differently than LDKF.
type Noise ¶
type Noise interface { Process(k int) *mat64.Vector // Returns the process noise w at step k Measurement(k int) *mat64.Vector // Returns the measurement noise w at step k ProcessMatrix() mat64.Symmetric // Returns the process noise matrix Q MeasurementMatrix() mat64.Symmetric // Returns the measurement noise matrix R Reset() // Reinitializes the noise String() string // Stringer interface implementation }
Noise allows to handle the noise for a KF.
type Noiseless ¶
Noiseless is noiseless and implements the Noise interface.
func NewNoiseless ¶
NewNoiseless creates new AWGN noise from the provided Q and R.
func (Noiseless) Measurement ¶
Measurement returns a vector of the correct size.
func (Noiseless) MeasurementMatrix ¶
MeasurementMatrix implements the Noise interface.
func (Noiseless) ProcessMatrix ¶
ProcessMatrix implements the Noise interface.
type SRIF ¶
SRIF defines a square root information filter for non-linear dynamical systems. Use NewSquareRootInformation to initialize.
func (*SRIF) Predict ¶
Predict computes only the time update (or prediction). Will return an error if the KF is locked (call Prepare to unlock).
func (*SRIF) SmoothAll ¶
func (kf *SRIF) SmoothAll(estimates []*SRIFEstimate) (err error)
SmoothAll will smooth all the previous estimates using the provided data. Returns the smoothed estimates. Will return an error if there are more estimates than there should be. WARNING: overwrites the provided array of estimates.
type SRIFEstimate ¶
type SRIFEstimate struct { Φ *mat64.Dense // Used for smoothing Δobs *mat64.Vector R *mat64.Dense // contains filtered or unexported fields }
SRIFEstimate is the output of each update state of the Vanilla KF. It implements the Estimate interface.
func NewSRIFEstimate ¶
func NewSRIFEstimate(Φ *mat64.Dense, sqinfoState, meas, Δobs *mat64.Vector, R0, predR0 *mat64.Dense) SRIFEstimate
NewSRIFEstimate initializes a new SRIFEstimate. NOTE: R0 and predR0 are mat64.Dense for simplicity of implementation, but they should be symmetric.
func (SRIFEstimate) Covariance ¶
func (e SRIFEstimate) Covariance() mat64.Symmetric
Covariance implements the Estimate interface.
func (SRIFEstimate) Innovation ¶
func (e SRIFEstimate) Innovation() *mat64.Vector
Innovation implements the Estimate interface.
func (SRIFEstimate) IsWithin2σ ¶
func (e SRIFEstimate) IsWithin2σ() bool
IsWithin2σ returns whether the estimation is within the 2σ bounds.
func (SRIFEstimate) IsWithinNσ ¶
func (e SRIFEstimate) IsWithinNσ(N float64) bool
IsWithinNσ returns whether the estimation is within the 2σ bounds.
func (SRIFEstimate) Measurement ¶
func (e SRIFEstimate) Measurement() *mat64.Vector
Measurement implements the Estimate interface.
func (SRIFEstimate) ObservationDev ¶
func (e SRIFEstimate) ObservationDev() *mat64.Vector
ObservationDev returns the observation deviation.
func (SRIFEstimate) PredCovariance ¶
func (e SRIFEstimate) PredCovariance() mat64.Symmetric
PredCovariance implements the Estimate interface.
func (SRIFEstimate) State ¶
func (e SRIFEstimate) State() *mat64.Vector
State implements the Estimate interface.
func (SRIFEstimate) String ¶
func (e SRIFEstimate) String() string
type SquareRoot ¶
type SquareRoot struct { F mat64.Matrix G mat64.Matrix H mat64.Matrix Noise Noise // contains filtered or unexported fields }
SquareRoot defines a square root kalman filter. Use NewSqrt to initialize.
func (*SquareRoot) GetInputControl ¶
func (kf *SquareRoot) GetInputControl() mat64.Matrix
GetInputControl returns the G matrix.
func (*SquareRoot) GetMeasurementMatrix ¶
func (kf *SquareRoot) GetMeasurementMatrix() mat64.Matrix
GetMeasurementMatrix returns the H matrix.
func (*SquareRoot) GetStateTransition ¶
func (kf *SquareRoot) GetStateTransition() mat64.Matrix
GetStateTransition returns the F matrix.
func (*SquareRoot) Reset ¶
func (kf *SquareRoot) Reset()
Reset reinitializes the KF with its initial estimate.
func (*SquareRoot) SetInputControl ¶
func (kf *SquareRoot) SetInputControl(G mat64.Matrix)
SetInputControl updates the G matrix.
func (*SquareRoot) SetMeasurementMatrix ¶
func (kf *SquareRoot) SetMeasurementMatrix(H mat64.Matrix)
SetMeasurementMatrix updates the H matrix.
func (*SquareRoot) SetStateTransition ¶
func (kf *SquareRoot) SetStateTransition(F mat64.Matrix)
SetStateTransition updates the F matrix.
type SquareRootEstimate ¶
type SquareRootEstimate struct {
// contains filtered or unexported fields
}
SquareRootEstimate is the output of each update state of the SquareRoot KF. It implements the Estimate interface.
func NewSqrtEstimate ¶
func NewSqrtEstimate(state, meas, innovation *mat64.Vector, stddev, predStddev, gain *mat64.Dense) SquareRootEstimate
NewSqrtEstimate initializes a new InformationEstimate.
func (SquareRootEstimate) Covariance ¶
func (e SquareRootEstimate) Covariance() mat64.Symmetric
Covariance implements the Estimate interface.
func (SquareRootEstimate) Gain ¶
func (e SquareRootEstimate) Gain() mat64.Matrix
Gain the Estimate interface.
func (SquareRootEstimate) Innovation ¶
func (e SquareRootEstimate) Innovation() *mat64.Vector
Innovation implements the Estimate interface.
func (SquareRootEstimate) IsWithin2σ ¶
func (e SquareRootEstimate) IsWithin2σ() bool
IsWithin2σ returns whether the estimation is within the 2σ bounds.
func (SquareRootEstimate) IsWithinNσ ¶
func (e SquareRootEstimate) IsWithinNσ(N float64) bool
IsWithinNσ returns whether the estimation is within the 2σ bounds.
func (SquareRootEstimate) Measurement ¶
func (e SquareRootEstimate) Measurement() *mat64.Vector
Measurement implements the Estimate interface.
func (SquareRootEstimate) PredCovariance ¶
func (e SquareRootEstimate) PredCovariance() mat64.Symmetric
PredCovariance implements the Estimate interface.
func (SquareRootEstimate) State ¶
func (e SquareRootEstimate) State() *mat64.Vector
State implements the Estimate interface.
func (SquareRootEstimate) String ¶
func (e SquareRootEstimate) String() string
type Vanilla ¶
type Vanilla struct { F mat64.Matrix G mat64.Matrix H mat64.Matrix Noise Noise // contains filtered or unexported fields }
Vanilla defines a vanilla kalman filter. Use NewVanilla to initialize.
func (*Vanilla) GetInputControl ¶
GetInputControl returns the G matrix.
func (*Vanilla) GetMeasurementMatrix ¶
GetMeasurementMatrix returns the H matrix.
func (*Vanilla) GetStateTransition ¶
GetStateTransition returns the F matrix.
func (*Vanilla) Reset ¶
func (kf *Vanilla) Reset()
Reset reinitializes the KF with its initial estimate.
func (*Vanilla) SetInputControl ¶
SetInputControl updates the F matrix.
func (*Vanilla) SetMeasurementMatrix ¶
SetMeasurementMatrix updates the H matrix.
func (*Vanilla) SetStateTransition ¶
SetStateTransition updates the F matrix.
type VanillaEstimate ¶
type VanillaEstimate struct {
// contains filtered or unexported fields
}
VanillaEstimate is the output of each update state of the Vanilla KF. It implements the Estimate interface.
func (VanillaEstimate) Covariance ¶
func (e VanillaEstimate) Covariance() mat64.Symmetric
Covariance implements the Estimate interface.
func (VanillaEstimate) Gain ¶
func (e VanillaEstimate) Gain() mat64.Matrix
Gain the Estimate interface.
func (VanillaEstimate) Innovation ¶
func (e VanillaEstimate) Innovation() *mat64.Vector
Innovation implements the Estimate interface.
func (VanillaEstimate) IsWithin2σ ¶
func (e VanillaEstimate) IsWithin2σ() bool
IsWithin2σ returns whether the estimation is within the 2σ bounds.
func (VanillaEstimate) IsWithinNσ ¶
func (e VanillaEstimate) IsWithinNσ(N float64) bool
IsWithinNσ returns whether the estimation is within the 2σ bounds.
func (VanillaEstimate) Measurement ¶
func (e VanillaEstimate) Measurement() *mat64.Vector
Measurement implements the Estimate interface.
func (VanillaEstimate) PredCovariance ¶
func (e VanillaEstimate) PredCovariance() mat64.Symmetric
PredCovariance implements the Estimate interface.
func (VanillaEstimate) State ¶
func (e VanillaEstimate) State() *mat64.Vector
State implements the Estimate interface.
func (VanillaEstimate) String ¶
func (e VanillaEstimate) String() string