Skip to content

Commit

Permalink
Sacado: Fix scalar assignment for hierarchical_dfad.
Browse files Browse the repository at this point in the history
As pointed out by Roger, there was a bug in scalar assignment when
hierarchical_dfad was enabled.  By an unfortunate set of coincidences,
this case fell through the cracks for testing.  I made sure this was
now properly tested and also tested Fad assignment.  I also made
partition_scalar<> work correctly when hierarchical_dfad was enabled.

Build/Test Cases Summary
Enabled Packages: Sacado
Disabled Packages: PyTrilinos,Claps,TriKota
Enabled all Forward Packages
0) MPI_RELEASE_DEBUG_SHARED_PT => passed: passed=1182,notpassed=0 (12.55 min)
  • Loading branch information
etphipp committed Jan 17, 2018
1 parent 0c34a26 commit c6c840d
Show file tree
Hide file tree
Showing 3 changed files with 109 additions and 15 deletions.
4 changes: 4 additions & 0 deletions packages/sacado/src/KokkosExp_View_Fad_Contiguous.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -183,6 +183,7 @@ namespace Sacado {
template <typename S> class GeneralFad;
}
}
#ifndef SACADO_VIEW_CUDA_HIERARCHICAL_DFAD
template <unsigned Stride, typename T, typename U>
KOKKOS_INLINE_FUNCTION
typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::DynamicStorage<T,U> >, Stride >::type
Expand All @@ -202,6 +203,7 @@ namespace Sacado {

return xp;
}
#endif
template <unsigned Stride, typename T, int N>
KOKKOS_INLINE_FUNCTION
typename LocalScalarType< Fad::Exp::GeneralFad< Fad::Exp::StaticStorage<T,N> >, Stride >::type
Expand Down Expand Up @@ -233,6 +235,7 @@ namespace Sacado {
template <typename T, int N> class SLFad;
template <typename T, int N> class SFad;
}
#ifndef SACADO_VIEW_CUDA_HIERARCHICAL_DFAD
template <unsigned Stride, typename T>
KOKKOS_INLINE_FUNCTION
typename LocalScalarType< Fad::DFad<T>, Stride >::type
Expand All @@ -252,6 +255,7 @@ namespace Sacado {

return xp;
}
#endif
template <unsigned Stride, typename T, int N>
KOKKOS_INLINE_FUNCTION
typename LocalScalarType< Fad::SLFad<T,N>, Stride >::type
Expand Down
9 changes: 3 additions & 6 deletions packages/sacado/src/Sacado_DynamicArrayTraits.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -495,11 +495,9 @@ namespace Sacado {
//! Copy array from \c src to \c dest of length \c sz
KOKKOS_INLINE_FUNCTION
static void strided_copy(const T* src, int src_stride,
T* dest, int dest_stride, int sz) {
T* dest, int dest_stride, int sz) {
for (int i=threadIdx.x; i<sz; i+=blockDim.x) {
*(dest) = *(src);
dest += dest_stride;
src += src_stride;
dest[i*dest_stride] = src[i*src_stride];
}
}

Expand All @@ -515,8 +513,7 @@ namespace Sacado {
KOKKOS_INLINE_FUNCTION
static void strided_zero(T* dest, int stride, int sz) {
for (int i=threadIdx.x; i<sz; i+=blockDim.x) {
*(dest) = T(0.);
dest += stride;
dest[i*stride] = T(0.);
}
}

Expand Down
111 changes: 102 additions & 9 deletions packages/sacado/test/UnitTests/Fad_KokkosTests.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -176,19 +176,15 @@ struct MultiplyKernel {
};

// Kernel to assign a constant to a view
template <typename ViewType, typename ScalarType>
template <typename ViewType>
struct ScalarAssignKernel {
typedef typename ViewType::execution_space execution_space;
typedef typename ViewType::size_type size_type;
typedef typename ViewType::value_type::value_type ScalarType;
typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
typedef Kokkos::RangePolicy< execution_space> range_policy_type;
typedef typename team_policy_type::member_type team_handle;
typedef typename Kokkos::ThreadLocalScalarType<ViewType>::type local_scalar_type;
#if defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;
#else
static const size_type stride = 32;
#endif

const ViewType m_v;
const ScalarType m_s;
Expand All @@ -199,8 +195,7 @@ struct ScalarAssignKernel {
// Multiply entries for row 'i' with a value
KOKKOS_INLINE_FUNCTION
void operator() (const size_type i) const {
local_scalar_type s = Sacado::partition_scalar<stride>(m_s);
m_v(i) = s;
m_v(i) = m_s;
}

KOKKOS_INLINE_FUNCTION
Expand Down Expand Up @@ -243,6 +238,71 @@ struct ScalarAssignKernel {
}
};

// Kernel to assign a constant to a view
template <typename ViewType>
struct ValueAssignKernel {
typedef typename ViewType::execution_space execution_space;
typedef typename ViewType::size_type size_type;
typedef typename ViewType::value_type ValueType;
typedef Kokkos::TeamPolicy< execution_space> team_policy_type;
typedef Kokkos::RangePolicy< execution_space> range_policy_type;
typedef typename team_policy_type::member_type team_handle;
typedef typename Kokkos::ThreadLocalScalarType<ViewType>::type local_scalar_type;
static const size_type stride = Kokkos::ViewScalarStride<ViewType>::stride;

const ViewType m_v;
const ValueType m_s;

ValueAssignKernel(const ViewType& v, const ValueType& s) :
m_v(v), m_s(s) {};

// Multiply entries for row 'i' with a value
KOKKOS_INLINE_FUNCTION
void operator() (const size_type i) const {
local_scalar_type s = Sacado::partition_scalar<stride>(m_s);
m_v(i) = s;
}

KOKKOS_INLINE_FUNCTION
void operator()( const team_handle& team ) const
{
const size_type i = team.league_rank()*team.team_size() + team.team_rank();
if (i < m_v.dimension_0())
(*this)(i);
}

// Kernel launch
static void apply(const ViewType& v, const ValueType& s) {
const size_type nrow = v.dimension_0();

#if defined (KOKKOS_HAVE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
const bool use_team =
std::is_same<execution_space, Kokkos::Cuda>::value &&
( Kokkos::is_view_fad_contiguous<ViewType>::value ||
Kokkos::is_dynrankview_fad_contiguous<ViewType>::value ) &&
( stride > 1 );
#elif defined (KOKKOS_HAVE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL_DFAD)
const bool use_team =
std::is_same<execution_space, Kokkos::Cuda>::value &&
( Kokkos::is_view_fad_contiguous<ViewType>::value ||
Kokkos::is_dynrankview_fad_contiguous<ViewType>::value ) &&
is_dfad<typename ViewType::non_const_value_type>::value;
#else
const bool use_team = false;
#endif

if (use_team) {
const size_type team_size = 256 / stride;
team_policy_type policy( (nrow+team_size-1)/team_size, team_size, stride );
Kokkos::parallel_for( policy, ValueAssignKernel(v,s) );
}
else {
range_policy_type policy( 0, nrow );
Kokkos::parallel_for( policy, ValueAssignKernel(v,s) );
}
}
};

// Kernel to assign a column of a rank-2 to a rank-1
template <typename InputViewType,
typename OutputViewType,
Expand Down Expand Up @@ -533,7 +593,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(

// Deep copy a constant scalar
value_type a = 2.3456;
ScalarAssignKernel<ViewType,value_type>::apply( v, a );
ScalarAssignKernel<ViewType>::apply( v, a );

// Copy to host
host_view_type hv = Kokkos::create_mirror_view(v);
Expand All @@ -551,6 +611,38 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
}
}

TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
Kokkos_View_Fad, ValueAssign, FadType, Layout, Device )
{
typedef Kokkos::View<FadType*,Layout,Device> ViewType;
typedef typename ViewType::size_type size_type;
typedef typename ViewType::HostMirror host_view_type;

const size_type num_rows = global_num_rows;
const size_type fad_size = global_fad_size;

// Create and fill view
ViewType v("view", num_rows, fad_size+1);
typename ViewType::array_type va = v;
Kokkos::deep_copy( va, 1.0 );

// Deep copy a constant scalar
FadType a(fad_size, 2.3456);
for (size_type i=0; i<fad_size; ++i)
a.fastAccessDx(i) = 7.89+i;
ValueAssignKernel<ViewType>::apply( v, a );

// Copy to host
host_view_type hv = Kokkos::create_mirror_view(v);
Kokkos::deep_copy(hv, v);

// Check
success = true;
for (size_type i=0; i<num_rows; ++i) {
success = success && checkFads(a, hv(i), out);
}
}

TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
Kokkos_View_Fad, Multiply, FadType, Layout, Device )
{
Expand Down Expand Up @@ -1603,6 +1695,7 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy_ConstantFad, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DeepCopy_ConstantFadFull, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ScalarAssign, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, ValueAssign, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Unmanaged, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Unmanaged2, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, UnmanagedConst, F, L, D ) \
Expand Down

0 comments on commit c6c840d

Please sign in to comment.