}
#endif
+// There is a similar issue with CUDA: The destructor of static objects might
+// run after the CUDA driver is unloaded. Hence, also release all memory
+// related to CUDA vectors.
+#ifdef DEAL_II_WITH_CUDA
+ GrowingVectorMemory<
+ LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>::
+ release_unused_memory();
+ GrowingVectorMemory<
+ LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>::
+ release_unused_memory();
+#endif
+
#ifdef DEAL_II_WITH_P4EST
// now end p4est and libsc
// Note: p4est has no finalize function
sparse_matrix_inst2.cc
)
+# Add CUDA wrapper files
+IF(DEAL_II_WITH_CUDA)
+ SET(_separate_src
+ ${_separate_src}
+ vector_memory.cu
+ )
+ENDIF()
+
SET(_inst
affine_constraints.inst.in
block_sparse_matrix.inst.in
// ---------------------------------------------------------------------
#include <deal.II/lac/block_vector.h>
-#include <deal.II/lac/cuda_vector.h>
#include <deal.II/lac/la_parallel_block_vector.h>
#include <deal.II/lac/la_parallel_vector.h>
#include <deal.II/lac/la_vector.h>
#include "vector_memory.inst"
-#ifdef DEAL_II_WITH_CUDA
-template class VectorMemory<LinearAlgebra::CUDAWrappers::Vector<float>>;
-template class VectorMemory<LinearAlgebra::CUDAWrappers::Vector<double>>;
-template class GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<float>>;
-template class GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<double>>;
-# ifdef DEAL_II_WITH_MPI
-template class VectorMemory<
- LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>;
-template class VectorMemory<
- LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>;
-template class GrowingVectorMemory<
- LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>;
-template class GrowingVectorMemory<
- LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>;
-# endif
-#endif
-
namespace internal
{
namespace GrowingVectorMemoryImplementation
{
+#ifdef DEAL_II_WITH_CUDA
+ void
+ release_all_unused_cuda_memory();
+#endif
+
void
release_all_unused_memory()
{
#include "vector_memory_release.inst"
-
#ifdef DEAL_II_WITH_CUDA
- dealii::GrowingVectorMemory<dealii::LinearAlgebra::CUDAWrappers::Vector<
- float>>::release_unused_memory();
- dealii::GrowingVectorMemory<dealii::LinearAlgebra::CUDAWrappers::Vector<
- double>>::release_unused_memory();
-# ifdef DEAL_II_WITH_MPI
- dealii::GrowingVectorMemory<
- dealii::LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>::
- release_unused_memory();
- dealii::GrowingVectorMemory<
- dealii::LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>::
- release_unused_memory();
-# endif
+ release_all_unused_cuda_memory();
#endif
}
} // namespace GrowingVectorMemoryImplementation
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2019 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.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/lac/cuda_vector.h>
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/vector_memory.templates.h>
+
+
+DEAL_II_NAMESPACE_OPEN
+
+template class VectorMemory<LinearAlgebra::CUDAWrappers::Vector<float>>;
+template class VectorMemory<LinearAlgebra::CUDAWrappers::Vector<double>>;
+template class GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<float>>;
+template class GrowingVectorMemory<LinearAlgebra::CUDAWrappers::Vector<double>>;
+template class VectorMemory<
+ LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>;
+template class VectorMemory<
+ LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>;
+template class GrowingVectorMemory<
+ LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>;
+template class GrowingVectorMemory<
+ LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>;
+
+namespace internal
+{
+ namespace GrowingVectorMemoryImplementation
+ {
+ void
+ release_all_unused_cuda_memory()
+ {
+ dealii::GrowingVectorMemory<dealii::LinearAlgebra::CUDAWrappers::Vector<
+ float>>::release_unused_memory();
+ dealii::GrowingVectorMemory<dealii::LinearAlgebra::CUDAWrappers::Vector<
+ double>>::release_unused_memory();
+ dealii::GrowingVectorMemory<
+ dealii::LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>::
+ release_unused_memory();
+ dealii::GrowingVectorMemory<
+ dealii::LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>::
+ release_unused_memory();
+ }
+ } // namespace GrowingVectorMemoryImplementation
+} // namespace internal
+
+DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2007 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test that we can successfully fill a GrowingVectorMemory pool
+// with LinearAlgebra::distributed::Vector<Number, MemorySpace::CUDA> objects.
+// Partially copied from lac/vector_memory.cc
+
+
+#include <deal.II/base/exceptions.h>
+
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include "../tests.h"
+
+
+using namespace dealii;
+
+template <typename VectorType>
+void
+test_stat()
+{
+ GrowingVectorMemory<VectorType> mem(1, true);
+ VectorType * v1 = mem.alloc();
+ VectorType * v2 = mem.alloc();
+ VectorType * v3 = mem.alloc();
+ VectorType * v4 = mem.alloc();
+ VectorType * v5 = mem.alloc();
+ VectorType * v6 = mem.alloc();
+ v1->reinit(5);
+ v2->reinit(5);
+ v3->reinit(5);
+ v4->reinit(5);
+ v5->reinit(5);
+ v6->reinit(5);
+ mem.free(v1);
+ mem.free(v2);
+ mem.free(v3);
+ mem.free(v4);
+ mem.free(v5);
+ mem.free(v6);
+ v1 = mem.alloc();
+ mem.free(v1);
+ v1 = mem.alloc();
+ mem.free(v1);
+ v1 = mem.alloc();
+ mem.free(v1);
+ v1 = mem.alloc();
+ mem.free(v1);
+}
+
+
+int
+main(int argc, char *argv[])
+{
+ initlog();
+ Utilities::MPI::MPI_InitFinalize mpi_init(argc, argv, 1);
+
+ test_stat<LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>();
+ test_stat<LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>();
+
+ return 0;
+}
--- /dev/null
+
+DEAL::GrowingVectorMemory:Overall allocated vectors: 10
+DEAL::GrowingVectorMemory:Maximum allocated vectors: 6
+DEAL::GrowingVectorMemory:Overall allocated vectors: 10
+DEAL::GrowingVectorMemory:Maximum allocated vectors: 6
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2007 - 2017 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.
+//
+// ---------------------------------------------------------------------
+
+
+// Test that memory leaks are detected correctly for a GrowingVectorMemory pool
+// with LinearAlgebra::distributed<Number, MemorySpace::CUDA> objects.
+// Partially copied from lac/vector_memory.cc
+
+
+#include <deal.II/base/exceptions.h>
+
+#include <deal.II/lac/la_parallel_vector.h>
+#include <deal.II/lac/vector_memory.h>
+
+#include "../tests.h"
+
+
+using namespace dealii;
+
+template <typename VectorType>
+void
+test_leak()
+{
+ GrowingVectorMemory<VectorType> mem;
+ VectorType * v = mem.alloc();
+ v->reinit(5);
+}
+
+int
+main(int argc, char *argv[])
+{
+ initlog();
+ Utilities::MPI::MPI_InitFinalize mpi_init(argc, argv, 1);
+ deal_II_exceptions::disable_abort_on_exception();
+
+ try
+ {
+ test_leak<
+ LinearAlgebra::distributed::Vector<double, MemorySpace::CUDA>>();
+ test_leak<LinearAlgebra::distributed::Vector<float, MemorySpace::CUDA>>();
+ }
+ catch (const StandardExceptions::ExcMemoryLeak &e)
+ {
+ deallog << "Exception: " << e.get_exc_name() << std::endl;
+ }
+
+ return 0;
+}
--- /dev/null
+
+DEAL::Exception: StandardExceptions::ExcMemoryLeak(current_alloc)
+DEAL::
+--------------------------------------------------------
+An error occurred in file <vector_memory.templates.h> in function
+ dealii::GrowingVectorMemory<VectorType>::~GrowingVectorMemory() [with VectorType = dealii::LinearAlgebra::distributed::Vector<double, dealii::MemorySpace::CUDA>]
+The violated condition was:
+ current_alloc == 0
+Additional information:
+ Destroying memory handler while 1 objects are still allocated.
+--------------------------------------------------------
+
+DEAL::Exception: StandardExceptions::ExcMemoryLeak(current_alloc)
+DEAL::
+--------------------------------------------------------
+An error occurred in file <vector_memory.templates.h> in function
+ dealii::GrowingVectorMemory<VectorType>::~GrowingVectorMemory() [with VectorType = dealii::LinearAlgebra::distributed::Vector<float, dealii::MemorySpace::CUDA>]
+The violated condition was:
+ current_alloc == 0
+Additional information:
+ Destroying memory handler while 1 objects are still allocated.
+--------------------------------------------------------
+