28#include "gsl/gsl_linalg.h"
29#include "gsl/gsl_eigen.h"
72 if(
_matrix !=
nullptr) gsl_matrix_complex_free( (gsl_matrix_complex*)
_matrix );
80 if (&mm ==
this)
return *
this;
84 gsl_matrix_memcpy((gsl_matrix*)
_matrix, (
const gsl_matrix*)mm.
_matrix);
91 if (&mm ==
this)
return *
this;
95 gsl_matrix_complex_memcpy((gsl_matrix_complex*)
_matrix, (
const gsl_matrix_complex*)mm.
_matrix);
104 _matrix = gsl_matrix_alloc(mm.num_row(), mm.num_col() );
105 gsl_matrix_memcpy( (gsl_matrix*)_matrix, (gsl_matrix*)mm._matrix );
114 _matrix = gsl_matrix_complex_alloc( mm.num_row(), mm.num_col() );
115 gsl_matrix_complex_memcpy( (gsl_matrix_complex*)_matrix, (gsl_matrix_complex*)mm._matrix );
119template <
class Tmplt>
125template <
class Tmplt>
129 for(
size_t a=1; a<=i; a++)
130 for(
size_t b=1; b<=j; b++)
131 operator()(a,b) = value;
134template <
class Tmplt>
140template <
class Tmplt>
144 for(
size_t a=1; a<=i; a++)
146 for(
size_t b=1; b< a; b++) mm(a, b) = off_diag_value;
147 for(
size_t b=a+1; b<=i; b++) mm(a, b) = off_diag_value;
148 mm(a,a) = diag_value;
154template <
class Tmplt>
165 _matrix = gsl_matrix_alloc(i, j);
171 _matrix = gsl_matrix_complex_alloc(i, j);
177 _matrix = gsl_matrix_alloc(i, j);
178 for(
size_t a=0; a<i; a++)
179 for(
size_t b=0; b<j; b++)
180 operator()(a+1,b+1) = data[b*i + a];
186 _matrix = gsl_matrix_complex_alloc(i, j);
187 for(
size_t a=0; a<i; a++)
188 for(
size_t b=0; b<j; b++)
189 operator()(a+1,b+1) = data[b*i + a];
201 gsl_permutation * p = gsl_permutation_alloc (
num_row());
203 gsl_linalg_complex_LU_decomp ((gsl_matrix_complex*)copy.
_matrix, p, &signum);
204 gsl_permutation_free(p);
205 return gsl_linalg_complex_LU_det((gsl_matrix_complex*)copy.
_matrix, signum);
213 gsl_permutation * p = gsl_permutation_alloc (
num_row());
215 gsl_linalg_LU_decomp ((gsl_matrix*)copy.
_matrix, p, &signum);
216 double d = gsl_linalg_LU_det((gsl_matrix*)copy.
_matrix, signum);
217 gsl_permutation_free(p);
221template <
class Tmplt>
234 gsl_permutation* perm = gsl_permutation_alloc(
num_row());
237 gsl_linalg_complex_LU_decomp ( (gsl_matrix_complex*)lu.
_matrix, perm, &s);
238 gsl_linalg_complex_LU_invert ( (gsl_matrix_complex*)lu.
_matrix, perm, (gsl_matrix_complex*)
_matrix);
239 gsl_permutation_free( perm );
240 for(
size_t i=1; i<=
num_row(); i++)
241 for(
size_t j=1; j<=
num_row(); j++)
242 if(
operator()(i,j) !=
operator()(i,j))
throw(
GeneralClassicException(
"MMatrix::invert()",
"Failed to invert matrix - singular?"));
249 gsl_permutation* perm = gsl_permutation_alloc(
num_row());
252 gsl_linalg_LU_decomp ( (gsl_matrix*)lu.
_matrix, perm, &s);
253 gsl_linalg_LU_invert ( (gsl_matrix*)lu.
_matrix, perm, (gsl_matrix*)
_matrix);
254 gsl_permutation_free( perm );
255 for(
size_t i=1; i<=
num_row(); i++)
256 for(
size_t j=1; j<=
num_row(); j++)
257 if(
operator()(i,j) !=
operator()(i,j))
throw(
GeneralClassicException(
"MMatrix::invert()",
"Failed to invert matrix - singular?"));
264 gsl_matrix_transpose_memcpy ((gsl_matrix*)in.
_matrix, (gsl_matrix*)
_matrix);
272 gsl_matrix_complex_transpose_memcpy ((gsl_matrix_complex*)in.
_matrix, (gsl_matrix_complex*)
_matrix);
276template <
class Tmplt>
280 for(
size_t i=min_row; i<=max_row; i++)
281 for(
size_t j=min_col; j<=max_col; j++)
282 sub_matrix(i-min_row+1, j-min_col+1) =
operator()(i,j);
289template <
class Tmplt>
293 for(
size_t i=2; i<=
num_row() && i<=
num_col(); i++) t+=
operator()(i,i);
300template <
class Tmplt>
306 gsl_eigen_nonsymm_workspace * w = gsl_eigen_nonsymm_alloc(
num_row());
307 gsl_eigen_nonsymm_params(0, 1, w);
309 gsl_eigen_nonsymm_free(w);
315template <
class Tmplt>
322 gsl_eigen_nonsymmv_workspace * w = gsl_eigen_nonsymmv_alloc(
num_row());
324 gsl_eigen_nonsymmv_free(w);
325 if(ierr != 0)
throw(
GeneralClassicException(
"MMatrix<Tmplt>::eigenvectors",
"Failed to calculate eigenvalue"));
335 gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1., (gsl_matrix*)m1.
_matrix, (gsl_matrix*)m2.
_matrix, 0., (gsl_matrix*)out.
_matrix);
352 for(
size_t i=1; i<=mat.
num_row(); i++)
355 for(
size_t j=1; j<mat.
num_col(); j++)
356 out << mat(i,j) <<
" ";
357 out << mat(i,mat.
num_col()) <<
"\n";
369 for(
size_t i=1; i<=
nr; i++)
370 for(
size_t j=1; j<=nc; j++)
381 for(
size_t i=1; i<=mc.
num_row(); i++)
382 for(
size_t j=1; j<=mc.
num_col(); j++)
383 md(i,j) =
re(mc(i,j));
390 for(
size_t i=1; i<=mc.
num_row(); i++)
391 for(
size_t j=1; j<=mc.
num_col(); j++)
392 md(i,j) =
im(mc(i,j));
399 for(
size_t i=1; i<=mc.
num_row(); i++)
400 for(
size_t j=1; j<=mc.
num_col(); j++)
408 throw(
GeneralClassicException(
"MMatrix<m_complex>::complex",
"Attempt to build complex matrix when real and imaginary matrix don't have the same size"));
410 for(
size_t i=1; i<=mc.
num_row(); i++)
411 for(
size_t j=1; j<=mc.
num_col(); j++)
416template <
class Tmplt>
m_complex m_complex_build(double r, double i)
MMatrix< m_complex > complex(MMatrix< double > real)
template std::istream & operator>>(std::istream &in, MMatrix< double > &mat)
MMatrix< double > re(MMatrix< m_complex > mc)
std::ostream & operator<<(std::ostream &out, const Mesh::Iterator &it)
MMatrix< double > & operator*=(MMatrix< double > &m1, MMatrix< double > m2)
MMatrix< double > im(MMatrix< m_complex > mc)
void invert()
turns this matrix into its inverse
static gsl_matrix * get_matrix(const MMatrix< double > &m)
MVector< Tmplt > get_mvector(size_t column) const
size_t num_col() const
returns number of columns in the matrix
MVector< m_complex > eigenvalues() const
returns a vector of eigenvalues. Throws an exception if either this matrix is not square or the eigen...
void build_matrix(size_t i, size_t j)
MMatrix(const MMatrix< double > &mm)
MMatrix< Tmplt > & operator=(const MMatrix< Tmplt > &mm)
MMatrix< Tmplt > T() const
returns matrix transpose T (such that M(i,j) = T(j,i))
size_t num_row() const
returns number of rows in the matrix
std::pair< MVector< m_complex >, MMatrix< m_complex > > eigenvectors() const
MMatrix< Tmplt > inverse() const
returns matrix inverse leaving this matrix unchanged
MMatrix< Tmplt > sub(size_t min_row, size_t max_row, size_t min_col, size_t max_col) const
static MMatrix< Tmplt > Diagonal(size_t i, Tmplt diag_value, Tmplt off_diag_value)
Construct a square matrix filling on and off diagonal values.
Tmplt determinant() const
returns matrix determinant; throws an exception if matrix is not square
MMatrix()
default constructor makes an empty MMatrix of size (0,0)
const Tmplt & operator()(size_t row, size_t column) const
Tmplt trace() const
returns sum of terms with row == column, even if matrix is not square
static gsl_vector * get_vector(const MVector< double > &m)