Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_ClenshawCurtisLegendreBasisImp.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 #ifdef HAVE_STOKHOS_DAKOTA
11 #include "sandia_rules.hpp"
12 #endif
14 
15 template <typename ordinal_type, typename value_type>
17 ClenshawCurtisLegendreBasis(ordinal_type p, bool normalize, bool isotropic_) :
18  LegendreBasis<ordinal_type, value_type>(p, normalize),
19  isotropic(isotropic_)
20 {
21 #ifdef HAVE_STOKHOS_DAKOTA
22  this->setSparseGridGrowthRule(webbur::level_to_order_exp_cc);
23 #endif
24 }
25 
26 template <typename ordinal_type, typename value_type>
29  const ClenshawCurtisLegendreBasis& basis) :
31  isotropic(basis.isotropic)
32 {
33 }
34 
35 template <typename ordinal_type, typename value_type>
38 {
39 }
40 
41 template <typename ordinal_type, typename value_type>
42 void
45  Teuchos::Array<value_type>& quad_points,
46  Teuchos::Array<value_type>& quad_weights,
47  Teuchos::Array< Teuchos::Array<value_type> >& quad_values) const
48 {
49 #ifdef HAVE_STOKHOS_DAKOTA
50  ordinal_type num_points;
51  if (quad_order % ordinal_type(2) == ordinal_type(1))
52  num_points = quad_order;
53  else
54  num_points = quad_order+1;
55  quad_points.resize(num_points);
56  quad_weights.resize(num_points);
57  quad_values.resize(num_points);
58 
59  webbur::clenshaw_curtis_compute(
60  num_points, &quad_points[0], &quad_weights[0]);
61 
62  for (ordinal_type i=0; i<num_points; i++) {
63  quad_weights[i] *= 0.5; // scale to unit measure
64  quad_values[i].resize(this->p+1);
65  this->evaluateBases(quad_points[i], quad_values[i]);
66  }
67 
68 #else
70  true, std::logic_error, "Clenshaw-Curtis requires TriKota to be enabled!");
71 #endif
72 }
73 
74 template <typename ordinal_type, typename value_type>
78 {
79  if (n % ordinal_type(2) == ordinal_type(1))
80  return n;
81  return n-ordinal_type(1);
82 }
83 
84 template <typename ordinal_type, typename value_type>
88 {
89  return
91 }
92 
93 template <typename ordinal_type, typename value_type>
97 {
98  if (n == ordinal_type(0))
99  return 0;
100  return (1 << (n-1)); // std::pow(2,n-1);
101 }
102 
103 template <typename ordinal_type, typename value_type>
107 {
108  return n;
109 }
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual ordinal_type quadDegreeOfExactness(ordinal_type n) const
virtual ordinal_type coefficientGrowth(ordinal_type n) const
Evaluate coefficient growth rule for Smolyak-type bases.
virtual Teuchos::RCP< OneDOrthogPolyBasis< ordinal_type, value_type > > cloneWithOrder(ordinal_type p) const
Clone this object with the option of building a higher order basis.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual ordinal_type pointGrowth(ordinal_type n) const
Evaluate point growth rule for Smolyak-type bases.
void resize(size_type new_size, const value_type &x=value_type())
Legendre polynomial basis.
Legendre polynomial basis using Clenshaw-Curtis quadrature points.
virtual void setSparseGridGrowthRule(LevelToOrderFnPtr ptr)
Set sparse grid rule.
ClenshawCurtisLegendreBasis(ordinal_type p, bool normalize=false, bool isotropic=false)
Constructor.
virtual void getQuadPoints(ordinal_type quad_order, Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &weights, Teuchos::Array< Teuchos::Array< value_type > > &values) const
Compute quadrature points, weights, and values of basis polynomials at given set of points points...