template<typename number> class Vector;
template<typename number> class BlockVector;
template<typename number> class FullMatrix;
+template<typename number> class SparseMatrix;
/**
*/
explicit LAPACKFullMatrix (const size_type n = 0);
+
/**
* Constructor. Initialize the matrix as a rectangular matrix.
*/
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>
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.
#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>
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)
}
-
template <typename number>
template <typename number2>
LAPACKFullMatrix<number> &
}
+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> &
}
-
template <typename number>
void
LAPACKFullMatrix<number>::vmult (
}
-
template <typename number>
void
LAPACKFullMatrix<number>::Tvmult (
}
-
template <typename number>
void
LAPACKFullMatrix<number>::mmult(LAPACKFullMatrix<number> &C,
}
-
template <typename number>
void
LAPACKFullMatrix<number>::mmult(FullMatrix<number> &C,
}
-
template <typename number>
void
LAPACKFullMatrix<number>::Tmmult(FullMatrix<number> &C,
}
-
template <typename number>
void
LAPACKFullMatrix<number>::mTmult(LAPACKFullMatrix<number> &C,
}
-
template <typename number>
void
LAPACKFullMatrix<number>::TmTmult(LAPACKFullMatrix<number> &C,
}
-
template <typename number>
void
LAPACKFullMatrix<number>::TmTmult(FullMatrix<number> &C,
}
-
template <typename number>
void
LAPACKFullMatrix<number>::compute_lu_factorization()
}
-
template <typename number>
void
LAPACKFullMatrix<number>::compute_svd()
}
-
template <typename number>
void
LAPACKFullMatrix<number>::compute_inverse_svd(const double threshold)
}
-
template <typename number>
void
LAPACKFullMatrix<number>::invert()
}
-
template <typename number>
void
LAPACKFullMatrix<number>::apply_lu_factorization(Vector<number> &v,
}
-
template <typename number>
void
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();
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);
Tvmult(w, v, true);
}
+
template <typename number>
void
LAPACKFullMatrix<number>::print_formatted (