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 // $Id$
2 // $Source$
3 // @HEADER
4 // ***********************************************************************
5 //
6 // Stokhos Package
7 // Copyright (2009) Sandia Corporation
8 //
9 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
10 // license for use of this work by or on behalf of the U.S. Government.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
40 //
41 // ***********************************************************************
42 // @HEADER
43 
44 #ifndef STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
45 #define STOKHOS_QUADORTHOGPOLYEXPANSION_HPP
46 
48 #include "Stokhos_Quadrature.hpp"
49 
50 #include "Teuchos_RCP.hpp"
51 #include "Teuchos_Array.hpp"
54 #include "Teuchos_BLAS.hpp"
55 
56 namespace Stokhos {
57 
59  template <typename ordinal_type, typename value_type,
61  value_type> >
63  public OrthogPolyExpansionBase<ordinal_type, value_type, node_type> {
64  public:
65 
72 
75 
76  void timesEqual(
78  const value_type& x);
79  void divideEqual(
81  const value_type& x);
82 
83  void timesEqual(
86  void divideEqual(
89 
90 
95  const value_type& a,
99  const value_type& b);
104  const value_type& a,
108  const value_type& b);
109 
124  const value_type& a,
128  const value_type& b);
151  const value_type& a,
155  const value_type& b);
162 
163  template <typename FuncT>
164  void nary_op(const FuncT& func,
167 
168  template <typename ExprT1, typename ExprT2>
169  value_type compute_times_coeff(ordinal_type k, const ExprT1& a,
170  const ExprT2& b) const;
171 
172  template <typename ExprT1, typename ExprT2>
174  const ExprT2& b) const;
175 
176  private:
177 
178  // Prohibit copying
180 
181  // Prohibit Assignment
183 
184  protected:
185 
188 
191 
194 
197 
200 
203 
206 
209 
212 
215 
218 
221 
224 
227 
230 
233 
236 
237  public:
238 
240  template <typename FuncT>
241  void unary_op(
242  const FuncT& func,
245 
247  template <typename FuncT>
248  void binary_op(
249  const FuncT& func,
253 
255  template <typename FuncT>
256  void binary_op(
257  const FuncT& func,
259  const value_type& a,
261 
263  template <typename FuncT>
264  void binary_op(
265  const FuncT& func,
268  const value_type& b);
269 
270  protected:
271 
272  struct times_quad_func {
273  value_type operator() (const value_type& a, const value_type& b) const {
274  return a * b;
275  }
276  };
277 
278  struct div_quad_func {
279  value_type operator() (const value_type& a, const value_type& b) const {
280  return a / b;
281  }
282  };
283 
284  struct exp_quad_func {
285  value_type operator() (const value_type& a) const {
286  return std::exp(a);
287  }
288  };
289 
290  struct log_quad_func {
291  value_type operator() (const value_type& a) const {
292  return std::log(a);
293  }
294  };
295 
296  struct log10_quad_func {
297  value_type operator() (const value_type& a) const {
298  return std::log10(a);
299  }
300  };
301 
302  struct sqrt_quad_func {
303  value_type operator() (const value_type& a) const {
304  return std::sqrt(a);
305  }
306  };
307 
308  struct cbrt_quad_func {
309  value_type operator() (const value_type& a) const {
310  return std::cbrt(a);
311  }
312  };
313 
314  struct pow_quad_func {
315  value_type operator() (const value_type& a, const value_type& b) const {
316  return std::pow(a,b);
317  }
318  };
319 
320  struct cos_quad_func {
321  value_type operator() (const value_type& a) const {
322  return std::cos(a);
323  }
324  };
325 
326  struct sin_quad_func {
327  value_type operator() (const value_type& a) const {
328  return std::sin(a);
329  }
330  };
331 
332  struct tan_quad_func {
333  value_type operator() (const value_type& a) const {
334  return std::tan(a);
335  }
336  };
337 
338  struct cosh_quad_func {
339  value_type operator() (const value_type& a) const {
340  return std::cosh(a);
341  }
342  };
343 
344  struct sinh_quad_func {
345  value_type operator() (const value_type& a) const {
346  return std::sinh(a);
347  }
348  };
349 
350  struct tanh_quad_func {
351  value_type operator() (const value_type& a) const {
352  return std::tanh(a);
353  }
354  };
355 
356  struct acos_quad_func {
357  value_type operator() (const value_type& a) const {
358  return std::acos(a);
359  }
360  };
361 
362  struct asin_quad_func {
363  value_type operator() (const value_type& a) const {
364  return std::asin(a);
365  }
366  };
367 
368  struct atan_quad_func {
369  value_type operator() (const value_type& a) const {
370  return std::atan(a);
371  }
372  };
373 
374  struct atan2_quad_func {
375  value_type operator() (const value_type& a, const value_type& b) const {
376  return std::atan2(a,b);
377  }
378  };
379 
380  struct acosh_quad_func {
381  value_type operator() (const value_type & a) const {
382  return std::log(a+std::sqrt(a*a-value_type(1.0)));
383  }
384  };
385 
386  struct asinh_quad_func {
387  value_type operator() (const value_type& a) const {
388  return std::log(a+std::sqrt(a*a+value_type(1.0)));
389  }
390  };
391 
392  struct atanh_quad_func {
393  value_type operator() (const value_type& a) const {
394  return 0.5*std::log((value_type(1.0)+a)/(value_type(1.0)-a));
395  }
396  };
397 
398  }; // class QuadOrthogPolyExpansion
399 
400 } // namespace Stokhos
401 
403 
404 #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)