//---------------------------------------------------------------------------
+namespace internal
+{
+ namespace ArrayViewHelper
+ {
+ template <typename MemorySpaceType>
+ inline bool
+ correct_memory_space(const void *const)
+ {
+ return std::is_same<MemorySpaceType, MemorySpace::Host>::value;
+ }
+
+#ifdef DEAL_II_COMPILER_CUDA_AWARE
+ template <>
+ inline bool
+ correct_memory_space<MemorySpace::Host>(const void *ptr)
+ {
+ cudaPointerAttributes attributes;
+ const cudaError_t cuda_error = cudaPointerGetAttributes(&attributes, ptr);
+ if (cuda_error != cudaErrorInvalidValue)
+ {
+ AssertCuda(cuda_error);
+ return attributes.memoryType == cudaMemoryTypeHost;
+ }
+ else
+ {
+ // ignore and reset the error since host pointers produce an error
+ cudaGetLastError();
+ return true;
+ }
+ }
+
+ template <>
+ inline bool
+ correct_memory_space<MemorySpace::CUDA>(const void *ptr)
+ {
+ cudaPointerAttributes attributes;
+ const cudaError_t cuda_error = cudaPointerGetAttributes(&attributes, ptr);
+ if (cuda_error != cudaErrorInvalidValue)
+ {
+ AssertCuda(cuda_error);
+ return attributes.memoryType == cudaMemoryTypeDevice;
+ }
+ else
+ {
+ // ignore and reset the error since host pointers produce an error
+ cudaGetLastError();
+ return false;
+ }
+ }
+#endif
+ } // namespace ArrayViewHelper
+} // namespace internal
+
+
template <typename ElementType, typename MemorySpaceType>
inline ArrayView<ElementType, MemorySpaceType>::ArrayView(
const std::size_t n_elements)
: starting_element(starting_element)
, n_elements(n_elements)
-{}
+{
+ Assert(internal::ArrayViewHelper::correct_memory_space<MemorySpaceType>(
+ starting_element),
+ ExcMessage(
+ "The memory space indicated by the template parameter "
+ "and the one derived from the pointer value do not match!"));
+}
"object has a const value_type. In other words, you can "
"only create an ArrayView to const values from a const "
"std::vector.");
+ Assert(internal::ArrayViewHelper::correct_memory_space<MemorySpaceType>(
+ vector.data()),
+ ExcMessage(
+ "The memory space indicated by the template parameter "
+ "and the one derived from the pointer value do not match!"));
}
return n_elements;
}
+
+
template <typename ElementType, typename MemorySpaceType>
inline typename ArrayView<ElementType, MemorySpaceType>::iterator
ArrayView<ElementType, MemorySpaceType>::begin() const
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// check that we detect creating ArrayView objects wit the wrong memory space.
+
+#include <deal.II/base/array_view.h>
+
+#include "../tests.h"
+
+int
+main(int argc, char **argv)
+{
+ deal_II_exceptions::disable_abort_on_exception();
+
+ initlog();
+
+ std::vector<unsigned int> dummy_host(2);
+ std::unique_ptr<unsigned int[], void (*)(unsigned int *)> dummy_cuda(
+ Utilities::CUDA::allocate_device_data<unsigned int>(2),
+ Utilities::CUDA::delete_device_data<unsigned int>);
+
+ deallog << "Testing host ArrayView with host memory" << std::endl;
+ ArrayView<unsigned int, MemorySpace::Host> view_1(dummy_host);
+
+ deallog << "Testing device ArrayView with host memory" << std::endl;
+ try
+ {
+ ArrayView<unsigned int, MemorySpace::CUDA> view_2(dummy_host);
+ }
+ catch (const ExceptionBase &exc)
+ {
+ deallog << exc.what() << std::endl;
+ }
+
+ deallog << "Testing host ArrayView with device memory" << std::endl;
+ try
+ {
+ ArrayView<unsigned int, MemorySpace::Host> view_3(dummy_cuda.get(), 2);
+ }
+ catch (const ExceptionBase &exc)
+ {
+ deallog << exc.what() << std::endl;
+ }
+
+ deallog << "Testing device ArrayView with device memory" << std::endl;
+ ArrayView<unsigned int, MemorySpace::CUDA> view_4(dummy_cuda.get(), 2);
+
+ return 0;
+}
--- /dev/null
+
+DEAL::Testing host ArrayView with host memory
+DEAL::Testing device ArrayView with host memory
+DEAL::
+--------------------------------------------------------
+An error occurred in file <array_view.h> in function
+ dealii::ArrayView<ElementType, MemorySpaceType>::ArrayView(dealii::ArrayView<ElementType, MemorySpaceType>::value_type*, std::size_t) [with ElementType = unsigned int; MemorySpaceType = dealii::MemorySpace::CUDA; dealii::ArrayView<ElementType, MemorySpaceType>::value_type = unsigned int; std::size_t = long unsigned int]
+The violated condition was:
+ internal::ArrayViewHelper::correct_memory_space<MemorySpaceType>( starting_element)
+Additional information:
+ The memory space indicated by the template parameter and the one derived from the pointer value do not match!
+--------------------------------------------------------
+
+DEAL::Testing host ArrayView with device memory
+DEAL::
+--------------------------------------------------------
+An error occurred in file <array_view.h> in function
+ dealii::ArrayView<ElementType, MemorySpaceType>::ArrayView(dealii::ArrayView<ElementType, MemorySpaceType>::value_type*, std::size_t) [with ElementType = unsigned int; MemorySpaceType = dealii::MemorySpace::Host; dealii::ArrayView<ElementType, MemorySpaceType>::value_type = unsigned int; std::size_t = long unsigned int]
+The violated condition was:
+ internal::ArrayViewHelper::correct_memory_space<MemorySpaceType>( starting_element)
+Additional information:
+ The memory space indicated by the template parameter and the one derived from the pointer value do not match!
+--------------------------------------------------------
+
+DEAL::Testing device ArrayView with device memory