]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove a few deprecated constructors and functions of PETSc parallel vectors.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 1 Feb 2015 22:03:53 +0000 (16:03 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 1 Feb 2015 22:04:50 +0000 (16:04 -0600)
doc/news/changes.h
include/deal.II/lac/constraint_matrix.templates.h
include/deal.II/lac/petsc_parallel_vector.h
source/lac/petsc_parallel_vector.cc

index f18722adc742fe33bf4cdc543aead27ca02259bd..39277d8a952ed3eeed56d17f539f0932045e5172 100644 (file)
@@ -210,6 +210,7 @@ inconvenience this causes.
   - Multigrid::vmult and friends.
   - Classes FEEvaluationDGP, FEEvaluationGeneral and FEEvaluationGL.
   - Classes PreconditionedMatrix and PreconditionLACSolver.
+  - PETScVectors::MPI::Vector constructors and reinit() variants.
   <br>
   This release also removes the deprecated class MGDoFHandler. The
   functionality of this class had previously been incorporated into
index 258923dea375028020d314a4f2715a01f50938f9..3e37eb9e384fe9910268d5f192ac1093e9c2087b 100644 (file)
@@ -979,7 +979,7 @@ namespace internal
                                        PETScWrappers::MPI::Vector       &output,
                                        const internal::bool2type<false>  /*is_block_vector*/)
     {
-      output.reinit (vec.get_mpi_communicator(), locally_owned_elements, needed_elements);
+      output.reinit (locally_owned_elements, needed_elements, vec.get_mpi_communicator());
       output = vec;
     }
 #endif
index 2e3ff3ea056f84ba3e6c543a27dd2eba5df00a7f..e2912a1abb2705cad45c41a8ac7db3cf05d16502 100644 (file)
@@ -235,30 +235,6 @@ namespace PETScWrappers
                        const VectorBase   &v,
                        const size_type     local_size);
 
-
-      /**
-       * Constructs a new parallel ghosted PETSc vector from an Indexset. Note
-       * that @p local must be contiguous and the global size of the vector is
-       * determined by local.size(). The global indices in @p ghost are
-       * supplied as ghost indices that can also be read locally.
-       *
-       * Note that the @p ghost IndexSet may be empty and that any indices
-       * already contained in @p local are ignored during construction. That
-       * way, the ghost parameter can equal the set of locally relevant
-       * degrees of freedom, see step-32.
-       *
-       * @deprecated Use Vector::Vector(const IndexSet &,const IndexSet
-       * &,const MPI_Comm &) instead.
-       *
-       * @note This operation always creates a ghosted vector.
-       *
-       * @see
-       * @ref GlossGhostedVector "vectors with ghost elements"
-       */
-      explicit Vector (const MPI_Comm     &communicator,
-                       const IndexSet   &local,
-                       const IndexSet &ghost) DEAL_II_DEPRECATED;
-
       /**
        * Constructs a new parallel ghosted PETSc vector from an Indexset. Note
        * that @p local must be contiguous and the global size of the vector is
@@ -279,16 +255,6 @@ namespace PETScWrappers
               const IndexSet &ghost,
               const MPI_Comm &communicator);
 
-      /**
-       * Constructs a new parallel PETSc vector from an IndexSet. This creates
-       * a non ghosted vector.
-       *
-       * @deprecated Use Vector::Vector(const IndexSet &,const MPI_Comm &)
-       * instead.
-       */
-      explicit Vector (const MPI_Comm &communicator,
-                       const IndexSet &local) DEAL_II_DEPRECATED;
-
       /**
        * Constructs a new parallel PETSc vector from an IndexSet. This creates
        * a non ghosted vector.
@@ -372,19 +338,6 @@ namespace PETScWrappers
       void reinit (const Vector &v,
                    const bool    fast = false);
 
-      /**
-       * Reinit as a ghosted vector. See constructor with same signature for
-       * more details.
-       *
-       * @deprecated Use Vector::reinit(const IndexSet &, const IndexSet &,
-       * const MPI_Comm &) instead.
-       *
-       * @see
-       * @ref GlossGhostedVector "vectors with ghost elements"
-       */
-      void reinit (const MPI_Comm     &communicator,
-                   const IndexSet   &local,
-                   const IndexSet &ghost) DEAL_II_DEPRECATED;
       /**
        * Reinit as a vector without ghost elements. See the constructor with
        * same signature for more details.
@@ -396,16 +349,6 @@ namespace PETScWrappers
                    const IndexSet &ghost,
                    const MPI_Comm &communicator);
 
-      /**
-       * Reinit as a vector without ghost elements. See constructor with same
-       * signature for more details.
-       *
-       * @deprecated Use Vector::reinit(const IndexSet &, const MPI_Comm &)
-       * instead.
-       */
-      void reinit (const MPI_Comm     &communicator,
-                   const IndexSet   &local) DEAL_II_DEPRECATED;
-
       /**
        * Reinit as a vector without ghost elements. See constructor with same
        * signature for more details.
index 7675bff7a42c431a48a2de4f6c568c483db6a3d6..301f7c0e1663a3f41311d86933603fe8302a3bd6 100644 (file)
@@ -67,20 +67,6 @@ namespace PETScWrappers
 
 
 
-    Vector::Vector (const MPI_Comm     &communicator,
-                    const IndexSet   &local,
-                    const IndexSet &ghost)
-      :
-      communicator (communicator)
-    {
-      Assert(local.is_contiguous(), ExcNotImplemented());
-
-      IndexSet ghost_set = ghost;
-      ghost_set.subtract_set(local);
-
-      Vector::create_vector(local.size(), local.n_elements(), ghost_set);
-    }
-
     Vector::Vector (const IndexSet   &local,
                     const IndexSet &ghost,
                     const MPI_Comm     &communicator)
@@ -106,15 +92,6 @@ namespace PETScWrappers
     }
 
 
-    Vector::Vector (const MPI_Comm     &communicator,
-                    const IndexSet   &local)
-      :
-      communicator (communicator)
-    {
-      Assert(local.is_contiguous(), ExcNotImplemented());
-      Vector::create_vector(local.size(), local.n_elements());
-    }
-
     void
     Vector::reinit (const MPI_Comm  &comm,
                     const size_type  n,
@@ -180,14 +157,6 @@ namespace PETScWrappers
 
 
 
-    void
-    Vector::reinit (const MPI_Comm     &comm,
-                    const IndexSet   &local,
-                    const IndexSet &ghost)
-    {
-      reinit(local, ghost, comm);
-    }
-
     void
     Vector::reinit (const IndexSet   &local,
                     const IndexSet &ghost,
@@ -211,13 +180,6 @@ namespace PETScWrappers
       create_vector(local.size(), local.n_elements(), ghost_set);
     }
 
-    void
-    Vector::reinit (const MPI_Comm     &comm,
-                    const IndexSet   &local)
-    {
-      reinit(local, comm);
-    }
-
     void
     Vector::reinit (const IndexSet &local,
                     const MPI_Comm &comm)

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.