Sacado Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
view_factory_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 #include "Sacado.hpp"
11 
12 #if !defined(__CUDACC__)
14 #include "Kokkos_ViewFactory.hpp"
15 #include "Teuchos_Assert.hpp"
16 
17 // Example to demonstrate the use of Kokkos::ViewFactory for constructing
18 // view's of Fad's without explicitly referencing the sacado dimension
19 
20 // An example function that takes two views of various ranks and scalar types
21 // and produces a third allocated using the ViewFactory
22 template <class View1, class View2>
23 Kokkos::View< typename Kokkos::ViewFactory<View1,View2>::value_type*,
24  typename View1::device_type >
25 my_func(const View1& v1, const View2& v2)
26 {
27  typedef Kokkos::ViewFactory<View1,View2> ViewFac;
28  typedef typename ViewFac::value_type value_type;
29  typedef typename View1::device_type device_type;
30  typedef typename View1::size_type size_type;
31  typedef Kokkos::View<value_type*,device_type> ViewTmp;
32 
33  static_assert( unsigned(View1::rank) == 2, "" );
34  static_assert( unsigned(View2::rank) == 1, "" );
35 
36  const size_type m = v1.extent(0);
37  const size_type n = v1.extent(1);
38  ViewTmp vtmp = ViewFac::template create_view<ViewTmp>(v1,v2,"tmp",m);
39 
40  Kokkos::parallel_for(m, KOKKOS_LAMBDA (const size_type i) {
41  value_type t = 0;
42  for (size_type j=0; j<n; ++j)
43  t += v1(i,j)*v2(j);
44  vtmp(i) = t;
45  });
46 
47  return vtmp;
48 }
49 
50 #if defined(HAVE_SACADO_KOKKOS)
51 // An example function that takes two dynamic-rank views of various ranks and
52 // scalar types and produces a third allocated using the ViewFactory
53 template <class View1, class View2>
54 Kokkos::DynRankView< typename Kokkos::ViewFactory<View1,View2>::value_type,
55  typename View1::device_type >
56 my_func_dynamic(const View1& v1, const View2& v2)
57 {
58  typedef Kokkos::ViewFactory<View1,View2> ViewFac;
59  typedef typename ViewFac::value_type value_type;
60  typedef typename View1::device_type device_type;
61  typedef typename View1::size_type size_type;
62  typedef Kokkos::DynRankView<value_type,device_type> ViewTmp;
63 
64  TEUCHOS_ASSERT( v1.rank() == 2 );
65  TEUCHOS_ASSERT( v2.rank() == 1 );
66 
67  const size_type m = v1.extent(0);
68  const size_type n = v1.extent(1);
69  ViewTmp vtmp = ViewFac::template create_view<ViewTmp>(v1,v2,"tmp",m);
70 
71  Kokkos::parallel_for(m, KOKKOS_LAMBDA (const size_type i) {
72  value_type t = 0;
73  for (size_type j=0; j<n; ++j)
74  t += v1(i,j)*v2(j);
75  vtmp(i) = t;
76  });
77 
78  return vtmp;
79 }
80 #endif
81 
82 #endif
83 
84 int main(int argc, char* argv[]) {
85 
86 #if !defined(__CUDACC__)
88  typedef Kokkos::DefaultExecutionSpace execution_space;
89 
90  Kokkos::initialize(argc, argv);
91 
92  const size_t m = 10;
93  const size_t n = 2;
94  const size_t deriv_dim = 1;
95 
96  // First use the static-rank view
97  {
98 
99  // Calculation with double's
100  Kokkos::View<double**,execution_space> v1("v1",m,n);
101  Kokkos::View<double* ,execution_space> v2("v2",n);
102 
103  Kokkos::deep_copy(v1, 2.0);
104  Kokkos::deep_copy(v2, 3.0);
105 
106  auto v3 = my_func(v1,v2);
107 
108  std::cout << "v3 = " << std::endl;
109  for (size_t i=0; i<m; ++i) {
110  std::cout << "\t" << v3(i) << std::endl;
111  }
112 
113  // Calculation with Fad's (initialize each component of v2_fad with a
114  // Fad object with value == 3.0, one derivative == 1
115  Kokkos::View<FadType*,execution_space> v2_fad("v2_fad",n,deriv_dim+1);
116  Kokkos::deep_copy(v2_fad, FadType(deriv_dim, 0, 3.0));
117 
118  auto v3_fad = my_func(v1,v2_fad);
119 
120  std::cout << "v3_fad = " << std::endl;
121  for (size_t i=0; i<m; ++i) {
122  std::cout << "\t" << v3_fad(i) << std::endl;
123  }
124 
125  }
126 
127 #if defined(HAVE_SACADO_KOKKOS)
128  // Now use the dynamic-rank view
129  {
130 
131  // Calculation with double's
132  Kokkos::DynRankView<double,execution_space> v1("v1",m,n);
133  Kokkos::DynRankView<double,execution_space> v2("v2",n);
134 
135  Kokkos::deep_copy(v1, 2.0);
136  Kokkos::deep_copy(v2, 3.0);
137 
138  auto v3 = my_func_dynamic(v1,v2);
139 
140  std::cout << "v3 = " << std::endl;
141  for (size_t i=0; i<m; ++i) {
142  std::cout << "\t" << v3(i) << std::endl;
143  }
144 
145  // Calculation with Fad's (initialize each component of v2_fad with a
146  // Fad object with value == 3.0, one derivative == 1
147  Kokkos::DynRankView<FadType,execution_space> v2_fad("v2_fad",n,deriv_dim+1);
148  Kokkos::deep_copy(v2_fad, FadType(deriv_dim, 0, 3.0));
149 
150  auto v3_fad = my_func_dynamic(v1,v2_fad);
151 
152  std::cout << "v3_fad = " << std::endl;
153  for (size_t i=0; i<m; ++i) {
154  std::cout << "\t" << v3_fad(i) << std::endl;
155  }
156 
157  }
158 #endif
159 
160  Kokkos::finalize();
161 #endif
162 
163  return 0;
164 }
Kokkos::View< typename Kokkos::ViewFactory< View1, View2 >::value_type *, typename View1::device_type > my_func(const View1 &v1, const View2 &v2)
Sacado::Fad::DFad< double > FadType
int main()
Definition: ad_example.cpp:171
#define TEUCHOS_ASSERT(assertion_test)
int n