n2p2 - A neural network potential package
nnp::KalmanFilter Class Reference

Implementation of the Kalman filter method. More...

#include <KalmanFilter.h>

Inheritance diagram for nnp::KalmanFilter:
Collaboration diagram for nnp::KalmanFilter:

Public Types

enum  KalmanType { KT_STANDARD , KT_FADINGMEMORY }
 Enumerate different Kalman filter types. More...
 

Public Member Functions

 KalmanFilter (std::size_t const sizeState, KalmanType const type)
 Kalman filter class constructor. More...
 
virtual ~KalmanFilter ()
 Destructor. More...
 
void setSizeObservation (std::size_t const sizeObservation)
 Set observation vector size. More...
 
void setState (double *state)
 Set pointer to current state. More...
 
void setError (double const *const error)
 Set pointer to current error vector. More...
 
void setError (double const *const error, std::size_t const size)
 Set pointer to current error vector. More...
 
void setJacobian (double const *const jacobian)
 Set pointer to current Jacobi matrix. More...
 
void setJacobian (double const *const jacobian, std::size_t const columns)
 Set pointer to current Jacobi matrix. More...
 
void update ()
 Update error covariance matrix and state vector. More...
 
void update (std::size_t const sizeObservation)
 Update error covariance matrix and state vector. More...
 
void setParametersStandard (double const epsilon, double const q0, double const qtau, double const qmin, double const eta0, double const etatau, double const etamax)
 Set parameters for standard Kalman filter. More...
 
void setParametersFadingMemory (double const epsilon, double const q0, double const qtau, double const qmin, double const lambda, double const nu)
 Set parameters for fading memory Kalman filter. More...
 
std::string status (std::size_t epoch) const
 Status report. More...
 
std::vector< std::string > statusHeader () const
 Header for status report file. More...
 
std::vector< std::string > info () const
 Information about Kalman filter settings. More...
 
KalmanType getType () const
 Getter for type. More...
 
std::size_t getSizeObservation () const
 Getter for sizeObservation. More...
 
std::size_t getNumUpdates () const
 Getter for numUpdates. More...
 
double getEta () const
 Getter for eta. More...
 
double getEpsilon () const
 Getter for epsilon. More...
 
double getQ0 () const
 Getter for q0. More...
 
double getQtau () const
 Getter for qtau. More...
 
double getQmin () const
 Getter for qmin. More...
 
double getLambda () const
 Getter for lambda. More...
 
double getNu () const
 Getter for nu. More...
 
double getGamma () const
 Getter for gamma. More...
 
- Public Member Functions inherited from nnp::Updater
virtual void setState (double *state)=0
 Set pointer to current state. More...
 
virtual void setError (double const *const error, std::size_t const size=1)=0
 Set pointer to current error vector. More...
 
virtual void setJacobian (double const *const jacobian, std::size_t const columns=1)=0
 Set pointer to current Jacobi matrix. More...
 
virtual void update ()=0
 Perform single update of state vector. More...
 
virtual std::string status (std::size_t epoch) const =0
 Status report. More...
 
virtual std::vector< std::string > statusHeader () const =0
 Header for status report file. More...
 
virtual std::vector< std::string > info () const =0
 Information about this updater. More...
 
virtual void setupTiming (std::string const &prefix="upd")
 Activate detailed timing. More...
 
virtual void resetTimingLoop ()
 Start a new timing loop (e.g. More...
 
virtual std::map< std::string, StopwatchgetTiming () const
 Return timings gathered in stopwatch map. More...
 

Private Attributes

KalmanType type
 Kalman filter type. More...
 
std::size_t sizeObservation
 Size of observation (measurement) vector. More...
 
std::size_t numUpdates
 Total number of updates performed. More...
 
double epsilon
 Error covariance initialization parameter \(\epsilon\). More...
 
double q
 Process noise \(q\). More...
 
double q0
 Process noise initial value \(q_0\). More...
 
double qtau
 Process noise exponential decay parameter \(q_{\tau}\). More...
 
double qmin
 Process noise minimum value \(q_{\text{min}}\). More...
 
double eta
 Learning rate \(\eta\). More...
 
double eta0
 Learning rate initial value \(\eta_0\). More...
 
double etatau
 Learning rate exponential increase parameter \(\eta_{\tau}\). More...
 
double etamax
 Learning rate maximum value \(\eta_{\text{max}}\). More...
 
double lambda
 Forgetting factor for fading memory Kalman filter. More...
 
double nu
 Parameter for fading memory Kalman filter. More...
 
double gamma
 Forgetting gain factor gamma for fading memory Kalman filter. More...
 
Eigen::Map< Eigen::VectorXd > * w
 State vector. More...
 
Eigen::Map< Eigen::VectorXd const > * xi
 Error vector. More...
 
Eigen::Map< Eigen::MatrixXd const > * H
 Derivative matrix. More...
 
Eigen::MatrixXd P
 Error covariance matrix. More...
 
Eigen::MatrixXd K
 Kalman gain matrix. More...
 
Eigen::MatrixXd X
 Intermediate result X = P . H. More...
 

Additional Inherited Members

- Protected Member Functions inherited from nnp::Updater
 Updater (std::size_t const sizeState)
 Constructor. More...
 
- Protected Attributes inherited from nnp::Updater
bool timing
 Whether detailed timing is enabled. More...
 
bool timingReset
 Internal loop timer reset switch. More...
 
std::size_t sizeState
 Number of neural network connections (weights + biases). More...
 
std::string prefix
 Prefix for timing stopwatches. More...
 
std::map< std::string, Stopwatchsw
 Stopwatch map for timing. More...
 

Detailed Description

Implementation of the Kalman filter method.

Definition at line 31 of file KalmanFilter.h.

Member Enumeration Documentation

◆ KalmanType

Enumerate different Kalman filter types.

Enumerator
KT_STANDARD 

Regular Kalman filter.

KT_FADINGMEMORY 

Kalman filtering with fading memory modification.

Definition at line 35 of file KalmanFilter.h.

36 {
41 };
@ KT_STANDARD
Regular Kalman filter.
Definition: KalmanFilter.h:38
@ KT_FADINGMEMORY
Kalman filtering with fading memory modification.
Definition: KalmanFilter.h:40

Constructor & Destructor Documentation

◆ KalmanFilter()

KalmanFilter::KalmanFilter ( std::size_t const  sizeState,
KalmanType const  type 
)

Kalman filter class constructor.

Parameters
[in]sizeStateSize of the state vector w.
[in]typeKalman filter type (see KalmanType).

Definition at line 27 of file KalmanFilter.cpp.

28 :
31 numUpdates (0 ),
32 epsilon (0.0 ),
33 q (0.0 ),
34 q0 (0.0 ),
35 qtau (0.0 ),
36 qmin (0.0 ),
37 eta (0.0 ),
38 eta0 (0.0 ),
39 etatau (0.0 ),
40 etamax (0.0 ),
41 lambda (0.0 ),
42 nu (0.0 ),
43 gamma (0.0 ),
44 w (NULL),
45 xi (NULL),
46 H (NULL)
47{
48 if (!(type == KT_STANDARD ||
50 {
51 throw runtime_error("ERROR: Unknown Kalman filter type.\n");
52 }
53
54 if (sizeState < 1)
55 {
56 throw runtime_error("ERROR: Wrong Kalman filter dimensions.\n");
57 }
58
59 this->type = type;
61
62 w = new Map<VectorXd >(0, sizeState);
63 xi = new Map<VectorXd const>(0, sizeObservation);
64 H = new Map<MatrixXd const>(0, sizeState, sizeObservation);
65 P.resize(sizeState, sizeState);
66 P.setIdentity();
67 // Prevent problems with unallocated K when log starts.
69 K.setZero();
70}
double eta0
Learning rate initial value .
Definition: KalmanFilter.h:214
Eigen::Map< Eigen::MatrixXd const > * H
Derivative matrix.
Definition: KalmanFilter.h:230
double nu
Parameter for fading memory Kalman filter.
Definition: KalmanFilter.h:222
double qtau
Process noise exponential decay parameter .
Definition: KalmanFilter.h:208
std::size_t sizeObservation
Size of observation (measurement) vector.
Definition: KalmanFilter.h:198
KalmanType type
Kalman filter type.
Definition: KalmanFilter.h:196
Eigen::Map< Eigen::VectorXd > * w
State vector.
Definition: KalmanFilter.h:226
double q0
Process noise initial value .
Definition: KalmanFilter.h:206
double lambda
Forgetting factor for fading memory Kalman filter.
Definition: KalmanFilter.h:220
double etatau
Learning rate exponential increase parameter .
Definition: KalmanFilter.h:216
double etamax
Learning rate maximum value .
Definition: KalmanFilter.h:218
double gamma
Forgetting gain factor gamma for fading memory Kalman filter.
Definition: KalmanFilter.h:224
std::size_t numUpdates
Total number of updates performed.
Definition: KalmanFilter.h:200
double q
Process noise .
Definition: KalmanFilter.h:204
double qmin
Process noise minimum value .
Definition: KalmanFilter.h:210
Eigen::MatrixXd K
Kalman gain matrix.
Definition: KalmanFilter.h:234
double eta
Learning rate .
Definition: KalmanFilter.h:212
double epsilon
Error covariance initialization parameter .
Definition: KalmanFilter.h:202
Eigen::Map< Eigen::VectorXd const > * xi
Error vector.
Definition: KalmanFilter.h:228
Eigen::MatrixXd P
Error covariance matrix.
Definition: KalmanFilter.h:232
Updater(std::size_t const sizeState)
Constructor.
Definition: Updater.cpp:22
std::size_t sizeState
Number of neural network connections (weights + biases).
Definition: Updater.h:110

References H, K, KT_FADINGMEMORY, KT_STANDARD, P, sizeObservation, nnp::Updater::sizeState, type, w, and xi.

◆ ~KalmanFilter()

KalmanFilter::~KalmanFilter ( )
virtual

Destructor.

Definition at line 72 of file KalmanFilter.cpp.

73{
74}

Member Function Documentation

◆ setSizeObservation()

void KalmanFilter::setSizeObservation ( std::size_t const  sizeObservation)

Set observation vector size.

Parameters
[in]sizeSize of the observation vector.

If the size of the observation vector is known in advance it can be set here.

Definition at line 76 of file KalmanFilter.cpp.

77{
78 sizeObservation = size;
79
80 return;
81}

References sizeObservation.

Referenced by nnp::Training::update().

Here is the caller graph for this function:

◆ setState()

void KalmanFilter::setState ( double *  state)
virtual

Set pointer to current state.

Parameters
[in,out]statePointer to state vector (weights vector), will be changed in-place upon calling update().

Implements nnp::Updater.

Definition at line 83 of file KalmanFilter.cpp.

84{
85 new (w) Map<VectorXd>(state, sizeState);
86
87 return;
88}

References nnp::Updater::sizeState, and w.

◆ setError() [1/2]

void KalmanFilter::setError ( double const *const  error)

Set pointer to current error vector.

Parameters
[in]errorPointer to error (difference between reference and neural network potential output).

Definition at line 90 of file KalmanFilter.cpp.

91{
93
94 return;
95}
void setError(double const *const error)
Set pointer to current error vector.

References setError(), and sizeObservation.

Referenced by setError().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setError() [2/2]

void nnp::KalmanFilter::setError ( double const *const  error,
std::size_t const  size 
)
virtual

Set pointer to current error vector.

Parameters
[in]errorPointer to error (difference between reference and neural network potential output).
[in]sizeNumber of error vector entries.

Implements nnp::Updater.

◆ setJacobian() [1/2]

void KalmanFilter::setJacobian ( double const *const  jacobian)

Set pointer to current Jacobi matrix.

Parameters
[in]jacobianDerivatives of error with respect to weights.

Definition at line 104 of file KalmanFilter.cpp.

105{
106 setJacobian(jacobian, sizeObservation);
107
108 return;
109}
void setJacobian(double const *const jacobian)
Set pointer to current Jacobi matrix.

References setJacobian(), and sizeObservation.

Referenced by setJacobian().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ setJacobian() [2/2]

void nnp::KalmanFilter::setJacobian ( double const *const  jacobian,
std::size_t const  columns 
)
virtual

Set pointer to current Jacobi matrix.

Parameters
[in]jacobianDerivatives of error with respect to weights.
[in]columnsNumber of gradients provided.
Note
If there are \(m\) errors and \(n\) weights, the Jacobi matrix is a \(n \times m\) matrix stored in column-major order.

Implements nnp::Updater.

◆ update() [1/2]

void KalmanFilter::update ( )
virtual

Update error covariance matrix and state vector.

Implements nnp::Updater.

Definition at line 119 of file KalmanFilter.cpp.

120{
122
123 return;
124}
void update()
Update error covariance matrix and state vector.

References sizeObservation, and update().

Referenced by update().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ update() [2/2]

void nnp::KalmanFilter::update ( std::size_t const  sizeObservation)

Update error covariance matrix and state vector.

Parameters
[in]sizeObservationSize of the observation vector.

◆ setParametersStandard()

void KalmanFilter::setParametersStandard ( double const  epsilon,
double const  q0,
double const  qtau,
double const  qmin,
double const  eta0,
double const  etatau,
double const  etamax 
)

Set parameters for standard Kalman filter.

Parameters
[in]epsilonError covariance initialization parameter \(\epsilon\).
[in]q0Process noise initial value \(q_0\).
[in]qtauProcess noise exponential decay parameter \(q_{\tau}\).
[in]qminProcess noise minimum value \(q_{\text{min}}\).
[in]eta0Initial learning rate \(\eta_0\).
[in]etatauLearning rate exponential increase parameter \(\eta_{\tau}\).
[in]etamaxLearning rate maximum value \(\eta_{\text{max}}\).

Definition at line 199 of file KalmanFilter.cpp.

206{
207 this->epsilon = epsilon;
208 this->q0 = q0 ;
209 this->qtau = qtau ;
210 this->qmin = qmin ;
211 this->eta0 = eta0 ;
212 this->etatau = etatau ;
213 this->etamax = etamax ;
214
215 q = q0;
216 eta = eta0;
217 P /= epsilon;
218
219 return;
220}

References epsilon, eta, eta0, etamax, etatau, P, q, q0, qmin, and qtau.

Referenced by nnp::Training::setupTraining().

Here is the caller graph for this function:

◆ setParametersFadingMemory()

void KalmanFilter::setParametersFadingMemory ( double const  epsilon,
double const  q0,
double const  qtau,
double const  qmin,
double const  lambda,
double const  nu 
)

Set parameters for fading memory Kalman filter.

Parameters
[in]epsilonError covariance initialization parameter \(\epsilon\).
[in]q0Process noise initial value \(q_0\).
[in]qtauProcess noise exponential decay parameter \(q_{\tau}\).
[in]qminProcess noise minimum value \(q_{\text{min}}\).
[in]lambdaForgetting factor \(\lambda\).
[in]nuFading memory parameter \(\nu\).

Definition at line 222 of file KalmanFilter.cpp.

228{
229 this->epsilon = epsilon;
230 this->q0 = q0 ;
231 this->qtau = qtau ;
232 this->qmin = qmin ;
233 this->lambda = lambda ;
234 this->nu = nu ;
235
236 q = q0;
237 P /= epsilon;
238 gamma = 1.0;
239
240 return;
241}

References epsilon, gamma, lambda, nu, P, q, q0, qmin, and qtau.

Referenced by nnp::Training::setupTraining().

Here is the caller graph for this function:

◆ status()

string KalmanFilter::status ( std::size_t  epoch) const
virtual

Status report.

Parameters
[in]epochCurrent epoch.
Returns
Line with current status information.

Implements nnp::Updater.

Definition at line 243 of file KalmanFilter.cpp.

244{
245
246 double Pasym = 0.5 * (P - P.transpose()).array().abs().mean();
247 double Pdiag = P.diagonal().array().abs().sum();
248 double Poffdiag = (P.array().abs().sum() - Pdiag)
249 / (sizeState * (sizeState - 1));
250 Pdiag /= sizeState;
251 double Kmean = K.array().abs().mean();
252
253 string s = strpr("%10zu %10zu %16.8E %16.8E %16.8E %16.8E %16.8E",
254 epoch, numUpdates, Pdiag, Poffdiag, Pasym, Kmean, q);
255 if (type == KT_STANDARD)
256 {
257 s += strpr(" %16.8E", eta);
258 }
259 else if (type == KT_FADINGMEMORY)
260 {
261 s += strpr(" %16.8E %16.8E", lambda, numUpdates * gamma);
262 }
263 s += '\n';
264
265 return s;
266}
string strpr(const char *format,...)
String version of printf function.
Definition: utility.cpp:90

References eta, gamma, K, KT_FADINGMEMORY, KT_STANDARD, lambda, numUpdates, P, q, nnp::Updater::sizeState, nnp::strpr(), and type.

Here is the call graph for this function:

◆ statusHeader()

vector< string > KalmanFilter::statusHeader ( ) const
virtual

Header for status report file.

Returns
Vector with header lines.

Implements nnp::Updater.

Definition at line 268 of file KalmanFilter.cpp.

269{
270 vector<string> header;
271
272 vector<string> title;
273 vector<string> colName;
274 vector<string> colInfo;
275 vector<size_t> colSize;
276 title.push_back("Kalman filter status report.");
277 colSize.push_back(10);
278 colName.push_back("epoch");
279 colInfo.push_back("Training epoch.");
280 colSize.push_back(10);
281 colName.push_back("nupdates");
282 colInfo.push_back("Number of updates performed.");
283 colSize.push_back(16);
284 colName.push_back("Pdiag");
285 colInfo.push_back("Mean of absolute diagonal values of error covariance "
286 "matrix P.");
287 colSize.push_back(16);
288 colName.push_back("Poffdiag");
289 colInfo.push_back("Mean of absolute off-diagonal values of error "
290 "covariance matrix P.");
291 colSize.push_back(16);
292 colName.push_back("Pasym");
293 colInfo.push_back("Asymmetry of P, i.e. mean of absolute values of "
294 "asymmetric part 0.5*(P - P^T).");
295 colSize.push_back(16);
296 colName.push_back("Kmean");
297 colInfo.push_back("Mean of abolute compontents of Kalman gain matrix K.");
298 colSize.push_back(16);
299 colName.push_back("q");
300 colInfo.push_back("Magnitude of process noise (= diagonal entries of Q).");
301 if (type == KT_STANDARD)
302 {
303 colSize.push_back(16);
304 colName.push_back("eta");
305 colInfo.push_back("Learning rate.");
306 }
307 else if (type == KT_FADINGMEMORY)
308 {
309 colSize.push_back(16);
310 colName.push_back("lambda");
311 colInfo.push_back("Forgetting factor for fading memory KF.");
312 colSize.push_back(16);
313 colName.push_back("kgamma");
314 colInfo.push_back("Forgetting gain k * gamma(k).");
315 }
316 header = createFileHeader(title, colSize, colName, colInfo);
317
318 return header;
319}
vector< string > createFileHeader(vector< string > const &title, vector< size_t > const &colSize, vector< string > const &colName, vector< string > const &colInfo, char const &commentChar)
Definition: utility.cpp:111

References nnp::createFileHeader(), KT_FADINGMEMORY, KT_STANDARD, and type.

Here is the call graph for this function:

◆ info()

vector< string > KalmanFilter::info ( ) const
virtual

Information about Kalman filter settings.

Returns
Vector with info lines.

Implements nnp::Updater.

Definition at line 321 of file KalmanFilter.cpp.

322{
323 vector<string> v;
324
325 if (type == KT_STANDARD)
326 {
327 v.push_back(strpr("KalmanType::KT_STANDARD (%d)\n", type));
328 v.push_back(strpr("sizeState = %zu\n", sizeState));
329 v.push_back(strpr("sizeObservation = %zu\n", sizeObservation));
330 v.push_back(strpr("epsilon = %12.4E\n", epsilon));
331 v.push_back(strpr("q0 = %12.4E\n", q0 ));
332 v.push_back(strpr("qtau = %12.4E\n", qtau ));
333 v.push_back(strpr("qmin = %12.4E\n", qmin ));
334 v.push_back(strpr("eta0 = %12.4E\n", eta0 ));
335 v.push_back(strpr("etatau = %12.4E\n", etatau ));
336 v.push_back(strpr("etamax = %12.4E\n", etamax ));
337 }
338 else if (type == KT_FADINGMEMORY)
339 {
340 v.push_back(strpr("KalmanType::KT_FADINGMEMORY (%d)\n", type));
341 v.push_back(strpr("sizeState = %zu\n", sizeState));
342 v.push_back(strpr("sizeObservation = %zu\n", sizeObservation));
343 v.push_back(strpr("epsilon = %12.4E\n", epsilon));
344 v.push_back(strpr("q0 = %12.4E\n", q0 ));
345 v.push_back(strpr("qtau = %12.4E\n", qtau ));
346 v.push_back(strpr("qmin = %12.4E\n", qmin ));
347 v.push_back(strpr("lambda = %12.4E\n", lambda));
348 v.push_back(strpr("nu = %12.4E\n", nu ));
349 }
350 v.push_back(strpr("OpenMP threads used: %d\n", nbThreads()));
351
352 return v;
353}

References epsilon, eta0, etamax, etatau, KT_FADINGMEMORY, KT_STANDARD, lambda, nu, q0, qmin, qtau, sizeObservation, nnp::Updater::sizeState, nnp::strpr(), and type.

Here is the call graph for this function:

◆ getType()

KalmanFilter::KalmanType nnp::KalmanFilter::getType ( ) const
inline

Getter for type.

Definition at line 243 of file KalmanFilter.h.

244{
245 return type;
246}

References type.

◆ getSizeObservation()

std::size_t nnp::KalmanFilter::getSizeObservation ( ) const
inline

Getter for sizeObservation.

Definition at line 248 of file KalmanFilter.h.

249{
250 return sizeObservation;
251}

References sizeObservation.

◆ getNumUpdates()

std::size_t nnp::KalmanFilter::getNumUpdates ( ) const
inline

Getter for numUpdates.

Definition at line 253 of file KalmanFilter.h.

254{
255 return numUpdates;
256}

References numUpdates.

◆ getEta()

double nnp::KalmanFilter::getEta ( ) const
inline

Getter for eta.

Definition at line 258 of file KalmanFilter.h.

259{
260 return eta;
261}

References eta.

◆ getEpsilon()

double nnp::KalmanFilter::getEpsilon ( ) const
inline

Getter for epsilon.

Definition at line 263 of file KalmanFilter.h.

264{
265 return epsilon;
266}

References epsilon.

◆ getQ0()

double nnp::KalmanFilter::getQ0 ( ) const
inline

Getter for q0.

Definition at line 268 of file KalmanFilter.h.

269{
270 return q0;
271}

References q0.

◆ getQtau()

double nnp::KalmanFilter::getQtau ( ) const
inline

Getter for qtau.

Definition at line 273 of file KalmanFilter.h.

274{
275 return qtau;
276}

References qtau.

◆ getQmin()

double nnp::KalmanFilter::getQmin ( ) const
inline

Getter for qmin.

Definition at line 278 of file KalmanFilter.h.

279{
280 return qmin;
281}

References qmin.

◆ getLambda()

double nnp::KalmanFilter::getLambda ( ) const
inline

Getter for lambda.

Definition at line 283 of file KalmanFilter.h.

284{
285 return lambda;
286}

References lambda.

◆ getNu()

double nnp::KalmanFilter::getNu ( ) const
inline

Getter for nu.

Definition at line 288 of file KalmanFilter.h.

289{
290 return nu;
291}

References nu.

◆ getGamma()

double nnp::KalmanFilter::getGamma ( ) const

Getter for gamma.

Member Data Documentation

◆ type

KalmanType nnp::KalmanFilter::type
private

Kalman filter type.

Definition at line 196 of file KalmanFilter.h.

Referenced by getType(), info(), KalmanFilter(), status(), and statusHeader().

◆ sizeObservation

std::size_t nnp::KalmanFilter::sizeObservation
private

Size of observation (measurement) vector.

Definition at line 198 of file KalmanFilter.h.

Referenced by getSizeObservation(), info(), KalmanFilter(), setError(), setJacobian(), setSizeObservation(), and update().

◆ numUpdates

std::size_t nnp::KalmanFilter::numUpdates
private

Total number of updates performed.

Definition at line 200 of file KalmanFilter.h.

Referenced by getNumUpdates(), and status().

◆ epsilon

double nnp::KalmanFilter::epsilon
private

Error covariance initialization parameter \(\epsilon\).

Definition at line 202 of file KalmanFilter.h.

Referenced by getEpsilon(), info(), setParametersFadingMemory(), and setParametersStandard().

◆ q

double nnp::KalmanFilter::q
private

Process noise \(q\).

Definition at line 204 of file KalmanFilter.h.

Referenced by setParametersFadingMemory(), setParametersStandard(), and status().

◆ q0

double nnp::KalmanFilter::q0
private

Process noise initial value \(q_0\).

Definition at line 206 of file KalmanFilter.h.

Referenced by getQ0(), info(), setParametersFadingMemory(), and setParametersStandard().

◆ qtau

double nnp::KalmanFilter::qtau
private

Process noise exponential decay parameter \(q_{\tau}\).

Definition at line 208 of file KalmanFilter.h.

Referenced by getQtau(), info(), setParametersFadingMemory(), and setParametersStandard().

◆ qmin

double nnp::KalmanFilter::qmin
private

Process noise minimum value \(q_{\text{min}}\).

Definition at line 210 of file KalmanFilter.h.

Referenced by getQmin(), info(), setParametersFadingMemory(), and setParametersStandard().

◆ eta

double nnp::KalmanFilter::eta
private

Learning rate \(\eta\).

Definition at line 212 of file KalmanFilter.h.

Referenced by getEta(), setParametersStandard(), and status().

◆ eta0

double nnp::KalmanFilter::eta0
private

Learning rate initial value \(\eta_0\).

Definition at line 214 of file KalmanFilter.h.

Referenced by info(), and setParametersStandard().

◆ etatau

double nnp::KalmanFilter::etatau
private

Learning rate exponential increase parameter \(\eta_{\tau}\).

Definition at line 216 of file KalmanFilter.h.

Referenced by info(), and setParametersStandard().

◆ etamax

double nnp::KalmanFilter::etamax
private

Learning rate maximum value \(\eta_{\text{max}}\).

Definition at line 218 of file KalmanFilter.h.

Referenced by info(), and setParametersStandard().

◆ lambda

double nnp::KalmanFilter::lambda
private

Forgetting factor for fading memory Kalman filter.

Definition at line 220 of file KalmanFilter.h.

Referenced by getLambda(), info(), setParametersFadingMemory(), and status().

◆ nu

double nnp::KalmanFilter::nu
private

Parameter for fading memory Kalman filter.

Definition at line 222 of file KalmanFilter.h.

Referenced by getNu(), info(), and setParametersFadingMemory().

◆ gamma

double nnp::KalmanFilter::gamma
private

Forgetting gain factor gamma for fading memory Kalman filter.

Definition at line 224 of file KalmanFilter.h.

Referenced by setParametersFadingMemory(), and status().

◆ w

Eigen::Map<Eigen::VectorXd>* nnp::KalmanFilter::w
private

State vector.

Definition at line 226 of file KalmanFilter.h.

Referenced by KalmanFilter(), and setState().

◆ xi

Eigen::Map<Eigen::VectorXd const>* nnp::KalmanFilter::xi
private

Error vector.

Definition at line 228 of file KalmanFilter.h.

Referenced by KalmanFilter().

◆ H

Eigen::Map<Eigen::MatrixXd const>* nnp::KalmanFilter::H
private

Derivative matrix.

Definition at line 230 of file KalmanFilter.h.

Referenced by KalmanFilter().

◆ P

Eigen::MatrixXd nnp::KalmanFilter::P
private

Error covariance matrix.

Definition at line 232 of file KalmanFilter.h.

Referenced by KalmanFilter(), setParametersFadingMemory(), setParametersStandard(), and status().

◆ K

Eigen::MatrixXd nnp::KalmanFilter::K
private

Kalman gain matrix.

Definition at line 234 of file KalmanFilter.h.

Referenced by KalmanFilter(), and status().

◆ X

Eigen::MatrixXd nnp::KalmanFilter::X
private

Intermediate result X = P . H.

Definition at line 236 of file KalmanFilter.h.


The documentation for this class was generated from the following files: