Loading...
Searching...
No Matches
02enKF_1DinverseHeatTransfer.C
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/>.
24Description
25 Example of a state reconstruction in 1D heat transfer problem using EnKF
26SourceFiles
27 02enKF_1DinverseHeatTransfer.C
28\*---------------------------------------------------------------------------*/
29
30#include <iostream>
31#include "fvCFD.H"
32#include "fvOptions.H"
33#include "simpleControl.H"
34#include "IOmanip.H"
35#include "Time.H"
36#include "laplacianProblem.H"
37#include "ITHACAutilities.H"
38#include <Eigen/Dense>
39#define _USE_MATH_DEFINES
40#include <cmath>
41#include "Foam2Eigen.H"
42
43#include "MUQ/Modeling/Distributions/Gaussian.h"
44
45#include "muq2ithaca.H"
46
48
49using namespace SPLINTER;
50
51
52int main(int argc, char* argv[])
53{
54 solverPerformance::debug = 1; //No verbose output
55 // Number of elements in the ensemble
56 int Nseeds = 50;
57 // Create the example object of the 02enKF_1DinverseHeatTransfer type
58 EnKF_1DinverseHeatTransfer example(argc, argv, Nseeds);
59 // Reading parameters from file
61 example._runTime());
62 example.k = para->ITHACAdict->lookupOrDefault<double>("thermalConductivity", 0);
63 M_Assert(example.k > 0, "thermalConductivity, k, not specified");
64 example.rho = para->ITHACAdict->lookupOrDefault<double>("density", 0);
65 M_Assert(example.rho > 0, "Density, rho, not specified");
66 example.Cp = para->ITHACAdict->lookupOrDefault<double>("heatCapacity", 0);
67 M_Assert(example.Cp > 0, "heatCapacity, Cp, not specified");
68 // Mesh and temperature field setup
69 fvMesh& mesh = example._mesh();
70 volScalarField& T = example._T();
71 int stateSize = T.size();
72 scalar initialField = T.internalField()[0];
73 // Performes true solution
74 example.solveDirect();
75 int Ntimes = example.Ntimes;
76 // Setting up the densities
77 // Gaussian densities are assumed
78 example.priorSetup(initialField, 0.7);
79 example.modelErrorSetup(0.0, 0.7);
80 example.measNoiseSetup(0.0, 0.05);
81 Eigen::MatrixXd posteriorSamples(stateSize, Nseeds);
82 Eigen::MatrixXd priorSamples(stateSize, Nseeds);
83 // Sampling prior density
84 example.priorSampling();
85 Eigen::MatrixXd posteriorMean(stateSize, Ntimes);
86 Eigen::MatrixXd minConfidence = posteriorMean;
87 Eigen::MatrixXd maxConfidence = minConfidence;
88 posteriorMean.col(0) = posteriorSamples.rowwise().mean();
89 // Performs reconstruction
90 example.reconstruct();
91 return 0;
92}
93//--------
97
int main(int argc, char *argv[])
Header file of the Foam2Eigen class.
Foam::fvMesh & mesh
Definition createMesh.H:47
#define M_Assert(Expr, Msg)
Header file of the ITHACAutilities namespace.
Class where the UQ tutorial number 2 is implemented.
void priorSetup(double mean, double cov)
Setup of the prior distribution.
void measNoiseSetup(double mean, double cov)
Setup of the measurement noise distribution.
void modelErrorSetup(double mean, double cov)
Setup of the model error distribution.
void reconstruct()
Reconstruction phase.
void solveDirect()
Preforming a true solution.
void priorSampling()
Samples the prior density.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
IOdictionary * ITHACAdict
Dictionary for input objects from file.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
autoPtr< volScalarField > _T
Temperature field.
autoPtr< fvMesh > _mesh
Mesh.
autoPtr< Time > _runTime
Time.
volScalarField & T
Definition createT.H:46
Header file of the laplacianProblem class.
Header file of the muq2ithaca namespace.