Teuchos - Trilinos Tools Package  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions | Protected Attributes | List of all members
Teuchos::Polynomial< CoeffT > Class Template Reference

Lightweight container class to represent a simple polynomial. More...

#include <Teuchos_PolynomialDecl.hpp>

Inheritance diagram for Teuchos::Polynomial< CoeffT >:
Teuchos::Describable Teuchos::LabeledObject

Public Types

typedef CoeffT coeff_type
 Typename of coefficients. More...
 
typedef
Teuchos::PolynomialTraits
< coeff_type >::scalar_type 
scalar_type
 Typename of scalars. More...
 

Public Member Functions

 Polynomial (unsigned int deg, const CoeffT &cloneCoeff, unsigned int reserve=0)
 Create a polynomial of degree deg. More...
 
 Polynomial (unsigned int deg, unsigned int reserve=0)
 Create a polynomial of degree deg without cloning. DANGEROUS! More...
 
 ~Polynomial ()
 Destructor. More...
 
unsigned int degree () const
 Return degree of polynomial. More...
 
void setDegree (unsigned int deg)
 Set degree of polynomial to deg. More...
 
Teuchos::RCP< CoeffT > getCoefficient (unsigned int i)
 Return ref-count pointer to coefficient i. More...
 
Teuchos::RCP< const CoeffT > getCoefficient (unsigned int i) const
 Return ref-count pointer to constant coefficient i. More...
 
void setCoefficient (unsigned int i, const CoeffT &v)
 Set coefficient i to c. More...
 
void setCoefficientPtr (unsigned int i, const Teuchos::RCP< CoeffT > &c_ptr)
 Set pointer for coefficient i to c_ptr. DANGEROUS! More...
 
void evaluate (typename Teuchos::Polynomial< CoeffT >::scalar_type &t, CoeffT *x, CoeffT *xdot=NULL) const
 Evaluate polynomial and possibly its derivative at time t. More...
 
- Public Member Functions inherited from Teuchos::Describable
virtual std::string description () const
 Return a simple one-line description of this object. More...
 
virtual void describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 Print the object with some verbosity level to a FancyOStream. More...
 
void describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 Version of describe() that takes an std::ostream instead of a FancyOStream. More...
 
virtual ~Describable ()
 Destructor (marked virtual for memory safety of derived classes). More...
 
- Public Member Functions inherited from Teuchos::LabeledObject
 LabeledObject ()
 Construct with an empty label. More...
 
virtual ~LabeledObject ()
 
virtual void setObjectLabel (const std::string &objectLabel)
 Set the object label (see LabeledObject). More...
 
virtual std::string getObjectLabel () const
 Get the object label (see LabeledObject). More...
 

Protected Attributes

unsigned int d
 Degree of polynomial. More...
 
unsigned int sz
 Size of polynomial (may be > d) More...
 
std::vector< Teuchos::RCP
< CoeffT > > 
coeff
 Vector of polynomial coefficients. More...
 

Additional Inherited Members

- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default = VERB_DEFAULT
 Default value for the verbLevel argument of describe(). More...
 

Detailed Description

template<typename CoeffT>
class Teuchos::Polynomial< CoeffT >

Lightweight container class to represent a simple polynomial.

This class represents a simple polynomial of the form:

\[ x(t) = \sum_{i=0}^d x_i t^i \]

where $d$ is the degree of the polynomial and $t$ is a scalar. The coefficients $x_i$, $i=0,\dots,d$ can be scalars, vectors, operators, etc., any type that can implement Teuchos::PolynomialTraits. A template specialization of Teuchos::PolynomialTraits must be provided for the coefficient type, however Teuchos::PolynomailTraits does provide a default template definition for scalars.

This class provides methods for creating a polynomial of some degree, setting and retreiving coefficients, and evaluating the polynomial and its derivative at some value $t$.

Definition at line 69 of file Teuchos_PolynomialDecl.hpp.

Member Typedef Documentation

template<typename CoeffT>
typedef CoeffT Teuchos::Polynomial< CoeffT >::coeff_type

Typename of coefficients.

Definition at line 73 of file Teuchos_PolynomialDecl.hpp.

template<typename CoeffT>
typedef Teuchos::PolynomialTraits<coeff_type>::scalar_type Teuchos::Polynomial< CoeffT >::scalar_type

Typename of scalars.

Definition at line 76 of file Teuchos_PolynomialDecl.hpp.

Constructor & Destructor Documentation

template<typename CoeffT >
Teuchos::Polynomial< CoeffT >::Polynomial ( unsigned int  deg,
const CoeffT &  cloneCoeff,
unsigned int  reserve = 0 
)

Create a polynomial of degree deg.

If reserve > deg, a polynomial of degree deg will be created, but space for reserve + 1 coefficients will be created to allow future increases in the degree of the polynomial to be more efficient. A clone of cloneCoeff will be placed at each coefficient.

Definition at line 49 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
Teuchos::Polynomial< CoeffT >::Polynomial ( unsigned int  deg,
unsigned int  reserve = 0 
)

Create a polynomial of degree deg without cloning. DANGEROUS!

In this version of the constructor, no clone object is provided, and therefore no storage will be created for each coefficient. In this case, setCoefficientPtr() below should be used to set each coefficient pointer to a valid cofficient. This constructor exists to be able to create an efficient "view" of another polynomial.

Definition at line 65 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
Teuchos::Polynomial< CoeffT >::~Polynomial ( )

Destructor.

Definition at line 78 of file Teuchos_Polynomial.hpp.

Member Function Documentation

template<typename CoeffT>
unsigned int Teuchos::Polynomial< CoeffT >::degree ( ) const
inline

Return degree of polynomial.

Definition at line 102 of file Teuchos_PolynomialDecl.hpp.

template<typename CoeffT >
void Teuchos::Polynomial< CoeffT >::setDegree ( unsigned int  deg)

Set degree of polynomial to deg.

Definition at line 84 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
Teuchos::RCP< CoeffT > Teuchos::Polynomial< CoeffT >::getCoefficient ( unsigned int  i)

Return ref-count pointer to coefficient i.

Definition at line 99 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
Teuchos::RCP< const CoeffT > Teuchos::Polynomial< CoeffT >::getCoefficient ( unsigned int  i) const

Return ref-count pointer to constant coefficient i.

Definition at line 113 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
void Teuchos::Polynomial< CoeffT >::setCoefficient ( unsigned int  i,
const CoeffT &  v 
)

Set coefficient i to c.

Definition at line 127 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
void Teuchos::Polynomial< CoeffT >::setCoefficientPtr ( unsigned int  i,
const Teuchos::RCP< CoeffT > &  c_ptr 
)

Set pointer for coefficient i to c_ptr. DANGEROUS!

Directly set the coefficient pointer to c_ptr. This method should be used with care since future calls to setCoefficient(i,c) will also modify the coefficient pointed to c_ptr. However, in certain situations it is necessary to do this for efficiency.

Definition at line 145 of file Teuchos_Polynomial.hpp.

template<typename CoeffT >
void Teuchos::Polynomial< CoeffT >::evaluate ( typename Teuchos::Polynomial< CoeffT >::scalar_type t,
CoeffT *  x,
CoeffT *  xdot = NULL 
) const

Evaluate polynomial and possibly its derivative at time t.

The value of the polynomial at t is computed and stored in *x. If xdot is not NULL, the derivative with respect to t is evaluated and stored in *xdot.

Horner's method is used to efficiently evaluate the polynomial and its derivative.

Definition at line 161 of file Teuchos_Polynomial.hpp.

Member Data Documentation

template<typename CoeffT>
unsigned int Teuchos::Polynomial< CoeffT >::d
protected

Degree of polynomial.

Definition at line 151 of file Teuchos_PolynomialDecl.hpp.

template<typename CoeffT>
unsigned int Teuchos::Polynomial< CoeffT >::sz
protected

Size of polynomial (may be > d)

Definition at line 154 of file Teuchos_PolynomialDecl.hpp.

template<typename CoeffT>
std::vector< Teuchos::RCP<CoeffT> > Teuchos::Polynomial< CoeffT >::coeff
protected

Vector of polynomial coefficients.

coeff[i] corresponds to the degree i term, i=0,...,d

Definition at line 160 of file Teuchos_PolynomialDecl.hpp.


The documentation for this class was generated from the following files: