Loading...
Searching...
No Matches
muq2ithaca.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 muq2ithaca
26Description
27 Set of function to convert ITHACA objects into MUQ2 objects and vice versa
28SourceFiles
29 muq2ithaca.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef muq2ithaca_H
38#define muq2ithaca_H
39
40#include <iostream>
41#include <Eigen/Eigen>
42#include "ITHACAassert.H"
43#include "Foam2Eigen.H"
44
45#include "MUQ/Modeling/Distributions/Gaussian.h"
46#include "MUQ/Modeling/Distributions/Density.h"
47
48#include "MUQ/SamplingAlgorithms/SamplingProblem.h"
49#include "MUQ/SamplingAlgorithms/SingleChainMCMC.h"
50#include "MUQ/SamplingAlgorithms/MCMCFactory.h"
51
52#include <boost/property_tree/ptree.hpp>
53
54
55namespace ITHACAmuq
56{
57namespace muq2ithaca
58{
59//--------------------------------------------------------------------------
69Eigen::MatrixXd EnsembleKalmanFilter(Eigen::MatrixXd prior,
70 Eigen::VectorXd measurements,
71 Eigen::MatrixXd measurementsCov,
72 Eigen::MatrixXd observedState);
73
74//--------------------------------------------------------------------------
84Eigen::MatrixXd EnsembleKalmanFilter(PtrList<volScalarField>& prior,
85 Eigen::VectorXd measurements,
86 Eigen::MatrixXd measurementsCov,
87 Eigen::MatrixXd observedState);
88
89//--------------------------------------------------------------------------
98double quantile(Eigen::VectorXd samps, double p, int method = 1);
99
100//--------------------------------------------------------------------------
109Eigen::VectorXd quantile(Eigen::MatrixXd samps, double p, int method = 1);
110}
111}
112
113#endif
Header file of the Foam2Eigen class.
Implementation of the assert function for ITHACA-FV.
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
volScalarField & p