Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_DerivOrthogPolyExpansion.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_DERIVORTHOGPOLYEXPANSION_HPP
45 #define STOKHOS_DERIVORTHOGPOLYEXPANSION_HPP
46 
48 #include "Stokhos_DerivBasis.hpp"
49 
50 #include "Teuchos_RCP.hpp"
51 #include "Teuchos_Array.hpp"
54 #include "Teuchos_LAPACK.hpp"
55 
56 namespace Stokhos {
57 
59  template <typename ordinal_type, typename value_type>
60  class DerivOrthogPolyExpansion : public OrthogPolyExpansion<ordinal_type, value_type> {
61  public:
62 
64 
71 
74 
76  ordinal_type size() const { return sz; }
77 
80  getBasis() const {return basis; }
81 
84  getTripleProduct() const { return Cijk; }
85 
86  // Operations
87  void unaryMinus(
90 
92  const value_type& x);
94  const value_type& x);
96  const value_type& x);
98  const value_type& x);
99 
100  void plusEqual(
103  void minusEqual(
106  void timesEqual(
109  void divideEqual(
112 
117  const value_type& a,
121  const value_type& b);
126  const value_type& a,
130  const value_type& b);
135  const value_type& a,
139  const value_type& b);
144  const value_type& a,
148  const value_type& b);
149 
164  const value_type& a,
168  const value_type& b);
187  template <typename OpT>
188  void quad(const OpT& quad_func,
198 // void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
199 // const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
200 // const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
201 // void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
202 // const T& a,
203 // const OrthogPolyApprox<ordinal_type, value_type, node_type>& b);
204 // void atan2(OrthogPolyApprox<ordinal_type, value_type, node_type>& c,
205 // const OrthogPolyApprox<ordinal_type, value_type, node_type>& a,
206 // const T& b);
221  const value_type& a,
225  const value_type& b);
230  const value_type& a,
234  const value_type& b);
237 
238  private:
239 
240  // Prohibit copying
242 
243  // Prohibit Assignment
245 
246  protected:
247 
250 
253 
256 
259 
262 
265 
268 
271 
274 
275  protected:
276 
279 
280  struct acos_quad_func {
281  value_type operator() (const value_type& a) const {
282  return std::acos(a);
283  }
284  };
285 
286  struct asin_quad_func {
287  value_type operator() (const value_type& a) const {
288  return std::asin(a);
289  }
290  };
291 
292  struct atan_quad_func {
293  value_type operator() (const value_type& a) const {
294  return std::atan(a);
295  }
296  };
297 
298  struct acosh_quad_func {
299  value_type operator() (const value_type & a) const {
300  return std::log(a+std::sqrt(a*a-value_type(1.0)));
301  }
302  };
303 
304  struct asinh_quad_func {
305  value_type operator() (const value_type& a) const {
306  return std::log(a+std::sqrt(a*a+value_type(1.0)));
307  }
308  };
309 
310  struct atanh_quad_func {
311  value_type operator() (const value_type& a) const {
312  return 0.5*std::log((value_type(1.0)+a)/(value_type(1.0)-a));
313  }
314  };
315 
316  }; // class DerivOrthogPolyExpansion
317 
318 } // namespace Stokhos
319 
321 
322 #endif // STOKHOS_DERIVORTHOGPOLYEXPANSION_HPP
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
void cbrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
DerivOrthogPolyExpansion(const Teuchos::RCP< const DerivBasis< ordinal_type, value_type > > &basis, const Teuchos::RCP< const Teuchos::SerialDenseMatrix< ordinal_type, value_type > > &Bij, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::RCP< const Stokhos::Dense3Tensor< ordinal_type, value_type > > &Dijk)
Constructor.
void acos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Abstract base class for multivariate orthogonal polynomials that support computing double and triple ...
Teuchos::SerialDenseMatrix< ordinal_type, value_type > A
Matrix.
virtual Teuchos::RCP< const Sparse3Tensor< ordinal_type, value_type > > getTripleProduct() const
Get triple product.
void minusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > getBasis() const
Get basis.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
void quad(const OpT &quad_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)
ordinal_type size() const
Get expansion size.
void cos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void plus(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 timesEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void asin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void cosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::RCP< const Teuchos::SerialDenseMatrix< ordinal_type, value_type > > Bij
Derivative double-product tensor.
KOKKOS_INLINE_FUNCTION PCE< Storage > acos(const PCE< Storage > &a)
void unaryMinus(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Abstract base class for orthogonal polynomial-based expansions.
Teuchos::RCP< const Stokhos::DerivBasis< ordinal_type, value_type > > basis
Basis.
void acosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void derivative(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void tan(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 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)
void min(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)
DerivOrthogPolyExpansion & operator=(const DerivOrthogPolyExpansion &b)
void atan(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::RCP< const Stokhos::Dense3Tensor< ordinal_type, value_type > > Dijk
Derivative Triple-product tensor.
void sinhcosh(OrthogPolyApprox< ordinal_type, value_type, node_type > &s, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void sqrt(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Othogonal polynomial expansions based on derivative calculations.
void minus(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 Stokhos::Sparse3Tensor< ordinal_type, value_type > > Cijk
Triple-product tensor.
KOKKOS_INLINE_FUNCTION PCE< Storage > atan(const PCE< Storage > &a)
void tanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void asinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void abs(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Class to store coefficients of a projection onto an orthogonal polynomial basis.
void sincos(OrthogPolyApprox< ordinal_type, value_type, node_type > &s, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void fabs(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::LAPACK< ordinal_type, value_type > lapack
LAPACK wrappers.
Teuchos::Array< ordinal_type > piv
Pivot array.
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)
void max(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 log(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
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)
Data structure storing a dense 3-tensor C(i,j,k).
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void atanh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Teuchos::SerialDenseMatrix< ordinal_type, value_type > B
RHS.
Stokhos::StandardStorage< ordinal_type, value_type > node_type
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
void sinh(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void plusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
KOKKOS_INLINE_FUNCTION PCE< Storage > asin(const PCE< Storage > &a)
ordinal_type solve(ordinal_type s, ordinal_type nrhs)
Solve linear system.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)