12 #if !defined(__CUDACC__)
22 template <
class View1,
class View2>
23 Kokkos::View< typename Kokkos::ViewFactory<View1,View2>::value_type*,
24 typename View1::device_type >
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;
33 static_assert(
unsigned(View1::rank) == 2,
"" );
34 static_assert(
unsigned(View2::rank) == 1,
"" );
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);
40 Kokkos::parallel_for(m, KOKKOS_LAMBDA (
const size_type
i) {
42 for (size_type j=0; j<n; ++j)
50 #if defined(HAVE_SACADO_KOKKOS)
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)
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;
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);
71 Kokkos::parallel_for(m, KOKKOS_LAMBDA (
const size_type
i) {
73 for (size_type j=0; j<n; ++j)
84 int main(
int argc,
char* argv[]) {
86 #if !defined(__CUDACC__)
88 typedef Kokkos::DefaultExecutionSpace execution_space;
90 Kokkos::initialize(argc, argv);
94 const size_t deriv_dim = 1;
100 Kokkos::View<double**,execution_space> v1(
"v1",m,n);
101 Kokkos::View<double* ,execution_space> v2(
"v2",n);
103 Kokkos::deep_copy(v1, 2.0);
104 Kokkos::deep_copy(v2, 3.0);
108 std::cout <<
"v3 = " << std::endl;
109 for (
size_t i=0; i<m; ++
i) {
110 std::cout <<
"\t" << v3(i) << std::endl;
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));
118 auto v3_fad =
my_func(v1,v2_fad);
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;
127 #if defined(HAVE_SACADO_KOKKOS)
132 Kokkos::DynRankView<double,execution_space> v1(
"v1",m,n);
133 Kokkos::DynRankView<double,execution_space> v2(
"v2",n);
135 Kokkos::deep_copy(v1, 2.0);
136 Kokkos::deep_copy(v2, 3.0);
138 auto v3 = my_func_dynamic(v1,v2);
140 std::cout <<
"v3 = " << std::endl;
141 for (
size_t i=0; i<m; ++
i) {
142 std::cout <<
"\t" << v3(i) << std::endl;
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));
150 auto v3_fad = my_func_dynamic(v1,v2_fad);
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;
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
#define TEUCHOS_ASSERT(assertion_test)