Loading...
Searching...
No Matches
Fang2017filter_wDF.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Class
25 Fang2017filter_wDF
26Description
27 A class that implements the Ensamble Kalman Filter in the version by Fang et al. 2017
28SourceFiles
29 Fang2017filter_wDF.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef fang2017filter_wDF_H
38#define fang2017filter_wDF_H
39
40#include "ITHACAutilities.H"
41#include "muq2ithaca.H"
42#include "ensembleClass.H"
43
44namespace ITHACAmuq
45{
46//--------------------------------------------------------------------------
49class Fang2017filter_wDF
50{
51 int Nsamples = 0;
52
53 int stateSize = 0;
54 int parameterSize = 0;
55 int observationSize = 0;
56
57 int timeStepI = 0;
58 double startTime;
59 double deltaTime;
60 double endTime;
61 int Ntimes;
62
63 bool initialStateFlag = 0;
64 bool parameterPriorFlag = 0;
65 bool modelErrorFlag = 0;
66 bool measurementNoiseFlag = 0;
67
69 Eigen::VectorXd initialStateMean;
70
72 Eigen::MatrixXd initialStateCov;
73
75 Eigen::VectorXd parameterPriorMean;
76
78 Eigen::MatrixXd parameterPriorCov;
79
81 Eigen::MatrixXd parameterMean;
82
84 Eigen::MatrixXd stateMean;
85
87 Eigen::VectorXd timeVector;
88
90 int observationStart;
91
93 int observationDelta;
94
96 Eigen::VectorXi observationBoolVec;
97
99 Eigen::MatrixXd observations;
100
102 std::shared_ptr<muq::Modeling::Gaussian> initialStateDensity;
103
105 std::shared_ptr<muq::Modeling::Gaussian> parameterPriorDensity;
106
107
108 int innerLoopI = 0;
109
112 bool univariateInitStateDensFlag = 0;
113
114
115 public:
116 // Constructors
117 Fang2017filter_wDF();
118
119 Fang2017filter_wDF(int _Nsamples);
120
121 virtual ~Fang2017filter_wDF() = default;
122
124 std::shared_ptr<muq::Modeling::Gaussian> modelErrorDensity;
125
127 std::shared_ptr<muq::Modeling::Gaussian> measNoiseDensity;
128
131
134
137
140
142 Eigen::MatrixXd state_maxConf;
143
145 Eigen::MatrixXd state_minConf;
146
148 Eigen::MatrixXd parameter_maxConf;
149
151 Eigen::MatrixXd parameter_minConf;
152
153 // Functions
154
155 //--------------------------------------------------------------------------
157 int getNumberOfSamples();
158
159 //--------------------------------------------------------------------------
161 double getTime();
162
163 //--------------------------------------------------------------------------
165 double getTime(int _timeStepI);
166
167 //--------------------------------------------------------------------------
169 int getTimeStep();
170
171 //--------------------------------------------------------------------------
173 Eigen::VectorXd getTimeVector();
174
175 //--------------------------------------------------------------------------
177 Eigen::MatrixXd getStateMean();
178
179 //--------------------------------------------------------------------------
181 Eigen::MatrixXd getParameterMean();
182
183 //--------------------------------------------------------------------------
185 Eigen::MatrixXd getParameterMaxConf(); // Kabir
186
187 //--------------------------------------------------------------------------
189 Eigen::MatrixXd getParameterMinConf(); // Kabir
190
191 //--------------------------------------------------------------------------
193 void setObservations(Eigen::MatrixXd _observations);
194
195 //--------------------------------------------------------------------------
197 void setTime(double _startTime, double _deltaTime, double _endTime);
198
199 //--------------------------------------------------------------------------
201 void setObservationTime(int _observationStart, int _observationDelta);
202
203 //--------------------------------------------------------------------------
205 void setModelError(double cov, bool univariate = 0);
206
207 //--------------------------------------------------------------------------
209 void setMeasNoise(double cov);
210
211 //--------------------------------------------------------------------------
213 void setInitialStateDensity(Eigen::VectorXd _mean, Eigen::MatrixXd _cov,
214 bool _univariateFlag = 0);
215
216 //--------------------------------------------------------------------------
218 void sampleInitialState();
219
220 //--------------------------------------------------------------------------
222 void setParameterPriorDensity(Eigen::VectorXd _mean, Eigen::MatrixXd _cov);
223
224 //--------------------------------------------------------------------------
226 void sampleParameterDist();
227
228 //--------------------------------------------------------------------------
230 Eigen::MatrixXd ensembleFromDensity(
231 std::shared_ptr<muq::Modeling::Gaussian> _density);
232
233 //--------------------------------------------------------------------------
235 void buildJointEns();
236
237 //--------------------------------------------------------------------------
239 virtual void stateProjection() = 0;
240
241 //--------------------------------------------------------------------------
243 virtual void observeState() = 0;
244
245 //--------------------------------------------------------------------------
247 void setObservationSize(int _size);
248
249 //--------------------------------------------------------------------------
251 void setStateSize(int _size);
252
253 //--------------------------------------------------------------------------
255 int getStateSize();
256
257 //--------------------------------------------------------------------------
259 void setParameterSize(int _size);
260
261 //--------------------------------------------------------------------------
263 int getParameterSize();
264
265 //--------------------------------------------------------------------------
267 void updateJointEns(Eigen::VectorXd _observation);
268
269 //--------------------------------------------------------------------------
271 void run(int innerLoopMax, word outputFolder);
272};
273}
274#endif
Header file of the ITHACAutilities namespace.
Eigen::MatrixXd parameter_minConf
parameter minimum confidence level
void setObservations(Eigen::MatrixXd _observations)
Set the observations matrix.
Eigen::MatrixXd parameter_maxConf
parameter maximum confidence level
void setMeasNoise(double cov)
Setup of the measurement noise distribution.
void run(int innerLoopMax, word outputFolder)
Run the filtering.
ensemble stateEns
State ensemble.
Eigen::MatrixXd getStateMean()
Return state mean.
Eigen::MatrixXd getParameterMaxConf()
Return parameter Max Confidence.
void setModelError(double cov, bool univariate=0)
Setup of the model error distribution.
void sampleInitialState()
Create initial state ensemble.
void setObservationTime(int _observationStart, int _observationDelta)
Setup the observation vector.
Eigen::MatrixXd getParameterMinConf()
Return parameter Min Confidence.
std::shared_ptr< muq::Modeling::Gaussian > measNoiseDensity
Measurement noise density.
ensemble jointEns
Joint state and parameter ensemble.
int getTimeStep()
Return timestep.
std::shared_ptr< muq::Modeling::Gaussian > modelErrorDensity
Model noise density.
void setParameterPriorDensity(Eigen::VectorXd _mean, Eigen::MatrixXd _cov)
Create parameter ensemble.
Eigen::MatrixXd ensembleFromDensity(std::shared_ptr< muq::Modeling::Gaussian > _density)
General class to sample from an input density.
void buildJointEns()
Concatenate state and parameter ensambles to create the joint ensamble.
ensemble parameterEns
Parameter ensemble.
Eigen::MatrixXd state_minConf
State minimum confidence level.
int getNumberOfSamples()
Return number of samples per ensamble.
void setInitialStateDensity(Eigen::VectorXd _mean, Eigen::MatrixXd _cov, bool _univariateFlag=0)
Create initial state ensemble.
void sampleParameterDist()
Create parameter ensemble.
Eigen::MatrixXd getParameterMean()
Return parameter mean.
ensemble observationEns
Observation ensemble.
void setTime(double _startTime, double _deltaTime, double _endTime)
Setup the time vector.
Eigen::MatrixXd state_maxConf
State maximum confidence level.
Eigen::VectorXd getTimeVector()
Return time vector.
Class for ensembles, each column is a sample.
Header file of the muq2ithaca namespace.