Stokhos Package Browser (Single Doxygen Collection)  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
Stokhos::KL::ExponentialRandomField< value_type, execution_space > Class Template Reference

Class representing a KL expansion of an exponential random field. More...

#include <Stokhos_KL_ExponentialRandomField.hpp>

Inheritance diagram for Stokhos::KL::ExponentialRandomField< value_type, execution_space >:
Inheritance graph
[legend]

Public Types

typedef
ExponentialOneDEigenFunction
< value_type > 
one_d_eigen_func_type
 
typedef OneDEigenPair
< one_d_eigen_func_type
one_d_eigen_pair_type
 
typedef ProductEigenPair
< one_d_eigen_func_type,
execution_space
product_eigen_pair_type
 
typedef Kokkos::View
< one_d_eigen_func_type
**, execution_space
eigen_func_array_type
 
typedef Kokkos::View
< value_type
*, execution_space
eigen_value_array_type
 

Public Member Functions

 ExponentialRandomField ()
 Default constructor. More...
 
 ExponentialRandomField (Teuchos::ParameterList &solverParams)
 Constructor. More...
 
KOKKOS_INLINE_FUNCTION ~ExponentialRandomField ()
 Destructor. More...
 
KOKKOS_INLINE_FUNCTION int spatialDimension () const
 Return spatial dimension of the field. More...
 
KOKKOS_INLINE_FUNCTION int stochasticDimension () const
 Return stochastic dimension of the field. More...
 
template<typename point_type , typename rv_type >
KOKKOS_INLINE_FUNCTION
Teuchos::PromotionTraits
< typename rv_type::value_type,
value_type >::promote 
evaluate (const point_type &point, const rv_type &random_variables) const
 Evaluate random field at a point. More...
 
template<typename point_type >
KOKKOS_INLINE_FUNCTION value_type evaluate_mean (const point_type &point) const
 Evaluate mean of random field at a point. More...
 
template<typename point_type >
KOKKOS_INLINE_FUNCTION
Teuchos::PromotionTraits
< typename
point_type::value_type,
value_type >::promote 
evaluate_standard_deviation (const point_type &point) const
 Evaluate standard deviation of random field at a point. More...
 
template<typename point_type >
KOKKOS_INLINE_FUNCTION
Teuchos::PromotionTraits
< typename
point_type::value_type,
value_type >::promote 
evaluate_eigenfunction (const point_type &point, int i) const
 Evaluate given eigenfunction at a point. More...
 
value_type KOKKOS_INLINE_FUNCTION eigenvalue (int i) const
 Return eigenvalue. More...
 
void print (std::ostream &os) const
 Print KL expansion. More...
 

Protected Attributes

int num_KL
 Number of KL terms. More...
 
int dim
 Dimension of expansion. More...
 
value_type mean
 Mean of random field. More...
 
value_type std_dev
 Standard deviation of random field. More...
 
eigen_func_array_type product_eigen_funcs
 Product eigenfunctions. More...
 
eigen_value_array_type product_eigen_values
 Product eigenvalues. More...
 

Detailed Description

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
class Stokhos::KL::ExponentialRandomField< value_type, execution_space >

Class representing a KL expansion of an exponential random field.

This class provides a means for evaluating a random field $a(x,\theta)$, $x\in D$, $\theta\in\Omega$ through the KL expansion

\[ a(x,\theta) \approx a_0 + \sigma\sum_{k=1}^M \sqrt{\lambda_k}b_k(x)\xi_k(\theta) \]

where $D$ is a $d$-dimensional hyper-rectangle, for the case when the covariance function of $a$ is exponential:

\[ \mbox{cov}(x,x') = \sigma\exp(-|x_1-x_1'|/L_1-...-|x_d-x_d'|/L_d). \]

In this case, the covariance function and domain factor into a product 1-dimensional covariance functions over 1-dimensional domains, and thus the eigenfunctions $b_k$ and eigenvalues $\lambda_k$ factor into a product of corresponding 1-dimensional eigenfunctions and values. These are computed by the OneDExponentialCovarianceFunction class For a given value of $M$, the code works by computing the $M$ eigenfunctions for each direction using this class. Then all possible tensor products of these one-dimensional eigenfunctions are formed, with corresponding eigenvalue given by the product of the one-dimensional eigenvalues. These eigenvalues are then sorted in increasing order, and the first $M$ are chosen as the $M$ KL eigenpairs. Then given values for the random variables $\xi_k$, the class provides a routine for evaluating a realization of the random field.

The idea for this approach was provided by Chris Miller.

All data is passed into this class through a Teuchos::ParameterList, which accepts the following parameters:

Additionally parameters for the OneDExponentialCovarianceFunction are also accepted.

Definition at line 76 of file Stokhos_KL_ExponentialRandomField.hpp.

Member Typedef Documentation

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
typedef ExponentialOneDEigenFunction<value_type> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::one_d_eigen_func_type

Definition at line 79 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
typedef OneDEigenPair<one_d_eigen_func_type> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::one_d_eigen_pair_type

Definition at line 80 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
typedef ProductEigenPair<one_d_eigen_func_type,execution_space> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::product_eigen_pair_type

Definition at line 81 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
typedef Kokkos::View<one_d_eigen_func_type**,execution_space> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::eigen_func_array_type

Definition at line 82 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
typedef Kokkos::View<value_type*,execution_space> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::eigen_value_array_type

Definition at line 83 of file Stokhos_KL_ExponentialRandomField.hpp.

Constructor & Destructor Documentation

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
Stokhos::KL::ExponentialRandomField< value_type, execution_space >::ExponentialRandomField ( )
inline

Default constructor.

Definition at line 86 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type , typename execution_space >
Stokhos::KL::ExponentialRandomField< value_type, execution_space >::ExponentialRandomField ( Teuchos::ParameterList solverParams)

Constructor.

Definition at line 18 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
KOKKOS_INLINE_FUNCTION Stokhos::KL::ExponentialRandomField< value_type, execution_space >::~ExponentialRandomField ( )
inline

Destructor.

Definition at line 93 of file Stokhos_KL_ExponentialRandomField.hpp.

Member Function Documentation

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
KOKKOS_INLINE_FUNCTION int Stokhos::KL::ExponentialRandomField< value_type, execution_space >::spatialDimension ( ) const
inline

Return spatial dimension of the field.

Definition at line 97 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
KOKKOS_INLINE_FUNCTION int Stokhos::KL::ExponentialRandomField< value_type, execution_space >::stochasticDimension ( ) const
inline

Return stochastic dimension of the field.

Definition at line 101 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type , typename execution_space >
template<typename point_type , typename rv_type >
KOKKOS_INLINE_FUNCTION Teuchos::PromotionTraits< typename rv_type::value_type, value_type >::promote Stokhos::KL::ExponentialRandomField< value_type, execution_space >::evaluate ( const point_type &  point,
const rv_type &  random_variables 
) const

Evaluate random field at a point.

Definition at line 135 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
template<typename point_type >
KOKKOS_INLINE_FUNCTION value_type Stokhos::KL::ExponentialRandomField< value_type, execution_space >::evaluate_mean ( const point_type &  point) const
inline

Evaluate mean of random field at a point.

Definition at line 114 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type , typename execution_space >
template<typename point_type >
KOKKOS_INLINE_FUNCTION Teuchos::PromotionTraits< typename point_type::value_type, value_type >::promote Stokhos::KL::ExponentialRandomField< value_type, execution_space >::evaluate_standard_deviation ( const point_type &  point) const

Evaluate standard deviation of random field at a point.

Definition at line 153 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.

template<typename value_type , typename execution_space >
template<typename point_type >
KOKKOS_INLINE_FUNCTION Teuchos::PromotionTraits< typename point_type::value_type, value_type >::promote Stokhos::KL::ExponentialRandomField< value_type, execution_space >::evaluate_eigenfunction ( const point_type &  point,
int  i 
) const

Evaluate given eigenfunction at a point.

Definition at line 171 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
value_type KOKKOS_INLINE_FUNCTION Stokhos::KL::ExponentialRandomField< value_type, execution_space >::eigenvalue ( int  i) const
inline

Return eigenvalue.

Definition at line 132 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type , typename execution_space >
void Stokhos::KL::ExponentialRandomField< value_type, execution_space >::print ( std::ostream &  os) const

Print KL expansion.

Definition at line 183 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.

Member Data Documentation

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
int Stokhos::KL::ExponentialRandomField< value_type, execution_space >::num_KL
protected

Number of KL terms.

Definition at line 140 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
int Stokhos::KL::ExponentialRandomField< value_type, execution_space >::dim
protected

Dimension of expansion.

Definition at line 143 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
value_type Stokhos::KL::ExponentialRandomField< value_type, execution_space >::mean
protected

Mean of random field.

Definition at line 146 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
value_type Stokhos::KL::ExponentialRandomField< value_type, execution_space >::std_dev
protected

Standard deviation of random field.

Definition at line 149 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
eigen_func_array_type Stokhos::KL::ExponentialRandomField< value_type, execution_space >::product_eigen_funcs
protected

Product eigenfunctions.

Definition at line 152 of file Stokhos_KL_ExponentialRandomField.hpp.

template<typename value_type, typename execution_space = Kokkos::DefaultExecutionSpace>
eigen_value_array_type Stokhos::KL::ExponentialRandomField< value_type, execution_space >::product_eigen_values
protected

Product eigenvalues.

Definition at line 155 of file Stokhos_KL_ExponentialRandomField.hpp.


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