Stokhos Package Browser (Single Doxygen Collection)
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
src
Stokhos_SmolyakSparseGridQuadratureImp.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
#include "Teuchos_TimeMonitor.hpp"
43
44
template
<
typename
ordinal_type,
typename
value_type,
typename
po
int
_compare_type>
45
template
<
typename
index_set_type>
46
Stokhos::SmolyakSparseGridQuadrature<ordinal_type, value_type, point_compare_type>::
47
SmolyakSparseGridQuadrature
(
48
const
Teuchos::RCP<
const
ProductBasis<ordinal_type,value_type>
>& product_basis,
49
const
index_set_type& index_set,
50
const
value_type
duplicate_tol,
51
const
point_compare_type& point_compare)
52
{
53
#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
54
TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos: Sparse Grid Generation"
);
55
#endif
56
57
typedef
SmolyakBasis<ordinal_type,value_type>
smolyak_basis_type;
58
smolyak_basis_type smolyak_basis(
59
product_basis->getCoordinateBases(), index_set, duplicate_tol);
60
61
typedef
SmolyakPseudoSpectralOperator<ordinal_type,value_type,point_compare_type>
smolyak_operator_type;
62
smolyak_operator_type smolyak_operator(smolyak_basis,
true
,
true
,
63
point_compare);
64
ordinal_type
nqp = smolyak_operator.point_size();
65
ordinal_type
npc = product_basis->
size
();
66
67
// Compute quad points, weights, values
68
quad_points
.resize(nqp);
69
quad_weights
.resize(nqp);
70
quad_values
.resize(nqp);
71
typedef
typename
smolyak_operator_type::const_set_iterator const_iterator;
72
ordinal_type
i = 0;
73
for
(const_iterator it = smolyak_operator.set_begin();
74
it != smolyak_operator.set_end(); ++it, ++i) {
75
quad_points
[i] = it->first.getTerm();
76
quad_weights
[i] = it->second.first;
77
quad_values
[i].resize(npc);
78
product_basis->evaluateBases(
quad_points
[i],
quad_values
[i]);
79
}
80
}
81
82
template
<
typename
ordinal_type,
typename
value_type,
typename
po
int
_compare_type>
83
const
Teuchos::Array< Teuchos::Array<value_type> >&
84
Stokhos::SmolyakSparseGridQuadrature<ordinal_type, value_type, point_compare_type>::
85
getQuadPoints
()
const
86
{
87
return
quad_points
;
88
}
89
90
template
<
typename
ordinal_type,
typename
value_type,
typename
po
int
_compare_type>
91
const
Teuchos::Array<value_type>&
92
Stokhos::SmolyakSparseGridQuadrature<ordinal_type, value_type, point_compare_type>::
93
getQuadWeights
()
const
94
{
95
return
quad_weights
;
96
}
97
98
template
<
typename
ordinal_type,
typename
value_type,
typename
po
int
_compare_type>
99
const
Teuchos::Array< Teuchos::Array<value_type> >&
100
Stokhos::SmolyakSparseGridQuadrature<ordinal_type, value_type, point_compare_type>::
101
getBasisAtQuadPoints
()
const
102
{
103
return
quad_values
;
104
}
105
106
template
<
typename
ordinal_type,
typename
value_type,
typename
po
int
_compare_type>
107
std::ostream&
108
Stokhos::SmolyakSparseGridQuadrature<ordinal_type, value_type, point_compare_type>::
109
print
(std::ostream& os)
const
110
{
111
ordinal_type
nqp =
quad_weights
.size();
112
os <<
"Smolyak Sparse Grid Quadrature with "
<< nqp <<
" points:"
113
<< std::endl <<
"Weight : Points"
<< std::endl;
114
for
(
ordinal_type
i=0; i<nqp; i++) {
115
os << i <<
": "
<<
quad_weights
[i] <<
" : "
;
116
for
(
ordinal_type
j
=0;
j<static_cast<ordinal_type>
(
quad_points
[i].
size
());
117
j
++)
118
os <<
quad_points
[i][
j
] <<
" "
;
119
os << std::endl;
120
}
121
os <<
"Basis values at quadrature points:"
<< std::endl;
122
for
(
ordinal_type
i=0; i<nqp; i++) {
123
os << i <<
" "
<<
": "
;
124
for
(
ordinal_type
j
=0;
j<static_cast<ordinal_type>
(
quad_values
[i].
size
());
125
j
++)
126
os <<
quad_values
[i][
j
] <<
" "
;
127
os << std::endl;
128
}
129
130
return
os;
131
}
j
j
Definition
Sacado_Fad_Exp_MP_Vector.hpp:527
Stokhos::ProductBasis
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
Definition
Stokhos_ProductBasis.hpp:66
Stokhos::value_type
Stokhos::ProductContainer::size
ordinal_type size() const
Return size.
Definition
Stokhos_ProductContainerImp.hpp:139
Stokhos::SmolyakBasis
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
Definition
Stokhos_SmolyakBasis.hpp:63
Stokhos::SmolyakPseudoSpectralOperator
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
Definition
Stokhos_SmolyakPseudoSpectralOperator.hpp:62
Stokhos::SmolyakSparseGridQuadrature::size
virtual ordinal_type size() const
Get number of quadrature points.
Definition
Stokhos_SmolyakSparseGridQuadrature.hpp:84
Stokhos::SmolyakSparseGridQuadrature::SmolyakSparseGridQuadrature
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.
Definition
Stokhos_SmolyakSparseGridQuadratureImp.hpp:47
Stokhos::SmolyakSparseGridQuadrature::getBasisAtQuadPoints
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const
Get values of basis at quadrature points.
Definition
Stokhos_SmolyakSparseGridQuadratureImp.hpp:101
Stokhos::SmolyakSparseGridQuadrature::quad_points
Teuchos::Array< Teuchos::Array< value_type > > quad_points
Quadrature points.
Definition
Stokhos_SmolyakSparseGridQuadrature.hpp:123
Stokhos::SmolyakSparseGridQuadrature::getQuadPoints
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getQuadPoints() const
Get quadrature points.
Definition
Stokhos_SmolyakSparseGridQuadratureImp.hpp:85
Stokhos::SmolyakSparseGridQuadrature::getQuadWeights
virtual const Teuchos::Array< value_type > & getQuadWeights() const
Get quadrature weights.
Definition
Stokhos_SmolyakSparseGridQuadratureImp.hpp:93
Stokhos::SmolyakSparseGridQuadrature::quad_weights
Teuchos::Array< value_type > quad_weights
Quadrature weights.
Definition
Stokhos_SmolyakSparseGridQuadrature.hpp:126
Stokhos::SmolyakSparseGridQuadrature::quad_values
Teuchos::Array< Teuchos::Array< value_type > > quad_values
Quadrature values.
Definition
Stokhos_SmolyakSparseGridQuadrature.hpp:129
Stokhos::SmolyakSparseGridQuadrature::print
virtual std::ostream & print(std::ostream &os) const
Print quadrature data.
Definition
Stokhos_SmolyakSparseGridQuadratureImp.hpp:109
Generated by
1.17.0