Loading...
Searching...
No Matches
LRSensitivity.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-------------------------------------------------------------------------------
12
13 License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30#ifndef LRSensitivity_H
31#define LRSensitivity_H
32#include <iostream>
33#include <vector>
34#include <string>
35#include <Eigen/Eigen>
36#include "FofM.H"
37#include "ITHACAstream.H"
38#include "ITHACAutilities.H"
39#include "ITHACAsampling.H"
40
41
43{
44 public:
45 //Constructors
47 explicit LRSensitivity(label Npara, label Np);
49
50 //Members
52 Eigen::MatrixXd trainingRange;
54 Eigen::MatrixXd MatX;
56 Eigen::MatrixXd MatXn;
58 autoPtr<FofM> M;
62 label Npoints;
64 Eigen::VectorXd y;
66 Eigen::VectorXd yn;
68 Eigen::VectorXd ylin;
70 double Ey;
72 double Vy;
74 Eigen::VectorXd EX;
76 Eigen::VectorXd VX;
78 Eigen::VectorXd betas;
81 double QI;
83 bool Xdone = false;
85 bool Ydone = false;
87 bool bdone = false;
88
89 //--------------------------------------------------------------------------------------------------------------------//
90 //Methods
91
95 void buildSamplingSet(std::vector<std::string>& pdflist, Eigen::MatrixXd plist);
97
99 void load_output();
101
103 void getXstats();
105
107 void getYstat();
109
111 void getBetas();
113
115 void assessQuality();
117
118
119
120
121
122 private:
123
125 void setAll();
126
127};
128
129#endif
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the ITHACAutilities namespace.
void getYstat()
Method to compute Ey and Vy, it sets Ydone=true.
void getXstats()
Method to compute EX and VX, it sets Xdone=true.
double Vy
Variance of the output.
Eigen::VectorXd y
Vectors to store the output of the model.
Eigen::VectorXd yn
Vector to store the normalized output of the model.
Eigen::MatrixXd MatXn
Matrices storing the independent variables' normalized values.
double QI
double to quantify the goodness of linear regression from 0 (no linearity) to 1 (perfect linearity) a...
Eigen::MatrixXd MatX
Matrices storing the independent variables' values.
Eigen::VectorXd VX
Vector to store the variances of the parameters.
double Ey
Mean value of the output.
label No_parameters
Number of parameters used for the analysis.
void getBetas()
Method to compute SRCs, it sets bdone=true.
void load_output()
Method to load the output of FOM/ROM inside member LRsensitivity member y.
label Npoints
Number of points used for Monte Carlo sampling.
bool bdone
boolean variable to check if SRCs are computed
Eigen::VectorXd ylin
Vector to store the result of linear approximation.
Eigen::MatrixXd trainingRange
Matrix to store the range used in training/offline stage.
Eigen::VectorXd betas
Vector to store the standardized regression coefficient (SRC)
void assessQuality()
Method to compute the quality of linear regression approximation.
void buildSamplingSet(std::vector< std::string > &pdflist, Eigen::MatrixXd plist)
Method to build MatX, it internally calls ITHACAsampling::samplingMC pdflist is the list with the nam...
Eigen::VectorXd EX
Vector to store the meanvalues of the parameters.
autoPtr< FofM > M
Figure of merit object.
bool Ydone
boolean variable to check if Ey and Vy are computed
bool Xdone
boolean variable to check if EX and VX are computed