10 #ifndef STOKHOS_FORUQTKORTHOGPOLYEXPANSION_HPP
11 #define STOKHOS_FORUQTKORTHOGPOLYEXPANSION_HPP
14 #ifdef HAVE_STOKHOS_FORUQTK
25 template <
typename ordinal_type,
typename value_type>
26 class ForUQTKOrthogPolyExpansion :
27 public OrthogPolyExpansionBase<ordinal_type, value_type,
28 Stokhos::StandardStorage<ordinal_type, value_type> > {
33 enum EXPANSION_METHOD {
39 ForUQTKOrthogPolyExpansion(
40 const Teuchos::RCP<
const OrthogPolyBasis<ordinal_type,value_type> >& basis,
42 EXPANSION_METHOD method = TAYLOR,
46 virtual ~ForUQTKOrthogPolyExpansion() {}
49 void timesEqual(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
51 void divideEqual(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
55 OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
56 const OrthogPolyApprox<ordinal_type, value_type, node_type>& x);
58 OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
59 const OrthogPolyApprox<ordinal_type, value_type, node_type>& x);
62 void times(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
63 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
64 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
65 void times(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
67 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
68 void times(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
69 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
71 void divide(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
72 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
73 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
74 void divide(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
76 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
77 void divide(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
78 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
81 void exp(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
82 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
83 void log(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
84 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
85 void log10(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
86 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
87 void sqrt(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
88 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
89 void cbrt(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
90 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
91 void pow(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
92 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
93 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
94 void pow(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
96 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
97 void pow(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
98 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
100 void cos(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
101 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
102 void sin(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
103 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
104 void tan(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
105 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
106 void cosh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
107 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
108 void sinh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
109 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
110 void tanh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
111 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
112 void acos(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
113 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
114 void asin(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
115 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
116 void atan(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
117 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
118 void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
119 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
120 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
121 void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
123 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
124 void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
125 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
127 void acosh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
128 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
129 void asinh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
130 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
131 void atanh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
132 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
137 ForUQTKOrthogPolyExpansion(
const ForUQTKOrthogPolyExpansion&);
140 ForUQTKOrthogPolyExpansion& operator=(
const ForUQTKOrthogPolyExpansion& b);
157 EXPANSION_METHOD method;
165 #endif // HAVE_STOKHOS_FORUQTK
167 #endif // STOKHOS_FORUQTKORTHOGPOLYEXPANSION_HPP
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > tan(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > sinh(const PCE< Storage > &a)
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
KOKKOS_INLINE_FUNCTION PCE< Storage > tanh(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > cbrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > acos(const PCE< Storage > &a)
atan2(expr1.val(), expr2.val())
KOKKOS_INLINE_FUNCTION PCE< Storage > cosh(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > atan(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > log10(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > asin(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > cos(const PCE< Storage > &a)