|
Stokhos Package Browser (Single Doxygen Collection)
Version of the Day
|
Class representing a KL expansion of an exponential random field. More...
#include <Stokhos_KL_ExponentialRandomField.hpp>

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... | |
Class representing a KL expansion of an exponential random field.
This class provides a means for evaluating a random field
,
,
through the KL expansion
where
is a
-dimensional hyper-rectangle, for the case when the covariance function of
is exponential:
In this case, the covariance function and domain factor into a product 1-dimensional covariance functions over 1-dimensional domains, and thus the eigenfunctions
and eigenvalues
factor into a product of corresponding 1-dimensional eigenfunctions and values. These are computed by the OneDExponentialCovarianceFunction class For a given value of
, the code works by computing the
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
are chosen as the
KL eigenpairs. Then given values for the random variables
, 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:
of KL terms
for each dimension
for each dimension
for each dimension
of the random field
of the random field Additionally parameters for the OneDExponentialCovarianceFunction are also accepted.
Definition at line 109 of file Stokhos_KL_ExponentialRandomField.hpp.
| typedef ExponentialOneDEigenFunction<value_type> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::one_d_eigen_func_type |
Definition at line 112 of file Stokhos_KL_ExponentialRandomField.hpp.
| typedef OneDEigenPair<one_d_eigen_func_type> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::one_d_eigen_pair_type |
Definition at line 113 of file Stokhos_KL_ExponentialRandomField.hpp.
| typedef ProductEigenPair<one_d_eigen_func_type,execution_space> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::product_eigen_pair_type |
Definition at line 114 of file Stokhos_KL_ExponentialRandomField.hpp.
| typedef Kokkos::View<one_d_eigen_func_type**,execution_space> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::eigen_func_array_type |
Definition at line 115 of file Stokhos_KL_ExponentialRandomField.hpp.
| typedef Kokkos::View<value_type*,execution_space> Stokhos::KL::ExponentialRandomField< value_type, execution_space >::eigen_value_array_type |
Definition at line 116 of file Stokhos_KL_ExponentialRandomField.hpp.
|
inline |
Default constructor.
Definition at line 119 of file Stokhos_KL_ExponentialRandomField.hpp.
| Stokhos::KL::ExponentialRandomField< value_type, execution_space >::ExponentialRandomField | ( | Teuchos::ParameterList & | solverParams | ) |
Constructor.
Definition at line 50 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.
|
inline |
Destructor.
Definition at line 126 of file Stokhos_KL_ExponentialRandomField.hpp.
|
inline |
Return spatial dimension of the field.
Definition at line 130 of file Stokhos_KL_ExponentialRandomField.hpp.
|
inline |
Return stochastic dimension of the field.
Definition at line 134 of file Stokhos_KL_ExponentialRandomField.hpp.
| 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 167 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.
|
inline |
Evaluate mean of random field at a point.
Definition at line 147 of file Stokhos_KL_ExponentialRandomField.hpp.
| 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 185 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.
| 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 203 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.
|
inline |
Return eigenvalue.
Definition at line 165 of file Stokhos_KL_ExponentialRandomField.hpp.
| void Stokhos::KL::ExponentialRandomField< value_type, execution_space >::print | ( | std::ostream & | os | ) | const |
Print KL expansion.
Definition at line 215 of file Stokhos_KL_ExponentialRandomFieldImp.hpp.
|
protected |
Number of KL terms.
Definition at line 173 of file Stokhos_KL_ExponentialRandomField.hpp.
|
protected |
Dimension of expansion.
Definition at line 176 of file Stokhos_KL_ExponentialRandomField.hpp.
|
protected |
Mean of random field.
Definition at line 179 of file Stokhos_KL_ExponentialRandomField.hpp.
|
protected |
Standard deviation of random field.
Definition at line 182 of file Stokhos_KL_ExponentialRandomField.hpp.
|
protected |
Product eigenfunctions.
Definition at line 185 of file Stokhos_KL_ExponentialRandomField.hpp.
|
protected |
Product eigenvalues.
Definition at line 188 of file Stokhos_KL_ExponentialRandomField.hpp.
1.8.5