inline
SparseMatrixEZ<Number> &
BlockSparseMatrixEZ<Number>::block (const unsigned int row,
- const unsigned int column)
+ const unsigned int column)
{
- Assert (row<n_rows(), ExcIndexRange (row, 0, n_rows()));
- Assert (column<n_cols(), ExcIndexRange (column, 0, n_cols()));
+ Assert (row<n_block_rows(), ExcIndexRange (row, 0, n_block_rows()));
+ Assert (column<n_block_cols(), ExcIndexRange (column, 0, n_block_cols()));
return blocks[row][column];
};
inline
const SparseMatrixEZ<Number> &
BlockSparseMatrixEZ<Number>::block (const unsigned int row,
- const unsigned int column) const
+ const unsigned int column) const
{
- Assert (row<n_rows(), ExcIndexRange (row, 0, n_rows()));
- Assert (column<n_cols(), ExcIndexRange (column, 0, n_cols()));
+ Assert (row<n_block_rows(), ExcIndexRange (row, 0, n_block_rows()));
+ Assert (column<n_block_cols(), ExcIndexRange (column, 0, n_block_cols()));
return blocks[row][column];
};
*/
DeclException2 (ExcInvalidEntry,
int, int,
- << "The entry with index <" << arg1 << ',' << arg2
- << "> does not exist.");
+ << "The entry with index (" << arg1 << ',' << arg2
+ << ") does not exist.");
+ DeclException2(ExcEntryAllocationFailure,
+ int, int,
+ << "An entry with index (" << arg1 << ',' << arg2
+ << ") cannot be allocated.");
private:
/**
* Find an entry. Return a
{
if (end >= row_info[row+1].start)
{
+ // Failure if increment 0
+ Assert(increment!=0,ExcEntryAllocationFailure(row,col));
+
// Insert new entries
data.insert(data.begin()+end, increment, Entry());
// Update starts of
}
-template <typename number>
-template <class STREAM>
-void
-SparseMatrixEZ<number>::print_statistics(STREAM& out, bool full)
-{
- typename std::vector<RowInfo>::const_iterator row = row_info.begin();
- const typename std::vector<RowInfo>::const_iterator endrow = row_info.end();
-
- // Add up entries actually used
- unsigned int entries_used = 0;
- unsigned int max_length = 0;
- for (; row != endrow ; ++ row)
- {
- entries_used += row->length;
- if (max_length < row->length)
- max_length = row->length;
- }
-
- // Number of entries allocated is
- // position of last entry used
- --row;
- unsigned int entries_alloc = row->start + row->length;
-
- out << "SparseMatrixEZ:used entries:" << entries_used << std::endl
- << "SparseMatrixEZ:alloc entries:" << entries_alloc << std::endl
- << "SparseMatrixEZ:reserved entries:" << data.capacity() << std::endl;
-
- if (full)
- {
- std::vector<unsigned int> length_used (max_length+1);
-
- for (row = row_info.begin() ; row != endrow; ++row)
- {
- ++length_used[row->length];
- }
- for (unsigned int i=0; i< length_used.size();++i)
- if (length_used[i] != 0)
- out << "SparseMatrixEZ:entries\t" << i
- << "\trows\t" << length_used[i]
- << std::endl;
- }
-}
-
#endif
/*---------------------------- sparse_matrix.h ---------------------------*/
template <typename number>
SparseMatrixEZ<number>&
-SparseMatrixEZ<number>::operator= (const SparseMatrixEZ<number>&)
+SparseMatrixEZ<number>::operator= (const SparseMatrixEZ<number>& m)
{
- Assert (false, ExcNotImplemented());
+ Assert (m.empty(), ExcInvalidConstructorCall());
return *this;
}
default_row_length = 5;
if (default_increment == Entry::invalid)
default_increment = 4;
- if (default_increment == 0)
- default_increment = 4;
increment = default_increment;
n_columns = n_cols;
}
+template <typename number>
+template <class STREAM>
+void
+SparseMatrixEZ<number>::print_statistics(STREAM& out, bool full)
+{
+ typename std::vector<RowInfo>::const_iterator row = row_info.begin();
+ const typename std::vector<RowInfo>::const_iterator endrow = row_info.end();
+
+ // Add up entries actually used
+ unsigned int entries_used = 0;
+ unsigned int max_length = 0;
+ for (; row != endrow ; ++ row)
+ {
+ entries_used += row->length;
+ if (max_length < row->length)
+ max_length = row->length;
+ }
+
+ // Number of entries allocated is
+ // position of last entry used
+ --row;
+ unsigned int entries_alloc = row->start + row->length;
+
+ out << "SparseMatrixEZ:used entries:" << entries_used << std::endl
+ << "SparseMatrixEZ:alloc entries:" << entries_alloc << std::endl
+ << "SparseMatrixEZ:reserved entries:" << data.capacity() << std::endl;
+
+ if (full)
+ {
+ std::vector<unsigned int> length_used (max_length+1);
+
+ for (row = row_info.begin() ; row != endrow; ++row)
+ {
+ ++length_used[row->length];
+ }
+ for (unsigned int i=0; i< length_used.size();++i)
+ if (length_used[i] != 0)
+ out << "SparseMatrixEZ:entries\t" << i
+ << "\trows\t" << length_used[i]
+ << std::endl;
+ }
+}
+
+
#endif // __deal2__sparse_matrix_ez_templates_h
-//---------------------------- sparse_matrix_ez.double.cc ---------------------------
+//------------------------------------------------------------------------
// $Id$
// Version: $Name$
//
// to the file deal.II/doc/license.html for the text and
// further information on this license.
//
-//---------------------------- sparse_matrix_ez.double.cc ---------------------------
+//------------------------------------------------------------------------
+#include <base/logstream.h>
#include <lac/sparse_matrix_ez.templates.h>
+#include <iostream>
#define TYPEMAT double
template class SparseMatrixEZ<TYPEMAT>;
+template void SparseMatrixEZ<TYPEMAT>::print_statistics(ostream&, bool);
+template void SparseMatrixEZ<TYPEMAT>::print_statistics(LogStream&, bool);
+
#define TYPE2 float
-//---------------------------- sparse_matrix_ez.float.cc ---------------------------
+//-------------------------------------------------------------------------
// $Id$
// Version: $Name$
//
// to the file deal.II/doc/license.html for the text and
// further information on this license.
//
-//---------------------------- sparse_matrix_ez.float.cc ---------------------------
-
+//-------------------------------------------------------------------------
+#include <base/logstream.h>
#include <lac/sparse_matrix_ez.templates.h>
+#include <iostream>
#define TYPEMAT float
template class SparseMatrixEZ<TYPEMAT>;
+template void SparseMatrixEZ<TYPEMAT>::print_statistics(ostream&, bool);
+template void SparseMatrixEZ<TYPEMAT>::print_statistics(LogStream&, bool);
#define TYPE2 float