<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.
}
}
+//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
* 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:
* 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
* 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
* 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
* 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;
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());
// 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,
// 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++);
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
{
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++);
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());
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());
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++);
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());
// $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
#include <lac/sparse_matrix.templates.h>
+#include <lac/block_vector.h>
#define TYPEMAT double
// $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
#include <lac/sparse_matrix.templates.h>
+#include <lac/block_vector.h>
#define TYPEMAT float
// $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>::