n2p2 - A neural network potential package
KalmanFilter.h
Go to the documentation of this file.
1// n2p2 - A neural network potential package
2// Copyright (C) 2018 Andreas Singraber (University of Vienna)
3//
4// This program is free software: you can redistribute it and/or modify
5// it under the terms of the GNU General Public License as published by
6// the Free Software Foundation, either version 3 of the License, or
7// (at your option) any later version.
8//
9// This program is distributed in the hope that it will be useful,
10// but WITHOUT ANY WARRANTY; without even the implied warranty of
11// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12// GNU General Public License for more details.
13//
14// You should have received a copy of the GNU General Public License
15// along with this program. If not, see <https://www.gnu.org/licenses/>.
16
17#ifndef KALMANFILTER_H
18#define KALMANFILTER_H
19
20#include <mpi.h>
21#include "Updater.h"
22#include <Eigen/Core>
23#include <cstddef>
24#include <string>
25#include <vector>
26
27namespace nnp
28{
29
31class KalmanFilter : public Updater
32{
33public:
36 {
41 };
42
48 KalmanFilter(std::size_t const sizeState,
49 KalmanType const type);
52 virtual ~KalmanFilter();
61 std::size_t const sizeObservation);
67 void setState(double* state);
73 void setError(double const* const error);
80 void setError(double const* const error,
81 std::size_t const size);
86 void setJacobian(double const* const jacobian);
96 void setJacobian(double const* const jacobian,
97 std::size_t const columns);
100 void update();
105 void update(std::size_t const sizeObservation);
119 void setParametersStandard(double const epsilon,
120 double const q0,
121 double const qtau,
122 double const qmin,
123 double const eta0,
124 double const etatau,
125 double const etamax);
137 void setParametersFadingMemory(double const epsilon,
138 double const q0,
139 double const qtau,
140 double const qmin,
141 double const lambda,
142 double const nu);
149 std::string status(std::size_t epoch) const;
154 std::vector<std::string> statusHeader() const;
159 std::vector<std::string> info() const;
162 KalmanType getType() const;
165 std::size_t getSizeObservation() const;
168 std::size_t getNumUpdates() const;
171 double getEta() const;
174 double getEpsilon() const;
177 double getQ0() const;
180 double getQtau() const;
183 double getQmin() const;
186 double getLambda() const;
189 double getNu() const;
192 double getGamma() const;
193
194private:
198 std::size_t sizeObservation;
200 std::size_t numUpdates;
202 double epsilon;
204 double q;
206 double q0;
208 double qtau;
210 double qmin;
212 double eta;
214 double eta0;
216 double etatau;
218 double etamax;
220 double lambda;
222 double nu;
224 double gamma;
226 Eigen::Map<Eigen::VectorXd>* w;
228 Eigen::Map<Eigen::VectorXd const>* xi;
230 Eigen::Map<Eigen::MatrixXd const>* H;
232 Eigen::MatrixXd P;
234 Eigen::MatrixXd K;
236 Eigen::MatrixXd X;
237};
238
240// Inlined function definitions //
242
244{
245 return type;
246}
247
248inline std::size_t KalmanFilter::getSizeObservation() const
249{
250 return sizeObservation;
251}
252
253inline std::size_t KalmanFilter::getNumUpdates() const
254{
255 return numUpdates;
256}
257
258inline double KalmanFilter::getEta() const
259{
260 return eta;
261}
262
263inline double KalmanFilter::getEpsilon() const
264{
265 return epsilon;
266}
267
268inline double KalmanFilter::getQ0() const
269{
270 return q0;
271}
272
273inline double KalmanFilter::getQtau() const
274{
275 return qtau;
276}
277
278inline double KalmanFilter::getQmin() const
279{
280 return qmin;
281}
282
283inline double KalmanFilter::getLambda() const
284{
285 return lambda;
286}
287
288inline double KalmanFilter::getNu() const
289{
290 return nu;
291}
292
293}
294
295#endif
Implementation of the Kalman filter method.
Definition: KalmanFilter.h:32
void update()
Update error covariance matrix and state vector.
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.
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::string status(std::size_t epoch) const
Status report.
std::size_t sizeObservation
Size of observation (measurement) vector.
Definition: KalmanFilter.h:198
double getQ0() const
Getter for q0.
Definition: KalmanFilter.h:268
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
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.
double getEpsilon() const
Getter for epsilon.
Definition: KalmanFilter.h:263
Eigen::MatrixXd X
Intermediate result X = P . H.
Definition: KalmanFilter.h:236
void setState(double *state)
Set pointer to current state.
void setJacobian(double const *const jacobian, std::size_t const columns)
Set pointer to current Jacobi matrix.
std::vector< std::string > statusHeader() const
Header for status report file.
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::vector< std::string > info() const
Information about Kalman filter settings.
void setError(double const *const error)
Set pointer to current error vector.
std::size_t numUpdates
Total number of updates performed.
Definition: KalmanFilter.h:200
KalmanType
Enumerate different Kalman filter types.
Definition: KalmanFilter.h:36
@ KT_STANDARD
Regular Kalman filter.
Definition: KalmanFilter.h:38
@ KT_FADINGMEMORY
Kalman filtering with fading memory modification.
Definition: KalmanFilter.h:40
void setError(double const *const error, std::size_t const size)
Set pointer to current error vector.
double q
Process noise .
Definition: KalmanFilter.h:204
std::size_t getNumUpdates() const
Getter for numUpdates.
Definition: KalmanFilter.h:253
KalmanFilter(std::size_t const sizeState, KalmanType const type)
Kalman filter class constructor.
double qmin
Process noise minimum value .
Definition: KalmanFilter.h:210
void setSizeObservation(std::size_t const sizeObservation)
Set observation vector size.
double getEta() const
Getter for eta.
Definition: KalmanFilter.h:258
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
double getNu() const
Getter for nu.
Definition: KalmanFilter.h:288
KalmanType getType() const
Getter for type.
Definition: KalmanFilter.h:243
std::size_t getSizeObservation() const
Getter for sizeObservation.
Definition: KalmanFilter.h:248
double getQtau() const
Getter for qtau.
Definition: KalmanFilter.h:273
Eigen::Map< Eigen::VectorXd const > * xi
Error vector.
Definition: KalmanFilter.h:228
double getGamma() const
Getter for gamma.
double getLambda() const
Getter for lambda.
Definition: KalmanFilter.h:283
Eigen::MatrixXd P
Error covariance matrix.
Definition: KalmanFilter.h:232
double getQmin() const
Getter for qmin.
Definition: KalmanFilter.h:278
virtual ~KalmanFilter()
Destructor.
void update(std::size_t const sizeObservation)
Update error covariance matrix and state vector.
void setJacobian(double const *const jacobian)
Set pointer to current Jacobi matrix.
Base class for different weight update methods.
Definition: Updater.h:32
std::size_t sizeState
Number of neural network connections (weights + biases).
Definition: Updater.h:110
Definition: Atom.h:28