]> https://gitweb.dealii.org/ - dealii.git/commitdiff
SparseMatrix operates on BlockVector
authorGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 1 Apr 2004 06:59:54 +0000 (06:59 +0000)
committerGuido Kanschat <dr.guido.kanschat@gmail.com>
Thu, 1 Apr 2004 06:59:54 +0000 (06:59 +0000)
git-svn-id: https://svn.dealii.org/trunk@8934 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/doc/news/c-4-0.html
deal.II/lac/include/lac/sparse_matrix.h
deal.II/lac/include/lac/sparse_matrix.templates.h
deal.II/lac/source/sparse_matrix.double.cc
deal.II/lac/source/sparse_matrix.float.cc
deal.II/lac/source/sparse_matrix_vector.in.h

index 0e2b98131e201122b9b071b5aed306b49eb75a90..5693a2465b5fab55a18571a520379deacddba501 100644 (file)
@@ -368,6 +368,13 @@ contributor's names are abbreviated by WB (Wolfgang Bangerth), GK
 <h3>lac</h3>
 
 <ol>
+  <li> <p> Improved: The matrix-vector operations of <code
+  class="class">SparseMatrix</code> accept arguments of type <code
+  class="class">BlockVector</code>.
+  <br>
+  (GK/2004/03/31)
+  </p>
+
   <li> <p> Fixed: The <code class="class">SparseMatrix</code> iterator classes
        had various interesting bugs when rows of the matrix were completely
        empty. These should now be fixed.
index 13bd6003f51d681cd48b21b341bc31c400116210..fe253ae39363a9676d509f6a1c3fb55cbcda40da 100644 (file)
@@ -370,18 +370,20 @@ namespace internals
   }
 }
 
+//TODO: Add multithreading to the other vmult functions.
 
 /**
- * Sparse matrix. This class implements the function to store values in the
- * locations of a sparse matrix denoted by a @ref SparsityPattern. The
- * separation of sparsity pattern and values is done since one can store data
- * elements of different type in these locations without the SparsityPattern
- * having to know this, and more importantly one can associate more than one
- * matrix with the same sparsity pattern.
+ * Sparse matrix. This class implements the function to store values
+ * in the locations of a sparse matrix denoted by a
+ * SparsityPattern. The separation of sparsity pattern and values is
+ * done since one can store data elements of different type in these
+ * locations without the SparsityPattern having to know this, and more
+ * importantly one can associate more than one matrix with the same
+ * sparsity pattern.
  *
  * @ref Instantiations some
  
- * @author several, 1994-2003
+ * @author several, 1994-2004
  */
 template <typename number>
 class SparseMatrix : public virtual Subscriptor
@@ -907,9 +909,9 @@ class SparseMatrix : public virtual Subscriptor
                                       * Source and destination must
                                       * not be the same vector.
                                      */
-    template <typename somenumber>
-    void vmult (Vector<somenumber>       &dst,
-               const Vector<somenumber> &src) const;
+    template <class OutVector, class InVector>
+    void vmult (OutVector& dst,
+               const InVector& src) const;
     
                                     /**
                                      * Matrix-vector multiplication:
@@ -922,9 +924,9 @@ class SparseMatrix : public virtual Subscriptor
                                       * Source and destination must
                                       * not be the same vector.
                                      */
-    template <typename somenumber>
-    void Tvmult (Vector<somenumber>       &dst,
-                const Vector<somenumber> &src) const;
+    template <class OutVector, class InVector>
+    void Tvmult (OutVector& dst,
+                const InVector& src) const;
   
                                     /**
                                      * Adding Matrix-vector
@@ -936,9 +938,9 @@ class SparseMatrix : public virtual Subscriptor
                                       * Source and destination must
                                       * not be the same vector.
                                      */
-    template <typename somenumber>
-    void vmult_add (Vector<somenumber>       &dst,
-                   const Vector<somenumber> &src) const;
+    template <class OutVector, class InVector>
+    void vmult_add (OutVector& dst,
+                   const InVector& src) const;
     
                                     /**
                                      * Adding Matrix-vector
@@ -953,9 +955,9 @@ class SparseMatrix : public virtual Subscriptor
                                       * Source and destination must
                                       * not be the same vector.
                                      */
-    template <typename somenumber>
-    void Tvmult_add (Vector<somenumber>       &dst,
-                    const Vector<somenumber> &src) const;
+    template <class OutVector, class InVector>
+    void Tvmult_add (OutVector& dst,
+                    const InVector& src) const;
   
                                     /**
                                      * Return the square of the norm
@@ -1508,9 +1510,9 @@ class SparseMatrix : public virtual Subscriptor
                                      * in the case of enabled
                                      * multithreading.
                                      */
-    template <typename somenumber>
-    void threaded_vmult (Vector<somenumber>       &dst,
-                        const Vector<somenumber> &src,
+    template <class OutVector, class InVector>
+    void threaded_vmult (OutVector& dst,
+                        const InVector& src,
                         const unsigned int        begin_row,
                         const unsigned int        end_row) const;
 
index a7d25e43261dc6e4214eb1ef17d47ef65e05dd93..0a3547f79b035e6e9184fe1ee4d367da837959f5 100644 (file)
@@ -277,10 +277,10 @@ SparseMatrix<number>::add_scaled (const number factor,
 
 
 template <typename number>
-template <typename somenumber>
+template <class OutVector, class InVector>
 void
-SparseMatrix<number>::vmult (Vector<somenumber>& dst,
-                            const Vector<somenumber>& src) const
+SparseMatrix<number>::vmult (OutVector& dst,
+                            const InVector& src) const
 {
   Assert (cols != 0, ExcMatrixNotInitialized());
   Assert (val != 0, ExcMatrixNotInitialized());
@@ -314,13 +314,13 @@ SparseMatrix<number>::vmult (Vector<somenumber>& dst,
                                       // correct type right away:
       typedef
        void (SparseMatrix<number>::*mem_fun_p)
-       (Vector<somenumber> &,
-        const Vector<somenumber> &,
+       (OutVector &,
+        const InVector &,
         const unsigned int        ,
         const unsigned int) const;
       const mem_fun_p comp
        = (&SparseMatrix<number>::
-           template threaded_vmult<somenumber>);
+           template threaded_vmult<InVector,OutVector>);
       Threads::ThreadGroup<> threads;
       for (unsigned int i=0; i<n_threads; ++i)
        threads += Threads::spawn (*this, comp) (dst, src,
@@ -336,10 +336,10 @@ SparseMatrix<number>::vmult (Vector<somenumber>& dst,
                                       // do it in an oldfashioned way
       const number       *val_ptr    = &val[cols->rowstart[0]];
       const unsigned int *colnum_ptr = &cols->colnums[cols->rowstart[0]];
-      somenumber         *dst_ptr    = &dst(0);
+      typename OutVector::iterator dst_ptr = dst.begin();
       for (unsigned int row=0; row<n_rows; ++row)
        {
-         somenumber s = 0.;
+         typename OutVector::value_type s = 0.;
          const number *const val_end_of_row = &val[cols->rowstart[row+1]];
          while (val_ptr != val_end_of_row)
            s += *val_ptr++ * src(*colnum_ptr++);
@@ -350,10 +350,10 @@ SparseMatrix<number>::vmult (Vector<somenumber>& dst,
 
 
 template <typename number>
-template <typename somenumber>
+template <class OutVector, class InVector>
 void
-SparseMatrix<number>::threaded_vmult (Vector<somenumber>       &dst,
-                                     const Vector<somenumber> &src,
+SparseMatrix<number>::threaded_vmult (OutVector       &dst,
+                                     const InVector &src,
                                      const unsigned int        begin_row,
                                      const unsigned int        end_row) const
 {
@@ -363,10 +363,10 @@ SparseMatrix<number>::threaded_vmult (Vector<somenumber>       &dst,
 
   const number       *val_ptr    = &val[cols->rowstart[begin_row]];
   const unsigned int *colnum_ptr = &cols->colnums[cols->rowstart[begin_row]];
-  somenumber         *dst_ptr    = &dst(begin_row);
+  typename OutVector::iterator dst_ptr = dst.begin() + begin_row;
   for (unsigned int row=begin_row; row<end_row; ++row)
     {
-      somenumber s = 0.;
+      typename OutVector::value_type s = 0.;
       const number *const val_end_of_row = &val[cols->rowstart[row+1]];
       while (val_ptr != val_end_of_row)
        s += *val_ptr++ * src(*colnum_ptr++);
@@ -376,10 +376,10 @@ SparseMatrix<number>::threaded_vmult (Vector<somenumber>       &dst,
 
 
 template <typename number>
-template <typename somenumber>
+template <class OutVector, class InVector>
 void
-SparseMatrix<number>::Tvmult (Vector<somenumber>& dst,
-                              const Vector<somenumber>& src) const
+SparseMatrix<number>::Tvmult (OutVector& dst,
+                              const InVector& src) const
 {
   Assert (val != 0, ExcMatrixNotInitialized());
   Assert (cols != 0, ExcMatrixNotInitialized());
@@ -402,10 +402,10 @@ SparseMatrix<number>::Tvmult (Vector<somenumber>& dst,
 
 
 template <typename number>
-template <typename somenumber>
+template <class OutVector, class InVector>
 void
-SparseMatrix<number>::vmult_add (Vector<somenumber>& dst,
-                                 const Vector<somenumber>& src) const
+SparseMatrix<number>::vmult_add (OutVector& dst,
+                                 const InVector& src) const
 {
   Assert (cols != 0, ExcMatrixNotInitialized());
   Assert (val != 0, ExcMatrixNotInitialized());
@@ -417,10 +417,10 @@ SparseMatrix<number>::vmult_add (Vector<somenumber>& dst,
   const unsigned int  n_rows     = m();
   const number       *val_ptr    = &val[cols->rowstart[0]];
   const unsigned int *colnum_ptr = &cols->colnums[cols->rowstart[0]];
-  somenumber         *dst_ptr    = &dst(0);
+  typename OutVector::iterator dst_ptr    = dst.begin();
   for (unsigned int row=0; row<n_rows; ++row)
     {
-      somenumber s = 0.;
+      typename OutVector::value_type s = 0.;
       const number *const val_end_of_row = &val[cols->rowstart[row+1]];
       while (val_ptr != val_end_of_row)
        s += *val_ptr++ * src(*colnum_ptr++);
@@ -430,10 +430,10 @@ SparseMatrix<number>::vmult_add (Vector<somenumber>& dst,
 
 
 template <typename number>
-template <typename somenumber>
+template <class OutVector, class InVector>
 void
-SparseMatrix<number>::Tvmult_add (Vector<somenumber>& dst,
-                                  const Vector<somenumber>& src) const
+SparseMatrix<number>::Tvmult_add (OutVector& dst,
+                                  const InVector& src) const
 {
   Assert (val != 0, ExcMatrixNotInitialized());
   Assert (cols != 0, ExcMatrixNotInitialized());
index 67a3b2b1923f0827be7631f9c010dbf24c2b181c..62a5533d775352052f0e5cb96f2dbd350ae27f4c 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -13,6 +13,7 @@
 
 
 #include <lac/sparse_matrix.templates.h>
+#include <lac/block_vector.h>
 
 #define TYPEMAT double
 
index 7aae94aac1db0190f93ad107725f0365f8b195bb..05c6bd3a74a0d1b4dc921ce0baac1ef202578132 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -13,6 +13,7 @@
 
 
 #include <lac/sparse_matrix.templates.h>
+#include <lac/block_vector.h>
 
 #define TYPEMAT float
 
index c30794a372789b9476b6574b797dc5168111b5ca..cdd09856dada4125f0f3e842a0e986dbfb765d7b 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
 // TYPEMAT and TYPEVEC are defined in sparsematrix?.cc
 
 template void SparseMatrix<TYPEMAT>::
-vmult<TYPEVEC> (Vector<TYPEVEC> &,
-                const Vector<TYPEVEC> &) const;
+vmult (Vector<TYPEVEC> &, const Vector<TYPEVEC> &) const;
 template void SparseMatrix<TYPEMAT>::
-Tvmult<TYPEVEC> (Vector<TYPEVEC> &,
-                 const Vector<TYPEVEC> &) const;
+Tvmult (Vector<TYPEVEC> &, const Vector<TYPEVEC> &) const;
 template void SparseMatrix<TYPEMAT>::
-vmult_add<TYPEVEC> (Vector<TYPEVEC> &,
-                    const Vector<TYPEVEC> &) const;
+vmult_add (Vector<TYPEVEC> &, const Vector<TYPEVEC> &) const;
 template void SparseMatrix<TYPEMAT>::
-Tvmult_add<TYPEVEC> (Vector<TYPEVEC> &,
-                     const Vector<TYPEVEC> &) const;
+Tvmult_add (Vector<TYPEVEC> &, const Vector<TYPEVEC> &) const;
+
+template void SparseMatrix<TYPEMAT>::
+vmult (BlockVector<TYPEVEC> &, const BlockVector<TYPEVEC> &) const;
+template void SparseMatrix<TYPEMAT>::
+Tvmult (BlockVector<TYPEVEC> &, const BlockVector<TYPEVEC> &) const;
+template void SparseMatrix<TYPEMAT>::
+vmult_add (BlockVector<TYPEVEC> &, const BlockVector<TYPEVEC> &) const;
+template void SparseMatrix<TYPEMAT>::
+Tvmult_add (BlockVector<TYPEVEC> &, const BlockVector<TYPEVEC> &) const;
 
 template TYPEVEC
 SparseMatrix<TYPEMAT>::

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.