Loading...
Searching...
No Matches
02enKF_1DinverseHeatTransfer.H
Go to the documentation of this file.
1#ifndef EnKF_1DinverseHeatTransfer_H
2#define EnKF_1DinverseHeatTransfer_H
3#include "simpleControl.H"
4#include "fvOptions.H"
5
6
11{
12 public:
13 explicit EnKF_1DinverseHeatTransfer(int argc, char* argv[], int Nseeds_)
14 :
15 laplacianProblem(argc, argv)
16 {
17 fvMesh& mesh = _mesh();
18 _simple = autoPtr<simpleControl>
19 (
20 new simpleControl
21 (
22 mesh
23 )
24 );
25 simpleControl& simple = _simple();
26 Time& runTime = _runTime();
27#include "createFvOptions.H"
28 left_ind = mesh.boundaryMesh().findPatchID("left");
29 Nseeds = Nseeds_;
30 Tensemble = PtrList<volScalarField>(Nseeds);
31 volScalarField& T(_T());
32 stateSize = T.size();
33 priorSamples = Eigen::MatrixXd::Zero(stateSize, Nseeds);
35 startTime = runTime.startTime().value();
36 deltaTime = runTime.deltaTValue();
37 endTime = runTime.endTime().value();
39 Info << "startTime = " << startTime << endl;
40 Info << "endTime = " << endTime << endl;
41 Info << "deltaTime = " << deltaTime << endl;
42 Info << "Ntimes = " << Ntimes << endl;
43 timeVector.resize(Ntimes);
44 probe_mean.resize(Nprobes, Ntimes);
48 }
49
51 autoPtr<simpleControl> _simple;
52 autoPtr<fv::options> _fvOptions;
54 double k;
56 double rho;
58 double Cp;
60 label left_ind;
61 Eigen::MatrixXd measurements;
65 int Nseeds;
66
68 PtrList<volScalarField> Tensemble;
70 PtrList<volScalarField> Tmean;
72
73 std::shared_ptr<muq::Modeling::Gaussian> priorDensity;
74 std::shared_ptr<muq::Modeling::Gaussian> modelErrorDensity;
75 std::shared_ptr<muq::Modeling::Gaussian> measNoiseDensity;
76 double meas_cov;
77
78 Eigen::MatrixXd priorSamples;
79 Eigen::MatrixXd posteriorSamples;
80
81 double startTime;
82 double deltaTime;
83 double endTime;
84 int Ntimes;
85
87 Eigen::VectorXd timeVector;
89 int Nprobes = 1;
91 Foam::vector probePosition = Foam::vector(0.5, 0.1, 0.005);
92 Eigen::MatrixXd probe_true;
93 Eigen::MatrixXd probe_mean;
94 Eigen::MatrixXd probe_MaxConfidence;
95 Eigen::MatrixXd probe_minConfidence;
96
98 int sampleI = 0;
99 int sampleFlag = 0;
100
103 double updateBC(bool reconstruction = 0)
104 {
105 Time& runTime = _runTime();
106 scalar time = runTime.value();
107 double BC = time * 5 + 2;
108
109 if (reconstruction)
110 {
111 BC = time * 5 + 0;
112 }
113
114 return BC;
115 }
116
119 {
120 samplingTimeI = 0;
121 sampleI = 0;
122 sampleFlag = 0;
123 }
124
126 Eigen::VectorXd observe(volScalarField field)
127 {
128 fvMesh& mesh = _mesh();
129 Time& runTime = _runTime();
130 IOdictionary measurementsDict
131 (
132 IOobject
133 (
134 "measurementsDict",
135 runTime.constant(),
136 mesh,
137 IOobject::MUST_READ,
138 IOobject::NO_WRITE
139 )
140 );
141 List<vector> measurementPoints(measurementsDict.lookup("positions"));
142 Eigen::VectorXd measures(measurementPoints.size());
143 forAll(measurementPoints, pntI)
144 {
145 //TODO
146 //add check if I put two measurements on the same cell
147 measures(pntI) = field[mesh.findCell(measurementPoints[pntI])];
148 }
149 return measures;
150 }
151
153 Eigen::MatrixXd solve(volScalarField& T, word outputFolder,
154 word outputFieldName, bool reconstruction = 0)
155 {
156 Eigen::MatrixXd obsMat;
157 fvMesh& mesh = _mesh();
158 Foam::Time& runTime = _runTime();
159 simpleControl& simple = _simple();
160 fv::options& fvOptions(_fvOptions());
161 dimensionedScalar diffusivity("diffusivity", dimensionSet(0, 2, -1, 0, 0, 0, 0),
162 k / (rho * Cp));
163 Info << "Thermal diffusivity = " << diffusivity << " m2/s" << endl;
164 int timeI = 0;
165
166 while (runTime.loop())
167 {
168 scalar time = runTime.value();
169
170 if (!reconstruction)
171 {
172 sampleFlag++;
173 timeVector(timeI) = time;
174 timeI++;
175 }
176
177 // Update BC
178 ITHACAutilities::assignBC(T, mesh.boundaryMesh().findPatchID("left"),
179 updateBC(reconstruction));
180
181 // Solve
182 while (simple.correctNonOrthogonal())
183 {
184 fvScalarMatrix TEqn
185 (
186 fvm::ddt(T) - fvm::laplacian(diffusivity, T)
187 );
188 fvOptions.constrain(TEqn);
189 TEqn.solve();
190 fvOptions.correct(T);
191 }
192
193 if (reconstruction)
194 {
195 //Adding model error
196 forAll(T.internalField(), cellI)
197 {
198 T.ref()[cellI] += modelErrorDensity->Sample()(cellI);
199 }
200 }
201 else
202 {
203 probe_true.col(timeI - 1) = sampleField(T, probePosition);
204 }
205
207 outputFolder,
208 outputFieldName);
209
211 {
212 Info << "Sampling at time = " << runTime.timeName() << nl << endl;
213 obsMat.conservativeResize(observe(T).size(), obsMat.cols() + 1);
214 obsMat.col(sampleI) = observe(T);
215
216 if (!reconstruction)
217 {
218 sampleFlag = 0;
219 sampleI++;
220 }
221 }
222
223 runTime.write();
224 }
225
226 if (!reconstruction)
227 {
228 ITHACAstream::exportMatrix(timeVector, "timeVector", "eigen", outputFolder);
229 ITHACAstream::exportMatrix(probe_true, "probe_true", "eigen", outputFolder);
230 }
231
232 return obsMat;
233 }
234
237 {
238 Info << "\n****************************************\n" << endl;
239 Info << "\nPerforming true solution\n" << endl;
240 word outputFolder = "./ITHACAoutput/direct/";
241 restart();
242 volScalarField& T(_T());
244 measurements = solve(T, outputFolder, "Tdirect");
245 ITHACAstream::exportMatrix(measurements, "trueMeasurements", "eigen",
246 outputFolder);
247 std::cout << "Number of samples in time = " << measurements.cols() << std::endl;
248 std::cout << "Number of sample in space = " << measurements.rows() << std::endl;
250 Info << "\nEND true solution\n" << endl;
251 Info << "\n****************************************\n" << endl;
252 }
253
254
257 {
258 restart();
259 Time& runTime = _runTime();
261 std::vector<Eigen::MatrixXd> observationTensor;
262 word outputFolder = "ITHACAoutput/reconstuction";
263 Tmean.resize(Ntimes);
264 Info << "\n****************************************\n" << endl;
265 Info << "\nStarting reconstruction\n" << endl;
266 int samplesCounter = 0;
267 // Read true field
268 PtrList<volScalarField> Ttrue;
269 ITHACAstream::read_fields(Ttrue, "Tdirect",
270 "ITHACAoutput/direct/");
271
272 for (int timeI = 0; timeI < Ntimes; timeI++)
273 {
274 double initialTime = timeI * deltaTime + startTime;
275 Info << "\n\nTime " << initialTime + deltaTime << endl;
276 Eigen::MatrixXd forecastOutput = forecastStep(initialTime, timeI,
277 initialTime + deltaTime);
278
279 if (forecastOutput.cols() > 0)
280 {
281 // This is a sampling step
282 samplesCounter++;
283 observationTensor.push_back(forecastOutput);
284 sampleFlag = 0;
285 Eigen::VectorXd meas = measurements.col(samplesCounter - 1);
286 // Kalman filter
288 Eigen::MatrixXd::Identity(meas.size(), meas.size()) * meas_cov, forecastOutput);
289
290 for (int i = 0; i < Nseeds; i++)
291 {
292 restart();
293 volScalarField& T = _T();
294 Eigen::VectorXd internalField = posteriorSamples.col(i);
295 assignIF(T, internalField);
296 Tensemble.Foam::PtrList<volScalarField>::set(i, T.clone());
297 }
298 }
299 else
300 {
301 Info << "debug : NOT a sampling step" << endl;
302 }
303
304 // Fill the mean value field
305 Eigen::VectorXd internalField = Foam2Eigen::PtrList2Eigen(
306 Tensemble).rowwise().mean();
307 volScalarField& Tmean = _T();
308 assignIF(Tmean, internalField);
309 ITHACAstream::exportSolution(Tmean, std::to_string(initialTime + deltaTime),
310 outputFolder, "Tmean");
311 ITHACAstream::exportSolution(Ttrue[timeI],
312 std::to_string(initialTime + deltaTime),
313 outputFolder, "Ttrue");
314 // Save values at measurements points
315 Eigen::MatrixXd ensambleProbe(Nprobes, Nseeds);
316
317 for (int i = 0; i < Nseeds; i++)
318 {
319 ensambleProbe.col(i) = sampleField(Tensemble[i], probePosition);
320 }
321
322 probe_mean.col(timeI) = ensambleProbe.rowwise().mean();
323 probe_MaxConfidence.col(timeI) = ITHACAmuq::muq2ithaca::quantile(ensambleProbe,
324 0.95);
325 probe_minConfidence.col(timeI) = ITHACAmuq::muq2ithaca::quantile(ensambleProbe,
326 0.05);
327 }
328
329 ITHACAstream::exportMatrix(probe_mean, "probe_mean", "eigen", outputFolder);
330 ITHACAstream::exportMatrix(probe_MaxConfidence, "probe_MaxConfidence", "eigen",
331 outputFolder);
332 ITHACAstream::exportMatrix(probe_minConfidence, "probe_minConfidence", "eigen",
333 outputFolder);
334 Info << "\n****************************************\n" << endl;
335 }
336
337
339 Eigen::MatrixXd forecastStep(double startTime_, int startIndex_,
340 double endTime_)
341 {
342 word outputFolder = "ITHACAoutput/forwardSamples";
343 bool reconstructionFlag = 1;
344 Eigen::MatrixXd observationMat;
345 sampleFlag++;
346
347 for (int i = 0; i < Nseeds; i++)
348 {
349 word fieldName = "Tsample" + std::to_string(i);
350 volScalarField T(_T());
351
352 if (startIndex_ == 0)
353 {
354 restart();
355 T = _T();
356 Eigen::VectorXd internalField = priorSamples.col(i);
357 assignIF(T, internalField);
358 }
359 else
360 {
361 T = Tensemble.set(i, nullptr)();
362 }
363
364 resetRunTime(startTime_, startIndex_, endTime_);
365
366 if (observationMat.cols() == 0)
367 {
368 observationMat = solve(T, outputFolder, fieldName, reconstructionFlag);
369 }
370 else
371 {
372 observationMat.conservativeResize(observationMat.rows(),
373 observationMat.cols() + 1);
374 Eigen::MatrixXd observationForCheck = solve(T, outputFolder, fieldName,
375 reconstructionFlag);
376 M_Assert(observationForCheck.cols() == 1,
377 "The forecast step passed trough more than a sampling step");
378 observationMat.col(i) = observationForCheck;
379 }
380
381 Tensemble.set(i, T.clone());
382 }
383
384 return observationMat;
385 }
386
388 void assignIF(volScalarField& field_, Eigen::VectorXd internalField_)
389 {
390 for (int i = 0; i < internalField_.size(); i++)
391 {
392 field_.ref()[i] = internalField_(i);
393 }
394 }
395
397 Eigen::VectorXd sampleField(volScalarField field_, Foam::vector probeLocation_)
398 {
399 Foam::fvMesh& mesh = _mesh();
400 Eigen::VectorXd output(1);
401 output(0) = field_[mesh.findCell(probeLocation_)];
402 return output;
403 }
404
406 void restart()
407 {
408 Time& runTime = _runTime();
409 instantList Times = runTime.times();
410 runTime.setTime(Times[1], 0);
411 _simple.clear();
412 _T.clear();
413 Foam::fvMesh& mesh = _mesh();
414 _simple = autoPtr<simpleControl>
415 (
416 new simpleControl
417 (
418 mesh
419 )
420 );
421 _T = autoPtr<volScalarField>
422 (
423 new volScalarField
424 (
425 IOobject
426 (
427 "T",
428 runTime.timeName(),
429 mesh,
430 IOobject::MUST_READ,
431 IOobject::AUTO_WRITE
432 ),
433 mesh
434 )
435 );
436 }
437
439 void resetRunTime(double startTime_, int startIndex_, double endTime_)
440 {
441 Time& runTime = _runTime();
442 instantList Times = runTime.times();
443 runTime.setTime(startTime_, startIndex_);
444 runTime.setEndTime(endTime_);
445 _simple.clear();
446 Foam::fvMesh& mesh = _mesh();
447 _simple = autoPtr<simpleControl>
448 (
449 new simpleControl
450 (
451 mesh
452 )
453 );
454 }
455
457 void priorSetup(double mean, double cov)
458 {
459 volScalarField& T(_T());
460 int stateSize = T.size();
461 Eigen::VectorXd prior_mu = Eigen::MatrixXd::Ones(stateSize, 1) * mean;
462 Eigen::MatrixXd prior_cov = Eigen::MatrixXd::Identity(stateSize,
463 stateSize) * cov;
464 priorDensity = std::make_shared<muq::Modeling::Gaussian>(prior_mu, prior_cov);
465 }
466
469 {
470 for (int i = 0; i < Nseeds; i++)
471 {
472 priorSamples.col(i) = priorDensity->Sample();
473 }
474
476 }
477
479 void modelErrorSetup(double mean, double cov)
480 {
481 volScalarField& T(_T());
482 int stateSize = T.size();
483 Eigen::VectorXd modelError_mu = Eigen::MatrixXd::Ones(stateSize, 1) * mean;
484 Eigen::MatrixXd modelError_cov = Eigen::MatrixXd::Identity(stateSize,
485 stateSize) * cov;
486 modelErrorDensity = std::make_shared<muq::Modeling::Gaussian>(modelError_mu,
487 modelError_cov);
488 }
489
491 void measNoiseSetup(double mean, double cov)
492 {
493 meas_cov = cov;
494 int obsSize = measurements.rows();
495 M_Assert(obsSize > 0, "Read measurements before setting up the noise");
496 Eigen::VectorXd measNoise_mu = Eigen::MatrixXd::Ones(obsSize, 1) * mean;
497 Eigen::MatrixXd measNoise_cov = Eigen::MatrixXd::Identity(obsSize,
498 obsSize) * cov;
499 measNoiseDensity = std::make_shared<muq::Modeling::Gaussian>(measNoise_mu,
500 measNoise_cov);
501 }
502};
503
504#endif
forAll(example_CG.gList, solutionI)
Definition CGtest.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
TEqn solve()
Class where the UQ tutorial number 2 is implemented.
Eigen::VectorXd observe(volScalarField field)
Returns the input field values at the measurement points defined in the measurementsDict.
PtrList< volScalarField > Tmean
Mean field of the ensemble.
label left_ind
Index of the left patch.
void resetRunTime(double startTime_, int startIndex_, double endTime_)
Reset runTime to given values.
std::shared_ptr< muq::Modeling::Gaussian > priorDensity
std::shared_ptr< muq::Modeling::Gaussian > measNoiseDensity
void priorSetup(double mean, double cov)
Setup of the prior distribution.
int Nseeds
Number of seeds in the ensemble.
Eigen::VectorXd timeVector
Timesteps vector.
void assignIF(volScalarField &field_, Eigen::VectorXd internalField_)
Assign internalField given an Eigen Vector.
EnKF_1DinverseHeatTransfer(int argc, char *argv[], int Nseeds_)
autoPtr< simpleControl > _simple
[tutorial02]
void measNoiseSetup(double mean, double cov)
Setup of the measurement noise distribution.
std::shared_ptr< muq::Modeling::Gaussian > modelErrorDensity
void modelErrorSetup(double mean, double cov)
Setup of the model error distribution.
double updateBC(bool reconstruction=0)
Updates the boundary conditions according to runTime.
Foam::vector probePosition
Probles position, probes are used only for visualization and postprocessing.
void reconstruct()
Reconstruction phase.
void resetSamplingCounters()
Set all the sampling counters to 0.
Eigen::MatrixXd forecastStep(double startTime_, int startIndex_, double endTime_)
Forecasting step.
PtrList< volScalarField > Tensemble
Ensemble.
void restart()
Restart before a now solution.
void solveDirect()
Preforming a true solution.
int NtimeStepsBetweenSamps
Number of timesteps between each sampling.
Eigen::MatrixXd solve(volScalarField &T, word outputFolder, word outputFieldName, bool reconstruction=0)
Solve the heat transfer problem.
void priorSampling()
Samples the prior density.
Eigen::VectorXd sampleField(volScalarField field_, Foam::vector probeLocation_)
Samples a field in a given point.
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
Class to implement a full order laplacian parametrized problem.
autoPtr< volScalarField > _T
Temperature field.
autoPtr< fvMesh > _mesh
Mesh.
autoPtr< Time > _runTime
Time.
volScalarField & T
Definition createT.H:46
fvScalarMatrix & TEqn
Definition TEqn.H:15
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
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
void assignBC(GeometricField< scalar, fvPatchField, volMesh > &s, label BC_ind, double value)
Assign uniform Boundary Condition to a volScalarField.
simpleControl simple(mesh)
label i
Definition pEqn.H:46