Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SmolyakSparseGridQuadratureImp.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 #include "Teuchos_TimeMonitor.hpp"
11 
12 template <typename ordinal_type, typename value_type, typename point_compare_type>
13 template <typename index_set_type>
16  const Teuchos::RCP<const ProductBasis<ordinal_type,value_type> >& product_basis,
17  const index_set_type& index_set,
18  const value_type duplicate_tol,
19  const point_compare_type& point_compare)
20 {
21 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
22  TEUCHOS_FUNC_TIME_MONITOR("Stokhos: Sparse Grid Generation");
23 #endif
24 
25  typedef SmolyakBasis<ordinal_type,value_type> smolyak_basis_type;
26  smolyak_basis_type smolyak_basis(
27  product_basis->getCoordinateBases(), index_set, duplicate_tol);
28 
30  smolyak_operator_type smolyak_operator(smolyak_basis, true, true,
31  point_compare);
32  ordinal_type nqp = smolyak_operator.point_size();
33  ordinal_type npc = product_basis->size();
34 
35  // Compute quad points, weights, values
36  quad_points.resize(nqp);
37  quad_weights.resize(nqp);
38  quad_values.resize(nqp);
39  typedef typename smolyak_operator_type::const_set_iterator const_iterator;
40  ordinal_type i = 0;
41  for (const_iterator it = smolyak_operator.set_begin();
42  it != smolyak_operator.set_end(); ++it, ++i) {
43  quad_points[i] = it->first.getTerm();
44  quad_weights[i] = it->second.first;
45  quad_values[i].resize(npc);
46  product_basis->evaluateBases(quad_points[i], quad_values[i]);
47  }
48 }
49 
50 template <typename ordinal_type, typename value_type, typename point_compare_type>
54 {
55  return quad_points;
56 }
57 
58 template <typename ordinal_type, typename value_type, typename point_compare_type>
62 {
63  return quad_weights;
64 }
65 
66 template <typename ordinal_type, typename value_type, typename point_compare_type>
70 {
71  return quad_values;
72 }
73 
74 template <typename ordinal_type, typename value_type, typename point_compare_type>
75 std::ostream&
77 print(std::ostream& os) const
78 {
79  ordinal_type nqp = quad_weights.size();
80  os << "Smolyak Sparse Grid Quadrature with " << nqp << " points:"
81  << std::endl << "Weight : Points" << std::endl;
82  for (ordinal_type i=0; i<nqp; i++) {
83  os << i << ": " << quad_weights[i] << " : ";
84  for (ordinal_type j=0; j<static_cast<ordinal_type>(quad_points[i].size());
85  j++)
86  os << quad_points[i][j] << " ";
87  os << std::endl;
88  }
89  os << "Basis values at quadrature points:" << std::endl;
90  for (ordinal_type i=0; i<nqp; i++) {
91  os << i << " " << ": ";
92  for (ordinal_type j=0; j<static_cast<ordinal_type>(quad_values[i].size());
93  j++)
94  os << quad_values[i][j] << " ";
95  os << std::endl;
96  }
97 
98  return os;
99 }
SmolyakSparseGridQuadrature(const Teuchos::RCP< const ProductBasis< ordinal_type, value_type > > &product_basis, const index_set_type &index_set, const value_type duplicate_tol=1.0e-12, const point_compare_type &point_compare=point_compare_type())
Constructor.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
virtual const Teuchos::Array< value_type > & getQuadWeights() const
Get quadrature weights.
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getQuadPoints() const
Get quadrature points.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const
Get values of basis at quadrature points.
virtual std::ostream & print(std::ostream &os) const
Print quadrature data.