From e79ee10b45c28613e955f76f1213be3ff5392099 Mon Sep 17 00:00:00 2001 From: Daniel Arndt Date: Wed, 7 Jun 2017 17:53:11 +0200 Subject: [PATCH] Remove deprecated member functions in vector classes --- contrib/utilities/print_deprecated.py | 0 .../incompatibilities/20170612DanielArndt | 4 + include/deal.II/lac/block_vector_base.h | 15 ---- include/deal.II/lac/matrix_lib.templates.h | 2 +- include/deal.II/lac/petsc_vector_base.h | 55 ------------ include/deal.II/lac/solver_qmrs.h | 2 +- include/deal.II/lac/vector.h | 10 --- include/deal.II/lac/vector.templates.h | 16 +--- source/lac/block_matrix_array.cc | 2 +- source/lac/petsc_vector_base.cc | 79 ---------------- tests/vector/complex_vector_38.cc | 90 ------------------- tests/vector/complex_vector_38.output | 2 - tests/vector/vector_38.cc | 88 ------------------ tests/vector/vector_38.output | 2 - 14 files changed, 10 insertions(+), 357 deletions(-) mode change 100644 => 100755 contrib/utilities/print_deprecated.py create mode 100644 doc/news/changes/incompatibilities/20170612DanielArndt delete mode 100644 tests/vector/complex_vector_38.cc delete mode 100644 tests/vector/complex_vector_38.output delete mode 100644 tests/vector/vector_38.cc delete mode 100644 tests/vector/vector_38.output diff --git a/contrib/utilities/print_deprecated.py b/contrib/utilities/print_deprecated.py old mode 100644 new mode 100755 diff --git a/doc/news/changes/incompatibilities/20170612DanielArndt b/doc/news/changes/incompatibilities/20170612DanielArndt new file mode 100644 index 0000000000..de9dea0042 --- /dev/null +++ b/doc/news/changes/incompatibilities/20170612DanielArndt @@ -0,0 +1,4 @@ +Changed: The deprecated member functions add(), normalize(), conjugate(), abs() +and mult() in the vector classes have been removed. +
+(Daniel Arndt, 2017/06/12) diff --git a/include/deal.II/lac/block_vector_base.h b/include/deal.II/lac/block_vector_base.h index 40f5072797..2e61d86c03 100644 --- a/include/deal.II/lac/block_vector_base.h +++ b/include/deal.II/lac/block_vector_base.h @@ -858,13 +858,6 @@ public: */ void add (const value_type s); - /** - * U+=V. Simple vector addition, equal to the operator +=. - * - * This function is deprecated use the operator += instead. - */ - void add (const BlockVectorBase &V) DEAL_II_DEPRECATED; - /** * U+=a*V. Simple addition of a scaled vector. */ @@ -1820,14 +1813,6 @@ void BlockVectorBase::add (const value_type a) -template -void BlockVectorBase::add (const BlockVectorBase &v) -{ - *this += v; -} - - - template void BlockVectorBase::add (const value_type a, const BlockVectorBase &v) diff --git a/include/deal.II/lac/matrix_lib.templates.h b/include/deal.II/lac/matrix_lib.templates.h index 247d23234a..af9abfe165 100644 --- a/include/deal.II/lac/matrix_lib.templates.h +++ b/include/deal.II/lac/matrix_lib.templates.h @@ -115,7 +115,7 @@ MeanValueFilter::vmult_add(BlockVector &dst, if (i == component) vmult_add(dst.block(i), src.block(i)); else - dst.block(i).add(src.block(i)); + dst.block(i)+=src.block(i); } diff --git a/include/deal.II/lac/petsc_vector_base.h b/include/deal.II/lac/petsc_vector_base.h index 0e3d4e624b..0306d34ef7 100644 --- a/include/deal.II/lac/petsc_vector_base.h +++ b/include/deal.II/lac/petsc_vector_base.h @@ -519,14 +519,6 @@ namespace PETScWrappers const VectorBase &V, const VectorBase &W); - /** - * Normalize vector by dividing by the $l_2$-norm of the vector. Return - * the vector norm before normalization. - * - * This function is deprecated. - */ - real_type normalize () const DEAL_II_DEPRECATED; - /** * Return the value of the vector element with the largest negative value. */ @@ -537,46 +529,6 @@ namespace PETScWrappers */ real_type max () const; - /** - * Replace every element in a vector with its absolute value. - * - * This function is deprecated. - */ - VectorBase &abs () DEAL_II_DEPRECATED; - - /** - * Conjugate a vector. - * - * This function is deprecated. - */ - VectorBase &conjugate () DEAL_II_DEPRECATED; - - /** - * A collective piecewise multiply operation on this vector - * with itself. TODO: The model for this function should be similar to add - * (). - * - * This function is deprecated. - */ - VectorBase &mult () DEAL_II_DEPRECATED; - - /** - * Same as above, but a collective piecewise multiply operation of - * this vector with v. - * - * This function is deprecated. - */ - VectorBase &mult (const VectorBase &v) DEAL_II_DEPRECATED; - - /** - * Same as above, but a collective piecewise multiply operation of - * u with v. - * - * This function is deprecated. - */ - VectorBase &mult (const VectorBase &u, - const VectorBase &v) DEAL_II_DEPRECATED; - /** * Return whether the vector contains only elements with value zero. This * is a collective operation. This function is expensive, because @@ -617,13 +569,6 @@ namespace PETScWrappers */ void add (const PetscScalar s); - /** - * Simple vector addition, equal to the operator +=. - * - * @deprecated Use the operator += instead. - */ - void add (const VectorBase &V) DEAL_II_DEPRECATED; - /** * Simple addition of a multiple of a vector, i.e. *this += a*V. */ diff --git a/include/deal.II/lac/solver_qmrs.h b/include/deal.II/lac/solver_qmrs.h index bc53ce6eac..60e0073b4c 100644 --- a/include/deal.II/lac/solver_qmrs.h +++ b/include/deal.II/lac/solver_qmrs.h @@ -391,7 +391,7 @@ SolverQMRS::iterate(const MatrixType &A, tau *= theta*psi; d.sadd(psi*theta_old, psi*alpha, p); - x.add(d); + x += d; print_vectors(step,x,v,d); // Step 5 diff --git a/include/deal.II/lac/vector.h b/include/deal.II/lac/vector.h index 49181d00bb..527b1999c8 100644 --- a/include/deal.II/lac/vector.h +++ b/include/deal.II/lac/vector.h @@ -640,16 +640,6 @@ public: */ void add (const Number s); - /** - * Simple vector addition, equal to the operator +=. - * - * @deprecated Use the operator += instead. - * - * @dealiiOperationIsMultithreaded - */ - void add (const Vector &V) DEAL_II_DEPRECATED; - - /** * Multiple addition of scaled vectors, i.e. *this += a*V+b*W. * diff --git a/include/deal.II/lac/vector.templates.h b/include/deal.II/lac/vector.templates.h index cf4f399007..d4ab64aa72 100644 --- a/include/deal.II/lac/vector.templates.h +++ b/include/deal.II/lac/vector.templates.h @@ -612,8 +612,10 @@ template Vector &Vector::operator += (const Vector &v) { Assert (vec_size!=0, ExcEmptyObject()); + Assert (vec_size == v.vec_size, ExcDimensionMismatch(vec_size, v.vec_size)); - add (v); + internal::VectorOperations::Vectorization_add_v vector_add(val, v.val); + internal::VectorOperations::parallel_for(vector_add,0,vec_size,thread_loop_partitioner); return *this; } @@ -644,18 +646,6 @@ void Vector::add (const Number v) -template -void Vector::add (const Vector &v) -{ - Assert (vec_size!=0, ExcEmptyObject()); - Assert (vec_size == v.vec_size, ExcDimensionMismatch(vec_size, v.vec_size)); - - internal::VectorOperations::Vectorization_add_v vector_add(val, v.val); - internal::VectorOperations::parallel_for(vector_add,0,vec_size,thread_loop_partitioner); -} - - - template void Vector::add (const Number a, const Vector &v, const Number b, const Vector &w) diff --git a/source/lac/block_matrix_array.cc b/source/lac/block_matrix_array.cc index c389f9d37f..acaa3eff79 100644 --- a/source/lac/block_matrix_array.cc +++ b/source/lac/block_matrix_array.cc @@ -374,7 +374,7 @@ BlockTrianglePrecondition::vmult_add BlockVectorType aux; aux.reinit(dst); vmult(aux, src); - dst.add(aux); + dst += aux; } diff --git a/source/lac/petsc_vector_base.cc b/source/lac/petsc_vector_base.cc index 46cf913bef..4506e2edf9 100644 --- a/source/lac/petsc_vector_base.cc +++ b/source/lac/petsc_vector_base.cc @@ -560,18 +560,6 @@ namespace PETScWrappers - VectorBase::real_type - VectorBase::normalize () const - { - real_type d; - Assert (!has_ghost_elements(), ExcGhostsPresent()); - const PetscErrorCode ierr = VecNormalize (vector, &d); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - return d; - } - - VectorBase::real_type VectorBase::min () const { @@ -598,65 +586,6 @@ namespace PETScWrappers } - VectorBase & - VectorBase::abs () - { - Assert (!has_ghost_elements(), ExcGhostsPresent()); - - const PetscErrorCode ierr = VecAbs (vector); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - return *this; - } - - - - VectorBase & - VectorBase::conjugate () - { - Assert (!has_ghost_elements(), ExcGhostsPresent()); - - const PetscErrorCode ierr = VecConjugate (vector); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - return *this; - } - - - - VectorBase & - VectorBase::mult () - { - Assert (!has_ghost_elements(), ExcGhostsPresent()); - - const PetscErrorCode ierr = VecPointwiseMult (vector,vector,vector); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - return *this; - } - - - VectorBase & - VectorBase::mult (const VectorBase &v) - { - Assert (!has_ghost_elements(), ExcGhostsPresent()); - const PetscErrorCode ierr = VecPointwiseMult (vector,vector,v); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - return *this; - } - - - VectorBase & - VectorBase::mult (const VectorBase &u, - const VectorBase &v) - { - Assert (!has_ghost_elements(), ExcGhostsPresent()); - const PetscErrorCode ierr = VecPointwiseMult (vector,u,v); - AssertThrow (ierr == 0, ExcPETScError(ierr)); - - return *this; - } bool VectorBase::all_zero () const @@ -810,14 +739,6 @@ namespace PETScWrappers - void - VectorBase::add (const VectorBase &v) - { - *this += v; - } - - - void VectorBase::add (const PetscScalar a, const VectorBase &v) diff --git a/tests/vector/complex_vector_38.cc b/tests/vector/complex_vector_38.cc deleted file mode 100644 index acf68d2562..0000000000 --- a/tests/vector/complex_vector_38.cc +++ /dev/null @@ -1,90 +0,0 @@ -// --------------------------------------------------------------------- -// -// Copyright (C) 2004 - 2016 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 Vector >::add(Vector) - -#include "../tests.h" -#include -#include -#include -#include - - -void test (Vector > &v, - Vector > &w) -{ - for (unsigned int i=0; i (i+1., i+2.); - } - - v.compress (); - w.compress (); - - v.add (w); - - // make sure we get the expected result - for (unsigned int i=0; i (i+1., i+2.), - ExcInternalError()); - AssertThrow (v(i) == 1.*i+std::complex (i+1., i+2.), - ExcInternalError()); - } - - deallog << "OK" << std::endl; -} - - - -int main () -{ - initlog(); - deallog.threshold_double(1.e-10); - - try - { - Vector > v (100); - Vector > w (100); - test (v,w); - } - 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/vector/complex_vector_38.output b/tests/vector/complex_vector_38.output deleted file mode 100644 index 0fd8fc12f0..0000000000 --- a/tests/vector/complex_vector_38.output +++ /dev/null @@ -1,2 +0,0 @@ - -DEAL::OK diff --git a/tests/vector/vector_38.cc b/tests/vector/vector_38.cc deleted file mode 100644 index 28ad1c24e4..0000000000 --- a/tests/vector/vector_38.cc +++ /dev/null @@ -1,88 +0,0 @@ -// --------------------------------------------------------------------- -// -// Copyright (C) 2004 - 2016 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 Vector::add(Vector) - -#include "../tests.h" -#include -#include -#include -#include - - -void test (Vector &v, - Vector &w) -{ - for (unsigned int i=0; i v (100); - Vector w (100); - test (v,w); - } - 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/vector/vector_38.output b/tests/vector/vector_38.output deleted file mode 100644 index 0fd8fc12f0..0000000000 --- a/tests/vector/vector_38.output +++ /dev/null @@ -1,2 +0,0 @@ - -DEAL::OK -- 2.39.5