Loading...
Searching...
No Matches
muq2ithaca.C
Go to the documentation of this file.
1#include "muq2ithaca.H"
2
3namespace ITHACAmuq
4{
5namespace muq2ithaca
6{
7Eigen::MatrixXd EnsembleKalmanFilter(Eigen::MatrixXd prior,
8 Eigen::VectorXd measurements,
9 Eigen::MatrixXd measurementsCov,
10 Eigen::MatrixXd observedState)
11{
12 M_Assert(measurements.rows() == observedState.rows(),
13 "The observed state should have the same dimention of the measurements");
14 M_Assert(observedState.cols() == prior.cols(),
15 "The input matrices should all have the samples on the columns");
16 M_Assert(measurementsCov.rows() == measurements.rows()
17 && measurementsCov.cols() == measurements.rows(),
18 "Wrong measurements covariance matrix");
19 unsigned Nseeds = prior.cols();
20 unsigned stateDim = prior.rows();
21 unsigned measDim = measurements.size();
22 Eigen::MatrixXd HA(measDim, Nseeds);
23 Eigen::VectorXd observedStateMean = observedState.rowwise().mean();
24 Eigen::MatrixXd A(stateDim, Nseeds);
25 Eigen::VectorXd priorMean = prior.rowwise().mean();
26 auto measNoise = std::make_shared<muq::Modeling::Gaussian>
27 (Eigen::VectorXd::Zero(measDim), measurementsCov);
28 Eigen::MatrixXd D(measDim, Nseeds);
29
30 for (unsigned i = 0; i < Nseeds; i++)
31 {
32 A.col(i) = prior.col(i) - priorMean;
33 HA.col(i) = observedState.col(i) - observedStateMean;
34 D.col(i) = measurements + measNoise->Sample();
35 }
36
37 Eigen::MatrixXd Y = D -
38 observedState; //diff measurement data and simulated data
39 Eigen::MatrixXd P = HA * HA.transpose() / (Nseeds - 1.) + measurementsCov;
40 Eigen::FullPivLU<Eigen::MatrixXd> lu_decomp(P);
41 auto P_rank = lu_decomp.rank();
42
43 if (P_rank < P.cols() || P_rank < P.rows())
44 {
45 std::cout << "Pseudo inverse of P should be implemented in the EnKF, exiting" <<
46 std::endl;
47 std::cout << P << std::endl;
48 exit(10);
49 }
50 else
51 {
52 P = P.inverse();
53 }
54
55 //KF update
56 Eigen::MatrixXd M = P * Y;
57 Eigen::MatrixXd Z = (1. / (Nseeds - 1.)) * A * HA.transpose();
58 //Kalman Gain
59 Eigen::MatrixXd K = (1. / (Nseeds - 1.)) * A * HA.transpose() * P;
60
61 if ((K.array() < 0.0).any())
62 {
63 std::cout << "WARNING: Kalman Gain is negative.\nK = \n" << K << std::endl <<
64 std::endl;
65 }
66
67 return prior + Z * M;
68}
69
70Eigen::MatrixXd EnsembleKalmanFilter(PtrList<volScalarField>& prior,
71 Eigen::VectorXd measurements,
72 Eigen::MatrixXd measurementsCov,
73 Eigen::MatrixXd observedState)
74{
75 M_Assert(measurements.rows() == observedState.rows(),
76 "The observed state should have the same dimention of the measurements");
77 M_Assert(observedState.cols() == prior.size(),
78 "The input matrices should all have the samples on the columns");
79 M_Assert(measurementsCov.rows() == measurements.rows()
80 && measurementsCov.cols() == measurements.rows(),
81 "Wrong measurements covariance matrix");
82 unsigned Nseeds = prior.size();
83 unsigned stateDim = prior[0].size();
84 unsigned measDim = measurements.size();
85 Eigen::MatrixXd priorMatrix = Foam2Eigen::PtrList2Eigen(prior);
86 Eigen::MatrixXd HA(measDim, Nseeds);
87 Eigen::VectorXd observedStateMean = observedState.rowwise().mean();
88 Eigen::MatrixXd A(stateDim, Nseeds);
89 Eigen::VectorXd priorMean = priorMatrix.rowwise().mean();
90 auto measNoise = std::make_shared<muq::Modeling::Gaussian>
91 (Eigen::VectorXd::Zero(measDim), measurementsCov);
92 Eigen::MatrixXd D(measDim, Nseeds);
93
94 for (unsigned i = 0; i < Nseeds; i++)
95 {
96 A.col(i) = priorMatrix.col(i) - priorMean;
97 HA.col(i) = observedState.col(i) - observedStateMean;
98 D.col(i) = measurements + measNoise->Sample();
99 }
100
101 Eigen::MatrixXd Y = D -
102 observedState; //diff measurement data and simulated data
103 Eigen::MatrixXd P = HA * HA.transpose() / (Nseeds - 1.) + measurementsCov;
104 Eigen::FullPivLU<Eigen::MatrixXd> lu_decomp(P);
105 auto P_rank = lu_decomp.rank();
106
107 if (P_rank < P.cols() || P_rank < P.rows())
108 {
109 std::cout << "Pseudo inverse of P should be implemented in the EnKF, exiting" <<
110 std::endl;
111 std::cout << P << std::endl;
112 exit(10);
113 }
114 else
115 {
116 P = P.inverse();
117 }
118
119 //KF update
120 Eigen::MatrixXd M = P * Y;
121 Eigen::MatrixXd Z = (1. / (Nseeds - 1.)) * A * HA.transpose();
122 //Kalman Gain
123 Eigen::MatrixXd K = (1. / (Nseeds - 1.)) * A * HA.transpose() * P;
124
125 if ((K.array() < 0.0).any())
126 {
127 std::cout << "WARNING: Kalman Gain is negative.\nK = \n" << K << std::endl <<
128 std::endl;
129 }
130
131 return priorMatrix + Z * M;
132}
133
134double quantile(Eigen::VectorXd samps, double p, int method)
135{
136 double m;
137 int n = samps.size();
138 int j;
139 double gamma = 0.0;
140 std::sort(samps.data(), samps.data() + samps.size());
141
142 if (method == 1)
143 {
144 m = 0;
145 j = std::floor(p * n + m);
146 double g = p * n + m - j * 1.0;
147
148 if (g > 0)
149 {
150 gamma = 1.0;
151 }
152 }
153 else if (method == 2)
154 {
155 gamma = 0.5;
156 m = 0;
157 j = std::floor(p * n + m);
158 double g = p * n + m - j * 1.0;
159
160 if (g > 0)
161 {
162 gamma = 1.0;
163 }
164 }
165 else if (method == 3)
166 {
167 gamma = 1.0;
168 m = -0.5;
169 j = std::floor(p * n + m);
170 double g = p * n + m - j * 1.0;
171
172 if (g == 0 && j % 2 == 0)
173 {
174 gamma = 0.0;
175 }
176 }
177 else
178 {
179 std::cout << "Quantile method not implemented. Exiting." << std::endl;
180 exit(10);
181 }
182
183 if (j >= samps.size() - 1)
184 {
185 j = samps.size() - 2;
186 }
187
188 return (1 - gamma) * samps(j) + gamma * samps(j + 1);
189}
190
191Eigen::VectorXd quantile(Eigen::MatrixXd samps, double p, int method)
192{
193 Eigen::VectorXd output(samps.rows());
194
195 for (int i = 0; i < samps.rows(); i++)
196 {
197 Eigen::VectorXd sampsRow = samps.row(i);
198 output(i) = quantile(sampsRow, p, method);
199 }
200
201 return output;
202}
203
204}
205}
#define M_Assert(Expr, Msg)
static Eigen::MatrixXd PtrList2Eigen(PtrList< GeometricField< Type, PatchField, GeoMesh > > &fields, label Nfields=-1)
Convert a PtrList of snapshots to Eigen matrix (only internal field)
Definition Foam2Eigen.C:649
volScalarField & A
volScalarField & D
Header file of the muq2ithaca namespace.
Eigen::MatrixXd EnsembleKalmanFilter(Eigen::MatrixXd prior, Eigen::VectorXd measurements, Eigen::MatrixXd measurementsCov, Eigen::MatrixXd observedState)
Ensemble Kalman Filter.
Definition muq2ithaca.C:7
double quantile(Eigen::VectorXd samps, double p, int method)
Returns quantile for a vector of samples.
Definition muq2ithaca.C:134
dimensionedVector & g
volScalarField & p
label i
Definition pEqn.H:46