Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SmolyakPseudoSpectralOperator.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_SMOLYAK_PSEUDO_SPECTRAL_OPERATOR_HPP
11 #define STOKHOS_SMOLYAK_PSEUDO_SPECTRAL_OPERATOR_HPP
12 
14 #include "Stokhos_SmolyakBasis.hpp"
17 #include "Teuchos_BLAS.hpp"
18 
19 namespace Stokhos {
20 
25  template <typename ordinal_t,
26  typename value_t,
27  typename point_compare_type =
28  typename DefaultPointCompare<ordinal_t,value_t>::type>
30  public PseudoSpectralOperator<ordinal_t,value_t,point_compare_type> {
31  public:
32 
33  typedef ordinal_t ordinal_type;
34  typedef value_t value_type;
36  typedef typename base_type::point_type point_type;
39  typedef typename base_type::iterator iterator;
43 
47 
49  template <typename coeff_compare_type>
52  bool use_smolyak_apply = true,
53  bool use_pst = true,
54  const point_compare_type& point_compare = point_compare_type());
55 
58 
60  ordinal_type point_size() const { return points.size(); }
61 
63  ordinal_type coeff_size() const { return coeff_sz; }
64 
66  iterator begin() { return point_map.begin(); }
67 
69  iterator end() { return point_map.end(); }
70 
72  const_iterator begin() const { return point_map.begin(); }
73 
75  const_iterator end() const { return point_map.end(); }
76 
78  set_iterator set_begin() { return points.begin(); }
79 
81  set_iterator set_end() { return points.end(); }
82 
84  const_set_iterator set_begin() const { return points.begin(); }
85 
87  const_set_iterator set_end() const { return points.end(); }
88 
91  const_set_iterator it = points.find(point);
93  it == points.end(), std::logic_error, "Invalid term " << point);
94  return it->second.second;
95  }
96 
98  const point_type& point(ordinal_type n) const { return point_map[n]; }
99 
101 
108  void transformQP2PCE(
109  const value_type& alpha,
112  const value_type& beta,
113  bool trans = false) const;
114 
116 
123  virtual void transformPCE2QP(
124  const value_type& alpha,
127  const value_type& beta,
128  bool trans = false) const;
129 
130  protected:
131 
133  void apply_direct(
135  const value_type& alpha,
138  const value_type& beta,
139  bool trans) const;
140 
146  const value_type& alpha,
149  const value_type& beta,
150  bool trans) const;
151 
157  const value_type& alpha,
160  const value_type& beta,
161  bool trans) const;
162 
163  void gather(
166  bool trans,
168 
169  void scatter(
172  bool trans,
174 
175  protected:
176 
179 
182 
185 
188 
191 
194 
197 
200 
203 
206 
209 
210  };
211 
212 }
213 
215 
216 #endif
Teuchos::SerialDenseMatrix< ordinal_type, value_type > pce2qp
Matrix mapping coefficients to points.
void apply_direct(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &A, const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans) const
Apply transformation operator using direct method.
std::map< point_type, std::pair< value_type, ordinal_type >, point_compare_type > point_set_type
point_map_type point_map
Map index to sparse grid term.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > qp2pce
Matrix mapping points to coefficients.
set_iterator set_end()
Iterator to end of point set.
Container storing a term in a generalized tensor product.
ordinal_type index(const point_type &point) const
Get point index for given point.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
An operator interface for building pseudo-spectral approximations.
iterator begin()
Iterator to begining of point set.
Teuchos::Array< Teuchos::RCP< operator_type > > operator_set_type
const_iterator begin() const
Iterator to begining of point set.
Teuchos::Array< Teuchos::Array< ordinal_type > > gather_maps
Gather maps for each operator for Smolyak apply.
Teuchos::Array< Teuchos::Array< ordinal_type > > scatter_maps
Scatter maps for each operator for Smolyak apply.
ordinal_type coeff_size() const
Number of coefficients.
void transformQP2PCE(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans=false) const
Transform values at quadrature points to PCE coefficients.
void transformPCE2QP_smolyak(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans) const
Transform PCE coefficients to values at quadrature points using Smolyak formula.
void scatter(const Teuchos::Array< ordinal_type > &map, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, bool trans, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result) const
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
Teuchos::BLAS< ordinal_type, value_type > blas
BLAS wrappers.
void gather(const Teuchos::Array< ordinal_type > &map, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, bool trans, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result) const
const_set_iterator set_end() const
Iterator to end of point set.
Teuchos::Array< value_type > smolyak_coeffs
Smolyak coefficients.
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
iterator end()
const_set_iterator set_begin() const
Iterator to begining of point set.
const_iterator end() const
Iterator to end of point set.
void transformQP2PCE_smolyak(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans) const
Transform values at quadrature points to PCE coefficients using Smolyak formula.
SmolyakPseudoSpectralOperator(const SmolyakBasis< ordinal_type, value_type, coeff_compare_type > &smolyak_basis, bool use_smolyak_apply=true, bool use_pst=true, const point_compare_type &point_compare=point_compare_type())
Constructor.
const point_type & point(ordinal_type n) const
Get point for given index.
An operator for building pseudo-spectral coefficients using tensor-product quadrature.
TensorProductPseudoSpectralOperator< ordinal_type, value_type > operator_type
set_iterator set_begin()
Iterator to begining of point set.
iterator begin()
virtual void transformPCE2QP(const value_type &alpha, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &result, const value_type &beta, bool trans=false) const
Transform PCE coefficients to quadrature values.
PseudoSpectralOperator< ordinal_type, value_type, point_compare_type > base_type
operator_set_type operators
Tensor pseudospectral operators.