23 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
26 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
35 MatVec(
const ViewTypeA& A_,
const ViewTypeB&
b_,
const ViewTypeC&
c_) :
37 m(
A.extent(0)),
n(
A.extent(1))
39 typedef typename ViewTypeC::execution_space execution_space;
40 Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>(0,m), *
this);
43 KOKKOS_INLINE_FUNCTION
44 void operator() (
const size_t i)
const {
45 typedef typename ViewTypeC::value_type scalar_type;
48 for (
size_t j=0; j<
n; ++j) {
49 scalar_type bb = b(j);
57 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
58 void run_mat_vec(
const ViewTypeA&
A,
const ViewTypeB& b,
const ViewTypeC&
c)
60 MatVec<ViewTypeA,ViewTypeB,ViewTypeC>
f(A,b,c);
65 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
74 MatVecDeriv(
const ViewTypeA& A_,
const ViewTypeB& b_,
const ViewTypeC& c_) :
76 m(
A.extent(0)),
n(
A.extent(1)),
p(
A.extent(2)-1)
78 typedef typename ViewTypeC::execution_space execution_space;
79 Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>(0,m), *
this);
82 KOKKOS_INLINE_FUNCTION
83 void operator() (
const size_t i)
const {
84 typedef typename ViewTypeC::value_type scalar_type;
87 for (
size_t k=0; k<
p; ++k) {
89 for (
size_t j=0; j<
n; ++j)
90 t +=
A(i,j,k)*b(j,p) +
A(i,j,p)*b(j,k);
96 for (
size_t j=0; j<n; ++j)
103 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
107 MatVecDeriv<ViewTypeA,ViewTypeB,ViewTypeC>
f(A,b,c);
112 int main(
int argc,
char* argv[]) {
115 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
117 Kokkos::initialize(argc, argv);
127 Kokkos::View<FadType**>
A(
"A",m,n,p+1);
128 Kokkos::View<FadType*> b(
"b",n,p+1);
129 Kokkos::View<FadType*>
c(
"c",m,p+1);
132 Kokkos::deep_copy( A,
FadType(p, 0, 2.0) );
133 Kokkos::deep_copy( b,
FadType(p, 0, 3.0) );
134 Kokkos::deep_copy( c, 0.0 );
140 std::cout <<
"\nc = A*b: Differentiated using Sacado:" << std::endl;
141 auto h_c = Kokkos::create_mirror_view(c);
142 Kokkos::deep_copy(h_c, c);
143 for (
size_t i=0; i<m; ++
i)
144 std::cout <<
"\tc(" << i <<
") = " << h_c(i) << std::endl;
149 Kokkos::View<FadType*> c2(
"c",m,p+1);
150 #ifndef SACADO_HAS_NEW_KOKKOS_VIEW_IMPL
151 Kokkos::View<double***> A_flat =
A;
152 Kokkos::View<double**> b_flat = b;
153 Kokkos::View<double**> c_flat = c2;
155 Kokkos::View<double***> A_flat(A.data_handle(), A.extent(0), A.extent(1), A.accessor().fad_size() + 1);
156 Kokkos::View<double**> b_flat(b.data_handle(), b.extent(0), b.accessor().fad_size() + 1);
157 Kokkos::View<double**> c_flat(c2.data_handle(), c2.extent(0), c2.accessor().fad_size() + 1);
162 std::cout <<
"\nc = A*b: Differentiated analytically:" << std::endl;
163 auto h_c2 = Kokkos::create_mirror_view(c2);
164 Kokkos::deep_copy(h_c2, c2);
165 for (
size_t i=0; i<m; ++
i)
166 std::cout <<
"\tc(" << i <<
") = " << h_c2(i) << std::endl;
170 for (
size_t i=0; i<m; ++
i) {
171 for (
size_t k=0; k<
p; ++k)
176 double tol = 1.0e-14;
179 std::cout <<
"\nExample passed!" << std::endl;
183 std::cout <<
"\nSomething is wrong, example failed!" << std::endl;
void run_mat_vec(const ViewTypeA &A, const ViewTypeB &b, const ViewTypeC &c)
Sacado::Fad::DFad< double > FadType
expr expr1 expr1 expr1 c expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr1 c expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr1 c *expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr1 c expr2 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr2 expr1 expr2 expr1 expr1 expr1 expr2 expr1 expr2 expr1 expr1 expr1 c
void run_mat_vec_deriv(const ViewTypeA &A, const ViewTypeB &b, const ViewTypeC &c)
expr expr expr fastAccessDx(i)) FAD_UNARYOP_MACRO(exp