42#ifndef STOKHOS_ORTHOGPOLYEXPANSIONBASE_HPP
43#define STOKHOS_ORTHOGPOLYEXPANSIONBASE_HPP
48#include "Teuchos_RCP.hpp"
49#include "Teuchos_ParameterList.hpp"
58 template <
typename ordinal_type,
typename value_type,
typename node_type>
67 const Teuchos::RCP<Teuchos::ParameterList>&
params = Teuchos::null);
76 Teuchos::RCP< const OrthogPolyBasis<ordinal_type, value_type> >
80 virtual Teuchos::RCP<const Sparse3Tensor<ordinal_type, value_type> >
186 Teuchos::RCP<const OrthogPolyBasis<ordinal_type, value_type> >
basis;
192 Teuchos::RCP<const Stokhos::Sparse3Tensor<ordinal_type, value_type> >
Cijk;
195 Teuchos::RCP<Teuchos::ParameterList>
params;
198 Teuchos::RCP<Stokhos::DivisionExpansionStrategy<ordinal_type,value_type,node_type> >
division_strategy;
Class to store coefficients of a projection onto an orthogonal polynomial basis.
Abstract base class for multivariate orthogonal polynomials.
void divide(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > Cijk
void minusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
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 max(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const value_type &b)
void plusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &x)
virtual Teuchos::RCP< const Sparse3Tensor< ordinal_type, value_type > > getTripleProduct() const
Get triple product.
void timesEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &x)
OrthogPolyExpansionBase(const Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > &basis, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor.
void times(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const value_type &b)
void fabs(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void minusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &x)
void divide(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const value_type &b)
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 minus(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void divideEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &x)
Teuchos::RCP< Teuchos::ParameterList > params
void unaryMinus(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 min(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void plus(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const value_type &b)
OrthogPolyExpansionBase(const OrthogPolyExpansionBase &)
void times(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
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 plus(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
OrthogPolyExpansionBase & operator=(const OrthogPolyExpansionBase &b)
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 OrthogPolyBasis< ordinal_type, value_type > > basis
Teuchos::RCP< Stokhos::DivisionExpansionStrategy< ordinal_type, value_type, node_type > > division_strategy
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)
ordinal_type size() const
Get expansion size.
Stokhos::Sparse3Tensor< ordinal_type, value_type > Cijk_type
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 plusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
void minus(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const value_type &b)
void max(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_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 value_type &b)
Teuchos::RCP< const OrthogPolyBasis< ordinal_type, value_type > > getBasis() const
Get basis.
virtual ~OrthogPolyExpansionBase()
Destructor.
void abs(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Top-level namespace for Stokhos classes and functions.