Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
random_field_example.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 // random_field_example
11 //
12 // usage:
13 // random_field_example
14 //
15 // output:
16 // prints out KL expansion for an exponential random field
17 
18 #include "Stokhos_ConfigDefs.h"
20 
21 // Functor for evaluating random field on the device
22 template <typename Value, typename Device = Kokkos::DefaultExecutionSpace>
23 struct RF {
24  typedef Value value_type;
25  typedef Device execution_space;
26 
28  typedef Kokkos::View<value_type*,execution_space> point_type;
29  typedef Kokkos::View<value_type*,execution_space> rv_type;
30  typedef Kokkos::View<value_type*,execution_space> result_type;
31  const rf_type rf;
32  const point_type x;
33  const rv_type rv;
34  const result_type y;
35 
36  RF(const rf_type& rf_, const point_type& x_, const rv_type& rv_) :
37  rf(rf_), x(x_), rv(rv_), y("y",1)
38  {
39  Kokkos::parallel_for(1, *this);
40  }
41 
42  KOKKOS_INLINE_FUNCTION
43  void operator() (const unsigned i) const {
44  y(0) = rf.evaluate(x, rv);
45  }
46 };
47 
48 int main(int argc, char **argv)
49 {
50  Kokkos::initialize();
51 
52  try {
53 
54  // Create random field
55  int M = 10;
56  Teuchos::ParameterList solverParams;
57  solverParams.set("Number of KL Terms", M);
58  solverParams.set("Mean", 1.0);
59  solverParams.set("Standard Deviation", 0.1);
60  int ndim = 3;
61  Teuchos::Array<double> domain_upper(ndim), domain_lower(ndim),
62  correlation_length(ndim);
63  for (int i=0; i<ndim; i++) {
64  domain_upper[i] = 1.0;
65  domain_lower[i] = 0.0;
66  correlation_length[i] = 10.0;
67  }
68  solverParams.set("Domain Upper Bounds", domain_upper);
69  solverParams.set("Domain Lower Bounds", domain_lower);
70  solverParams.set("Correlation Lengths", correlation_length);
72  rf.print(std::cout);
73 
74  // Evaluate random field at a point
75  Teuchos::Array<double> x(ndim);
76  for (int i=0; i<ndim; i++)
77  x[i] = (domain_upper[i] + domain_lower[i])/2.0 +
78  0.1*(domain_upper[i] - domain_lower[i])/2.0;
79  Teuchos::Array<double> rvar(M);
80  for (int i=0; i<M; i++)
81  rvar[i] = 1.5;
82  double result = rf.evaluate(x, rvar);
83  std::cout << "result (host) = " << result << std::endl;
84 
85  // Evaluate random field in a functor on device
86  typedef Kokkos::View<double*> view_type;
87  typedef view_type::HostMirror host_view_type;
88  view_type x_view("x", ndim);
89  host_view_type host_x = Kokkos::create_mirror_view(x_view);
90  for (int i=0; i<ndim; i++)
91  host_x(i) = x[i];
92  Kokkos::deep_copy(x_view, host_x);
93  view_type rvar_view("rvar", M);
94  host_view_type host_rvar = Kokkos::create_mirror_view(rvar_view);
95  for (int i=0; i<M; i++)
96  host_rvar(i) = rvar[i];
97  Kokkos::deep_copy(rvar_view, host_rvar);
98  RF<double> rf_func(rf, x_view, rvar_view);
99  host_view_type host_y = Kokkos::create_mirror_view(rf_func.y);
100  Kokkos::deep_copy(host_y, rf_func.y);
101  double result2 = host_y(0);
102  std::cout << "result (device)= " << result2 << std::endl;
103  }
104  catch (std::exception& e) {
105  std::cout << e.what() << std::endl;
106  }
107 
108  Kokkos::finalize();
109 }
KOKKOS_INLINE_FUNCTION void operator()(const unsigned i) const
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.
const rv_type rv
void print(std::ostream &os) const
Print KL expansion.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Kokkos::View< value_type *, execution_space > point_type
const rf_type rf
void deep_copy(const Stokhos::CrsMatrix< ValueType, DstDevice, Layout > &dst, const Stokhos::CrsMatrix< ValueType, SrcDevice, Layout > &src)
Kokkos::View< value_type *, execution_space > result_type
int main(int argc, char **argv)
const point_type x
Stokhos::KL::ExponentialRandomField< value_type, execution_space > rf_type
Value value_type
Class representing a KL expansion of an exponential random field.
const result_type y
RF(const rf_type &rf_, const point_type &x_, const rv_type &rv_)
Device execution_space
Stokhos::CrsMatrix< ValueType, Device, Layout >::HostMirror create_mirror_view(const Stokhos::CrsMatrix< ValueType, Device, Layout > &A)
Kokkos::View< value_type *, execution_space > rv_type