42 #ifndef STOKHOS_FORUQTKORTHOGPOLYEXPANSION_HPP
43 #define STOKHOS_FORUQTKORTHOGPOLYEXPANSION_HPP
46 #ifdef HAVE_STOKHOS_FORUQTK
57 template <
typename ordinal_type,
typename value_type>
58 class ForUQTKOrthogPolyExpansion :
59 public OrthogPolyExpansionBase<ordinal_type, value_type,
60 Stokhos::StandardStorage<ordinal_type, value_type> > {
65 enum EXPANSION_METHOD {
71 ForUQTKOrthogPolyExpansion(
72 const Teuchos::RCP<
const OrthogPolyBasis<ordinal_type,value_type> >& basis,
74 EXPANSION_METHOD method = TAYLOR,
78 virtual ~ForUQTKOrthogPolyExpansion() {}
81 void timesEqual(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
83 void divideEqual(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
87 OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
88 const OrthogPolyApprox<ordinal_type, value_type, node_type>& x);
90 OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
91 const OrthogPolyApprox<ordinal_type, value_type, node_type>& x);
94 void times(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
95 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
96 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
97 void times(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
99 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
100 void times(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
101 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
103 void divide(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
104 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
105 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
106 void divide(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
108 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
109 void divide(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
110 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
113 void exp(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
114 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
115 void log(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
116 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
117 void log10(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
118 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
119 void sqrt(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
120 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
121 void cbrt(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
122 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
123 void pow(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
124 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
125 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
126 void pow(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
128 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
129 void pow(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
130 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
132 void cos(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
133 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
134 void sin(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
135 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
136 void tan(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
137 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
138 void cosh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
139 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
140 void sinh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
141 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
142 void tanh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
143 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
144 void acos(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
145 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
146 void asin(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
147 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
148 void atan(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
149 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
150 void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
151 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
152 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
153 void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
155 const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
156 void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
157 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
159 void acosh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
160 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
161 void asinh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
162 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
163 void atanh(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
164 const OrthogPolyApprox<ordinal_type, value_type, node_type>& a);
169 ForUQTKOrthogPolyExpansion(
const ForUQTKOrthogPolyExpansion&);
172 ForUQTKOrthogPolyExpansion& operator=(
const ForUQTKOrthogPolyExpansion& b);
189 EXPANSION_METHOD method;
197 #endif // HAVE_STOKHOS_FORUQTK
199 #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)