18 namespace ForUQTKExpansionUnitTest {
21 template <
typename OrdinalType,
typename ValueType>
30 Stokhos::OrthogPolyApprox<OrdinalType,ValueType> x,
y,
u,
u2,
cx,
cu,
cu2,
sx,
su,
su2;
39 const OrdinalType d = 2;
40 const OrdinalType p = 7;
44 for (OrdinalType i=0; i<d; i++)
59 Teuchos::rcp(
new Stokhos::ForUQTKOrthogPolyExpansion<OrdinalType,ValueType>(
basis,
Cijk, Stokhos::ForUQTKOrthogPolyExpansion<OrdinalType,ValueType>::INTEGRATION));
75 for (OrdinalType i=0; i<d; i++) {
80 for (OrdinalType i=0; i<d; i++)
94 OrdinalType nqp = weights.
size();
97 for (OrdinalType i=0; i<c.
size(); i++)
102 for (OrdinalType k=0; k<nqp; k++) {
105 for (
int i=0; i<c.
size(); i++)
110 template <
class Func>
121 OrdinalType nqp = weights.
size();
124 for (OrdinalType i=0; i<c.
size(); i++)
129 for (OrdinalType k=0; k<nqp; k++) {
130 ValueType val1 = a.
evaluate(points[k], values[k]);
131 ValueType val2 = b.
evaluate(points[k], values[k]);
132 ValueType
val = func(val1, val2);
133 for (
int i=0; i<c.
size(); i++)
138 template <
class Func>
150 OrdinalType nqp = weights.
size();
153 for (OrdinalType i=0; i<c.
size(); i++)
158 for (OrdinalType k=0; k<nqp; k++) {
159 ValueType val2 = b.
evaluate(points[k], values[k]);
160 ValueType
val = func(a, val2);
161 for (
int i=0; i<c.
size(); i++)
166 template <
class Func>
178 OrdinalType nqp = weights.
size();
181 for (OrdinalType i=0; i<c.
size(); i++)
186 for (OrdinalType k=0; k<nqp; k++) {
187 ValueType val1 = a.
evaluate(points[k], values[k]);
188 ValueType
val = func(val1, b);
189 for (
int i=0; i<c.
size(); i++)
255 return 0.5*
std::log((1.0+a)/(1.0-a));
260 double operator() (
double a,
double b)
const {
return a + b; }
263 double operator() (
double a,
double b)
const {
return a - b; }
266 double operator() (
double a,
double b)
const {
return a * b; }
269 double operator() (
double a,
double b)
const {
return a / b; }
787 setup.exp->plus(ru, v, w);
903 setup.exp->minus(ru, v, w);
1019 setup.exp->times(ru, v, w);
1087 setup.exp->divide(ru, v, w);
1196 setup.exp->plusEqual(ru, v);
1249 setup.exp->minusEqual(ru, v);
1250 setup.exp->unaryMinus(v, v);
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > u
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > tan(const PCE< Storage > &a)
double operator()(double a, double b) const
double operator()(double a) const
double operator()(double a, double b) const
double operator()(double a) const
KOKKOS_INLINE_FUNCTION PCE< Storage > sinh(const PCE< Storage > &a)
value_type evaluate(const Teuchos::Array< value_type > &point) const
Evaluate polynomial approximation at a point.
Teuchos::RCP< Stokhos::Sparse3Tensor< int, double > > Cijk
UnitTestSetup< int, double > setup
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > cu2
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
double operator()(double a) const
bool comparePCEs(const PCEType &a1, const std::string &a1_name, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &a2, const std::string &a2_name, const ValueType &rel_tol, const ValueType &abs_tol, Teuchos::FancyOStream &out)
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const =0
Get values of basis at quadrature points.
double operator()(double a) const
double operator()(double a, double b) const
double operator()(double a, double b) const
KOKKOS_INLINE_FUNCTION PCE< Storage > tanh(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > cbrt(const PCE< Storage > &a)
double operator()(double a) const
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > su2
double operator()(double a) const
TEUCHOS_UNIT_TEST(Stokhos_ForUQTKExpansion, UMinus)
KOKKOS_INLINE_FUNCTION PCE< Storage > acos(const PCE< Storage > &a)
double operator()(double a) const
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > cu
virtual const Teuchos::Array< value_type > & getQuadWeights() const =0
Get quadrature weights.
double operator()(double a) const
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getQuadPoints() const =0
Get quadrature points.
static int runUnitTestsFromMain(int argc, char *argv[])
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > sx
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > su
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void reset(const Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > &new_basis, ordinal_type sz=0)
Reset to a new basis.
double operator()(double a) const
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
KOKKOS_INLINE_FUNCTION PCE< Storage > cosh(const PCE< Storage > &a)
Teuchos::RCP< Stokhos::ForUQTKOrthogPolyExpansion< OrdinalType, ValueType > > exp
double operator()(double a) const
double operator()(double a) const
double operator()(double a, double b) const
void computePCE1(Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &c, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > atan(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
double operator()(double a) const
void computePCE2LC(Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &c, ValueType a, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &b)
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > cx
Teuchos::RCP< const Stokhos::CompletePolynomialBasis< OrdinalType, ValueType > > basis
int main(int argc, char **argv)
double operator()(double a) const
Teuchos::RCP< const Stokhos::Quadrature< OrdinalType, ValueType > > quad
void computePCE2RC(Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &c, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &a, ValueType b)
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
double operator()(double a) const
void computePCE2(Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &c, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &a, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &b)
virtual ordinal_type size() const
Return total size of basis.
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > x
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
ordinal_type size() const
Return size.
double operator()(double a) const
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > y
double operator()(double a) const
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > log10(const PCE< Storage > &a)
double operator()(double a) const
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > u2
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
KOKKOS_INLINE_FUNCTION PCE< Storage > asin(const PCE< Storage > &a)
double operator()(double a) const
KOKKOS_INLINE_FUNCTION PCE< Storage > cos(const PCE< Storage > &a)
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...
reference term(ordinal_type dimension, ordinal_type order)
Get coefficient term for given dimension and order.