]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added second template argument to BlockMatrixArray representing generic BlockVector
authorMaien Hamed <maien@maien.uct.ac.za>
Thu, 25 Jun 2015 12:35:56 +0000 (14:35 +0200)
committerMaien Hamed <maien@maien.uct.ac.za>
Thu, 25 Jun 2015 12:35:56 +0000 (14:35 +0200)
include/deal.II/lac/block_matrix_array.h
include/deal.II/lac/pointer_matrix.h
source/lac/block_matrix_array.cc

index 9cdc19c4974ede5e2aa32e1e87d65ec586da2a5c..640d0e9921d8b878706d17fb1f432b2c24248467 100644 (file)
@@ -31,9 +31,6 @@
 
 DEAL_II_NAMESPACE_OPEN
 
-template <typename> class BlockVector;
-template <typename> class Vector;
-
 /*! @addtogroup Matrix2
  *@{
  */
@@ -117,7 +114,7 @@ template <typename> class Vector;
  * @author Guido Kanschat
  * @date 2000-2005, 2010
  */
-template <typename number = double>
+template <typename number = double, typename BLOCK_VECTOR=BlockVector<number> >
 class BlockMatrixArray : public Subscriptor
 {
 public:
@@ -165,7 +162,7 @@ public:
   void enter (const MATRIX       &matrix,
               const unsigned int  row,
               const unsigned int  col,
-              const double        prefix = 1.,
+              const number        prefix = 1.,
               const bool          transpose = false);
 
   /**
@@ -186,39 +183,39 @@ public:
   /**
    * Matrix-vector multiplication.
    */
-  void vmult (BlockVector<number> &dst,
-              const BlockVector<number> &src) const;
+  void vmult (BLOCK_VECTOR &dst,
+              const BLOCK_VECTOR &src) const;
 
   /**
    * Matrix-vector multiplication adding to <tt>dst</tt>.
    */
-  void vmult_add (BlockVector<number> &dst,
-                  const BlockVector<number> &src) const;
+       void vmult_add(BLOCK_VECTOR &dst,
+                  const BLOCK_VECTOR &src) const;
 
   /**
    * Transposed matrix-vector multiplication.
    */
-  void Tvmult (BlockVector<number> &dst,
-               const BlockVector<number> &src) const;
+  void Tvmult (BLOCK_VECTOR &dst,
+               const BLOCK_VECTOR &src) const;
 
   /**
    * Transposed matrix-vector multiplication adding to <tt>dst</tt>.
    */
-  void Tvmult_add (BlockVector<number> &dst,
-                   const BlockVector<number> &src) const;
+  void Tvmult_add (BLOCK_VECTOR &dst,
+                   const BLOCK_VECTOR &src) const;
 
   /**
    * Matrix scalar product between two vectors (at least for a symmetric
    * matrix).
    */
-  number matrix_scalar_product (const BlockVector<number> &u,
-                                const BlockVector<number> &v) const;
+  number matrix_scalar_product (const BLOCK_VECTOR &u,
+                                const BLOCK_VECTOR &v) const;
 
   /**
    * Compute $u^T M u$. This is the square of the norm induced by the matrix
    * assuming the matrix is symmetric positive definitive.
    */
-  number matrix_norm_square (const BlockVector<number> &u) const;
+  number matrix_norm_square (const BLOCK_VECTOR &u) const;
 
   /**
    * Print the block structure as a LaTeX-array. This output will not be very
@@ -282,7 +279,7 @@ protected:
     template<class MATRIX>
     Entry (const MATRIX &matrix,
            size_type row, size_type col,
-           double prefix, bool transpose);
+           number prefix, bool transpose);
 
     /**
      * Copy constructor invalidating the old object. Since it is only used for
@@ -312,7 +309,7 @@ protected:
     /**
      * Factor in front of the matrix block.
      */
-    double prefix;
+    number prefix;
 
     /**
      * Indicates that matrix block must be transposed for multiplication.
@@ -322,7 +319,7 @@ protected:
     /**
      * The matrix block itself.
      */
-    PointerMatrixBase<Vector<number> > *matrix;
+    PointerMatrixBase<typename BLOCK_VECTOR::BlockType > *matrix;
   };
 
   /**
@@ -398,9 +395,9 @@ private:
  * @ingroup Preconditioners
  * @author Guido Kanschat, 2001, 2005
  */
-template <typename number = double>
+template <typename number = double, typename BLOCK_VECTOR = BlockVector<number> >
 class BlockTrianglePrecondition
-  : private BlockMatrixArray<number>
+  : private BlockMatrixArray<number,BLOCK_VECTOR>
 {
 public:
   /**
@@ -434,32 +431,32 @@ public:
   void enter (const MATRIX   &matrix,
               const size_type row,
               const size_type col,
-              const double    prefix = 1.,
+              const number    prefix = 1.,
               const bool      transpose = false);
 
   /**
    * Preconditioning.
    */
-  void vmult (BlockVector<number> &dst,
-              const BlockVector<number> &src) const;
+  void vmult (BLOCK_VECTOR &dst,
+              const BLOCK_VECTOR &src) const;
 
   /**
    * Preconditioning adding to <tt>dst</tt>.
    */
-  void vmult_add (BlockVector<number> &dst,
-                  const BlockVector<number> &src) const;
+  void vmult_add (BLOCK_VECTOR &dst,
+                  const BLOCK_VECTOR &src) const;
 
   /**
    * Transposed preconditioning
    */
-  void Tvmult (BlockVector<number> &dst,
-               const BlockVector<number> &src) const;
+  void Tvmult (BLOCK_VECTOR &dst,
+               const BLOCK_VECTOR &src) const;
 
   /**
    * Transposed preconditioning adding to <tt>dst</tt>.
    */
-  void Tvmult_add (BlockVector<number> &dst,
-                   const BlockVector<number> &src) const;
+  void Tvmult_add (BLOCK_VECTOR &dst,
+                   const BLOCK_VECTOR &src) const;
 
   /**
    * Make function of base class available.
@@ -506,7 +503,7 @@ private:
    * Add all off-diagonal contributions and return the entry of the diagonal
    * element for one row.
    */
-  void do_row (BlockVector<number> &dst,
+  void do_row (BLOCK_VECTOR &dst,
                size_type row_num) const;
 
   /**
@@ -519,34 +516,34 @@ private:
 #ifndef DOXYGEN
 //---------------------------------------------------------------------------
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 template <class MATRIX>
 inline
-BlockMatrixArray<number>::Entry::Entry (
+BlockMatrixArray<number, BLOCK_VECTOR>::Entry::Entry (
   const MATRIX &m,
   size_type row,
   size_type col,
-  double prefix,
+  number prefix,
   bool transpose)
   :
   row (row),
   col (col),
   prefix (prefix),
   transpose (transpose),
-  matrix (new_pointer_matrix_base(m, Vector<number>(), typeid(*this).name()))
+  matrix (new_pointer_matrix_base(m, typename BLOCK_VECTOR::BlockType(), typeid(*this).name()))
 {}
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 template <class MATRIX>
 inline
 void
-BlockMatrixArray<number>::enter (
+BlockMatrixArray<number, BLOCK_VECTOR>::enter (
   const MATRIX &matrix,
   unsigned int row,
   unsigned int col,
-  double prefix,
+  number prefix,
   bool transpose)
 {
   Assert(row<n_block_rows(), ExcIndexRange(row, 0, n_block_rows()));
@@ -555,11 +552,11 @@ BlockMatrixArray<number>::enter (
 }
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 template <class STREAM>
 inline
 void
-BlockMatrixArray<number>::print_latex (STREAM &out) const
+BlockMatrixArray<number, BLOCK_VECTOR>::print_latex (STREAM &out) const
 {
   out << "\\begin{array}{"
       << std::string(n_block_cols(), 'c')
@@ -567,7 +564,7 @@ BlockMatrixArray<number>::print_latex (STREAM &out) const
 
   Table<2,std::string> array(n_block_rows(), n_block_cols());
 
-  typedef std::map<const PointerMatrixBase<Vector<number> > *, std::string> NameMap;
+  typedef std::map<const PointerMatrixBase<typename BLOCK_VECTOR::BlockType > *, std::string> NameMap;
   NameMap matrix_names;
 
   typename std::vector<Entry>::const_iterator m = entries.begin();
@@ -580,7 +577,7 @@ BlockMatrixArray<number>::print_latex (STREAM &out) const
         {
           std::pair<typename NameMap::iterator, bool> x =
             matrix_names.insert(
-              std::pair<const PointerMatrixBase<Vector<number> >*, std::string> (m->matrix,
+              std::pair<const PointerMatrixBase<typename BLOCK_VECTOR::BlockType >*, std::string> (m->matrix,
                   std::string("M")));
           std::ostringstream stream;
           stream << matrix_number++;
@@ -618,15 +615,15 @@ BlockMatrixArray<number>::print_latex (STREAM &out) const
   out << "\\end{array}" << std::endl;
 }
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 template <class MATRIX>
 inline
 void
-BlockTrianglePrecondition<number>::enter (const MATRIX &matrix,
+BlockTrianglePrecondition<number, BLOCK_VECTOR>::enter (const MATRIX &matrix,
                                           size_type row, size_type col,
-                                          double prefix, bool transpose)
+                                          number prefix, bool transpose)
 {
-  BlockMatrixArray<number>::enter(matrix, row, col, prefix, transpose);
+  BlockMatrixArray<number, BLOCK_VECTOR>::enter(matrix, row, col, prefix, transpose);
 }
 
 
index f4a5bac95c843067b89e0be63aa1f81b4f3c5377..6181ca95e9052979a318a4572f902c03cca840f4 100644 (file)
@@ -32,9 +32,8 @@ template <typename number> class SparseMatrix;
 template <typename number> class BlockSparseMatrix;
 template <typename number> class SparseMatrixEZ;
 template <typename number> class BlockSparseMatrixEZ;
-template <typename number> class BlockMatrixArray;
 template <typename number> class TridiagonalMatrix;
-
+template <typename number, typename BlockVectorType> class BlockMatrixArray;
 
 /*! @addtogroup Matrix2
  *@{
@@ -549,11 +548,11 @@ new_pointer_matrix_base(const BlockSparseMatrixEZ<numberm> &matrix, const VECTOR
  *
  * @relates PointerMatrixBase @relates PointerMatrix
  */
-template <typename numberv, typename numberm>
-PointerMatrixBase<BlockVector<numberv> > *
-new_pointer_matrix_base(const BlockMatrixArray<numberm> &matrix, const BlockVector<numberv> &, const char *name = "PointerMatrix")
+template <typename numberv, typename numberm, typename BLOCK_VECTOR=BlockVector<numberv> >
+PointerMatrixBase<BLOCK_VECTOR> *
+new_pointer_matrix_base(const BlockMatrixArray<numberm,BLOCK_VECTOR> &matrix, const BLOCK_VECTOR &, const char *name = "PointerMatrix")
 {
-  return new PointerMatrix<BlockMatrixArray<numberm>, BlockVector<numberv> >(&matrix, name);
+  return new PointerMatrix<BlockMatrixArray<numberm,BLOCK_VECTOR>, BlockVector<numberv> >(&matrix, name);
 }
 
 
index e5ac5df4db65aed959d5467ae0e62a1c5fce9d5e..edfe73058c0de15283050e7d82a070181a85f398 100644 (file)
 #include <deal.II/lac/block_matrix_array.h>
 #include <deal.II/lac/vector.h>
 #include <deal.II/lac/block_vector.h>
+#include <deal.II/lac/trilinos_block_vector.h>
 
 DEAL_II_NAMESPACE_OPEN
 
 
-template <typename number>
-BlockMatrixArray<number>::Entry::Entry (const Entry &e)
+template <typename number, typename BLOCK_VECTOR>
+BlockMatrixArray<number,BLOCK_VECTOR>::Entry::Entry (const Entry &e)
   :
   row(e.row),
   col(e.col),
@@ -36,25 +37,24 @@ BlockMatrixArray<number>::Entry::Entry (const Entry &e)
 
 
 
-template <typename number>
-BlockMatrixArray<number>::Entry::~Entry ()
+template <typename number, typename BLOCK_VECTOR>
+BlockMatrixArray<number,BLOCK_VECTOR>::Entry::~Entry ()
 {
   if (matrix)
     delete matrix;
 }
 
 
-
-template <typename number>
-BlockMatrixArray<number>::BlockMatrixArray ()
+template <typename number, typename BLOCK_VECTOR>
+BlockMatrixArray<number,BLOCK_VECTOR>::BlockMatrixArray ()
   : block_rows (0),
     block_cols (0)
 {}
 
 
 
-template <typename number>
-BlockMatrixArray<number>::BlockMatrixArray (
+template <typename number, typename BLOCK_VECTOR>
+BlockMatrixArray<number,BLOCK_VECTOR>::BlockMatrixArray (
   const unsigned int n_block_rows,
   const unsigned int n_block_cols)
   : block_rows (n_block_rows),
@@ -62,9 +62,9 @@ BlockMatrixArray<number>::BlockMatrixArray (
 {}
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::initialize (
+BlockMatrixArray<number,BLOCK_VECTOR>::initialize (
   const unsigned int n_block_rows,
   const unsigned int n_block_cols)
 {
@@ -74,9 +74,9 @@ BlockMatrixArray<number>::initialize (
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::reinit (
+BlockMatrixArray<number,BLOCK_VECTOR>::reinit (
   const unsigned int n_block_rows,
   const unsigned int n_block_cols)
 {
@@ -87,27 +87,27 @@ BlockMatrixArray<number>::reinit (
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::clear ()
+BlockMatrixArray<number,BLOCK_VECTOR>::clear ()
 {
   entries.clear();
 }
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::vmult_add (BlockVector<number> &dst,
-                                     const BlockVector<number> &src) const
+BlockMatrixArray<number,BLOCK_VECTOR>::vmult_add (BLOCK_VECTOR &dst,
+                                     const BLOCK_VECTOR &src) const
 {
-  GrowingVectorMemory<Vector<number> > mem;
+  GrowingVectorMemory<typename BLOCK_VECTOR::BlockType > mem;
   Assert (dst.n_blocks() == block_rows,
           ExcDimensionMismatch(dst.n_blocks(), block_rows));
   Assert (src.n_blocks() == block_cols,
           ExcDimensionMismatch(src.n_blocks(), block_cols));
 
-  typename VectorMemory<Vector<number> >::Pointer p_aux(mem);
-  Vector<number> &aux = *p_aux;
+  typename VectorMemory<typename BLOCK_VECTOR::BlockType >::Pointer p_aux(mem);
+  typename BLOCK_VECTOR::BlockType &aux = *p_aux;
 
   typename std::vector<Entry>::const_iterator m = entries.begin();
   typename std::vector<Entry>::const_iterator end = entries.end();
@@ -126,10 +126,10 @@ BlockMatrixArray<number>::vmult_add (BlockVector<number> &dst,
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::vmult (BlockVector<number> &dst,
-                                 const BlockVector<number> &src) const
+BlockMatrixArray<number,BLOCK_VECTOR>::vmult (BLOCK_VECTOR &dst,
+                                 const BLOCK_VECTOR &src) const
 {
   dst = 0.;
   vmult_add (dst, src);
@@ -138,12 +138,12 @@ BlockMatrixArray<number>::vmult (BlockVector<number> &dst,
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::Tvmult_add (BlockVector<number> &dst,
-                                      const BlockVector<number> &src) const
+BlockMatrixArray<number,BLOCK_VECTOR>::Tvmult_add (BLOCK_VECTOR &dst,
+                                      const BLOCK_VECTOR &src) const
 {
-  GrowingVectorMemory<Vector<number> > mem;
+  GrowingVectorMemory<typename BLOCK_VECTOR::BlockType > mem;
   Assert (dst.n_blocks() == block_cols,
           ExcDimensionMismatch(dst.n_blocks(), block_cols));
   Assert (src.n_blocks() == block_rows,
@@ -152,8 +152,8 @@ BlockMatrixArray<number>::Tvmult_add (BlockVector<number> &dst,
   typename std::vector<Entry>::const_iterator m = entries.begin();
   typename std::vector<Entry>::const_iterator end = entries.end();
 
-  typename VectorMemory<Vector<number> >::Pointer p_aux(mem);
-  Vector<number> &aux = *p_aux;
+  typename VectorMemory<typename BLOCK_VECTOR::BlockType >::Pointer p_aux(mem);
+  typename BLOCK_VECTOR::BlockType &aux = *p_aux;
 
   for (; m != end ; ++m)
     {
@@ -168,10 +168,10 @@ BlockMatrixArray<number>::Tvmult_add (BlockVector<number> &dst,
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockMatrixArray<number>::Tvmult (BlockVector<number> &dst,
-                                  const BlockVector<number> &src) const
+BlockMatrixArray<number,BLOCK_VECTOR>::Tvmult (BLOCK_VECTOR &dst,
+                                  const BLOCK_VECTOR &src) const
 {
   dst = 0.;
   Tvmult_add (dst, src);
@@ -180,20 +180,20 @@ BlockMatrixArray<number>::Tvmult (BlockVector<number> &dst,
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 number
-BlockMatrixArray<number>::matrix_scalar_product (
-  const BlockVector<number> &u,
-  const BlockVector<number> &v) const
+BlockMatrixArray<number,BLOCK_VECTOR>::matrix_scalar_product (
+  const BLOCK_VECTOR &u,
+  const BLOCK_VECTOR &v) const
 {
-  GrowingVectorMemory<Vector<number> > mem;
+  GrowingVectorMemory<typename BLOCK_VECTOR::BlockType > mem;
   Assert (u.n_blocks() == block_rows,
           ExcDimensionMismatch(u.n_blocks(), block_rows));
   Assert (v.n_blocks() == block_cols,
           ExcDimensionMismatch(v.n_blocks(), block_cols));
 
-  typename VectorMemory<Vector<number> >::Pointer p_aux(mem);
-  Vector<number> &aux = *p_aux;
+  typename VectorMemory<typename BLOCK_VECTOR::BlockType >::Pointer p_aux(mem);
+  typename BLOCK_VECTOR::BlockType &aux = *p_aux;
 
   typename std::vector<Entry>::const_iterator m;
   typename std::vector<Entry>::const_iterator end = entries.end();
@@ -220,28 +220,28 @@ BlockMatrixArray<number>::matrix_scalar_product (
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 number
-BlockMatrixArray<number>::matrix_norm_square (
-  const BlockVector<number> &u) const
+BlockMatrixArray<number,BLOCK_VECTOR>::matrix_norm_square (
+  const BLOCK_VECTOR &u) const
 {
   return matrix_scalar_product(u,u);
 }
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 unsigned int
-BlockMatrixArray<number>::n_block_rows () const
+BlockMatrixArray<number,BLOCK_VECTOR>::n_block_rows () const
 {
   return block_rows;
 }
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 unsigned int
-BlockMatrixArray<number>::n_block_cols () const
+BlockMatrixArray<number,BLOCK_VECTOR>::n_block_cols () const
 {
   return block_cols;
 }
@@ -250,47 +250,47 @@ BlockMatrixArray<number>::n_block_cols () const
 
 //---------------------------------------------------------------------------
 
-template <typename number>
-BlockTrianglePrecondition<number>::BlockTrianglePrecondition()
-  : BlockMatrixArray<number> (),
+template <typename number, typename BLOCK_VECTOR>
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::BlockTrianglePrecondition()
+  : BlockMatrixArray<number,BLOCK_VECTOR> (),
     backward(false)
 {}
 
 
-template <typename number>
-BlockTrianglePrecondition<number>::BlockTrianglePrecondition(
+template <typename number, typename BLOCK_VECTOR>
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::BlockTrianglePrecondition(
   const unsigned int block_rows)
   :
-  BlockMatrixArray<number> (block_rows, block_rows),
+  BlockMatrixArray<number,BLOCK_VECTOR> (block_rows, block_rows),
   backward(false)
 {}
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockTrianglePrecondition<number>::reinit (
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::reinit (
   const unsigned int n)
 {
-  BlockMatrixArray<number>::reinit(n,n);
+  BlockMatrixArray<number,BLOCK_VECTOR>::reinit(n,n);
 }
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockTrianglePrecondition<number>::do_row (
-  BlockVector<number> &dst,
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::do_row (
+  BLOCK_VECTOR &dst,
   size_type row_num) const
 {
-  GrowingVectorMemory<Vector<number> > mem;
-  typename std::vector<typename BlockMatrixArray<number>::Entry>::const_iterator
+  GrowingVectorMemory<typename BLOCK_VECTOR::BlockType > mem;
+  typename std::vector<typename BlockMatrixArray<number,BLOCK_VECTOR>::Entry>::const_iterator
   m = this->entries.begin();
-  typename std::vector<typename BlockMatrixArray<number>::Entry>::const_iterator
+  typename std::vector<typename BlockMatrixArray<number,BLOCK_VECTOR>::Entry>::const_iterator
   end = this->entries.end();
-  std::vector<typename std::vector<typename BlockMatrixArray<number>::Entry>::const_iterator>
+  std::vector<typename std::vector<typename BlockMatrixArray<number,BLOCK_VECTOR>::Entry>::const_iterator>
   diagonals;
 
-  typename VectorMemory<Vector<number> >::Pointer p_aux(mem);
-  Vector<number> &aux = *p_aux;
+  typename VectorMemory<typename BLOCK_VECTOR::BlockType >::Pointer p_aux(mem);
+  typename BLOCK_VECTOR::BlockType &aux = *p_aux;
 
   aux.reinit(dst.block(row_num), true);
 
@@ -357,18 +357,18 @@ BlockTrianglePrecondition<number>::do_row (
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockTrianglePrecondition<number>::vmult_add (
-  BlockVector<number> &dst,
-  const BlockVector<number> &src) const
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::vmult_add (
+  BLOCK_VECTOR &dst,
+  const BLOCK_VECTOR &src) const
 {
   Assert (dst.n_blocks() == n_block_rows(),
           ExcDimensionMismatch(dst.n_blocks(), n_block_rows()));
   Assert (src.n_blocks() == n_block_cols(),
           ExcDimensionMismatch(src.n_blocks(), n_block_cols()));
 
-  BlockVector<number> aux;
+  BLOCK_VECTOR aux;
   aux.reinit(dst);
   vmult(aux, src);
   dst.add(aux);
@@ -376,11 +376,11 @@ BlockTrianglePrecondition<number>::vmult_add (
 
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockTrianglePrecondition<number>::vmult (
-  BlockVector<number> &dst,
-  const BlockVector<number> &src) const
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::vmult (
+  BLOCK_VECTOR &dst,
+  const BLOCK_VECTOR &src) const
 {
   Assert (dst.n_blocks() == n_block_rows(),
           ExcDimensionMismatch(dst.n_blocks(), n_block_rows()));
@@ -402,21 +402,21 @@ BlockTrianglePrecondition<number>::vmult (
 
 }
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockTrianglePrecondition<number>::Tvmult (
-  BlockVector<number> &,
-  const BlockVector<number> &) const
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::Tvmult (
+  BLOCK_VECTOR &,
+  const BLOCK_VECTOR &) const
 {
   Assert (false, ExcNotImplemented());
 }
 
 
-template <typename number>
+template <typename number, typename BLOCK_VECTOR>
 void
-BlockTrianglePrecondition<number>::Tvmult_add (
-  BlockVector<number> &,
-  const BlockVector<number> &) const
+BlockTrianglePrecondition<number,BLOCK_VECTOR>::Tvmult_add (
+  BLOCK_VECTOR &,
+  const BLOCK_VECTOR &) const
 {
   Assert (false, ExcNotImplemented());
 }
@@ -426,4 +426,7 @@ template class BlockMatrixArray<double>;
 template class BlockTrianglePrecondition<float>;
 template class BlockTrianglePrecondition<double>;
 
+template class BlockMatrixArray<double, TrilinosWrappers::MPI::BlockVector>;
+template class BlockMatrixArray<float, TrilinosWrappers::MPI::BlockVector>;
+
 DEAL_II_NAMESPACE_CLOSE

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.