template <int rows, int columns> class BlockSparsityPattern;
template <typename number> class SparseMatrix;
template <typename number, int rows, int columns> class BlockSparseMatrix;
-template <int n_blocks> class BlockIndices;
+class BlockIndices;
/**
* This class represents the matrix denoting the distribution of the degrees
Assert (sparsity.get_column_indices() == sparsity.get_row_indices(),
ExcMatrixNotSquare());
- const BlockIndices<blocks> &
+ const BlockIndices &
index_mapping = sparsity.get_column_indices();
// store for each index whether it
Assert (sparsity.get_column_indices() == sparsity.get_row_indices(),
ExcMatrixNotSquare());
- const BlockIndices<blocks> &
+ const BlockIndices &
index_mapping = sparsity.get_column_indices();
// store for each index whether it
// the row and column mappings are
// equal, store a pointer on only
// one of them
- const BlockIndices<blocks> &
+ const BlockIndices &
index_mapping = sparsity_pattern.get_column_indices();
// now loop over all boundary dofs
* useful if a matrix is composed of several blocks, where you have to
* translate global row and column indices to local ones.
*
- * @author Wolfgang Bangerth, 2000
+ * @author Wolfgang Bangerth, Guido Kanschat, 2000
*/
-template <int n_blocks>
class BlockIndices
{
public:
/**
- * Defaulut constructor. Set all
- * indices denoting the start of
- * the blocks to zero.
+ * Default
+ * constructor. Initialize for
+ * @p{n_blocks} blocks and set
+ * all block sizes to zero.
*/
- BlockIndices ();
+ BlockIndices (const unsigned int n_blocks);
/**
* Constructor. Initialize the
- * number of indices within each
- * block from the given
- * argument. The size of the
- * vector shall be equal to
- * @p{n_blocks}.
+ * number of entries in each
+ * block @p{i} as @p{n[i]}. The
+ * number of blocks will be the
+ * size of the vector
*/
BlockIndices (const vector<unsigned int> &n);
/**
- * Reset the number of indices
- * within each block from the
- * given argument. The size of
- * the vector shall be equal to
- * @p{n_blocks}.
+ * Reinitialize the number of
+ * indices within each block from
+ * the given argument. The number
+ * of blocks will be adjusted to
+ * the size of @p{n} and the size
+ * of block @p{i} is set to
+ * @p{n[i]}.
*/
void reinit (const vector<unsigned int> &n);
* for the global index @p{i}. The
* first element of the pair is
* the block, the second the
- * index within that.
+ * index within it.
*/
pair<unsigned int,unsigned int>
global_to_local (const unsigned int i) const;
unsigned int local_to_global (const unsigned int block,
const unsigned int index) const;
+ /**
+ * Number of blocks in index field.
+ */
+ unsigned int size () const;
+
/**
* Return the total number of
* indices accumulated over all
- * blocks.
+ * blocks, that is, the dimension
+ * of the vector space of the
+ * block vector.
*/
unsigned int total_size () const;
/**
* Copy operator.
*/
- BlockIndices<n_blocks> & operator = (const BlockIndices<n_blocks> &b);
+ BlockIndices & operator = (const BlockIndices &b);
/**
* Compare whether two objects
* starting indices of all blocks
* are equal.
*/
- bool operator == (const BlockIndices<n_blocks> &b) const;
+ bool operator == (const BlockIndices &b) const;
/**
* Swap the contents of these two
* objects.
*/
- void swap (BlockIndices<n_blocks> &b);
+ void swap (BlockIndices &b);
private:
+ /**
+ * Number of blocks. This is made
+ * constant to avoid accidental
+ * changes during lifetime.
+ */
+ unsigned int n_blocks;
+
/**
* Global starting index of each
* vector. The last and redundant
* value is the total number of
* entries.
*/
- unsigned int start_indices[n_blocks+1];
+ vector<unsigned int> start_indices;
};
/* ---------------------- template and inline functions ------------------- */
-template <int n_blocks>
inline
-BlockIndices<n_blocks>::BlockIndices ()
+BlockIndices::BlockIndices (unsigned int n_blocks)
+ : n_blocks(n_blocks),
+ start_indices(n_blocks)
{
for (unsigned int i=0; i<=n_blocks; ++i)
start_indices[i] = 0;
-template <int n_blocks>
inline
-BlockIndices<n_blocks>::BlockIndices (const vector<unsigned int> &n)
+BlockIndices::BlockIndices (const vector<unsigned int> &n)
{
reinit (n);
};
-template <int n_blocks>
inline
void
-BlockIndices<n_blocks>::reinit (const vector<unsigned int> &n)
+BlockIndices::reinit (const vector<unsigned int> &n)
{
- Assert(n.size()==n_blocks,
- ExcDimensionMismatch(n.size(), n_blocks));
-
+ if (start_indices.size() != n.size()+1)
+ {
+ n_blocks = n.size();
+ start_indices.resize(n_blocks+1);
+ }
start_indices[0] = 0;
for (unsigned int i=1; i<=n_blocks; ++i)
start_indices[i] = start_indices[i-1] + n[i-1];
-template <int n_blocks>
inline
pair<unsigned int,unsigned int>
-BlockIndices<n_blocks>::global_to_local (const unsigned int i) const
+BlockIndices::global_to_local (const unsigned int i) const
{
Assert (i<total_size(), ExcIndexRange(i, 0, total_size()));
};
-template <int n_blocks>
inline
unsigned int
-BlockIndices<n_blocks>::local_to_global (const unsigned int block,
+BlockIndices::local_to_global (const unsigned int block,
const unsigned int index) const
{
Assert (block < n_blocks, ExcIndexRange(block, 0, n_blocks));
};
+inline
+unsigned int
+BlockIndices::size () const
+{
+ return n_blocks;
+}
+
+
-template <int n_blocks>
inline
unsigned int
-BlockIndices<n_blocks>::total_size () const
+BlockIndices::total_size () const
{
return start_indices[n_blocks];
};
-template <int n_blocks>
inline
-BlockIndices<n_blocks> &
-BlockIndices<n_blocks>::operator = (const BlockIndices<n_blocks> &b)
+BlockIndices &
+BlockIndices::operator = (const BlockIndices &b)
{
for (unsigned int i=0; i<=n_blocks; ++i)
start_indices[i] = b.start_indices[i];
-template <int n_blocks>
inline
bool
-BlockIndices<n_blocks>::operator == (const BlockIndices<n_blocks> &b) const
+BlockIndices::operator == (const BlockIndices &b) const
{
for (unsigned int i=0; i<=n_blocks; ++i)
if (start_indices[i] != b.start_indices[i])
-template <int n_blocks>
inline
void
-BlockIndices<n_blocks>::swap (BlockIndices<n_blocks> &b)
+BlockIndices::swap (BlockIndices &b)
{
for (unsigned int i=0; i<=n_blocks; ++i)
std::swap (start_indices[i], b.start_indices[i]);
*
* @author Wolfgang Bangerth, 2000
*/
-template <int n_blocks>
inline
-void swap (BlockIndices<n_blocks> &u,
- BlockIndices<n_blocks> &v)
+void swap (BlockIndices &u, BlockIndices &v)
{
u.swap (v);
};
*/
virtual void clear ();
+ /**
+ * Return the number of blocks in a
+ * column.
+ */
+ unsigned int n_block_rows () const;
+
+ /**
+ * Return the number of blocks in a
+ * row.
+ */
+ unsigned int n_block_cols () const;
+
/**
* Return whether the object is
* empty. It is empty if either
};
+template <typename number, int rows, int columns>
+inline
+unsigned int
+BlockSparseMatrix<number,rows,columns>::n_block_cols () const
+{
+ return columns;
+}
+
+
+
+template <typename number, int rows, int columns>
+inline
+unsigned int
+BlockSparseMatrix<number,rows,columns>::n_block_rows () const
+{
+ return rows;
+}
+
#endif // __deal2__block_sparse_matrix_h
* row indices to the individual
* blocks.
*/
- const BlockIndices<rows> &
+ const BlockIndices &
get_row_indices () const;
/**
* column indices to the individual
* blocks.
*/
- const BlockIndices<columns> &
+ const BlockIndices &
get_column_indices () const;
/**
*/
void compress ();
+ /**
+ * Return the number of blocks in a
+ * column.
+ */
+ unsigned int n_block_rows () const;
+
+ /**
+ * Return the number of blocks in a
+ * row.
+ */
+ unsigned int n_block_cols () const;
+
/**
* Return whether the object is
* empty. It is empty if no
* indices to indices of the
* sub-objects.
*/
- BlockIndices<rows> row_indices;
+ BlockIndices row_indices;
/**
* Object storing and managing
* indices to indices of the
* sub-objects.
*/
- BlockIndices<columns> column_indices;
+ BlockIndices column_indices;
/**
* Make the block sparse matrix a
template <int rows, int columns>
inline
-const BlockIndices<rows> &
+const BlockIndices &
BlockSparsityPattern<rows,columns>::get_row_indices () const
{
return row_indices;
template <int rows, int columns>
inline
-const BlockIndices<columns> &
+const BlockIndices &
BlockSparsityPattern<rows,columns>::get_column_indices () const
{
return column_indices;
+template <int rows, int columns>
+inline
+unsigned int
+BlockSparsityPattern<rows,columns>::n_block_cols () const
+{
+ return columns;
+}
+
+
+
+template <int rows, int columns>
+inline
+unsigned int
+BlockSparsityPattern<rows,columns>::n_block_rows () const
+{
+ return rows;
+}
+
#endif
template <int rows, int columns>
BlockSparsityPattern<rows,columns>::BlockSparsityPattern ()
+ : row_indices (rows), column_indices(columns)
{
Assert (rows>0, ExcInvalidSize (rows));
Assert (columns>0, ExcInvalidSize (columns));
* Return a reference on the
* object that describes the
* mapping between block and
- * global indices.
+ * global indices. The use of
+ * this function is highly
+ * deprecated and it should
+ * vanish in one of the next
+ * versions
*/
- const BlockIndices<n_blocks> &
+ const BlockIndices &
get_block_indices () const;
/**
void block_read (istream &in);
//@}
- /**
- * Exception
- */
- DeclException2 (ExcDimensionsDontMatch,
- int, int,
- << "The dimensions " << arg1 << " and " << arg2
- << " do not match here.");
- /**
- * Exception
- */
- DeclException2 (ExcInvalidIndex,
- int, int,
- << "The given index " << arg1
- << " should be less than " << arg2 << ".");
- /**
- * Exception
- */
- DeclException1 (ExcInvalidNumber,
- int,
- << "The provided number is invalid here: " << arg1);
- /**
- * Exception
- */
- DeclException0 (ExcOutOfMemory);
- /**
- * Exception
- */
- DeclException0 (ExcEmptyVector);
- /**
- * Exception
- */
- DeclException0 (ExcIO);
-
protected:
/**
* Pointer to the array of components.
* indices and indices within the
* different blocks.
*/
- BlockIndices<n_blocks> block_indices;
+ BlockIndices block_indices;
};
template <int n_blocks, typename Number>
inline
-const BlockIndices<n_blocks>&
+const BlockIndices&
BlockVector<n_blocks,Number>::get_block_indices () const
{
return block_indices;
template <int n_blocks, typename Number>
BlockVector<n_blocks,Number>::BlockVector ()
+ : block_indices(n_blocks)
{}
template <int n_blocks, typename Number>
BlockVector<n_blocks,Number>::BlockVector (const vector<unsigned int> &n)
+ : block_indices(n_blocks)
{
reinit (n, false);
}
template <int n_blocks, typename Number>
BlockVector<n_blocks,Number>::BlockVector (const BlockVector<n_blocks,Number>& v)
+ : block_indices(n_blocks)
{
block_indices = v.block_indices;
for (unsigned int i=0; i<n_blocks; ++i)
void block_read (istream &in);
//@}
- /**
- * Exception
- */
- DeclException2 (ExcDimensionsDontMatch,
- int, int,
- << "The dimensions " << arg1 << " and " << arg2
- << " do not match here.");
- /**
- * Exception
- */
- DeclException1 (ExcInvalidNumber,
- int,
- << "The provided number is invalid here: " << arg1);
- /**
- * Exception
- */
- DeclException0 (ExcOutOfMemory);
/**
* Exception
*/
DeclException0 (ExcEmptyVector);
- /**
- * Exception
- */
- DeclException0 (ExcIO);
protected:
if (&v == this)
return norm_sqr();
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
Number sum0 = 0,
sum1 = 0,
Vector<Number>& Vector<Number>::operator -= (const Vector<Number>& v)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
iterator i_ptr = begin(),
i_end = end();
void Vector<Number>::add (const Vector<Number>& v)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
iterator i_ptr = begin(),
i_end = end();
void Vector<Number>::add (const Number a, const Vector<Number>& v)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
iterator i_ptr = begin(),
i_end = end();
const Number b, const Vector<Number>& w)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
- Assert (dim == w.dim, ExcDimensionsDontMatch(dim, w.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
+ Assert (dim == w.dim, ExcDimensionMismatch(dim, w.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator v_ptr = v.begin(),
void Vector<Number>::sadd (const Number x, const Vector<Number>& v)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator v_ptr = v.begin();
void Vector<Number>::sadd (const Number x, const Number a, const Vector<Number>& v)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator v_ptr = v.begin();
const Vector<Number>& v, const Number b, const Vector<Number>& w)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
- Assert (dim == w.dim, ExcDimensionsDontMatch(dim, w.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
+ Assert (dim == w.dim, ExcDimensionMismatch(dim, w.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator v_ptr = v.begin(),
const Vector<Number>& w, const Number c, const Vector<Number>& y)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
- Assert (dim == w.dim, ExcDimensionsDontMatch(dim, w.dim));
- Assert (dim == y.dim, ExcDimensionsDontMatch(dim, y.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
+ Assert (dim == w.dim, ExcDimensionMismatch(dim, w.dim));
+ Assert (dim == y.dim, ExcDimensionMismatch(dim, y.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator v_ptr = v.begin(),
const Number b, const Vector<Number>& v)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == u.dim, ExcDimensionsDontMatch(dim, u.dim));
- Assert (dim == v.dim, ExcDimensionsDontMatch(dim, v.dim));
+ Assert (dim == u.dim, ExcDimensionMismatch(dim, u.dim));
+ Assert (dim == v.dim, ExcDimensionMismatch(dim, v.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator u_ptr = u.begin(),
void Vector<Number>::equ (const Number a, const Vector<Number>& u)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (dim == u.dim, ExcDimensionsDontMatch(dim, u.dim));
+ Assert (dim == u.dim, ExcDimensionMismatch(dim, u.dim));
iterator i_ptr = begin(),
i_end = end();
const_iterator u_ptr = u.begin();
void Vector<Number>::ratio (const Vector<Number> &a, const Vector<Number> &b)
{
Assert (dim!=0, ExcEmptyVector());
- Assert (a.dim == b.dim, ExcDimensionsDontMatch (a.dim, b.dim));
+ Assert (a.dim == b.dim, ExcDimensionMismatch (a.dim, b.dim));
// no need to reinit with zeros, since
// we overwrite them anyway
# mgbase.check mg.check are removed until the sutructure of the
# multigrid classes will be fixed.
-all: full_matrix.check solver.check
+all: vector-vector.check block_vector.check block_matrices.check full_matrix.check solver.check
exe: $(all:.check=.testcase) benchmark
run: $(all:.check=.output)
@$(CXX) $(LDFLAGS) -g -o $@ $^ $(LIBS)
+############################################################
+
+block_matrices-cc-files = block_matrices.cc
+
+ifeq ($(debug-mode),on)
+block_matrices-o-files = $(block_matrices-cc-files:.cc=.go)
+else
+block_matrices-o-files = $(block_matrices-cc-files:.cc=.o)
+endif
+
+block_matrices.testcase: $(block_matrices-o-files) $(libraries)
+ @echo =====linking======= $<
+ @$(CXX) $(LDFLAGS) -g -o $@ $^ $(LIBS)
+
+
+############################################################
+
+block_vector-cc-files = block_vector.cc
+
+ifeq ($(debug-mode),on)
+block_vector-o-files = $(block_vector-cc-files:.cc=.go)
+else
+block_vector-o-files = $(block_vector-cc-files:.cc=.o)
+endif
+
+block_vector.testcase: $(block_vector-o-files) $(libraries)
+ @echo =====linking======= $<
+ @$(CXX) $(LDFLAGS) -g -o $@ $^ $(LIBS)
+
+
############################################################
--- /dev/null
+//--------------------------------------------------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 2000 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//--------------------------------------------------------------------
+
+
+#include <base/logstream.h>
+#include <lac/block_sparsity_pattern.h>
+#include <lac/block_sparse_matrix.h>
+#include <lac/block_vector.h>
+#include <fstream>
+#include <vector>
+
+void abort ()
+{}
+
+
+void test ()
+{
+ deallog.push("BlockIndices");
+
+ vector<unsigned int> ivector(4);
+ ivector[0] = 3;
+ ivector[1] = 0;
+ ivector[2] = 1;
+ ivector[3] = 2;
+
+ BlockIndices i1(ivector);
+ BlockIndices i2 = i1;
+ BlockIndices i3(13);
+ i3.reinit(ivector);
+
+ deallog.push("global->local");
+
+ unsigned int n = i1.total_size();
+ for (unsigned int i=0;i<n;++i)
+ {
+ deallog << i
+ << '\t' << i1.global_to_local(i).first
+ << '\t' << i1.global_to_local(i).second
+ << '\t' << i2.global_to_local(i).first
+ << '\t' << i2.global_to_local(i).second
+ << '\t' << i3.global_to_local(i).first
+ << '\t' << i3.global_to_local(i).second
+ << endl;
+ }
+
+ deallog.pop();
+
+ deallog.push("local->global");
+ for (unsigned int i=0 ; i<i1.size() ; ++i)
+ for (unsigned int j=0 ; j<ivector[i] ; ++j)
+ deallog << i << '\t' << j << '\t'
+ << i1.local_to_global(i,j) << endl;
+
+ deallog.pop();
+
+ deallog.push("reinit");
+
+ ivector.insert(ivector.begin(), 5);
+ i1.reinit(ivector);
+ n = i1.total_size();
+ for (unsigned int i=0;i<n;++i)
+ {
+ deallog << i
+ << '\t' << i1.global_to_local(i).first
+ << '\t' << i1.global_to_local(i).second
+ << endl;
+ }
+ deallog << "---" << endl;
+
+ ivector.erase(ivector.begin());
+ ivector.erase(ivector.begin());
+ ivector.erase(ivector.begin());
+ i1.reinit(ivector);
+ n = i1.total_size();
+ for (unsigned int i=0;i<n;++i)
+ {
+ deallog << i
+ << '\t' << i1.global_to_local(i).first
+ << '\t' << i1.global_to_local(i).second
+ << endl;
+ }
+ deallog.pop();
+
+ deallog << "Range" << endl;
+ unsigned int i = i1.global_to_local(3).first;
+ i = i1.local_to_global(1,2);
+ i = i1.local_to_global(2,0);
+}
+
+
+
+
+int main ()
+{
+ ofstream logfile("block_vector.output");
+ logfile.setf(ios::fixed);
+ logfile.precision(3);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ cerr = logfile;
+
+ try
+ {
+ test ();
+ }
+ catch (exception &e)
+ {
+ cerr << endl << endl
+ << "----------------------------------------------------"
+ << endl;
+ cerr << "Exception on processing: " << e.what() << endl
+ << "Aborting!" << endl
+ << "----------------------------------------------------"
+ << endl;
+ // abort
+ return 2;
+ }
+ catch (...)
+ {
+ cerr << endl << endl
+ << "----------------------------------------------------"
+ << endl;
+ cerr << "Unknown exception!" << endl
+ << "Aborting!" << endl
+ << "----------------------------------------------------"
+ << endl;
+ // abort
+ return 3;
+ };
+
+
+ return 0;
+};
+
--- /dev/null
+
+DEAL:BlockIndices:global->local::0 0 0 0 0 0 0
+DEAL:BlockIndices:global->local::1 0 1 0 1 0 1
+DEAL:BlockIndices:global->local::2 0 2 0 2 0 2
+DEAL:BlockIndices:global->local::3 2 0 2 0 2 0
+DEAL:BlockIndices:global->local::4 3 0 3 0 3 0
+DEAL:BlockIndices:global->local::5 3 1 3 1 3 1
+DEAL:BlockIndices:local->global::0 0 0
+DEAL:BlockIndices:local->global::0 1 1
+DEAL:BlockIndices:local->global::0 2 2
+DEAL:BlockIndices:local->global::2 0 3
+DEAL:BlockIndices:local->global::3 0 4
+DEAL:BlockIndices:local->global::3 1 5
+DEAL:BlockIndices:reinit::0 0 0
+DEAL:BlockIndices:reinit::1 0 1
+DEAL:BlockIndices:reinit::2 0 2
+DEAL:BlockIndices:reinit::3 0 3
+DEAL:BlockIndices:reinit::4 0 4
+DEAL:BlockIndices:reinit::5 1 0
+DEAL:BlockIndices:reinit::6 1 1
+DEAL:BlockIndices:reinit::7 1 2
+DEAL:BlockIndices:reinit::8 3 0
+DEAL:BlockIndices:reinit::9 4 0
+DEAL:BlockIndices:reinit::10 4 1
+DEAL:BlockIndices:reinit::---
+DEAL:BlockIndices:reinit::0 0 0
+DEAL:BlockIndices:reinit::1 1 0
+DEAL:BlockIndices:reinit::2 1 1
+DEAL:BlockIndices::Range
+--------------------------------------------------------
+An error occurred in line <174> of file </home/kanschat/deal.II/lac/include/lac/block_indices.h> in function
+ struct pair<unsigned int,unsigned int> BlockIndices::global_to_local(unsigned int) const
+The violated condition was:
+ i<total_size()
+The name and call sequence of the exception was:
+ ExcIndexRange(i, 0, total_size())
+Additional Information:
+Index 3 is not in [0,3[
+--------------------------------------------------------
+--------------------------------------------------------
+An error occurred in line <191> of file </home/kanschat/deal.II/lac/include/lac/block_indices.h> in function
+ unsigned int BlockIndices::local_to_global(unsigned int, unsigned int) const
+The violated condition was:
+ index < start_indices[block+1]-start_indices[block]
+The name and call sequence of the exception was:
+ ExcIndexRange (index, 0, start_indices[block+1]-start_indices[block])
+Additional Information:
+Index 2 is not in [0,2[
+--------------------------------------------------------
+--------------------------------------------------------
+An error occurred in line <189> of file </home/kanschat/deal.II/lac/include/lac/block_indices.h> in function
+ unsigned int BlockIndices::local_to_global(unsigned int, unsigned int) const
+The violated condition was:
+ block < n_blocks
+The name and call sequence of the exception was:
+ ExcIndexRange(block, 0, n_blocks)
+Additional Information:
+Index 2 is not in [0,2[
+--------------------------------------------------------