Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_GaussPattersonLegendreBasisImp.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 GaussPattersonLegendreBasis(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_gp);
23 #endif
24 }
25 
26 template <typename ordinal_type, typename value_type>
29  const GaussPattersonLegendreBasis& 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  // Gauss-Patterson points have the following structure
51  // (cf. http://people.sc.fsu.edu/~jburkardt/f_src/patterson_rule/patterson_rule.html):
52  // Level l Num points n Precision p
53  // -----------------------------------
54  // 0 1 1
55  // 1 3 5
56  // 2 7 11
57  // 3 15 23
58  // 4 31 47
59  // 5 63 95
60  // 6 127 191
61  // 7 255 383
62  // Thus for l > 0, n = 2^{l+1}-1 and p = 3*2^l-1. So for a given quadrature
63  // order p, we find the smallest l s.t. 3*s^l-1 >= p and then compute the
64  // number of points n from the above. In this case, l = ceil(log2((p+1)/3))
65  ordinal_type num_points;
66  if (quad_order <= ordinal_type(1))
67  num_points = 1;
68  else {
69  ordinal_type l = std::ceil(std::log((quad_order+1.0)/3.0)/std::log(2.0));
70  num_points = (1 << (l+1)) - 1; // std::pow(2,l+1)-1;
71  }
72 
73  quad_points.resize(num_points);
74  quad_weights.resize(num_points);
75  quad_values.resize(num_points);
76 
77  webbur::patterson_lookup(num_points, &quad_points[0], &quad_weights[0]);
78 
79  for (ordinal_type i=0; i<num_points; i++) {
80  quad_weights[i] *= 0.5; // scale to unit measure
81  quad_values[i].resize(this->p+1);
82  this->evaluateBases(quad_points[i], quad_values[i]);
83  }
84 
85 #else
87  true, std::logic_error, "Clenshaw-Curtis requires TriKota to be enabled!");
88 #endif
89 }
90 
91 template <typename ordinal_type, typename value_type>
95 {
96  // Based on the above structure, we find the largest l s.t. 2^{l+1}-1 <= n,
97  // which is floor(log2(n+1)-1) and compute p = 3*2^l-1
98  if (n == ordinal_type(1))
99  return 1;
100  ordinal_type l = std::floor(std::log(n+1.0)/std::log(2.0)-1.0);
101  return (3 << l) - 1; // 3*std::pow(2,l)-1;
102 }
103 
104 template <typename ordinal_type, typename value_type>
108 {
109  return
111 }
112 
113 template <typename ordinal_type, typename value_type>
117 {
118  // Gauss-Patterson rules have precision 3*2^l-1, which is odd.
119  // Since discrete orthogonality requires integrating polynomials of
120  // order 2*p, setting p = 3*2^{l-1}-1 will yield the largest p such that
121  // 2*p <= 3*2^l-1
122  if (n == 0)
123  return 0;
124  return (3 << (n-1)) - 1; // 3*std::pow(2,n-1) - 1;
125 }
126 
127 template <typename ordinal_type, typename value_type>
131 {
132  return n;
133 }
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.
virtual ordinal_type pointGrowth(ordinal_type n) const
Evaluate point growth rule for Smolyak-type bases.
#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.
Legendre polynomial basis using Gauss-Patterson quadrature points.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
KOKKOS_INLINE_FUNCTION PCE< Storage > ceil(const PCE< Storage > &a)
void resize(size_type new_size, const value_type &x=value_type())
Legendre polynomial basis.
GaussPattersonLegendreBasis(ordinal_type p, bool normalize=false, bool isotropic=false)
Constructor.
virtual void setSparseGridGrowthRule(LevelToOrderFnPtr ptr)
Set sparse grid rule.
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
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...