| 
    Stokhos Package Browser (Single Doxygen Collection)
    Version of the Day
    
   | 
 
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_type > | eig_pair | 
| Eigenpairs.  More... | |
Private Member Functions | |
| OneDExponentialCovarianceFunction (const OneDExponentialCovarianceFunction &) | |
| Prohibit copying.  More... | |
| OneDExponentialCovarianceFunction & | operator= (const OneDExponentialCovarianceFunction &) | 
| Prohibit copying.  More... | |
Class representing an exponential covariance function and its KL eigevalues/eigenfunctions.
This class provides the exponential covariance function
 The corresponding eigenfunctions can be shown to be 
 and 
 for 
 where 
 and 
 satisfy 
and
 respectively, where 
 and 
. Then 
and the corresponding eigenvalues are given by
and
 It is straightforward to show that for each 
, 
, and thus 
. 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 
, the code works by computing the 
 eigenfunctions using a bisection root solver to compute the frequencies 
 and 
.
Data for the root solver is passed through a Teuchos::ParameterList, which accepts the following parameters:
 for bounding the frequencies 
 in the bisection algorithm Definition at line 87 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
| typedef ExponentialOneDEigenFunction<value_type> Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::eigen_function_type | 
Definition at line 90 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
| typedef OneDEigenPair<eigen_function_type> Stokhos::KL::OneDExponentialCovarianceFunction< value_type >::eigen_pair_type | 
Definition at line 91 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
      
  | 
  protected | 
Definition at line 125 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
| 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.
      
  | 
  inline | 
Destructor.
Definition at line 102 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
      
  | 
  private | 
Prohibit copying.
      
  | 
  inline | 
Evaluate covariance.
Definition at line 105 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
      
  | 
  inline | 
Get eigenpairs.
Definition at line 111 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
      
  | 
  private | 
Prohibit copying.
      
  | 
  protected | 
A basic root finder based on Newton's method.
Definition at line 68 of file Stokhos_KL_OneDExponentialCovarianceFunctionImp.hpp.
      
  | 
  protected | 
A basic root finder based on bisection.
Definition at line 92 of file Stokhos_KL_OneDExponentialCovarianceFunctionImp.hpp.
      
  | 
  protected | 
Correlation length.
Definition at line 128 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
      
  | 
  protected | 
Eigenpairs.
Definition at line 131 of file Stokhos_KL_OneDExponentialCovarianceFunction.hpp.
 1.8.5