From fdd71bdb2a030cf219c78de7937163075a78f4a2 Mon Sep 17 00:00:00 2001 From: turcksin Date: Mon, 9 Jun 2014 21:13:28 +0000 Subject: [PATCH] Add new wrappers to initialize Paralution and modify the current tests to use them. git-svn-id: https://svn.dealii.org/branches/branch_paralution@33028 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/include/deal.II/base/paralution.h | 175 ++++++++++++++++++ .../deal.II/lac/paralution_sparse_matrix.h | 7 +- .../source/lac/paralution_sparse_matrix.cc | 50 +++++ tests/lac/paralution_solver_01.cc | 5 +- tests/lac/paralution_solver_02.cc | 5 +- tests/lac/paralution_sparse_matrix_01.cc | 6 +- tests/lac/paralution_vector_01.cc | 6 +- tests/lac/paralution_vector_02.cc | 6 +- tests/lac/paralution_vector_03.cc | 6 +- tests/lac/paralution_vector_04.cc | 6 +- tests/lac/paralution_vector_05.cc | 6 +- tests/lac/paralution_vector_06.cc | 6 +- tests/tests.h | 5 +- 13 files changed, 252 insertions(+), 37 deletions(-) create mode 100644 deal.II/include/deal.II/base/paralution.h diff --git a/deal.II/include/deal.II/base/paralution.h b/deal.II/include/deal.II/base/paralution.h new file mode 100644 index 0000000000..315df97e77 --- /dev/null +++ b/deal.II/include/deal.II/base/paralution.h @@ -0,0 +1,175 @@ +// --------------------------------------------------------------------- +// $Id: paralution.h 32926 2014-05-16 17:12:16Z turcksin $ +// +// Copyright (C) 2014 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + +#ifndef __deal2__paralution_h +#define __deal2__paralution_h + +#include + +#ifdef DEAL_II_WITH_PARALUTION + +#include + +#include + +DEAL_II_NAMESPACE_OPEN + + +namespace Utilities +{ + /** + * A namespace for utility functions that abstract certain operations using + * the Paralution library. + */ + namespace Paralution + { + /** + * A class that is used to initialize and to stop Paralution. If a program + * uses Paralution one would typically just create an object of this type at + * the beginninh of main(). The constructor of this class then + * runs paralution::init_paralution(). At the end of the + * program, the compiler will invoke the destructor of this object which in + * turns calls paralution::stop_paralution(). + */ + class Paralution_InitFinalize + { + public : + /** + * Constructor. Initialize Paralution by calling paralution::init_paralution + * and set the number of threads used by Paralution to the given parameters. + */ + Paralution_InitFinalize(const unsigned int max_num_threads); + + /** + * Destructor. Calls paralution::stop_paralution. + */ + ~Paralution_InitFinalize(); + + /** + * Set a device. + */ + void set_device(const unsigned int device) const; + + /** + * Set a specific gpu. + */ + void set_gpu_cuda(const unsigned int gpu) const; + + /** + * Set OpenCL compute units. + */ + void set_ocl_compute_units(const unsigned int compute_unit) const; + + /** + * Set a specific OpenCL platform. + */ + void set_ocl_platform(const unsigned int platform) const; + + /** + * Set a specific OpenCL platform and device. + */ + void set_ocl(const unsigned int platform, const unsigned int device) const; + + /** + * Set the number of OpenMP thread. + */ + void set_omp_threads(const unsigned int n_threads) const; + + /** + * Print information about the platform. + */ + void info() const; + }; + + + +// ------------------- inline functions -------------- + inline Paralution_InitFinalize::Paralution_InitFinalize(const unsigned int max_num_threads) + { + // Initialize paralution + paralution::init_paralution(); + + // Set the number of OpenMP threads + set_omp_threads(max_num_threads); + + // Set the number of TBB threads + multithread_info.set_thread_limit(max_num_threads); + } + + + inline Paralution_InitFinalize::~Paralution_InitFinalize() + { + paralution::stop_paralution(); + } + + + + inline void Paralution_InitFinalize::set_device(const unsigned int device) const + { + paralution::set_device_paralution(device); + } + + + + inline void Paralution_InitFinalize::set_gpu_cuda(const unsigned int gpu) const + { + paralution::set_gpu_cuda_paralution(gpu); + } + + + + inline void Paralution_InitFinalize::set_ocl_compute_units(const unsigned int compute_unit) const + { + paralution::set_ocl_compute_units_paralution(compute_unit); + } + + + + inline void Paralution_InitFinalize::set_ocl_platform(const unsigned int platform) const + { + paralution::set_ocl_platform_paralution(platform); + } + + + + inline void Paralution_InitFinalize::set_ocl(const unsigned int platform, + const unsigned int device) const + { + paralution::set_ocl_paralution(platform, device); + } + + + + inline void Paralution_InitFinalize::set_omp_threads(const unsigned int n_threads) const + { + paralution::set_omp_threads_paralution(n_threads); + } + + + + inline void Paralution_InitFinalize::info() const + { + paralution::info_paralution(); + } + } +} + + +DEAL_II_NAMESPACE_CLOSE + +#endif // DEAL_II_WITH_PARALUTION + +#endif diff --git a/deal.II/include/deal.II/lac/paralution_sparse_matrix.h b/deal.II/include/deal.II/lac/paralution_sparse_matrix.h index c98e9d84f1..40f365597b 100644 --- a/deal.II/include/deal.II/lac/paralution_sparse_matrix.h +++ b/deal.II/include/deal.II/lac/paralution_sparse_matrix.h @@ -24,7 +24,6 @@ #include #include #include -#include #include @@ -371,6 +370,12 @@ namespace ParalutionWrappers */ void sync(); + /** + * Convert the paralution matrix to a different format. This function should + * only be called after convert_to_paralution_csr. + */ + void convert_format(matrix_format format); + /** * Return a constant reference to the underlying dealii::SparseMatrix. */ diff --git a/deal.II/source/lac/paralution_sparse_matrix.cc b/deal.II/source/lac/paralution_sparse_matrix.cc index 7f20e32150..530f34b1a8 100644 --- a/deal.II/source/lac/paralution_sparse_matrix.cc +++ b/deal.II/source/lac/paralution_sparse_matrix.cc @@ -137,6 +137,56 @@ namespace ParalutionWrappers else local_matrix.CopyFromAsync(sparse_matrix.paralution_matrix()); } + + template + void SparseMatrix::convert_format(matrix_format format) + { + if (is_local_matrix==false) + { + switch (format) + { + case DENSE : + { + local_matrix.ConvertToDense(); + break; + } + case CSR : + { + local_matrix.ConvertToCSR(); + break; + } + case MCSR : + { + local_matrix.ConvertToMCSR(); + break; + } + case BCSR : + { + local_matrix.ConvertToBCSR(); + break; + } + case COO : + { + local_matrix.ConvertToDIA(); + break; + } + case ELL : + { + local_matrix.ConvertToELL(); + break; + } + case HYB : + { + local_matrix.ConvertToHYB(); + break; + } + default : + { + AssertThrow(false,ExcMessage("Wrong format of Paralution matrix.")); + } + } + } + } } // Explicit instantiations diff --git a/tests/lac/paralution_solver_01.cc b/tests/lac/paralution_solver_01.cc index 720f95e074..8057394373 100644 --- a/tests/lac/paralution_solver_01.cc +++ b/tests/lac/paralution_solver_01.cc @@ -23,6 +23,7 @@ #include #include #include +#include #include #include #include @@ -30,7 +31,7 @@ int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -69,6 +70,4 @@ int main() } solver.solve(A,u,f,p); } - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_solver_02.cc b/tests/lac/paralution_solver_02.cc index 7bc614a68e..34b3cca2cd 100644 --- a/tests/lac/paralution_solver_02.cc +++ b/tests/lac/paralution_solver_02.cc @@ -23,6 +23,7 @@ #include #include #include +#include #include #include #include @@ -30,7 +31,7 @@ int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -69,6 +70,4 @@ int main() } solver.solve(A,u,f,p); } - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_sparse_matrix_01.cc b/tests/lac/paralution_sparse_matrix_01.cc index 44335a3de8..ca8b7661ea 100644 --- a/tests/lac/paralution_sparse_matrix_01.cc +++ b/tests/lac/paralution_sparse_matrix_01.cc @@ -18,9 +18,9 @@ // Test ParalutionWrappers::SparseMatrix #include "../tests.h" +#include #include #include -#include "paralution.hpp" template void check() @@ -55,7 +55,7 @@ void check() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -66,7 +66,5 @@ int main() check(); check(); - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_vector_01.cc b/tests/lac/paralution_vector_01.cc index 73311620d6..e7997c38d9 100644 --- a/tests/lac/paralution_vector_01.cc +++ b/tests/lac/paralution_vector_01.cc @@ -18,8 +18,8 @@ // Test ParalutionWrappers::Vector constructor #include "../tests.h" +#include #include -#include "paralution.hpp" template void check() @@ -35,7 +35,7 @@ void check() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -47,6 +47,4 @@ int main() check(); check(); check(); - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_vector_02.cc b/tests/lac/paralution_vector_02.cc index d3d3068465..05adabc943 100644 --- a/tests/lac/paralution_vector_02.cc +++ b/tests/lac/paralution_vector_02.cc @@ -18,8 +18,8 @@ // Test copy constructor and access elements of ParalutionWrappers::Vector #include "../tests.h" +#include #include -#include "paralution.hpp" void check() { @@ -39,7 +39,7 @@ void check() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -49,6 +49,4 @@ int main() deallog.threshold_double(1.e-10); check(); - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_vector_03.cc b/tests/lac/paralution_vector_03.cc index d4f298679d..51ff86d802 100644 --- a/tests/lac/paralution_vector_03.cc +++ b/tests/lac/paralution_vector_03.cc @@ -18,8 +18,8 @@ // Test ParalutionWrappers::Vector iterator #include "../tests.h" +#include #include -#include "paralution.hpp" void check() { @@ -40,7 +40,7 @@ void check() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -50,6 +50,4 @@ int main() deallog.threshold_double(1.e-10); check(); - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_vector_04.cc b/tests/lac/paralution_vector_04.cc index c21c6fd418..6a77a6db57 100644 --- a/tests/lac/paralution_vector_04.cc +++ b/tests/lac/paralution_vector_04.cc @@ -18,8 +18,8 @@ // Test reinit for ParalutionWrappers::Vector #include "../tests.h" +#include #include -#include "paralution.hpp" void check() { @@ -39,7 +39,7 @@ void check() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -49,6 +49,4 @@ int main() deallog.threshold_double(1.e-10); check(); - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_vector_05.cc b/tests/lac/paralution_vector_05.cc index 6e55260b08..79bd2df364 100644 --- a/tests/lac/paralution_vector_05.cc +++ b/tests/lac/paralution_vector_05.cc @@ -19,8 +19,8 @@ // operator-= #include "../tests.h" +#include #include -#include "paralution.hpp" void check_1() { @@ -58,7 +58,7 @@ void check_2() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -69,6 +69,4 @@ int main() check_1(); check_2(); - - paralution::stop_paralution(); } diff --git a/tests/lac/paralution_vector_06.cc b/tests/lac/paralution_vector_06.cc index fe00d96ed6..6ac90a3ab0 100644 --- a/tests/lac/paralution_vector_06.cc +++ b/tests/lac/paralution_vector_06.cc @@ -19,8 +19,8 @@ #include "../tests.h" +#include #include -#include "paralution.hpp" void check() { @@ -37,7 +37,7 @@ void check() int main() { - paralution::init_paralution(); + Utilities::Paralution::Paralution_InitFinalize paralution(1); std::ofstream logfile("output"); deallog << std::fixed; @@ -47,6 +47,4 @@ int main() deallog.threshold_double(1.e-10); check(); - - paralution::stop_paralution(); } diff --git a/tests/tests.h b/tests/tests.h index ef37ca2046..9d8680ed62 100644 --- a/tests/tests.h +++ b/tests/tests.h @@ -157,9 +157,10 @@ void unify_pretty_function (const std::string &filename) * Note that we can't do this if we run in MPI mode because then * MPI_InitFinalize already calls this function. Since every test * calls MPI_InitFinalize itself, we can't adjust the thread count - * for this here. + * for this here. The same is true for Paralution, Paralution_InitFinalize + * calls set_thread_limit. */ -#ifndef DEAL_II_WITH_MPI +#if !(defined DEAL_II_WITH_MPI) && !(defined DEAL_II_WITH_PARALUTION) struct LimitConcurrency { LimitConcurrency () -- 2.39.5