]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Get rid of PETScWrappers::Vector in the PETSc matrix classes.
authorDavid Wells <wellsd2@rpi.edu>
Sat, 29 Apr 2017 16:35:55 +0000 (12:35 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Thu, 4 May 2017 01:45:41 +0000 (21:45 -0400)
We can simply use PETScWrappers::VectorBase instead. In addition, get rid of
some methods in the derived class that are identical to those in MatrixBase.

include/deal.II/lac/petsc_sparse_matrix.h
source/lac/petsc_matrix_base.cc
source/lac/petsc_sparse_matrix.cc

index 41a1481f71c10a7b538a0b1d605b3b4886eb9aa3..ce88ecc7805670dbde81e521a0e92431176b85e6 100644 (file)
@@ -198,33 +198,6 @@ namespace PETScWrappers
      */
     virtual const MPI_Comm &get_mpi_communicator () const;
 
-    /**
-     * Return the square of the norm of the vector $v$ with respect to the
-     * norm induced by this matrix, i.e. $\left(v,Mv\right)$. This is useful,
-     * e.g. in the finite element context, where the $L_2$ norm of a function
-     * equals the matrix norm with respect to the mass matrix of the vector
-     * representing the nodal values of the finite element function.
-     *
-     * Obviously, the matrix needs to be quadratic for this operation.
-     *
-     * The implementation of this function is not as efficient as the one in
-     * the @p MatrixBase class used in deal.II (i.e. the original one, not the
-     * PETSc wrapper class) since PETSc doesn't support this operation and
-     * needs a temporary vector.
-     */
-    PetscScalar matrix_norm_square (const VectorBase &v) const;
-
-    /**
-     * Compute the matrix scalar product $\left(u,Mv\right)$.
-     *
-     * The implementation of this function is not as efficient as the one in
-     * the @p MatrixBase class used in deal.II (i.e. the original one, not the
-     * PETSc wrapper class) since PETSc doesn't support this operation and
-     * needs a temporary vector.
-     */
-    PetscScalar matrix_scalar_product (const VectorBase &u,
-                                       const VectorBase &v) const;
-
     /**
      * Return the number of rows of this matrix.
      */
index ac5c4b77db930cdaa7f902bb0a4b386aea2289d7..41dc140cd1094487a0cc4f8ae6a41afc8da5e2b4 100644 (file)
@@ -22,7 +22,7 @@
 #  include <deal.II/lac/petsc_full_matrix.h>
 #  include <deal.II/lac/petsc_sparse_matrix.h>
 #  include <deal.II/lac/petsc_parallel_sparse_matrix.h>
-#  include <deal.II/lac/petsc_vector.h>
+#  include <deal.II/lac/petsc_vector_base.h>
 
 DEAL_II_NAMESPACE_OPEN
 
@@ -384,8 +384,8 @@ namespace PETScWrappers
   PetscScalar
   MatrixBase::matrix_norm_square (const VectorBase &v) const
   {
-    Vector tmp(v.size());
-    vmult (tmp, v);
+    VectorBase tmp(v);
+    vmult(tmp, v);
     return tmp*v;
   }
 
@@ -394,8 +394,8 @@ namespace PETScWrappers
   MatrixBase::matrix_scalar_product (const VectorBase &u,
                                      const VectorBase &v) const
   {
-    Vector tmp(v.size());
-    vmult (tmp, v);
+    VectorBase tmp(u);
+    vmult(tmp, v);
     return u*tmp;
   }
 
index 913490b7e8d547c46ffed0de90f7e35d8e6e5d8f..e71099b7129f98b4592e92ca6e45806958073077 100644 (file)
@@ -19,7 +19,7 @@
 
 #  include <deal.II/lac/exceptions.h>
 #  include <deal.II/lac/petsc_compatibility.h>
-#  include <deal.II/lac/petsc_vector.h>
+#  include <deal.II/lac/petsc_vector_base.h>
 #  include <deal.II/lac/sparsity_pattern.h>
 #  include <deal.II/lac/dynamic_sparsity_pattern.h>
 
@@ -284,23 +284,6 @@ namespace PETScWrappers
   template void
   SparseMatrix::do_reinit (const DynamicSparsityPattern &,
                            const bool);
-
-  PetscScalar
-  SparseMatrix::matrix_norm_square (const VectorBase &v) const
-  {
-    Vector tmp (v.size());
-    vmult (tmp, v);
-    return tmp*v;
-  }
-
-  PetscScalar
-  SparseMatrix::matrix_scalar_product (const VectorBase &u,
-                                       const VectorBase &v) const
-  {
-    Vector tmp (v.size());
-    vmult (tmp, v);
-    return u*tmp;
-  }
 }
 
 

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.