Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_ProductLanczosPCEBasis.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stokhos Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef STOKHOS_PRODUCT_LANCZOS_PCE_BASIS_HPP
43 #define STOKHOS_PRODUCT_LANCZOS_PCE_BASIS_HPP
44 
45 #include "Teuchos_RCP.hpp"
46 #include "Teuchos_Array.hpp"
50 
51 #include "Stokhos_ProductBasis.hpp"
54 #include "Stokhos_Quadrature.hpp"
57 
58 namespace Stokhos {
59 
70  template <typename ordinal_type, typename value_type>
72  public ProductBasis<ordinal_type,value_type>,
73  public ReducedPCEBasis<ordinal_type,value_type> {
74  public:
75 
77 
83  ordinal_type p,
88 
90  virtual ~ProductLanczosPCEBasis();
91 
93 
94 
96  ordinal_type order() const;
97 
99  ordinal_type dimension() const;
100 
102  virtual ordinal_type size() const;
103 
105 
109  virtual const Teuchos::Array<value_type>& norm_squared() const;
110 
112  virtual const value_type& norm_squared(ordinal_type i) const;
113 
115 
121  virtual
124 
126  virtual
129 
131  virtual value_type evaluateZero(ordinal_type i) const;
132 
134 
138  virtual void evaluateBases(
140  Teuchos::Array<value_type>& basis_vals) const;
141 
143  virtual void print(std::ostream& os) const;
144 
146  virtual const std::string& getName() const;
147 
149 
151 
152 
154 
159  virtual const MultiIndex<ordinal_type>& term(ordinal_type i) const;
160 
162 
166  virtual ordinal_type index(const MultiIndex<ordinal_type>& term) const;
167 
169 
173  value_type> > >
174  getCoordinateBases() const;
175 
177  virtual MultiIndex<ordinal_type> getMaxOrders() const;
178 
180 
182 
183 
185  virtual void
187  value_type *out,
188  ordinal_type ncol = 1,
189  bool transpose = false) const;
190 
192  virtual void
194  value_type *out,
195  ordinal_type ncol = 1,
196  bool transpose = false) const;
197 
200  getReducedQuadrature() const;
201 
203 
204  protected:
205 
206  // Determine if a pce is linear, in that it has a total degree of at
207  // most 1. If the pce is nonlinear, return -2, or if it is constant,
208  // return -1, otherwise return the index of the variable the pce is
209  // linear in, ie, if return value is i, the pce = a_{i+1}*\xi_i
212 
213  private:
214 
215  // Prohibit copying
217 
218  // Prohibit Assignment
220 
221  protected:
222 
225 
227  std::string name;
228 
231 
234 
237 
240 
243 
246 
249 
250  }; // class ProductLanczosPCEBasis
251 
252 } // Namespace Stokhos
253 
254 // Include template definitions
256 
257 #endif
Teuchos::SerialDenseMatrix< ordinal_type, value_type > SDM
virtual void print(std::ostream &os) const
Print basis to stream os.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
SDM Phi
Values of transformed basis at quadrature points.
Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product me...
SDM A
Transition matrix from reduced basis to original.
SDM Ainv
Projection matrix from original matrix to reduced.
virtual Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > getReducedQuadrature() const
Get reduced quadrature object.
Stokhos::OrthogPolyApprox< ordinal_type, value_type > tmp_pce
Temporary pce used in invariant subspace calculations.
ProductLanczosPCEBasis(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::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::ParameterList &params=Teuchos::ParameterList())
Constructor.
virtual const std::string & getName() const
Return string name of basis.
virtual ordinal_type index(const MultiIndex< ordinal_type > &term) const
Get index of the multivariate polynomial given orders of each coordinate.
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.
ordinal_type dimension() const
Return dimension of basis.
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.
ordinal_type isInvariant(const Stokhos::OrthogPolyApprox< ordinal_type, value_type > &pce) const
Teuchos::ParameterList params
Algorithm parameters.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
Teuchos::RCP< Stokhos::CompletePolynomialBasis< ordinal_type, value_type > > tensor_lanczos_basis
Product Lanczos basis.
virtual MultiIndex< ordinal_type > getMaxOrders() const
Return maximum order allowable for each coordinate basis.
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 ordinal_type size() const
Return total size of basis.
Abstract base class for 1-D orthogonal polynomials.
virtual const MultiIndex< ordinal_type > & term(ordinal_type i) const
Get orders of each coordinate polynomial given an index i.
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
Teuchos::SerialDenseVector< ordinal_type, value_type > SDV
ordinal_type order() const
Return order of basis.
ProductLanczosPCEBasis & operator=(const ProductLanczosPCEBasis &)
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > getCoordinateBases() const
Return coordinate bases.
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.