MUQ  0.4.3
PolynomialChaosExpansion.h
Go to the documentation of this file.
1 #ifndef POLYNOMIALCHAOSEXPANSION_H_
2 #define POLYNOMIALCHAOSEXPANSION_H_
3 
6 
7 #include <vector>
8 #include <set>
9 
10 namespace muq {
11 namespace Approximation {
12 
13  class PCEFactory;
14 
29  friend class PCEFactory;
30 
31  public:
32 
33  PolynomialChaosExpansion(std::shared_ptr<OrthogonalPolynomial> const& basisCompsIn,
34  std::shared_ptr<muq::Utilities::MultiIndexSet> const& multisIn,
35  Eigen::MatrixXd const& coeffsIn);
36 
37  PolynomialChaosExpansion(std::shared_ptr<OrthogonalPolynomial> const& basisCompsIn,
38  std::shared_ptr<muq::Utilities::MultiIndexSet> const& multisIn,
39  unsigned int outputDim);
40 
41  PolynomialChaosExpansion(std::vector<std::shared_ptr<IndexedScalarBasis>> const& basisCompsIn,
42  std::shared_ptr<muq::Utilities::MultiIndexSet> const& multisIn,
43  Eigen::MatrixXd const& coeffsIn);
44 
45  PolynomialChaosExpansion(std::vector<std::shared_ptr<IndexedScalarBasis>> const& basisCompsIn,
46  std::shared_ptr<muq::Utilities::MultiIndexSet> const& multisIn,
47  unsigned int outputDim);
48 
49 
50  virtual ~PolynomialChaosExpansion() = default;
51 
53  Eigen::VectorXd Variance() const;
54  Eigen::MatrixXd Covariance() const;
55  Eigen::VectorXd Mean() const;
56 
58  Eigen::VectorXd Magnitude() const;
59 
66  static std::shared_ptr<PolynomialChaosExpansion> ComputeWeightedSum(std::vector<std::shared_ptr<PolynomialChaosExpansion>> expansions,
67  Eigen::VectorXd const& weights);
68 
69 
71  Eigen::VectorXd TotalSensitivity(unsigned int targetDim) const;
72 
74  Eigen::MatrixXd TotalSensitivity() const;
75 
77  Eigen::VectorXd SobolSensitivity(unsigned int targetDim) const;
78 
81  Eigen::VectorXd SobolSensitivity(std::vector<unsigned int> const& targetDims) const;
82 
86  Eigen::MatrixXd MainSensitivity() const;
87 
92  virtual void ToHDF5(muq::Utilities::H5Object &group) const override;
93 
106  static std::shared_ptr<PolynomialChaosExpansion> FromHDF5(std::string filename, std::string groupName="/");
107 
114  static std::shared_ptr<PolynomialChaosExpansion> FromHDF5(muq::Utilities::H5Object &group);
115 
116 
117 
118  private:
119 
125  Eigen::VectorXd GetNormalizationVec() const;
126 
132  static std::shared_ptr<PolynomialChaosExpansion> ComputeWeightedSum(std::vector<std::shared_ptr<PolynomialChaosExpansion>> expansions,
133  Eigen::VectorXd const& weights,
134  std::shared_ptr<muq::Utilities::MultiIndexSet> const& polynomials,
135  std::vector<std::vector<unsigned int>> const& locToGlob);
136 
137  };
138 
139 } // namespace muq
140 } // namespace Approximation
141 
142 
143 
144 #endif /* POLYNOMIALCHAOSEXPANSION_H_ */
Class for defining expansions of basis functions defined by a MultiIndexSet and collection of IndexSc...
Factory class for constructing a pseudo-spectral polynomial chaos approximation using a fixed quadrat...
Definition: PCEFactory.h:22
A class for representing and using expansions of orthogonal multivariate polynomials.
static std::shared_ptr< PolynomialChaosExpansion > FromHDF5(std::string filename, std::string groupName="/")
Loads an expansion from an HDF5 file.
Eigen::VectorXd Variance() const
compute the variance of the current expansion
Eigen::VectorXd Magnitude() const
Compute the L2 norm of each output.
static std::shared_ptr< PolynomialChaosExpansion > ComputeWeightedSum(std::vector< std::shared_ptr< PolynomialChaosExpansion >> expansions, Eigen::VectorXd const &weights)
Eigen::MatrixXd TotalSensitivity() const
Compute all Sobol total sensitivities. Rows are outputs, each column is an input.
PolynomialChaosExpansion(std::shared_ptr< OrthogonalPolynomial > const &basisCompsIn, std::shared_ptr< muq::Utilities::MultiIndexSet > const &multisIn, Eigen::MatrixXd const &coeffsIn)
virtual void ToHDF5(muq::Utilities::H5Object &group) const override
Write PCE to HDF5 file.
Eigen::VectorXd SobolSensitivity(unsigned int targetDim) const
Compute the main sensitivity index for the input dimension, for each output dimension.