From: Daniel Arndt Date: Thu, 1 Jun 2017 15:57:22 +0000 (+0200) Subject: Allow using VectorTools::subtract_mean_value for all distributed vector classes X-Git-Tag: v9.0.0-rc1~1544^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F4463%2Fhead;p=dealii.git Allow using VectorTools::subtract_mean_value for all distributed vector classes --- diff --git a/include/deal.II/numerics/vector_tools.h b/include/deal.II/numerics/vector_tools.h index b4c0699602..bdd53a4df1 100644 --- a/include/deal.II/numerics/vector_tools.h +++ b/include/deal.II/numerics/vector_tools.h @@ -2830,8 +2830,8 @@ namespace VectorTools * not equal to $(1,1,\ldots,1)^T$. For such elements, a different procedure * has to be used when subtracting the mean value. * - * @warning This function is only implemented for Vector and BlockVector. It - * is not implemented for any of the distributed vector classes. + * @warning This function can only be used for distributed vector classes + * provided the boolean mask is empty, i.e. selecting the whole vector. */ template void subtract_mean_value(VectorType &v, diff --git a/include/deal.II/numerics/vector_tools.templates.h b/include/deal.II/numerics/vector_tools.templates.h index 7ee26717ae..482940d769 100644 --- a/include/deal.II/numerics/vector_tools.templates.h +++ b/include/deal.II/numerics/vector_tools.templates.h @@ -7753,48 +7753,67 @@ namespace VectorTools return gradient[0]; } + namespace + { + template + typename std::enable_if::value == true>::type + internal_subtract_mean_value(VectorType &v, + const std::vector &p_select) + { + if (p_select.size() == 0) + { + // In case of an empty boolean mask operate on the whole vector: + v.add( - v.mean_value() ); + } + else + { + const unsigned int n = v.size(); + Assert(p_select.size() == n, + ExcDimensionMismatch(p_select.size(), n)); - template - void - subtract_mean_value(VectorType &v, - const std::vector &p_select) - { - if (p_select.size() == 0) - { - // In case of an empty boolean mask operate on the whole vector: - v.add( - v.mean_value() ); - } - else - { - // This function is not implemented for distributed vectors. - Assert((dealii::is_serial_vector< VectorType >::value == true), - ExcNotImplemented()); + typename VectorType::value_type s = 0.; + unsigned int counter = 0; + for (unsigned int i=0; i 0, ComponentMask::ExcNoComponentSelected()); - const unsigned int n = v.size(); + s /= counter; - Assert(p_select.size() == n, - ExcDimensionMismatch(p_select.size(), n)); + for (unsigned int i=0; i 0, ComponentMask::ExcNoComponentSelected()); - s /= counter; - for (unsigned int i=0; i + typename std::enable_if::value == false>::type + internal_subtract_mean_value(VectorType &v, + const std::vector &p_select) + { + (void) p_select; + Assert(p_select.size() == 0, ExcNotImplemented()); + // In case of an empty boolean mask operate on the whole vector: + v.add( - v.mean_value() ); + } + } + + + template + void + subtract_mean_value(VectorType &v, + const std::vector &p_select) + { + internal_subtract_mean_value(v,p_select); } namespace diff --git a/source/numerics/vector_tools_mean_value.cc b/source/numerics/vector_tools_mean_value.cc index 2861b76a98..43b8e51e98 100644 --- a/source/numerics/vector_tools_mean_value.cc +++ b/source/numerics/vector_tools_mean_value.cc @@ -18,19 +18,6 @@ DEAL_II_NAMESPACE_OPEN -#if defined(DEAL_II_WITH_TRILINOS) && defined(DEAL_II_WITH_MPI) -namespace VectorTools -{ - template <> - void - subtract_mean_value(LinearAlgebra::EpetraWrappers::Vector &, - const std::vector &) - { - Assert(false,ExcNotImplemented()); - } -} -#endif - // ---------------------------- explicit instantiations -------------------- #include "vector_tools_mean_value.inst" diff --git a/tests/numerics/subtract_mean_value_01.cc b/tests/numerics/subtract_mean_value_01.cc new file mode 100644 index 0000000000..3a938d9c35 --- /dev/null +++ b/tests/numerics/subtract_mean_value_01.cc @@ -0,0 +1,110 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// check VectorTools::subtract_mean_value() for deal.II serial vectors + +#include "../tests.h" +#include +#include +#include +#include +#include +#include +#include + +template +void test (VectorType &v) +{ + std::vector filter(v.size(), false); + // set some elements of the vector + for (unsigned int i=0; i v (10); + test (v); + } + + { + Vector v (10); + test (v); + } + + { + BlockVector v (std::vector(1,10)); + test (v); + } + + { + BlockVector v (std::vector< types::global_dof_index>(1,10)); + test (v); + } + + { + LinearAlgebra::Vector v(10); + test (v); + } + + { + LinearAlgebra::Vector v(10); + test (v); + } + } + catch (std::exception &exc) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Exception on processing: " << std::endl + << exc.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + + return 1; + } + catch (...) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + }; +} diff --git a/tests/numerics/subtract_mean_value_01.output b/tests/numerics/subtract_mean_value_01.output new file mode 100644 index 0000000000..0aa61ff573 --- /dev/null +++ b/tests/numerics/subtract_mean_value_01.output @@ -0,0 +1,7 @@ + +DEAL::OK +DEAL::OK +DEAL::OK +DEAL::OK +DEAL::OK +DEAL::OK diff --git a/tests/numerics/subtract_mean_value_02.cc b/tests/numerics/subtract_mean_value_02.cc new file mode 100644 index 0000000000..828181bd17 --- /dev/null +++ b/tests/numerics/subtract_mean_value_02.cc @@ -0,0 +1,110 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// check VectorTools::subtract_mean_value() for deal.II parallel vector + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + +template +void test (VectorType &v) +{ + // set some elements of the vector + unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + for (unsigned int i=5*my_id; i<5*(my_id+1); ++i) + { + v(i) = i; + } + v.compress (VectorOperation::insert); + + // then check the norm + VectorTools::subtract_mean_value(v); + AssertThrow (std::fabs(v.mean_value()) < 1e-10*v.l2_norm(), + ExcInternalError()); + + deallog << "OK" << std::endl; +} + +int main(int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1); + mpi_initlog(); + + unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + try + { + IndexSet local_range(10); + local_range.add_range(5*my_id,5*(my_id+1)); + IndexSet ghost_indices(10); + ghost_indices.add_range(3, 8); + { + LinearAlgebra::distributed::Vector v + (local_range, ghost_indices, MPI_COMM_WORLD); + test (v); + } + + { + LinearAlgebra::distributed::Vector v + (local_range, ghost_indices, MPI_COMM_WORLD); + test (v); + } + + { + LinearAlgebra::distributed::BlockVector v + (std::vector(1, local_range), + std::vector(1, ghost_indices), MPI_COMM_WORLD); + test (v); + } + + { + LinearAlgebra::distributed::BlockVector v + (std::vector(1, local_range), + std::vector(1, ghost_indices), MPI_COMM_WORLD); + test (v); + } + } + catch (std::exception &exc) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Exception on processing: " << std::endl + << exc.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + + return 1; + } + catch (...) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + }; +} diff --git a/tests/numerics/subtract_mean_value_02.with_mpi=true.mpirun=2.output b/tests/numerics/subtract_mean_value_02.with_mpi=true.mpirun=2.output new file mode 100644 index 0000000000..5f42bb230f --- /dev/null +++ b/tests/numerics/subtract_mean_value_02.with_mpi=true.mpirun=2.output @@ -0,0 +1,5 @@ + +DEAL::OK +DEAL::OK +DEAL::OK +DEAL::OK diff --git a/tests/petsc/subtract_mean_value_03.cc b/tests/petsc/subtract_mean_value_03.cc new file mode 100644 index 0000000000..5c76824e13 --- /dev/null +++ b/tests/petsc/subtract_mean_value_03.cc @@ -0,0 +1,94 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// check VectorTools::subtract_mean_value() for PETSc vectors + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + +template +void test (VectorType &v) +{ + // set some elements of the vector + unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + for (unsigned int i=5*my_id; i<5*(my_id+1); ++i) + { + v(i) = i; + } + v.compress (VectorOperation::insert); + + // then check the norm + VectorTools::subtract_mean_value(v); + AssertThrow (std::fabs(v.mean_value()) < 1e-10*v.l2_norm(), + ExcInternalError()); + + deallog << "OK" << std::endl; +} + +int main (int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1); + mpi_initlog(); + + unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + try + { + IndexSet local_range(10); + local_range.add_range(5*my_id,5*(my_id+1)); + { + PETScWrappers::MPI::Vector v + (local_range, MPI_COMM_WORLD); + test (v); + } + + { + PETScWrappers::MPI::BlockVector v + (std::vector(1, local_range), MPI_COMM_WORLD); + test (v); + } + } + catch (std::exception &exc) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Exception on processing: " << std::endl + << exc.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + + return 1; + } + catch (...) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + }; +} diff --git a/tests/petsc/subtract_mean_value_03.with_petsc=true.with_mpi=true.mpirun=2.output b/tests/petsc/subtract_mean_value_03.with_petsc=true.with_mpi=true.mpirun=2.output new file mode 100644 index 0000000000..8b3b075900 --- /dev/null +++ b/tests/petsc/subtract_mean_value_03.with_petsc=true.with_mpi=true.mpirun=2.output @@ -0,0 +1,3 @@ + +DEAL::OK +DEAL::OK diff --git a/tests/trilinos/subtract_mean_value_04.cc b/tests/trilinos/subtract_mean_value_04.cc new file mode 100644 index 0000000000..f52250b44c --- /dev/null +++ b/tests/trilinos/subtract_mean_value_04.cc @@ -0,0 +1,101 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + + +// check VectorTools::subtract_mean_value() for Trilinos vectors + +#include "../tests.h" +#include +#include +#include +#include +#include +#include + +template +void test (VectorType &v) +{ + // set some elements of the vector + unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + for (unsigned int i=5*my_id; i<5*(my_id+1); ++i) + { + v(i) = i; + } + v.compress (VectorOperation::insert); + + // then check the norm + VectorTools::subtract_mean_value(v); + AssertThrow (std::fabs(v.mean_value()) < 1e-10*v.l2_norm(), + ExcInternalError()); + + deallog << "OK" << std::endl; +} + +int main (int argc, char *argv[]) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization (argc, argv, 1); + mpi_initlog(); + + unsigned int my_id = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + + try + { + IndexSet local_range(10); + local_range.add_range(5*my_id,5*(my_id+1)); + { + TrilinosWrappers::MPI::Vector v1 (local_range, MPI_COMM_WORLD); + test(v1); + + LinearAlgebra::ReadWriteVector v_tmp (local_range); + LinearAlgebra::EpetraWrappers::Vector v2 (local_range, MPI_COMM_WORLD); + v_tmp.import(v1, VectorOperation::insert); + v2.import(v_tmp, VectorOperation::insert); + VectorTools::subtract_mean_value(v2); + AssertThrow (std::fabs(v2.mean_value()) < 1e-10*v2.l2_norm(), + ExcInternalError()); + deallog << "OK" << std::endl; + } + { + TrilinosWrappers::MPI::BlockVector v (std::vector(1, local_range), + MPI_COMM_WORLD); + test(v); + } + } + catch (std::exception &exc) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Exception on processing: " << std::endl + << exc.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + + return 1; + } + catch (...) + { + deallog << std::endl << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + }; +} diff --git a/tests/trilinos/subtract_mean_value_04.with_trilinos=true.with_mpi=true.mpirun=2.output b/tests/trilinos/subtract_mean_value_04.with_trilinos=true.with_mpi=true.mpirun=2.output new file mode 100644 index 0000000000..fb71de2867 --- /dev/null +++ b/tests/trilinos/subtract_mean_value_04.with_trilinos=true.with_mpi=true.mpirun=2.output @@ -0,0 +1,4 @@ + +DEAL::OK +DEAL::OK +DEAL::OK