Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
Stokhos::KL::OneDExponentialCovarianceFunction< value_type > Class Template Reference

Class representing an exponential covariance function and its KL eigevalues/eigenfunctions. More...

#include <Stokhos_KL_OneDExponentialCovarianceFunction.hpp>

Classes

struct  EigFuncCos
 Nonlinear function whose roots define eigenvalues for cos() eigenfunction. More...
 
struct  EigFuncSin
 Nonlinear function whose roots define eigenvalues for sin() eigenfunction. More...
 

Public Types

typedef
ExponentialOneDEigenFunction
< value_type > 
eigen_function_type
 
typedef OneDEigenPair
< eigen_function_type
eigen_pair_type
 

Public Member Functions

 OneDExponentialCovarianceFunction (int M, const value_type &a, const value_type &b, const value_type &L, const int dim_name, Teuchos::ParameterList &solverParams)
 Constructor. More...
 
 ~OneDExponentialCovarianceFunction ()
 Destructor. More...
 
value_type evaluateCovariance (const value_type &x, const value_type &xp) const
 Evaluate covariance. More...
 
const Teuchos::Array
< eigen_pair_type > & 
getEigenPairs () const
 Get eigenpairs. More...
 

Protected Types

typedef Teuchos::ScalarTraits
< value_type >::magnitudeType 
magnitude_type
 

Protected Member Functions

template<class Func >
value_type newton (const Func &func, const value_type &a, const value_type &b, magnitude_type tol, int max_num_its)
 A basic root finder based on Newton's method. More...
 
template<class Func >
value_type bisection (const Func &func, const value_type &a, const value_type &b, magnitude_type tol, int max_num_its)
 A basic root finder based on bisection. More...
 

Protected Attributes

value_type L
 Correlation length. More...
 
Teuchos::Array< eigen_pair_typeeig_pair
 Eigenpairs. More...
 

Private Member Functions

 OneDExponentialCovarianceFunction (const OneDExponentialCovarianceFunction &)
 Prohibit copying. More...
 
OneDExponentialCovarianceFunctionoperator= (const OneDExponentialCovarianceFunction &)
 Prohibit copying. More...
 

Detailed Description

template<typename value_type>
class Stokhos::KL::OneDExponentialCovarianceFunction< value_type >

Class representing an exponential covariance function and its KL eigevalues/eigenfunctions.

This class provides the exponential covariance function

\[ \mbox{cov}(x,x') = \exp(-|x-x'|/L). \]

The corresponding eigenfunctions can be shown to be $A_n\cos(\omega_n(x-\beta))$ and $B_n\sin(\omega^\ast_n(x-\beta))$ for $x\in[a,b]$ where $\omega_n$ and $\omega^\ast_n$ satisfy

\[ 1 - L\omega_n\tan(\omega_n\alpha) = 0 \]

and

\[ L\omega^\ast_n + \tan(\omega^\ast_n\alpha) = 0 \]

respectively, where $\alpha=(b-a)/2$ and $\beta=(b+a)/2$. Then

\[ A_n = \frac{1}{\left(\int_a^b\cos^2(\omega_n(x-\beta)) dx\right)^{1/2}} = \frac{1}{\sqrt{\alpha + \frac{\sin(2\omega_n\alpha)}{2\omega_n}}}, \]

\[ B_n = \frac{1}{\left(\int_a^b\sin^2(\omega^\ast_n(x-\beta)) dx\right)^{1/2}} = \frac{1}{\sqrt{\alpha - \frac{\sin(2\omega_n^\ast\alpha)}{2\omega^\ast_n}}} \]

and the corresponding eigenvalues are given by

\[ \lambda_n = \frac{2L}{(L\omega_n)^2 + 1} \]

and

\[ \lambda^\ast_n = \frac{2L}{(L\omega^\ast_n)^2 + 1}. \]

It is straightforward to show that for each $n$, $\omega^\ast_n < \omega_n$, and thus $\lambda_n < \lambda^\ast_n$. Hence when sorted on decreasing eigenvalue, the eigenfunctions alternate between cosine and sine. See "Stochastic Finite Elements" by Ghanem and Spanos for a complete description of how to derive these relationships.

For a given value of $M$, the code works by computing the $M$ eigenfunctions using a bisection root solver to compute the frequencies $\omega_n$ and $\omega^\ast_n$.

Data for the root solver is passed through a Teuchos::ParameterList, which accepts the following parameters:

Definition at line 87 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.

Member Typedef Documentation

template<typename value_type>
typedef ExponentialOneDEigenFunction<value_type> Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::eigen_function_type
template<typename value_type>
typedef Teuchos::ScalarTraits<value_type>::magnitudeType Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::magnitude_type
protected

Constructor & Destructor Documentation

template<typename value_type >
Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::OneDExponentialCovarianceFunction ( int  M,
const value_type &  a,
const value_type &  b,
const value_type &  L,
const int  dim_name,
Teuchos::ParameterList solverParams 
)

Constructor.

Definition at line 14 of file Stokhos_KL_OneDExponentialCovarianceFunctionImp.hpp.

template<typename value_type>
Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::~OneDExponentialCovarianceFunction ( )
inline

Destructor.

Definition at line 102 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.

template<typename value_type>
Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::OneDExponentialCovarianceFunction ( const OneDExponentialCovarianceFunction< value_type > &  )
private

Prohibit copying.

Member Function Documentation

template<typename value_type>
value_type Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::evaluateCovariance ( const value_type &  x,
const value_type &  xp 
) const
inline

Evaluate covariance.

Definition at line 105 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.

template<typename value_type>
const Teuchos::Array<eigen_pair_type>& Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::getEigenPairs ( ) const
inline

Get eigenpairs.

Definition at line 111 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.

template<typename value_type>
OneDExponentialCovarianceFunction& Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::operator= ( const OneDExponentialCovarianceFunction< value_type > &  )
private

Prohibit copying.

template<typename value_type >
template<class Func >
value_type Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::newton ( const Func &  func,
const value_type &  a,
const value_type &  b,
magnitude_type  tol,
int  max_num_its 
)
protected

A basic root finder based on Newton's method.

Definition at line 68 of file Stokhos_KL_OneDExponentialCovarianceFunctionImp.hpp.

template<typename value_type >
template<class Func >
value_type Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::bisection ( const Func &  func,
const value_type &  a,
const value_type &  b,
magnitude_type  tol,
int  max_num_its 
)
protected

A basic root finder based on bisection.

Definition at line 92 of file Stokhos_KL_OneDExponentialCovarianceFunctionImp.hpp.

Member Data Documentation

template<typename value_type>
value_type Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::L
protected

Correlation length.

Definition at line 128 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.

template<typename value_type>
Teuchos::Array<eigen_pair_type> Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::eig_pair
protected

Eigenpairs.

Definition at line 131 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.


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