]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Add an assignment operator from a regular SparseMatrix.
authoryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 31 Mar 2014 07:14:35 +0000 (07:14 +0000)
committeryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Mon, 31 Mar 2014 07:14:35 +0000 (07:14 +0000)
git-svn-id: https://svn.dealii.org/trunk@32706 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/lac/lapack_full_matrix.h
deal.II/source/lac/lapack_full_matrix.cc
deal.II/source/lac/lapack_full_matrix.inst.in

index 92e0dd09e8b5227c0364886f232ffb89e458289c..ee8339ba60147e065532d2a3a82ecdbd347c35e7 100644 (file)
@@ -34,6 +34,7 @@ DEAL_II_NAMESPACE_OPEN
 template<typename number> class Vector;
 template<typename number> class BlockVector;
 template<typename number> class FullMatrix;
+template<typename number> class SparseMatrix;
 
 
 /**
@@ -70,6 +71,7 @@ public:
    */
   explicit LAPACKFullMatrix (const size_type n = 0);
 
+
   /**
    * Constructor. Initialize the matrix as a rectangular matrix.
    */
@@ -95,14 +97,23 @@ public:
   operator = (const LAPACKFullMatrix<number> &);
 
   /**
-   * Assignment operator for a regular FullMatrix. Note that since LAPACK
-   * expects matrices in transposed order, this transposition is included
-   * here.
+   * Assignment operator from a regular FullMatrix. @note Since LAPACK
+   * expects matrices in transposed order, this transposition is
+   * included here.
    */
   template <typename number2>
   LAPACKFullMatrix<number> &
   operator = (const FullMatrix<number2> &);
 
+  /**
+   * Assignment operator from a regular SparseMatrix. @note Since
+   * LAPACK expects matrices in transposed order, this transposition
+   * is included here.
+   */
+  template <typename number2>
+  LAPACKFullMatrix<number> &
+  operator = (const SparseMatrix<number2> &);
+
   /**
    * This operator assigns a scalar to a matrix. To avoid confusion with
    * constructors, zero is the only value allowed for <tt>d</tt>
@@ -544,6 +555,10 @@ public:
                         const double        threshold   = 0.) const;
 
 private:
+  /**
+   * n_rows
+   */
+
   /**
    * Since LAPACK operations notoriously change the meaning of the matrix
    * entries, we record the current state after the last operation here.
index 2cf6cbee9a094895cc68ed801ce407e9ee37fbf0..98aa0a1d52dcc967b83f113c6126f07b3b56eefb 100644 (file)
@@ -19,6 +19,7 @@
 #include <deal.II/lac/lapack_templates.h>
 #include <deal.II/lac/lapack_support.h>
 #include <deal.II/lac/full_matrix.h>
+#include <deal.II/lac/sparse_matrix.h>
 #include <deal.II/lac/vector.h>
 #include <deal.II/lac/block_vector.h>
 
@@ -30,34 +31,30 @@ DEAL_II_NAMESPACE_OPEN
 using namespace LAPACKSupport;
 
 template <typename number>
-LAPACKFullMatrix<number>::LAPACKFullMatrix(const size_type n)
+LAPACKFullMatrix<number>::LAPACKFullMatrix (const size_type n)
   :
   TransposeTable<number> (n,n),
-  state(matrix)
+  state (matrix)
 {}
 
 
-
 template <typename number>
-LAPACKFullMatrix<number>::LAPACKFullMatrix(
-  const size_type m,
-  const size_type n)
+LAPACKFullMatrix<number>::LAPACKFullMatrix (const size_type m,
+                                           const size_type n)
   :
   TransposeTable<number> (m,n),
-  state(matrix)
+  state (matrix)
 {}
 
 
-
 template <typename number>
-LAPACKFullMatrix<number>::LAPACKFullMatrix(const LAPACKFullMatrix &M)
+LAPACKFullMatrix<number>::LAPACKFullMatrix (const LAPACKFullMatrix &M)
   :
   TransposeTable<number> (M),
-  state(matrix)
+  state (matrix)
 {}
 
 
-
 template <typename number>
 LAPACKFullMatrix<number> &
 LAPACKFullMatrix<number>::operator = (const LAPACKFullMatrix<number> &M)
@@ -68,7 +65,6 @@ LAPACKFullMatrix<number>::operator = (const LAPACKFullMatrix<number> &M)
 }
 
 
-
 template <typename number>
 template <typename number2>
 LAPACKFullMatrix<number> &
@@ -85,6 +81,21 @@ LAPACKFullMatrix<number>::operator = (const FullMatrix<number2> &M)
 }
 
 
+template <typename number>
+template <typename number2>
+LAPACKFullMatrix<number> &
+LAPACKFullMatrix<number>::operator = (const SparseMatrix<number2> &M)
+{
+  Assert (this->n_rows() == M.n(), ExcDimensionMismatch(this->n_rows(), M.n()));
+  Assert (this->n_cols() == M.m(), ExcDimensionMismatch(this->n_cols(), M.m()));
+  for (size_type i=0; i<this->n_rows(); ++i)
+    for (size_type j=0; j<this->n_cols(); ++j)
+      (*this)(i,j) = M.el(i,j);
+
+  state = LAPACKSupport::matrix;
+  return *this;
+}
+
 
 template <typename number>
 LAPACKFullMatrix<number> &
@@ -100,7 +111,6 @@ LAPACKFullMatrix<number>::operator = (const double d)
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::vmult (
@@ -159,7 +169,6 @@ LAPACKFullMatrix<number>::vmult (
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::Tvmult (
@@ -220,7 +229,6 @@ LAPACKFullMatrix<number>::Tvmult (
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::mmult(LAPACKFullMatrix<number>       &C,
@@ -244,7 +252,6 @@ LAPACKFullMatrix<number>::mmult(LAPACKFullMatrix<number>       &C,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::mmult(FullMatrix<number>             &C,
@@ -293,7 +300,6 @@ LAPACKFullMatrix<number>::Tmmult(LAPACKFullMatrix<number>       &C,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::Tmmult(FullMatrix<number>             &C,
@@ -318,7 +324,6 @@ LAPACKFullMatrix<number>::Tmmult(FullMatrix<number>             &C,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::mTmult(LAPACKFullMatrix<number>       &C,
@@ -367,7 +372,6 @@ LAPACKFullMatrix<number>::mTmult(FullMatrix<number>             &C,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::TmTmult(LAPACKFullMatrix<number>       &C,
@@ -391,7 +395,6 @@ LAPACKFullMatrix<number>::TmTmult(LAPACKFullMatrix<number>       &C,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::TmTmult(FullMatrix<number>             &C,
@@ -416,7 +419,6 @@ LAPACKFullMatrix<number>::TmTmult(FullMatrix<number>             &C,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::compute_lu_factorization()
@@ -436,7 +438,6 @@ LAPACKFullMatrix<number>::compute_lu_factorization()
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::compute_svd()
@@ -487,7 +488,6 @@ LAPACKFullMatrix<number>::compute_svd()
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::compute_inverse_svd(const double threshold)
@@ -509,7 +509,6 @@ LAPACKFullMatrix<number>::compute_inverse_svd(const double threshold)
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::invert()
@@ -542,7 +541,6 @@ LAPACKFullMatrix<number>::invert()
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::apply_lu_factorization(Vector<number> &v,
@@ -565,7 +563,6 @@ LAPACKFullMatrix<number>::apply_lu_factorization(Vector<number> &v,
 }
 
 
-
 template <typename number>
 void
 LAPACKFullMatrix<number>::apply_lu_factorization(LAPACKFullMatrix<number> &B,
@@ -588,12 +585,10 @@ LAPACKFullMatrix<number>::apply_lu_factorization(LAPACKFullMatrix<number> &B,
 }
 
 
-
 template <typename number>
 void
-LAPACKFullMatrix<number>::compute_eigenvalues(
-  const bool right,
-  const bool left)
+LAPACKFullMatrix<number>::compute_eigenvalues(const bool right,
+                                             const bool left)
 {
   Assert(state == matrix, ExcState(state));
   const int nn = this->n_cols();
@@ -662,12 +657,11 @@ LAPACKFullMatrix<number>::compute_eigenvalues(
 
 template <typename number>
 void
-LAPACKFullMatrix<number>::compute_eigenvalues_symmetric(
-  const number lower_bound,
-  const number upper_bound,
-  const number abs_accuracy,
-  Vector<number> &eigenvalues,
-  FullMatrix<number> &eigenvectors)
+LAPACKFullMatrix<number>::compute_eigenvalues_symmetric(const number        lower_bound,
+                                                       const number        upper_bound,
+                                                       const number        abs_accuracy,
+                                                       Vector<number>     &eigenvalues,
+                                                       FullMatrix<number> &eigenvectors)
 {
   Assert(state == matrix, ExcState(state));
   const int nn = (this->n_cols() > 0 ? this->n_cols() : 1);
@@ -968,6 +962,7 @@ LAPACKFullMatrix<number>::Tvmult_add (
   Tvmult(w, v, true);
 }
 
+
 template <typename number>
 void
 LAPACKFullMatrix<number>::print_formatted (
index 9db46db4a07de69bf2bf3c6cb005ad3203b1ef4d..23d91651bd5b666e88a3705b8ef058e00c389acd 100644 (file)
@@ -26,4 +26,7 @@ for (S1, S2 : REAL_SCALARS)
   {
     template LAPACKFullMatrix<S1> &
     LAPACKFullMatrix<S1>::operator = (const FullMatrix<S2> &M);
+
+    template LAPACKFullMatrix<S1> &
+    LAPACKFullMatrix<S1>::operator = (const SparseMatrix<S2> &M);
   }

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.