Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_GSReducedPCEBasisBase.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_GS_REDUCED_PCE_BASIS_BASE_HPP
11 #define STOKHOS_GS_REDUCED_PCE_BASIS_BASE_HPP
12 
13 #include "Teuchos_RCP.hpp"
14 #include "Teuchos_Array.hpp"
15 #include "Teuchos_BLAS.hpp"
19 
22 #include "Stokhos_Quadrature.hpp"
24 
25 namespace Stokhos {
26 
36  template <typename ordinal_type, typename value_type>
38  public ReducedPCEBasis<ordinal_type,value_type> {
39  public:
40 
42 
54 
56  virtual ~GSReducedPCEBasisBase();
57 
59 
60 
62  ordinal_type order() const;
63 
65  ordinal_type dimension() const;
66 
68  virtual ordinal_type size() const;
69 
71 
75  virtual const Teuchos::Array<value_type>& norm_squared() const;
76 
78  virtual const value_type& norm_squared(ordinal_type i) const;
79 
81 
87  virtual
90 
92  virtual
95 
97  virtual value_type evaluateZero(ordinal_type i) const;
98 
100 
104  virtual void evaluateBases(
106  Teuchos::Array<value_type>& basis_vals) const;
107 
109  virtual void print(std::ostream& os) const;
110 
112 
114 
115 
117  virtual void
119  value_type *out,
120  ordinal_type ncol = 1,
121  bool transpose = false) const;
122 
124  virtual void
126  value_type *out,
127  ordinal_type ncol = 1,
128  bool transpose = false) const;
129 
132  getReducedQuadrature() const;
133 
135 
136  protected:
137 
138  void setup(
139  ordinal_type p,
142 
144 
147  virtual ordinal_type
149  ordinal_type max_p,
150  value_type threshold,
153  const Teuchos::Array<value_type>& weights,
155  Teuchos::Array<ordinal_type>& num_terms_,
158 
159  private:
160 
161  // Prohibit copying
163 
164  // Prohibit Assignment
166 
167  protected:
168 
172 
174  std::string name;
175 
178 
181 
184 
187 
190 
193 
196 
199 
202 
205 
208 
211 
213  bool verbose;
214 
217 
220 
222 
223  }; // class GSReducedPCEBasisBase
224 
225 } // Namespace Stokhos
226 
227 // Include template definitions
229 
230 #endif
ordinal_type sz
Total size of basis.
SDM Qp
Coefficients of transformed basis in original basis.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
Utilities for indexing a multi-variate complete polynomial basis.
bool verbose
Whether to print a bunch of stuff out.
Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product me...
SDM Q
Values of transformed basis at quadrature points.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > pce_basis
Original pce basis.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > SDM
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > reduced_quad
Reduced quadrature object.
virtual void evaluateBases(const Teuchos::ArrayView< const value_type > &point, Teuchos::Array< value_type > &basis_vals) const
Evaluate basis polynomials at given point point.
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
void setup(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad)
virtual void print(std::ostream &os) const
Print basis to stream os.
Abstract base class for quadrature methods.
Abstract base class for reduced basis strategies built from polynomial chaos expansions in some other...
virtual void transformToOriginalBasis(const value_type *in, value_type *out, ordinal_type ncol=1, bool transpose=false) const
Transform coefficients to original basis from this basis.
Stokhos::CompletePolynomialBasisUtils< ordinal_type, value_type > CPBUtils
Teuchos::Array< value_type > norms
Norms.
GSReducedPCEBasisBase(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, const Teuchos::ParameterList &params=Teuchos::ParameterList())
Constructor.
virtual ordinal_type size() const
Return total size of basis.
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.
std::string orthogonalization_method
Orthogonalization method.
Teuchos::Array< Stokhos::MultiIndex< ordinal_type > > terms
2-D array of basis terms
ordinal_type order() const
Return order of basis.
ordinal_type p
Total order of basis.
Teuchos::ParameterList params
Algorithm parameters.
Teuchos::SerialDenseVector< ordinal_type, value_type > SDV
ordinal_type dimension() const
Return dimension of basis.
ordinal_type pce_sz
Size of original pce basis.
Teuchos::Array< ordinal_type > num_terms
Number of terms up to each order.
GSReducedPCEBasisBase & operator=(const GSReducedPCEBasisBase &b)
ordinal_type d
Total dimension of basis.
Teuchos::BLAS< ordinal_type, value_type > blas
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
virtual ordinal_type buildReducedBasis(ordinal_type max_p, value_type threshold, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &A, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &F, const Teuchos::Array< value_type > &weights, Teuchos::Array< Stokhos::MultiIndex< ordinal_type > > &terms_, Teuchos::Array< ordinal_type > &num_terms_, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Qp_, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Q_)=0
Build the reduced basis, parameterized by total order max_p.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
virtual void transformFromOriginalBasis(const value_type *in, value_type *out, ordinal_type ncol=1, bool transpose=false) const
Transform coefficients from original basis to this basis.