]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Copy operations from identity matrix
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 24 Apr 2006 23:06:59 +0000 (23:06 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 24 Apr 2006 23:06:59 +0000 (23:06 +0000)
git-svn-id: https://svn.dealii.org/trunk@12878 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/lac/include/lac/full_matrix.h
deal.II/lac/include/lac/full_matrix.templates.h
deal.II/lac/include/lac/identity_matrix.h
deal.II/lac/include/lac/sparse_matrix.h
deal.II/lac/include/lac/sparse_matrix.templates.h

index a5ea36872aeada24c55fa979fa84999d3ccde6ea..bbb148ed4381f18082a493017122517faa886c71 100644 (file)
@@ -17,6 +17,7 @@
 #include <base/config.h>
 #include <base/table.h>
 #include <lac/exceptions.h>
+#include <lac/identity_matrix.h>
 
 #include <vector>
 
@@ -245,6 +246,20 @@ class FullMatrix : public Table<2,number>
                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.
                                      */
@@ -271,6 +286,20 @@ class FullMatrix : public Table<2,number>
     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
index ae08654650121444e95a86058d9b28d1eb50cccd..c7b0e60b2add93dbee6d77c8b7f404fb87c3c603 100644 (file)
@@ -57,6 +57,19 @@ FullMatrix<number>::FullMatrix (const FullMatrix &m)
 {}
 
 
+
+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)
@@ -76,6 +89,20 @@ FullMatrix<number>::operator = (const FullMatrix<number2>& 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
index a244633cd467c3990fd48221e13993a75d67c456..e41cc413f6a719c6f764cce5cf53ae34c2f32f1e 100644 (file)
  @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>
index 4c375109817c80e1574284f67a799e87bbf050d3..1508743b22b4055fa2e0efa978d6b1aec932c039 100644 (file)
@@ -18,6 +18,7 @@
 #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;
@@ -564,6 +565,21 @@ class SparseMatrix : public virtual Subscriptor
                                      * 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
@@ -578,6 +594,21 @@ class SparseMatrix : public virtual Subscriptor
                                      */
     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
index 1bf79f2d9e5685575ef9d39bfc3baf838a46669e..dfa36e645176e9f6f2f7368aa567335a8fd70897 100644 (file)
@@ -93,6 +93,24 @@ SparseMatrix<number>::SparseMatrix (const SparsityPattern &c)
 
 
 
+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 ()
 {
@@ -121,6 +139,24 @@ SparseMatrix<number>::operator = (const double d)
 
 
 
+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)

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.