42#ifndef SACADO_PCE_ORTHOGPOLY_HPP
43#define SACADO_PCE_ORTHOGPOLY_HPP
47#ifdef HAVE_STOKHOS_SACADO
49#include "Teuchos_RCP.hpp"
51#include "Sacado_Traits.hpp"
52#include "Sacado_Handle.hpp"
53#include "Sacado_mpl_apply.hpp"
73 template <
typename T,
typename Storage >
90 typedef Stokhos::OrthogPolyBasis<ordinal_type,T>
basis_type;
93 typedef Stokhos::OrthogPolyExpansion<ordinal_type,T,Storage> expansion_type;
96 typedef Stokhos::OrthogPolyApprox<ordinal_type,T,Storage> approx_type;
98 typedef typename approx_type::pointer pointer;
99 typedef typename approx_type::const_pointer const_pointer;
100 typedef typename approx_type::reference reference;
101 typedef typename approx_type::const_reference const_reference;
104 template <
typename S>
106 typedef typename Sacado::mpl::apply<Storage,ordinal_type,S>::type
storage_type;
107 typedef OrthogPoly<S,storage_type> type;
120 OrthogPoly(
const value_type& x);
126 OrthogPoly(
const Teuchos::RCP<expansion_type>& expansion);
132 OrthogPoly(
const Teuchos::RCP<expansion_type>& expansion,
136 OrthogPoly(
const OrthogPoly& x);
142 void init(
const T& v) { th->init(v); }
145 void init(
const T* v) { th->init(v); }
148 template <
typename S>
149 void init(
const OrthogPoly<T,S>& v) { th->init(v.getOrthogPolyApprox()); }
152 void load(T* v) { th->load(v); }
155 template <
typename S>
156 void load(OrthogPoly<T,S>& v) { th->load(v.getOrthogPolyApprox()); }
162 void reset(
const Teuchos::RCP<expansion_type>& expansion);
168 void reset(
const Teuchos::RCP<expansion_type>& expansion,
181 void copyForWrite() { th.makeOwnCopy(); }
184 value_type evaluate(
const Teuchos::Array<value_type>& point)
const;
187 value_type evaluate(
const Teuchos::Array<value_type>& point,
188 const Teuchos::Array<value_type>& bvals)
const;
194 value_type standard_deviation()
const {
return th->standard_deviation(); }
197 value_type two_norm()
const {
return th->two_norm(); }
200 value_type two_norm_squared()
const {
return th->two_norm_squared(); }
203 value_type inner_product(
const OrthogPoly& b)
const {
204 return th->inner_product(b.getOrthogPolyApprox()); }
207 std::ostream& print(std::ostream& os)
const {
return th->print(os); }
210 bool isEqualTo(
const OrthogPoly& x)
const;
218 OrthogPoly<T,Storage>& operator=(
const value_type&
val);
221 OrthogPoly<T,Storage>& operator=(
const OrthogPoly<T,Storage>& x);
231 Teuchos::RCP<const basis_type> basis()
const {
return th->basis(); }
234 Teuchos::RCP<expansion_type> expansion()
const {
return expansion_; }
244 const_reference
val()
const {
return (*th)[0]; }
247 reference
val() {
return (*th)[0]; }
260 bool hasFastAccess(ordinal_type sz)
const {
return th->size()>=sz;}
263 const_pointer coeff()
const {
return th->coeff();}
266 pointer coeff() {
return th->coeff();}
279 reference term(ordinal_type dimension, ordinal_type order) {
280 return th->term(dimension, order); }
283 const_reference term(ordinal_type dimension, ordinal_type order)
const {
284 return th->term(dimension, order); }
287 Teuchos::Array<ordinal_type> order(ordinal_type term)
const {
288 return th->order(term); }
304 OrthogPoly<T,Storage>& operator += (
const value_type& x);
307 OrthogPoly<T,Storage>& operator -= (
const value_type& x);
310 OrthogPoly<T,Storage>& operator *= (
const value_type& x);
313 OrthogPoly<T,Storage>& operator /= (
const value_type& x);
316 OrthogPoly<T,Storage>& operator += (
const OrthogPoly<T,Storage>& x);
319 OrthogPoly<T,Storage>& operator -= (
const OrthogPoly<T,Storage>& x);
322 OrthogPoly<T,Storage>& operator *= (
const OrthogPoly<T,Storage>& x);
325 OrthogPoly<T,Storage>& operator /= (
const OrthogPoly<T,Storage>& x);
330 const approx_type& getOrthogPolyApprox()
const {
return *th; }
333 approx_type& getOrthogPolyApprox() {
return *th; }
338 Teuchos::RCP<expansion_type> expansion_;
341 Teuchos::RCP<expansion_type> const_expansion_;
343 Sacado::Handle< Stokhos::OrthogPolyApprox<int,value_type,Storage> > th;
398 template <
typename T,
typename Storage>
void
496 template <
typename T,
typename Storage>
bool
500 template <
typename T,
typename Storage>
bool
504 template <
typename T,
typename Storage>
bool
508 template <
typename T,
typename Storage>
bool
512 template <
typename T,
typename Storage>
bool
516 template <
typename T,
typename Storage>
bool
520 template <
typename T,
typename Storage>
bool
524 template <
typename T,
typename Storage>
bool
528 template <
typename T,
typename Storage>
bool
532 template <
typename T,
typename Storage>
bool
536 template <
typename T,
typename Storage>
bool
540 template <
typename T,
typename Storage>
bool
544 template <
typename T,
typename Storage>
bool
548 template <
typename T,
typename Storage>
bool
552 template <
typename T,
typename Storage>
bool
556 template <
typename T,
typename Storage>
bool
560 template <
typename T,
typename Storage>
bool
564 template <
typename T,
typename Storage>
bool
568 template <
typename T,
typename Storage> std::ostream&
571 template <
typename T,
typename Storage> std::istream&
expr1 expr1 expr1 expr2 expr1 expr1 c expr2 expr1 c fastAccessCoeff(j) - expr2.val(j)
atan2(expr1.val(), expr2.val())
Stokhos::StandardStorage< int, double > storage_type
Stokhos::LegendreBasis< int, double > basis_type
Stokhos::StandardStorage< int, double > Storage
std::istream & operator>>(std::istream &is, OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > cosh(const OrthogPoly< T, Storage > &a)
bool operator>=(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > operator+(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > operator/(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > tan(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > asin(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > sqrt(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > abs(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > exp(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > cos(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > fabs(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > tanh(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > cbrt(const OrthogPoly< T, Storage > &a)
bool operator!=(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > acos(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > operator-(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
bool operator==(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > log(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > log10(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > atan(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > acosh(const OrthogPoly< T, Storage > &a)
bool operator<=(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
std::ostream & operator<<(std::ostream &os, const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > min(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > max(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
bool operator<(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > sin(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > pow(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > sinh(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > atanh(const OrthogPoly< T, Storage > &a)
OrthogPoly< T, Storage > operator*(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
bool operator>(const OrthogPoly< T, Storage > &a, const OrthogPoly< T, Storage > &b)
OrthogPoly< T, Storage > asinh(const OrthogPoly< T, Storage > &a)