Skip to content

Commit

Permalink
Merge pull request #9857 from NexGenAnalytics/Sacado-UVM-off
Browse files Browse the repository at this point in the history
Sacado: Re-enable and fix Sacado tests for UVM-off PR testers
  • Loading branch information
etphipp authored Oct 29, 2021
2 parents 91e55fa + 6940295 commit b414ad3
Show file tree
Hide file tree
Showing 9 changed files with 52 additions and 53 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -148,7 +148,6 @@ set (Domi_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (Kokkos_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (KokkosKernels_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (ROL_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (Sacado_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (SEACAS_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (ShyLU_DD_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
set (STK_ENABLE_TESTS OFF CACHE BOOL "Turn off tests for non-UVM build")
Expand Down
4 changes: 0 additions & 4 deletions packages/sacado/test/UnitTests/Fad_Fad_KokkosTests_Cuda.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -35,11 +35,7 @@
#include "Kokkos_Core.hpp"

// Instantiate tests for Cuda device
#if defined(KOKKOS_ENABLE_CUDA_UVM)
#define SACADO_TEST_DFAD 1
#else
#define SACADO_TEST_DFAD 0
#endif
using Kokkos::Cuda;
VIEW_FAD_TESTS_D( Cuda )

Expand Down
20 changes: 13 additions & 7 deletions packages/sacado/test/UnitTests/Fad_KokkosAtomicTests.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -261,26 +261,30 @@ bool testAtomic(const TagType& tag, Teuchos::FancyOStream& out)

// Create and fill view
ViewType v;
ScalarViewType s0;
#if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
v = ViewType ("view", num_rows);
v = ViewType ("view", num_rows);
s0 = ScalarViewType ("");
#else
v = ViewType ("view", num_rows, fad_size+1);
v = ViewType ("view", num_rows, fad_size+1);
s0 = ScalarViewType ("", fad_size+1);
#endif
host_view_type h_v = Kokkos::create_mirror_view(v);
for (size_type i=0; i<num_rows; ++i)
h_v(i) =
generate_fad<FadType>(num_rows, size_type(1), fad_size, i, size_type(0));
Kokkos::deep_copy(v, h_v);

Kokkos::deep_copy(s0, tag.init());

// Create scalar view
ScalarViewType s;
FadType s0 = FadType(fad_size,tag.init());
#if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
s = ScalarViewType ("scalar view");
#else
s = ScalarViewType ("scalar view", fad_size+1);
#endif
Kokkos::deep_copy( s, s0 );
Kokkos::deep_copy( s, tag.init() );

// Call atomic_add kernel, which adds up entries in v
AtomicKernel<ViewType,ScalarViewType,OperFetch>::apply( tag, v, s );
Expand All @@ -290,12 +294,14 @@ bool testAtomic(const TagType& tag, Teuchos::FancyOStream& out)
Kokkos::deep_copy(hs, s);

// Compute correct result
FadType b = s0;
auto b = Kokkos::create_mirror_view(s0);
Kokkos::deep_copy(b, s0);

for (size_type i=0; i<num_rows; ++i)
b = tag.apply(b, h_v(i));
b() = tag.apply(b(), h_v(i));

// Check
bool success = checkFads(b, hs(), out);
bool success = checkFads(b(), hs(), out);

return success;
}
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -32,11 +32,8 @@

#include "Kokkos_Macros.hpp"

#if defined(KOKKOS_ENABLE_CUDA_UVM)
#define SACADO_TEST_DFAD 1
#else
#define SACADO_TEST_DFAD 0
#endif

#include "Fad_KokkosAtomicTests.hpp"

// Instantiate tests for Cuda device. We can only test DFad is UVM is enabled.
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -37,11 +37,8 @@

#include "Kokkos_Macros.hpp"

#if defined(KOKKOS_ENABLE_CUDA_UVM)
#define SACADO_TEST_DFAD 1
#else
#define SACADO_TEST_DFAD 0
#endif

#include "Fad_KokkosAtomicTests.hpp"

typedef Kokkos::LayoutContiguous<Kokkos::LayoutLeft,32> LeftContiguous32;
Expand Down
56 changes: 35 additions & 21 deletions packages/sacado/test/UnitTests/Fad_KokkosTests.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -239,7 +239,7 @@ struct ScalarAssignKernel {
};

// Kernel to assign a constant to a view
template <typename ViewType>
template <typename ViewType, typename ScalarViewType>
struct ValueAssignKernel {
typedef typename ViewType::execution_space execution_space;
typedef typename ViewType::size_type size_type;
Expand All @@ -250,16 +250,16 @@ struct ValueAssignKernel {
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;
const ViewType m_v;
const ScalarViewType m_s;

ValueAssignKernel(const ViewType& v, const ValueType& s) :
ValueAssignKernel(const ViewType& v, const ScalarViewType& 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);
local_scalar_type s = Sacado::partition_scalar<stride>(m_s());
m_v(i) = s;
}

Expand All @@ -272,7 +272,7 @@ struct ValueAssignKernel {
}

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

#if defined (KOKKOS_ENABLE_CUDA) && defined (SACADO_VIEW_CUDA_HIERARCHICAL)
Expand Down Expand Up @@ -654,11 +654,10 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
FadType a(fad_size, 2.3456);
for (size_type i=0; i<fad_size; ++i)
a.fastAccessDx(i) = 7.89 + (i+1);
Kokkos::deep_copy( v, a );

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

// Check
success = true;
Expand Down Expand Up @@ -714,36 +713,47 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
Kokkos_View_Fad, ValueAssign, FadType, Layout, Device )
{
typedef Kokkos::View<FadType*,Layout,Device> ViewType;
typedef Kokkos::View<FadType,Layout,Device> ScalarViewType;
typedef typename ViewType::size_type size_type;
typedef typename ViewType::HostMirror host_view_type;
typedef typename ScalarViewType::HostMirror host_scalar_view_type;

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

// Create and fill view
ViewType v;
ScalarViewType a;
#if defined (SACADO_DISABLE_FAD_VIEW_SPEC)
v = ViewType ("view", num_rows);
a = ScalarViewType ("fad");
#else
v = ViewType ("view", num_rows, fad_size+1);
a = ScalarViewType ("fad", fad_size+1);
#endif
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 );
Kokkos::deep_copy(a, 2.3456);

Kokkos::parallel_for(Kokkos::RangePolicy< Device>(0, fad_size), KOKKOS_LAMBDA(const int i) {
a().fastAccessDx(i) = 7.89 + i;
});

ValueAssignKernel<ViewType, ScalarViewType>::apply( v, a );

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

host_scalar_view_type ha = Kokkos::create_mirror_view(a);
Kokkos::deep_copy(ha, a);

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

Expand Down Expand Up @@ -1016,10 +1026,12 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
#else
v = ViewType ("view", num_rows, fad_size+1);
#endif
FadType a(fad_size, 2.3456);
for (size_type i=0; i<fad_size; ++i)
a.fastAccessDx(i) = 7.89+i;
Kokkos::deep_copy( v, a );
Kokkos::deep_copy(v, 2.3456);

Kokkos::parallel_for(Kokkos::RangePolicy<Device>(0, num_rows), KOKKOS_LAMBDA(const size_type i) {
for (size_type j = 0; j < fad_size; ++j)
v(i).fastAccessDx(j) = 7.89 + j;
});

// Create scalar view
ScalarViewType s;
Expand All @@ -1028,17 +1040,19 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
#else
s = ScalarViewType ("scalar view", fad_size+1);
#endif
Kokkos::deep_copy( s, FadType(fad_size,0.0) );

// Call atomic_add kernel, which adds up entries in v
// Call atomic_add kernel, which adds up entries in v
AtomicAddKernel<ViewType,ScalarViewType>::apply( v, s );

// Copy to host
host_scalar_view_type hs = Kokkos::create_mirror_view(s);
Kokkos::deep_copy(hs, s);

// Check
FadType b = num_rows*a;
auto hv = Kokkos::create_mirror_view(v);
Kokkos::deep_copy(hv, v);

FadType b = num_rows*hv(0);
success = checkFads(b, hs(), out);
}

Expand Down Expand Up @@ -2232,9 +2246,9 @@ TEUCHOS_UNIT_TEST_TEMPLATE_3_DECL(
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, MultiplyUpdate, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, MultiplyConst, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, MultiplyMixed, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AtomicAdd, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Rank8, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, Roger, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AtomicAdd, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, AssignDifferentStrides, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankDimensionScalar, F, L, D ) \
TEUCHOS_UNIT_TEST_TEMPLATE_3_INSTANT( Kokkos_View_Fad, DynRankAssignStatic0, F, L, D ) \
Expand Down
4 changes: 0 additions & 4 deletions packages/sacado/test/UnitTests/Fad_KokkosTests_Cuda.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -32,11 +32,7 @@

#include "Kokkos_Macros.hpp"

#if defined(KOKKOS_ENABLE_CUDA_UVM)
#define SACADO_TEST_DFAD 1
#else
#define SACADO_TEST_DFAD 0
#endif
#include "Fad_KokkosTests.hpp"

// Instantiate tests for Cuda device. We can only test DFad is UVM is enabled.
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -37,11 +37,8 @@

#include "Kokkos_Macros.hpp"

#if defined(KOKKOS_ENABLE_CUDA_UVM)
#define SACADO_TEST_DFAD 1
#else
#define SACADO_TEST_DFAD 0
#endif

#include "Fad_KokkosTests.hpp"

typedef Kokkos::LayoutContiguous<Kokkos::LayoutLeft,32> LeftContiguous32;
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -35,11 +35,8 @@

#include "Kokkos_Macros.hpp"

#if defined(KOKKOS_ENABLE_CUDA_UVM)
#define SACADO_TEST_DFAD 1
#else
#define SACADO_TEST_DFAD 0
#endif

#include "Fad_KokkosTests.hpp"

// Instantiate tests for Cuda device
Expand Down

0 comments on commit b414ad3

Please sign in to comment.