Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_QuadOrthogPolyExpansion.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 #ifndef STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
11 #define STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
12 
14 #include "Stokhos_Quadrature.hpp"
15 
16 #include "Teuchos_RCP.hpp"
17 #include "Teuchos_Array.hpp"
20 #include "Teuchos_BLAS.hpp"
21 
22 namespace Stokhos {
23 
25  template <typename ordinal_type, typename value_type,
27  value_type> >
29  public OrthogPolyExpansionBase<ordinal_type, value_type, node_type> {
30  public:
31 
38 
41 
42  void timesEqual(
44  const value_type& x);
45  void divideEqual(
47  const value_type& x);
48 
49  void timesEqual(
52  void divideEqual(
55 
56 
61  const value_type& a,
65  const value_type& b);
70  const value_type& a,
74  const value_type& b);
75 
90  const value_type& a,
94  const value_type& b);
117  const value_type& a,
121  const value_type& b);
128 
129  template <typename FuncT>
130  void nary_op(const FuncT& func,
133 
134  template <typename ExprT1, typename ExprT2>
135  value_type compute_times_coeff(ordinal_type k, const ExprT1& a,
136  const ExprT2& b) const;
137 
138  template <typename ExprT1, typename ExprT2>
140  const ExprT2& b) const;
141 
142  private:
143 
144  // Prohibit copying
146 
147  // Prohibit Assignment
149 
150  protected:
151 
154 
157 
160 
163 
166 
169 
172 
175 
178 
181 
184 
187 
190 
193 
196 
199 
202 
203  public:
204 
206  template <typename FuncT>
207  void unary_op(
208  const FuncT& func,
211 
213  template <typename FuncT>
214  void binary_op(
215  const FuncT& func,
219 
221  template <typename FuncT>
222  void binary_op(
223  const FuncT& func,
225  const value_type& a,
227 
229  template <typename FuncT>
230  void binary_op(
231  const FuncT& func,
234  const value_type& b);
235 
236  protected:
237 
238  struct times_quad_func {
239  value_type operator() (const value_type& a, const value_type& b) const {
240  return a * b;
241  }
242  };
243 
244  struct div_quad_func {
245  value_type operator() (const value_type& a, const value_type& b) const {
246  return a / b;
247  }
248  };
249 
250  struct exp_quad_func {
251  value_type operator() (const value_type& a) const {
252  return std::exp(a);
253  }
254  };
255 
256  struct log_quad_func {
257  value_type operator() (const value_type& a) const {
258  return std::log(a);
259  }
260  };
261 
262  struct log10_quad_func {
263  value_type operator() (const value_type& a) const {
264  return std::log10(a);
265  }
266  };
267 
268  struct sqrt_quad_func {
269  value_type operator() (const value_type& a) const {
270  return std::sqrt(a);
271  }
272  };
273 
274  struct cbrt_quad_func {
275  value_type operator() (const value_type& a) const {
276  return std::cbrt(a);
277  }
278  };
279 
280  struct pow_quad_func {
281  value_type operator() (const value_type& a, const value_type& b) const {
282  return std::pow(a,b);
283  }
284  };
285 
286  struct cos_quad_func {
287  value_type operator() (const value_type& a) const {
288  return std::cos(a);
289  }
290  };
291 
292  struct sin_quad_func {
293  value_type operator() (const value_type& a) const {
294  return std::sin(a);
295  }
296  };
297 
298  struct tan_quad_func {
299  value_type operator() (const value_type& a) const {
300  return std::tan(a);
301  }
302  };
303 
304  struct cosh_quad_func {
305  value_type operator() (const value_type& a) const {
306  return std::cosh(a);
307  }
308  };
309 
310  struct sinh_quad_func {
311  value_type operator() (const value_type& a) const {
312  return std::sinh(a);
313  }
314  };
315 
316  struct tanh_quad_func {
317  value_type operator() (const value_type& a) const {
318  return std::tanh(a);
319  }
320  };
321 
322  struct acos_quad_func {
323  value_type operator() (const value_type& a) const {
324  return std::acos(a);
325  }
326  };
327 
328  struct asin_quad_func {
329  value_type operator() (const value_type& a) const {
330  return std::asin(a);
331  }
332  };
333 
334  struct atan_quad_func {
335  value_type operator() (const value_type& a) const {
336  return std::atan(a);
337  }
338  };
339 
340  struct atan2_quad_func {
341  value_type operator() (const value_type& a, const value_type& b) const {
342  return std::atan2(a,b);
343  }
344  };
345 
346  struct acosh_quad_func {
347  value_type operator() (const value_type & a) const {
348  return std::log(a+std::sqrt(a*a-value_type(1.0)));
349  }
350  };
351 
352  struct asinh_quad_func {
353  value_type operator() (const value_type& a) const {
354  return std::log(a+std::sqrt(a*a+value_type(1.0)));
355  }
356  };
357 
358  struct atanh_quad_func {
359  value_type operator() (const value_type& a) const {
360  return 0.5*std::log((value_type(1.0)+a)/(value_type(1.0)-a));
361  }
362  };
363 
364  }; // class QuadOrthogPolyExpansion
365 
366 } // namespace Stokhos
367 
369 
370 #endif // STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > tan(const PCE< Storage > &a)
const Teuchos::Array< value_type > & norms
Norms of basis vectors.
const Teuchos::Array< Teuchos::Array< value_type > > & quad_points
Array of Quad points.
void binary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Nonlinear binary function.
Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > basis
Basis.
OrthogPolyExpansionBase< ordinal_type, value_type, node_type >::Cijk_type Cijk_type
Short-hand for Cijk.
Teuchos::RCP< Teuchos::ParameterList > params
Parameter list.
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void times(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
KOKKOS_INLINE_FUNCTION PCE< Storage > sinh(const PCE< Storage > &a)
void sqrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
ordinal_type nqp
Number of Quad points.
value_type operator()(const value_type &a, const value_type &b) const
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
void acos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
const Teuchos::Array< Teuchos::Array< value_type > > & quad_values
Values of basis at Quad points.
QuadOrthogPolyExpansion(const Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > &basis, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::RCP< const Quadrature< ordinal_type, value_type > > &quad, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
value_type operator()(const value_type &a, const value_type &b) const
Teuchos::Array< Teuchos::Array< Teuchos::Array< value_type > > > navals
Temporary array for values of n-ary arguments at quad points.
void atanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void log10(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void cos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Base class for consolidating common expansion implementations.
void asin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void atan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
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::Serial node_type
void divide(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Teuchos::RCP< const Quadrature< ordinal_type, value_type > > quad
Quadrature routine.
value_type compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
value_type operator()(const value_type &a, const value_type &b) const
Abstract base class for multivariate orthogonal polynomials.
void sinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
value_type fast_compute_times_coeff(ordinal_type k, const ExprT1 &a, const ExprT2 &b) const
Abstract base class for quadrature methods.
Teuchos::Array< value_type > fvals
Temporary array for values of operation at quad points.
void atan2(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void nary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > **a)
value_type operator()(const value_type &a, const value_type &b) const
bool use_quad_for_times
Use quadrature for times functions.
KOKKOS_INLINE_FUNCTION PCE< Storage > cosh(const PCE< Storage > &a)
Teuchos::Array< value_type > sqv
Quad values scaled by norms.
QuadOrthogPolyExpansion & operator=(const QuadOrthogPolyExpansion &b)
KOKKOS_INLINE_FUNCTION PCE< Storage > atan(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
void log(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
bool use_quad_for_division
Use quadrature for division functions.
Class to store coefficients of a projection onto an orthogonal polynomial basis.
void timesEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void cbrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
const Teuchos::Array< value_type > & quad_weights
Array of Quad weights.
Teuchos::Array< value_type > bvals
Temporary array for values of second argument at quad points.
void pow(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
void tanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void tan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::Array< value_type > avals
Temporary array for values of first argument at quad points.
void cosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > Cijk
Triple-product tensor.
Orthogonal polynomial expansions based on numerical quadrature.
void asinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &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)
Teuchos::Array< value_type > qv
Reshaped quad values into 1D array.
KOKKOS_INLINE_FUNCTION PCE< Storage > cos(const PCE< Storage > &a)
void unary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Nonlinear unary function.
Teuchos::BLAS< ordinal_type, value_type > blas
BLAS wrappers.
void acosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)