54 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
57 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
66 MatVec(
const ViewTypeA& A_,
const ViewTypeB& b_,
const ViewTypeC& c_) :
68 m(
A.extent(0)),
n(
A.extent(1))
70 typedef typename ViewTypeC::execution_space execution_space;
71 Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>(0,m), *
this);
75 void operator() (
const size_t i)
const {
76 typedef typename ViewTypeC::value_type scalar_type;
79 for (
size_t j=0; j<
n; ++j) {
80 scalar_type bb = b(j);
88 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
89 void run_mat_vec(
const ViewTypeA&
A,
const ViewTypeB& b,
const ViewTypeC&
c)
91 MatVec<ViewTypeA,ViewTypeB,ViewTypeC>
f(A,b,c);
96 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
101 const size_t m,
n, p;
105 MatVecDeriv(
const ViewTypeA& A_,
const ViewTypeB& b_,
const ViewTypeC& c_) :
107 m(
A.extent(0)),
n(
A.extent(1)), p(
A.extent(2)-1)
109 typedef typename ViewTypeC::execution_space execution_space;
110 Kokkos::parallel_for(Kokkos::RangePolicy<execution_space>(0,m), *
this);
114 void operator() (
const size_t i)
const {
115 typedef typename ViewTypeC::value_type scalar_type;
118 for (
size_t k=0; k<p; ++k) {
120 for (
size_t j=0; j<
n; ++j)
121 t +=
A(i,j,k)*b(j,p) +
A(i,j,p)*b(j,k);
127 for (
size_t j=0; j<n; ++j)
128 t +=
A(i,j,p)*b(j,p);
134 template <
typename ViewTypeA,
typename ViewTypeB,
typename ViewTypeC>
138 MatVecDeriv<ViewTypeA,ViewTypeB,ViewTypeC>
f(A,b,c);
143 int main(
int argc,
char* argv[]) {
146 #if defined(HAVE_SACADO_VIEW_SPEC) && !defined(SACADO_DISABLE_FAD_VIEW_SPEC)
148 Kokkos::initialize(argc, argv);
158 Kokkos::View<FadType**>
A(
"A",m,n,p+1);
159 Kokkos::View<FadType*> b(
"b",n,p+1);
160 Kokkos::View<FadType*>
c(
"c",m,p+1);
163 Kokkos::deep_copy( A,
FadType(p, 0, 2.0) );
164 Kokkos::deep_copy( b,
FadType(p, 0, 3.0) );
165 Kokkos::deep_copy( c, 0.0 );
171 std::cout <<
"\nc = A*b: Differentiated using Sacado:" << std::endl;
172 auto h_c = Kokkos::create_mirror_view(c);
173 Kokkos::deep_copy(h_c, c);
174 for (
size_t i=0; i<m; ++i)
175 std::cout <<
"\tc(" << i <<
") = " << h_c(i) << std::endl;
180 Kokkos::View<FadType*> c2(
"c",m,p+1);
181 Kokkos::View<double***> A_flat =
A;
182 Kokkos::View<double**> b_flat = b;
183 Kokkos::View<double**> c_flat = c2;
187 std::cout <<
"\nc = A*b: Differentiated analytically:" << std::endl;
188 auto h_c2 = Kokkos::create_mirror_view(c2);
189 Kokkos::deep_copy(h_c2, c2);
190 for (
size_t i=0; i<m; ++i)
191 std::cout <<
"\tc(" << i <<
") = " << h_c2(i) << std::endl;
195 for (
size_t i=0; i<m; ++i) {
196 for (
size_t k=0; k<p; ++k)
201 double tol = 1.0e-14;
204 std::cout <<
"\nExample passed!" << std::endl;
208 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
#define KOKKOS_INLINE_FUNCTION
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
expr expr expr fastAccessDx(i)) FAD_UNARYOP_MACRO(exp
void run_mat_vec_deriv(const ViewTypeA &A, const ViewTypeB &b, const ViewTypeC &c)