]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove deprecated member functions in vector classes
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Wed, 7 Jun 2017 15:53:11 +0000 (17:53 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Mon, 12 Jun 2017 20:08:33 +0000 (22:08 +0200)
14 files changed:
contrib/utilities/print_deprecated.py [changed mode: 0644->0755]
doc/news/changes/incompatibilities/20170612DanielArndt [new file with mode: 0644]
include/deal.II/lac/block_vector_base.h
include/deal.II/lac/matrix_lib.templates.h
include/deal.II/lac/petsc_vector_base.h
include/deal.II/lac/solver_qmrs.h
include/deal.II/lac/vector.h
include/deal.II/lac/vector.templates.h
source/lac/block_matrix_array.cc
source/lac/petsc_vector_base.cc
tests/vector/complex_vector_38.cc [deleted file]
tests/vector/complex_vector_38.output [deleted file]
tests/vector/vector_38.cc [deleted file]
tests/vector/vector_38.output [deleted file]

old mode 100644 (file)
new mode 100755 (executable)
diff --git a/doc/news/changes/incompatibilities/20170612DanielArndt b/doc/news/changes/incompatibilities/20170612DanielArndt
new file mode 100644 (file)
index 0000000..de9dea0
--- /dev/null
@@ -0,0 +1,4 @@
+Changed: The deprecated member functions add(), normalize(), conjugate(), abs()
+and mult() in the vector classes have been removed.
+<br>
+(Daniel Arndt, 2017/06/12)
index 40f50727972a5e63d53c05913445a1b00dab50bf..2e61d86c032fa2a57e7ab79ea76c9e70edbaf925 100644 (file)
@@ -858,13 +858,6 @@ public:
    */
   void add (const value_type s);
 
-  /**
-   * U+=V. Simple vector addition, equal to the <tt>operator +=</tt>.
-   *
-   * This function is deprecated use the <tt>operator +=</tt> instead.
-   */
-  void add (const BlockVectorBase &V) DEAL_II_DEPRECATED;
-
   /**
    * U+=a*V. Simple addition of a scaled vector.
    */
@@ -1820,14 +1813,6 @@ void BlockVectorBase<VectorType>::add (const value_type a)
 
 
 
-template <class VectorType>
-void BlockVectorBase<VectorType>::add (const BlockVectorBase<VectorType> &v)
-{
-  *this += v;
-}
-
-
-
 template <class VectorType>
 void BlockVectorBase<VectorType>::add (const value_type a,
                                        const BlockVectorBase<VectorType> &v)
index 247d23234a6230edc7a5e0cb1cd2f9cb83a35787..af9abfe1654fe9dc6ea9d5841c52f4aaad6eafa2 100644 (file)
@@ -115,7 +115,7 @@ MeanValueFilter::vmult_add(BlockVector<number> &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);
 }
 
 
index 0e3d4e624bbe6220c73c402b967f2202ef23df55..0306d34ef7bd3688e5fc0733919fb0b942c93af8 100644 (file)
@@ -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 <code>this</code> 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
-     * <code>this</code> vector with <b>v</b>.
-     *
-     * This function is deprecated.
-     */
-    VectorBase &mult (const VectorBase &v) DEAL_II_DEPRECATED;
-
-    /**
-     * Same as above, but a collective piecewise multiply operation of
-     * <b>u</b> with <b>v</b>.
-     *
-     * 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 <tt>operator +=</tt>.
-     *
-     * @deprecated Use the <tt>operator +=</tt> instead.
-     */
-    void add (const VectorBase &V) DEAL_II_DEPRECATED;
-
     /**
      * Simple addition of a multiple of a vector, i.e. <tt>*this += a*V</tt>.
      */
index bc53ce6eaceebde9c21f9f00eff32b4cf858629d..60e0073b4c8eedc4cfe3dfaaa0d61a839ef8f1cb 100644 (file)
@@ -391,7 +391,7 @@ SolverQMRS<VectorType>::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
index 49181d00bbccaf40563b76691f86b54168d7876e..527b1999c8ccc338c02686da0bd791ccb8fe0943 100644 (file)
@@ -640,16 +640,6 @@ public:
    */
   void add (const Number s);
 
-  /**
-   * Simple vector addition, equal to the <tt>operator +=</tt>.
-   *
-   * @deprecated Use the <tt>operator +=</tt> instead.
-   *
-   * @dealiiOperationIsMultithreaded
-   */
-  void add (const Vector<Number> &V) DEAL_II_DEPRECATED;
-
-
   /**
    * Multiple addition of scaled vectors, i.e. <tt>*this += a*V+b*W</tt>.
    *
index cf4f399007e4072ca691478b1d38718b93319af0..d4ab64aa724eae2613f9d42af5aa4cf5a83f99ac 100644 (file)
@@ -612,8 +612,10 @@ template <typename Number>
 Vector<Number> &Vector<Number>::operator += (const Vector<Number> &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<Number> 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<Number>::add (const Number v)
 
 
 
-template <typename Number>
-void Vector<Number>::add (const Vector<Number> &v)
-{
-  Assert (vec_size!=0, ExcEmptyObject());
-  Assert (vec_size == v.vec_size, ExcDimensionMismatch(vec_size, v.vec_size));
-
-  internal::VectorOperations::Vectorization_add_v<Number> vector_add(val, v.val);
-  internal::VectorOperations::parallel_for(vector_add,0,vec_size,thread_loop_partitioner);
-}
-
-
-
 template <typename Number>
 void Vector<Number>::add (const Number a, const Vector<Number> &v,
                           const Number b, const Vector<Number> &w)
index c389f9d37fc9975acbd4b0b1ab094219ad4a8897..acaa3eff7976de26c1db2970b90b690feafc170f 100644 (file)
@@ -374,7 +374,7 @@ BlockTrianglePrecondition<number,BlockVectorType>::vmult_add
   BlockVectorType aux;
   aux.reinit(dst);
   vmult(aux, src);
-  dst.add(aux);
+  dst += aux;
 }
 
 
index 46cf913bef8067ef092da37b10bb7119f169563b..4506e2edf9a292630e8e8f598938ed5ee54d96fc 100644 (file)
@@ -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 (file)
index acf68d2..0000000
+++ /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<std::complex<double> >::add(Vector)
-
-#include "../tests.h"
-#include <deal.II/lac/vector.h>
-#include <fstream>
-#include <iomanip>
-#include <vector>
-
-
-void test (Vector<std::complex<double> > &v,
-           Vector<std::complex<double> > &w)
-{
-  for (unsigned int i=0; i<v.size(); ++i)
-    {
-      v(i) = i;
-      w(i) = std::complex<double> (i+1., i+2.);
-    }
-
-  v.compress ();
-  w.compress ();
-
-  v.add (w);
-
-  // make sure we get the expected result
-  for (unsigned int i=0; i<v.size(); ++i)
-    {
-      AssertThrow (w(i) == std::complex<double> (i+1., i+2.),
-                   ExcInternalError());
-      AssertThrow (v(i) == 1.*i+std::complex<double> (i+1., i+2.),
-                   ExcInternalError());
-    }
-
-  deallog << "OK" << std::endl;
-}
-
-
-
-int main ()
-{
-  initlog();
-  deallog.threshold_double(1.e-10);
-
-  try
-    {
-      Vector<std::complex<double> > v (100);
-      Vector<std::complex<double> > 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 (file)
index 0fd8fc1..0000000
+++ /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 (file)
index 28ad1c2..0000000
+++ /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<double>::add(Vector)
-
-#include "../tests.h"
-#include <deal.II/lac/vector.h>
-#include <fstream>
-#include <iomanip>
-#include <vector>
-
-
-void test (Vector<double> &v,
-           Vector<double> &w)
-{
-  for (unsigned int i=0; i<v.size(); ++i)
-    {
-      v(i) = i;
-      w(i) = i+1.;
-    }
-
-  v.compress ();
-  w.compress ();
-
-  v.add (w);
-
-  // make sure we get the expected result
-  for (unsigned int i=0; i<v.size(); ++i)
-    {
-      Assert (w(i) == i+1., ExcInternalError());
-      Assert (v(i) == i+(i+1.), ExcInternalError());
-    }
-
-  deallog << "OK" << std::endl;
-}
-
-
-
-int main ()
-{
-  initlog();
-  deallog.threshold_double(1.e-10);
-
-  try
-    {
-      Vector<double> v (100);
-      Vector<double> 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 (file)
index 0fd8fc1..0000000
+++ /dev/null
@@ -1,2 +0,0 @@
-
-DEAL::OK

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.