#include <base/config.h>
#include <base/table.h>
#include <lac/exceptions.h>
+#include <lac/identity_matrix.h>
#include <vector>
const unsigned int cols,
const number* entries);
+ /**
+ * Construct a full matrix that
+ * equals the identity matrix of
+ * the size of the
+ * argument. Using this
+ * constructor, one can easily
+ * create an identity matrix of
+ * size <code>n</code> by saying
+ * @verbatim
+ * FullMatrix<double> M(IdentityMatrix(n));
+ * @endverbatim
+ */
+ explicit FullMatrix (const IdentityMatrix &id);
+
/**
* Assignment operator.
*/
FullMatrix<number> &
operator = (const double d);
+ /**
+ * Copy operator to create a full
+ * matrix that equals the
+ * identity matrix of the size of
+ * the argument. This way, one can easily
+ * create an identity matrix of
+ * size <code>n</code> by saying
+ * @verbatim
+ * M = IdentityMatrix(n);
+ * @endverbatim
+ */
+ FullMatrix<number> &
+ operator = (const IdentityMatrix &id);
+
/**
* Assignment from different
* matrix classes. This
{}
+
+template <typename number>
+FullMatrix<number>::FullMatrix (const IdentityMatrix &id)
+ :
+ Table<2,number> (id.m(), id.n())
+{
+ for (unsigned int i=0; i<id.m(); ++i)
+ (*this)(i,i) = 1;
+}
+
+
+
+
template <typename number>
FullMatrix<number>&
FullMatrix<number>::operator = (const FullMatrix<number>& M)
}
+
+template <typename number>
+FullMatrix<number>&
+FullMatrix<number>::operator = (const IdentityMatrix &id)
+{
+ reinit (id.m(), id.n());
+ for (unsigned int i=0; i<id.m(); ++i)
+ (*this)(i,i) = 1.;
+
+ return *this;
+}
+
+
+
template <typename number>
bool
FullMatrix<number>::all_zero () const
@endverbatim
*
* This creates a $10\times 10$ matrix with ones on the diagonal and
- * zeros everywhere else. All matrix types have conversion
- * constructors and can therefore be filled rather easily with
- * identity matrices.
+ * zeros everywhere else. Most matrix types, in particular FullMatrix
+ * and SparseMatrix, have conversion constructors and assignment
+ * operators for IdentityMatrix, and can therefore be filled rather
+ * easily with identity matrices.
*
*
* <h4>Preconditioning</h4>
#include <base/subscriptor.h>
#include <base/smartpointer.h>
#include <lac/sparsity_pattern.h>
+#include <lac/identity_matrix.h>
#include <lac/exceptions.h>
template<typename number> class Vector;
* generated then.
*/
explicit SparseMatrix (const SparsityPattern &sparsity);
+
+ /**
+ * Copy constructor: initialize
+ * the matrix with the identity
+ * matrix. This constructor will
+ * throw an exception if the
+ * sizes of the sparsity pattern
+ * and the identity matrix do not
+ * coincide, or if the sparsity
+ * pattern does not provide for
+ * nonzero entries on the entire
+ * diagonal.
+ */
+ SparseMatrix (const SparsityPattern &sparsity,
+ const IdentityMatrix &id);
/**
* Destructor. Free all memory, but do not
*/
SparseMatrix<number>& operator = (const SparseMatrix<number> &);
+ /**
+ * Copy operator: initialize
+ * the matrix with the identity
+ * matrix. This operator will
+ * throw an exception if the
+ * sizes of the sparsity pattern
+ * and the identity matrix do not
+ * coincide, or if the sparsity
+ * pattern does not provide for
+ * nonzero entries on the entire
+ * diagonal.
+ */
+ SparseMatrix<number> &
+ operator= (const IdentityMatrix &id);
+
/**
* This operator assigns a scalar to
* a matrix. Since this does usually
+template <typename number>
+SparseMatrix<number>::SparseMatrix (const SparsityPattern &c,
+ const IdentityMatrix &id)
+ :
+ cols(0),
+ val(0),
+ max_len(0)
+{
+ Assert (c.n_rows() == id.m(), ExcDimensionMismatch (c.n_rows(), id.m()));
+ Assert (c.n_cols() == id.n(), ExcDimensionMismatch (c.n_cols(), id.n()));
+
+ reinit (c);
+ for (unsigned int i=0; i<n(); ++i)
+ this->set(i,i,1.);
+}
+
+
+
template <typename number>
SparseMatrix<number>::~SparseMatrix ()
{
+template <typename number>
+SparseMatrix<number> &
+SparseMatrix<number>::operator= (const IdentityMatrix &id)
+{
+ Assert (cols->n_rows() == id.m(),
+ ExcDimensionMismatch (cols->n_rows(), id.m()));
+ Assert (cols->n_cols() == id.n(),
+ ExcDimensionMismatch (cols->n_cols(), id.n()));
+
+ *this = 0;
+ for (unsigned int i=0; i<n(); ++i)
+ this->set(i,i,1.);
+
+ return *this;
+}
+
+
+
template <typename number>
void
SparseMatrix<number>::reinit (const SparsityPattern &sparsity)