]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Get rid of VectorSpaceVector in the LinearAlgebra::EpetraWrappers::Vector class.
authorWolfgang Bangerth <bangerth@colostate.edu>
Mon, 3 Jul 2023 17:13:46 +0000 (11:13 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Mon, 3 Jul 2023 17:13:46 +0000 (11:13 -0600)
include/deal.II/lac/trilinos_epetra_vector.h
source/lac/trilinos_epetra_vector.cc

index b7a630548bfd27056179b681928be37df7e4cd75..e131be901966da207eb2159c880dfae021235096 100644 (file)
@@ -54,17 +54,13 @@ namespace LinearAlgebra
   {
     /**
      * This class implements a wrapper to the Trilinos distributed vector
-     * class Epetra_FEVector. This class is derived from the
-     * LinearAlgebra::VectorSpaceVector class. Note however that Epetra only
-     * works with Number = double. This class requires Trilinos to be compiled
+     * class Epetra_FEVector. This class requires Trilinos to be compiled
      * with MPI support.
      *
      * @ingroup TrilinosWrappers
      * @ingroup Vectors
      */
-    class Vector : public VectorSpaceVector<double>,
-                   public ReadVector<double>,
-                   public Subscriptor
+    class Vector : public ReadVector<double>, public Subscriptor
     {
     public:
       using value_type = double;
@@ -105,9 +101,8 @@ namespace LinearAlgebra
        * Change the dimension to that of the vector @p V. The elements of @p V are not
        * copied.
        */
-      virtual void
-      reinit(const VectorSpaceVector<double> &V,
-             const bool omit_zeroing_entries = false) override;
+      void
+      reinit(const Vector &V, const bool omit_zeroing_entries = false);
 
       /**
        * Extract a range of elements all at once.
@@ -129,8 +124,8 @@ namespace LinearAlgebra
        * Sets all elements of the vector to the scalar @p s. This operation is
        * only allowed if @p s is equal to zero.
        */
-      virtual Vector &
-      operator=(const double s) override;
+      Vector &
+      operator=(const double s);
 
       /**
        * Imports all the elements present in the vector's IndexSet from the
@@ -140,22 +135,22 @@ namespace LinearAlgebra
        * communication pattern is used multiple times. This can be used to
        * improve performance.
        */
-      virtual void
+      void
       import_elements(
         const ReadWriteVector<double> &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
+      void
       import(const ReadWriteVector<double> &V,
              VectorOperation::values        operation,
              std::shared_ptr<const Utilities::MPI::CommunicationPatternBase>
-               communication_pattern = {}) override
+               communication_pattern = {})
       {
         import_elements(V, operation, communication_pattern);
       }
@@ -163,65 +158,60 @@ namespace LinearAlgebra
       /**
        * Multiply the entire vector by a fixed factor.
        */
-      virtual Vector &
-      operator*=(const double factor) override;
+      Vector &
+      operator*=(const double factor);
 
       /**
        * Divide the entire vector by a fixed factor.
        */
-      virtual Vector &
-      operator/=(const double factor) override;
+      Vector &
+      operator/=(const double factor);
 
       /**
        * Add the vector @p V to the present one.
        */
-      virtual Vector &
-      operator+=(const VectorSpaceVector<double> &V) override;
+      Vector &
+      operator+=(const Vector &V);
 
       /**
        * Subtract the vector @p V from the present one.
        */
-      virtual Vector &
-      operator-=(const VectorSpaceVector<double> &V) override;
+      Vector &
+      operator-=(const Vector &V);
 
       /**
        * Return the scalar product of two vectors. The vectors need to have the
        * same layout.
        */
-      virtual double
-      operator*(const VectorSpaceVector<double> &V) const override;
+      double
+      operator*(const Vector &V) const;
 
       /**
        * Add @p a to all components. Note that @p is a scalar not a vector.
        */
-      virtual void
-      add(const double a) override;
+      void
+      add(const double a);
 
       /**
        * Simple addition of a multiple of a vector, i.e. <tt>*this +=
        * a*V</tt>. The vectors need to have the same layout.
        */
-      virtual void
-      add(const double a, const VectorSpaceVector<double> &V) override;
+      void
+      add(const double a, const Vector &V);
 
       /**
        * Multiple addition of multiple of a vector, i.e. <tt>*this> +=
        * a*V+b*W</tt>. The vectors need to have the same layout.
        */
-      virtual void
-      add(const double                     a,
-          const VectorSpaceVector<double> &V,
-          const double                     b,
-          const VectorSpaceVector<double> &W) override;
+      void
+      add(const double a, const Vector &V, const double b, const Vector &W);
 
       /**
        * Scaling and simple addition of a multiple of a vector, i.e. <tt>*this
        * = s*(*this)+a*V</tt>.
        */
-      virtual void
-      sadd(const double                     s,
-           const double                     a,
-           const VectorSpaceVector<double> &V) override;
+      void
+      sadd(const double s, const double a, const Vector &V);
 
       /**
        * Scale each element of this vector by the corresponding element in the
@@ -229,47 +219,47 @@ namespace LinearAlgebra
        * (and immediate re-assignment) by a diagonal scaling matrix. The
        * vectors need to have the same layout.
        */
-      virtual void
-      scale(const VectorSpaceVector<double> &scaling_factors) override;
+      void
+      scale(const Vector &scaling_factors);
 
       /**
        * Assignment <tt>*this = a*V</tt>.
        */
-      virtual void
-      equ(const double a, const VectorSpaceVector<double> &V) override;
+      void
+      equ(const double a, const Vector &V);
 
       /**
        * Return whether the vector contains only elements with value zero.
        */
-      virtual bool
-      all_zero() const override;
+      bool
+      all_zero() const;
 
       /**
        * Return the mean value of the element of this vector.
        */
-      virtual double
-      mean_value() const override;
+      double
+      mean_value() const;
 
       /**
        * 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 double
-      l1_norm() const override;
+      double
+      l1_norm() const;
 
       /**
        * Return the l<sub>2</sub> norm of the vector (i.e., the square root of
        * the sum of the square of all entries among all processors).
        */
-      virtual double
-      l2_norm() const override;
+      double
+      l2_norm() const;
 
       /**
        * Return the maximum norm of the vector (i.e., the maximum absolute value
        * among all entries and among all processors).
        */
-      virtual double
-      linfty_norm() const override;
+      double
+      linfty_norm() const;
 
       /**
        * Performs a combined operation of a vector addition and a subsequent
@@ -293,10 +283,8 @@ namespace LinearAlgebra
        * implemented as
        * $\left<v,w\right>=\sum_i v_i \bar{w_i}$.
        */
-      virtual double
-      add_and_dot(const double                     a,
-                  const VectorSpaceVector<double> &V,
-                  const VectorSpaceVector<double> &W) override;
+      double
+      add_and_dot(const double a, const Vector &V, const Vector &W);
       /**
        * This function always returns false and is present only for backward
        * compatibility.
@@ -335,8 +323,21 @@ 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;
+
+      /**
+       * Compress the underlying representation of the Trilinos object, i.e.
+       * flush the buffers of the vector object if it has any. This function is
+       * necessary after writing into a vector element-by-element and before
+       * anything else can be done on it.
+       *
+       * See
+       * @ref GlossCompress "Compressing distributed objects"
+       * for more information.
+       */
+      void
+      compress(const VectorOperation::values operation);
 
       /**
        * Return a const reference to the underlying Trilinos
@@ -355,17 +356,17 @@ namespace LinearAlgebra
       /**
        * Prints 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;
 
       /**
        * The vectors have different partitioning, i.e. their IndexSet objects
index 9a69a287c8dcf09ed31c1830d03030ba86bdaa8f..4cc1411abd669ca12910a5ab9efcdaf409aab088 100644 (file)
@@ -82,18 +82,10 @@ namespace LinearAlgebra
 
 
     void
-    Vector::reinit(const VectorSpaceVector<double> &V,
-                   const bool                       omit_zeroing_entries)
+    Vector::reinit(const Vector &V, const bool omit_zeroing_entries)
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast V. If fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
-
-      reinit(down_V.locally_owned_elements(),
-             down_V.get_mpi_communicator(),
+      reinit(V.locally_owned_elements(),
+             V.get_mpi_communicator(),
              omit_zeroing_entries);
     }
 
@@ -243,30 +235,23 @@ namespace LinearAlgebra
 
 
     Vector &
-    Vector::operator+=(const VectorSpaceVector<double> &V)
+    Vector::operator+=(const Vector &V)
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast V. If fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
       // If the maps are the same we can Update right away.
-      if (vector->Map().SameAs(down_V.trilinos_vector().Map()))
+      if (vector->Map().SameAs(V.trilinos_vector().Map()))
         {
-          const int ierr = vector->Update(1., down_V.trilinos_vector(), 1.);
+          const int ierr = vector->Update(1., V.trilinos_vector(), 1.);
           Assert(ierr == 0, ExcTrilinosError(ierr));
           (void)ierr;
         }
       else
         {
-          Assert(this->size() == down_V.size(),
-                 ExcDimensionMismatch(this->size(), down_V.size()));
+          Assert(this->size() == V.size(),
+                 ExcDimensionMismatch(this->size(), V.size()));
 
 #  if DEAL_II_TRILINOS_VERSION_GTE(11, 11, 0)
-          Epetra_Import data_exchange(vector->Map(),
-                                      down_V.trilinos_vector().Map());
-          const int     ierr = vector->Import(down_V.trilinos_vector(),
+          Epetra_Import data_exchange(vector->Map(), V.trilinos_vector().Map());
+          const int     ierr = vector->Import(V.trilinos_vector(),
                                           data_exchange,
                                           Epetra_AddLocalAlso);
           Assert(ierr == 0, ExcTrilinosError(ierr));
@@ -276,11 +261,9 @@ namespace LinearAlgebra
           // adding Hence, we provide a workaround in this case
 
           Epetra_MultiVector dummy(vector->Map(), 1, false);
-          Epetra_Import      data_exchange(dummy.Map(),
-                                      down_V.trilinos_vector().Map());
+          Epetra_Import data_exchange(dummy.Map(), V.trilinos_vector().Map());
 
-          int ierr =
-            dummy.Import(down_V.trilinos_vector(), data_exchange, Insert);
+          int ierr = dummy.Import(V.trilinos_vector(), data_exchange, Insert);
           Assert(ierr == 0, ExcTrilinosError(ierr));
 
           ierr = vector->Update(1.0, dummy, 1.0);
@@ -295,7 +278,7 @@ namespace LinearAlgebra
 
 
     Vector &
-    Vector::operator-=(const VectorSpaceVector<double> &V)
+    Vector::operator-=(const Vector &V)
     {
       this->add(-1., V);
 
@@ -305,21 +288,15 @@ namespace LinearAlgebra
 
 
     double
-    Vector::operator*(const VectorSpaceVector<double> &V) const
+    Vector::operator*(const Vector &V) const
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast V. If fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
-      Assert(this->size() == down_V.size(),
-             ExcDimensionMismatch(this->size(), down_V.size()));
-      Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
+      Assert(this->size() == V.size(),
+             ExcDimensionMismatch(this->size(), V.size()));
+      Assert(vector->Map().SameAs(V.trilinos_vector().Map()),
              ExcDifferentParallelPartitioning());
 
       double    result(0.);
-      const int ierr = vector->Dot(down_V.trilinos_vector(), &result);
+      const int ierr = vector->Dot(V.trilinos_vector(), &result);
       Assert(ierr == 0, ExcTrilinosError(ierr));
       (void)ierr;
 
@@ -340,19 +317,13 @@ namespace LinearAlgebra
 
 
     void
-    Vector::add(const double a, const VectorSpaceVector<double> &V)
+    Vector::add(const double a, const Vector &V)
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast V. If fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
       AssertIsFinite(a);
-      Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
+      Assert(vector->Map().SameAs(V.trilinos_vector().Map()),
              ExcDifferentParallelPartitioning());
 
-      const int ierr = vector->Update(a, down_V.trilinos_vector(), 1.);
+      const int ierr = vector->Update(a, V.trilinos_vector(), 1.);
       Assert(ierr == 0, ExcTrilinosError(ierr));
       (void)ierr;
     }
@@ -360,31 +331,20 @@ namespace LinearAlgebra
 
 
     void
-    Vector::add(const double                     a,
-                const VectorSpaceVector<double> &V,
-                const double                     b,
-                const VectorSpaceVector<double> &W)
-    {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&W) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast V. If fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
-      // Downcast W. If fails, throws an exception.
-      const Vector &down_W = dynamic_cast<const Vector &>(W);
-      Assert(vector->Map().SameAs(down_V.trilinos_vector().Map()),
+    Vector::add(const double  a,
+                const Vector &V,
+                const double  b,
+                const Vector &W)
+    {
+      Assert(vector->Map().SameAs(V.trilinos_vector().Map()),
              ExcDifferentParallelPartitioning());
-      Assert(vector->Map().SameAs(down_W.trilinos_vector().Map()),
+      Assert(vector->Map().SameAs(W.trilinos_vector().Map()),
              ExcDifferentParallelPartitioning());
       AssertIsFinite(a);
       AssertIsFinite(b);
 
-      const int ierr = vector->Update(
-        a, down_V.trilinos_vector(), b, down_W.trilinos_vector(), 1.);
+      const int ierr =
+        vector->Update(a, V.trilinos_vector(), b, W.trilinos_vector(), 1.);
       Assert(ierr == 0, ExcTrilinosError(ierr));
       (void)ierr;
     }
@@ -392,18 +352,10 @@ namespace LinearAlgebra
 
 
     void
-    Vector::sadd(const double                     s,
-                 const double                     a,
-                 const VectorSpaceVector<double> &V)
+    Vector::sadd(const double s, const double a, const Vector &V)
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-
       *this *= s;
-      // Downcast V. It fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
-      Vector        tmp(down_V);
+      Vector tmp(V);
       tmp *= a;
       *this += tmp;
     }
@@ -411,22 +363,13 @@ namespace LinearAlgebra
 
 
     void
-    Vector::scale(const VectorSpaceVector<double> &scaling_factors)
+    Vector::scale(const Vector &scaling_factors)
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&scaling_factors) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast scaling_factors. If fails, throws an exception.
-      const Vector &down_scaling_factors =
-        dynamic_cast<const Vector &>(scaling_factors);
-      Assert(vector->Map().SameAs(down_scaling_factors.trilinos_vector().Map()),
+      Assert(vector->Map().SameAs(scaling_factors.trilinos_vector().Map()),
              ExcDifferentParallelPartitioning());
 
-      const int ierr = vector->Multiply(1.0,
-                                        down_scaling_factors.trilinos_vector(),
-                                        *vector,
-                                        0.0);
+      const int ierr =
+        vector->Multiply(1.0, scaling_factors.trilinos_vector(), *vector, 0.0);
       Assert(ierr == 0, ExcTrilinosError(ierr));
       (void)ierr;
     }
@@ -434,21 +377,15 @@ namespace LinearAlgebra
 
 
     void
-    Vector::equ(const double a, const VectorSpaceVector<double> &V)
+    Vector::equ(const double a, const Vector &V)
     {
-      // Check that casting will work.
-      Assert(dynamic_cast<const Vector *>(&V) != nullptr,
-             ExcVectorTypeNotCompatible());
-
-      // Downcast V. If fails, throws an exception.
-      const Vector &down_V = dynamic_cast<const Vector &>(V);
       // If we don't have the same map, copy.
-      if (vector->Map().SameAs(down_V.trilinos_vector().Map()) == false)
+      if (vector->Map().SameAs(V.trilinos_vector().Map()) == false)
         this->sadd(0., a, V);
       else
         {
           // Otherwise, just update
-          int ierr = vector->Update(a, down_V.trilinos_vector(), 0.);
+          int ierr = vector->Update(a, V.trilinos_vector(), 0.);
           Assert(ierr == 0, ExcTrilinosError(ierr));
           (void)ierr;
         }
@@ -539,9 +476,7 @@ namespace LinearAlgebra
 
 
     double
-    Vector::add_and_dot(const double                     a,
-                        const VectorSpaceVector<double> &V,
-                        const VectorSpaceVector<double> &W)
+    Vector::add_and_dot(const double a, const Vector &V, const Vector &W)
     {
       this->add(a, V);
 
@@ -614,6 +549,11 @@ namespace LinearAlgebra
     }
 
 
+    void
+    Vector::compress(const VectorOperation::values /*operation*/)
+    {}
+
+
 
     const Epetra_FEVector &
     Vector::trilinos_vector() const

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.