]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove derivation of LinearAlgebra::distributed::Vector from VectorSpaceVector.
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 6 Jul 2023 13:02:54 +0000 (07:02 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Fri, 7 Jul 2023 05:16:01 +0000 (23:16 -0600)
include/deal.II/lac/la_parallel_vector.h
include/deal.II/lac/la_parallel_vector.templates.h

index 0d6dc96e31704579dc614813fd048b01ed30127b..6b117fd4a74e7691a9fdb1724878ef51c1f141ea 100644 (file)
@@ -248,9 +248,7 @@ namespace LinearAlgebra
      * @endcode
      */
     template <typename Number, typename MemorySpace = MemorySpace::Host>
-    class Vector : public ::dealii::LinearAlgebra::VectorSpaceVector<Number>,
-                   public ::dealii::ReadVector<Number>,
-                   public Subscriptor
+    class Vector : public ::dealii::ReadVector<Number>, public Subscriptor
     {
     public:
       using memory_space    = MemorySpace;
@@ -338,7 +336,7 @@ namespace LinearAlgebra
       /**
        * Destructor.
        */
-      virtual ~Vector() override;
+      ~Vector();
 
       /**
        * Set the global size of the vector to @p size without any actual
@@ -537,8 +535,8 @@ namespace LinearAlgebra
        * with VectorOperation::max two times subsequently, the maximal value
        * after the second calculation will be zero.
        */
-      virtual void
-      compress(VectorOperation::values operation) override;
+      void
+      compress(VectorOperation::values operation);
 
       /**
        * Fills the data field for ghost indices with the values stored in the
@@ -710,41 +708,33 @@ namespace LinearAlgebra
       /** @} */
 
       /**
-       * @name 3: Implementation of VectorSpaceVector
+       * @name 3: Implementation of vector space operations
        */
       /** @{ */
 
-      /**
-       * Change the dimension to that of the vector V. The elements of V are not
-       * copied.
-       */
-      virtual void
-      reinit(const VectorSpaceVector<Number> &V,
-             const bool omit_zeroing_entries = false) override;
-
       /**
        * Multiply the entire vector by a fixed factor.
        */
-      virtual Vector<Number, MemorySpace> &
-      operator*=(const Number factor) override;
+      Vector<Number, MemorySpace> &
+      operator*=(const Number factor);
 
       /**
        * Divide the entire vector by a fixed factor.
        */
-      virtual Vector<Number, MemorySpace> &
-      operator/=(const Number factor) override;
+      Vector<Number, MemorySpace> &
+      operator/=(const Number factor);
 
       /**
        * Add the vector @p V to the present one.
        */
-      virtual Vector<Number, MemorySpace> &
-      operator+=(const VectorSpaceVector<Number> &V) override;
+      Vector<Number, MemorySpace> &
+      operator+=(const Vector<Number, MemorySpace> &V);
 
       /**
        * Subtract the vector @p V from the present one.
        */
-      virtual Vector<Number, MemorySpace> &
-      operator-=(const VectorSpaceVector<Number> &V) override;
+      Vector<Number, MemorySpace> &
+      operator-=(const Vector<Number, MemorySpace> &V);
 
       /**
        * Import all the elements present in the vector's IndexSet from the input
@@ -757,21 +747,21 @@ namespace LinearAlgebra
        * @note If the MemorySpace is Default, the data in the ReadWriteVector will
        * be moved to the @ref GlossDevice "device".
        */
-      virtual void
+      void
       import_elements(
         const LinearAlgebra::ReadWriteVector<Number> &V,
         VectorOperation::values                       operation,
         std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
-          communication_pattern = {}) override;
+          communication_pattern = {});
 
       /**
        * @deprecated Use import_elements() instead.
        */
-      DEAL_II_DEPRECATED virtual void
+      DEAL_II_DEPRECATED void
       import(const LinearAlgebra::ReadWriteVector<Number> &V,
              VectorOperation::values                       operation,
              std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
-               communication_pattern = {}) override
+               communication_pattern = {})
       {
         import_elements(V, operation, communication_pattern);
       }
@@ -779,35 +769,35 @@ namespace LinearAlgebra
       /**
        * Return the scalar product of two vectors.
        */
-      virtual Number
-      operator*(const VectorSpaceVector<Number> &V) const override;
+      Number
+      operator*(const Vector<Number, MemorySpace> &V) const;
 
       /**
        * Add @p a to all components. Note that @p a is a scalar not a vector.
        */
-      virtual void
-      add(const Number a) override;
+      void
+      add(const Number a);
 
       /**
        * Simple addition of a multiple of a vector, i.e. <tt>*this += a*V</tt>.
        */
-      virtual void
-      add(const Number a, const VectorSpaceVector<Number> &V) override;
+      void
+      add(const Number a, const Vector<Number, MemorySpace> &V);
 
       /**
        * Multiple addition of scaled vectors, i.e. <tt>*this += a*V+b*W</tt>.
        */
-      virtual void
-      add(const Number                     a,
-          const VectorSpaceVector<Number> &V,
-          const Number                     b,
-          const VectorSpaceVector<Number> &W) override;
+      void
+      add(const Number                       a,
+          const Vector<Number, MemorySpace> &V,
+          const Number                       b,
+          const Vector<Number, MemorySpace> &W);
 
       /**
        * A collective add operation: This function adds a whole set of values
        * stored in @p values to the vector components specified by @p indices.
        */
-      virtual void
+      void
       add(const std::vector<size_type> &indices,
           const std::vector<Number> &   values);
 
@@ -815,38 +805,38 @@ namespace LinearAlgebra
        * Scaling and simple addition of a multiple of a vector, i.e. <tt>*this =
        * s*(*this)+a*V</tt>.
        */
-      virtual void
-      sadd(const Number                     s,
-           const Number                     a,
-           const VectorSpaceVector<Number> &V) override;
+      void
+      sadd(const Number                       s,
+           const Number                       a,
+           const Vector<Number, MemorySpace> &V);
 
       /**
        * Scale each element of this vector by the corresponding element in the
        * argument. This function is mostly meant to simulate multiplication (and
        * immediate re-assignment) by a diagonal scaling matrix.
        */
-      virtual void
-      scale(const VectorSpaceVector<Number> &scaling_factors) override;
+      void
+      scale(const Vector<Number, MemorySpace> &scaling_factors);
 
       /**
        * Assignment <tt>*this = a*V</tt>.
        */
-      virtual void
-      equ(const Number a, const VectorSpaceVector<Number> &V) override;
+      void
+      equ(const Number a, const Vector<Number, MemorySpace> &V);
 
       /**
        * Return the l<sub>1</sub> norm of the vector (i.e., the sum of the
        * absolute values of all entries among all processors).
        */
-      virtual real_type
-      l1_norm() const override;
+      real_type
+      l1_norm() const;
 
       /**
        * Return the $l_2$ norm of the vector (i.e., the square root of
        * the sum of the square of all entries among all processors).
        */
-      virtual real_type
-      l2_norm() const override;
+      real_type
+      l2_norm() const;
 
       /**
        * Return the square of the $l_2$ norm of the vector.
@@ -858,8 +848,8 @@ namespace LinearAlgebra
        * Return the maximum norm of the vector (i.e., the maximum absolute value
        * among all entries and among all processors).
        */
-      virtual real_type
-      linfty_norm() const override;
+      real_type
+      linfty_norm() const;
 
       /**
        * Perform a combined operation of a vector addition and a subsequent
@@ -881,10 +871,10 @@ namespace LinearAlgebra
        * implemented as
        * $\left<v,w\right>=\sum_i v_i \bar{w_i}$.
        */
-      virtual Number
-      add_and_dot(const Number                     a,
-                  const VectorSpaceVector<Number> &V,
-                  const VectorSpaceVector<Number> &W) override;
+      Number
+      add_and_dot(const Number                       a,
+                  const Vector<Number, MemorySpace> &V,
+                  const Vector<Number, MemorySpace> &W);
 
       /**
        * Return the global size of the vector, equal to the sum of the number of
@@ -904,27 +894,27 @@ namespace LinearAlgebra
        *  vec.locally_owned_elements() == complete_index_set(vec.size())
        * @endcode
        */
-      virtual dealii::IndexSet
-      locally_owned_elements() const override;
+      dealii::IndexSet
+      locally_owned_elements() const;
 
       /**
        * Print the vector to the output stream @p out.
        */
-      virtual void
+      void
       print(std::ostream &     out,
             const unsigned int precision  = 3,
             const bool         scientific = true,
-            const bool         across     = true) const override;
+            const bool         across     = true) const;
 
       /**
        * Return the memory consumption of this class in bytes.
        */
-      virtual std::size_t
-      memory_consumption() const override;
+      std::size_t
+      memory_consumption() const;
       /** @} */
 
       /**
-       * @name 4: Other vector operations not included in VectorSpaceVector
+       * @name 4: Other vector operations
        */
       /** @{ */
 
@@ -933,8 +923,8 @@ namespace LinearAlgebra
        * zero, also ghost elements are set to zero, otherwise they remain
        * unchanged.
        */
-      virtual Vector<Number, MemorySpace> &
-      operator=(const Number s) override;
+      Vector<Number, MemorySpace> &
+      operator=(const Number s);
 
       /**
        * This is a collective add operation that adds a whole set of values
@@ -1168,14 +1158,14 @@ namespace LinearAlgebra
        * This is a collective operation. This function is expensive, because
        * potentially all elements have to be checked.
        */
-      virtual bool
-      all_zero() const override;
+      bool
+      all_zero() const;
 
       /**
        * Compute the mean value of all the entries in the vector.
        */
-      virtual Number
-      mean_value() const override;
+      Number
+      mean_value() const;
 
       /**
        * $l_p$-norm of the vector. The pth root of the sum of the pth powers
@@ -1301,16 +1291,16 @@ namespace LinearAlgebra
        * without MPI communication.
        */
       void
-      add_local(const Number a, const VectorSpaceVector<Number> &V);
+      add_local(const Number a, const Vector<Number, MemorySpace> &V);
 
       /**
        * Scaling and simple addition of a multiple of a vector, i.e. <tt>*this =
        * s*(*this)+a*V</tt> without MPI communication.
        */
       void
-      sadd_local(const Number                     s,
-                 const Number                     a,
-                 const VectorSpaceVector<Number> &V);
+      sadd_local(const Number                       s,
+                 const Number                       a,
+                 const Vector<Number, MemorySpace> &V);
 
       /**
        * Local part of the inner product of two vectors.
index 967a6bbb25bcdf3d4e8e85addccb3edac46fca32..84bbd3b1f3f49d634f8bd3194657966374722ec9 100644 (file)
@@ -1423,33 +1423,11 @@ namespace LinearAlgebra
 
 
 
-    template <typename Number, typename MemorySpaceType>
-    void
-    Vector<Number, MemorySpaceType>::reinit(const VectorSpaceVector<Number> &V,
-                                            const bool omit_zeroing_entries)
-    {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &down_V = dynamic_cast<const VectorType &>(V);
-
-      reinit(down_V, omit_zeroing_entries);
-    }
-
-
-
     template <typename Number, typename MemorySpaceType>
     Vector<Number, MemorySpaceType> &
     Vector<Number, MemorySpaceType>::operator+=(
-      const VectorSpaceVector<Number> &vv)
+      const Vector<Number, MemorySpaceType> &v)
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&vv) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       AssertDimension(locally_owned_size(), v.locally_owned_size());
 
       dealii::internal::VectorOperations::
@@ -1470,14 +1448,8 @@ namespace LinearAlgebra
     template <typename Number, typename MemorySpaceType>
     Vector<Number, MemorySpaceType> &
     Vector<Number, MemorySpaceType>::operator-=(
-      const VectorSpaceVector<Number> &vv)
+      const Vector<Number, MemorySpaceType> &v)
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&vv) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       AssertDimension(locally_owned_size(), v.locally_owned_size());
 
       dealii::internal::VectorOperations::
@@ -1514,15 +1486,9 @@ namespace LinearAlgebra
     template <typename Number, typename MemorySpaceType>
     void
     Vector<Number, MemorySpaceType>::add_local(
-      const Number                     a,
-      const VectorSpaceVector<Number> &vv)
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &v)
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&vv) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       AssertIsFinite(a);
       AssertDimension(locally_owned_size(), v.locally_owned_size());
 
@@ -1543,8 +1509,9 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpaceType>
     void
-    Vector<Number, MemorySpaceType>::add(const Number                     a,
-                                         const VectorSpaceVector<Number> &vv)
+    Vector<Number, MemorySpaceType>::add(
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &vv)
     {
       add_local(a, vv);
 
@@ -1556,20 +1523,12 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpaceType>
     void
-    Vector<Number, MemorySpaceType>::add(const Number                     a,
-                                         const VectorSpaceVector<Number> &vv,
-                                         const Number                     b,
-                                         const VectorSpaceVector<Number> &ww)
-    {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&vv) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-      Assert(dynamic_cast<const VectorType *>(&ww) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &w = dynamic_cast<const VectorType &>(ww);
-
+    Vector<Number, MemorySpaceType>::add(
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &v,
+      const Number                           b,
+      const Vector<Number, MemorySpaceType> &w)
+    {
       AssertIsFinite(a);
       AssertIsFinite(b);
 
@@ -1631,16 +1590,10 @@ namespace LinearAlgebra
     template <typename Number, typename MemorySpaceType>
     void
     Vector<Number, MemorySpaceType>::sadd_local(
-      const Number                     x,
-      const Number                     a,
-      const VectorSpaceVector<Number> &vv)
+      const Number                           x,
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &v)
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert((dynamic_cast<const VectorType *>(&vv) != nullptr),
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       AssertIsFinite(x);
       AssertIsFinite(a);
       AssertDimension(locally_owned_size(), v.locally_owned_size());
@@ -1659,11 +1612,12 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpaceType>
     void
-    Vector<Number, MemorySpaceType>::sadd(const Number                     x,
-                                          const Number                     a,
-                                          const VectorSpaceVector<Number> &vv)
+    Vector<Number, MemorySpaceType>::sadd(
+      const Number                           x,
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &v)
     {
-      sadd_local(x, a, vv);
+      sadd_local(x, a, v);
 
       if (vector_is_ghosted)
         update_ghost_values();
@@ -1704,14 +1658,9 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpaceType>
     void
-    Vector<Number, MemorySpaceType>::scale(const VectorSpaceVector<Number> &vv)
+    Vector<Number, MemorySpaceType>::scale(
+      const Vector<Number, MemorySpaceType> &v)
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&vv) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       AssertDimension(locally_owned_size(), v.locally_owned_size());
 
       dealii::internal::VectorOperations::
@@ -1726,15 +1675,10 @@ namespace LinearAlgebra
 
     template <typename Number, typename MemorySpaceType>
     void
-    Vector<Number, MemorySpaceType>::equ(const Number                     a,
-                                         const VectorSpaceVector<Number> &vv)
+    Vector<Number, MemorySpaceType>::equ(
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &v)
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert(dynamic_cast<const VectorType *>(&vv) != nullptr,
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       AssertIsFinite(a);
       AssertDimension(locally_owned_size(), v.locally_owned_size());
 
@@ -1803,14 +1747,8 @@ namespace LinearAlgebra
     template <typename Number, typename MemorySpaceType>
     Number
     Vector<Number, MemorySpaceType>::operator*(
-      const VectorSpaceVector<Number> &vv) const
+      const Vector<Number, MemorySpaceType> &v) const
     {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert((dynamic_cast<const VectorType *>(&vv) != nullptr),
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-
       Number local_result = inner_product_local(v);
       if (partitioner->n_mpi_processes() > 1)
         return Utilities::MPI::sum(local_result,
@@ -2021,19 +1959,10 @@ namespace LinearAlgebra
     template <typename Number, typename MemorySpaceType>
     Number
     Vector<Number, MemorySpaceType>::add_and_dot(
-      const Number                     a,
-      const VectorSpaceVector<Number> &vv,
-      const VectorSpaceVector<Number> &ww)
-    {
-      // Downcast. Throws an exception if invalid.
-      using VectorType = Vector<Number, MemorySpaceType>;
-      Assert((dynamic_cast<const VectorType *>(&vv) != nullptr),
-             ExcVectorTypeNotCompatible());
-      const VectorType &v = dynamic_cast<const VectorType &>(vv);
-      Assert((dynamic_cast<const VectorType *>(&ww) != nullptr),
-             ExcVectorTypeNotCompatible());
-      const VectorType &w = dynamic_cast<const VectorType &>(ww);
-
+      const Number                           a,
+      const Vector<Number, MemorySpaceType> &v,
+      const Vector<Number, MemorySpaceType> &w)
+    {
       Number local_result = add_and_dot_local(a, v, w);
       if (partitioner->n_mpi_processes() > 1)
         return Utilities::MPI::sum(local_result,

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.