OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
MMatrix.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2015, Chris Rogers
3  * All rights reserved.
4  * Redistribution and use in source and binary forms, with or without
5  * modification, are permitted provided that the following conditions are met:
6  * 1. Redistributions of source code must retain the above copyright notice,
7  * this list of conditions and the following disclaimer.
8  * 2. Redistributions in binary form must reproduce the above copyright notice,
9  * this list of conditions and the following disclaimer in the documentation
10  * and/or other materials provided with the distribution.
11  * 3. Neither the name of STFC nor the names of its contributors may be used to
12  * endorse or promote products derived from this software without specific
13  * prior written permission.
14  *
15  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
16  * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
17  * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
18  * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
19  * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
20  * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
21  * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
22  * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
23  * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
24  * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
25  * POSSIBILITY OF SUCH DAMAGE.
26  */
27 
28 #ifndef _SRC_COMMON_CPP_MATHS_MMATRIX_HH_
29 #define _SRC_COMMON_CPP_MATHS_MMATRIX_HH_
30 
31 #include <iostream>
32 #include <vector>
33 
34 #include "gsl/gsl_matrix.h"
35 #include "gsl/gsl_blas.h"
36 
38 
40 
41 namespace interpolation {
42 
44 
61 template <class Tmplt>
62 class MMatrix
63 {
64 public:
67  MMatrix();
68 
71  MMatrix(const MMatrix<Tmplt>& mv );
72 
81  MMatrix(size_t nrows, size_t ncols, Tmplt* data_beginning );
82 
89  MMatrix(size_t nrows, size_t ncols, Tmplt value);
90 
96  MMatrix(size_t nrows, size_t ncols);
97 
106  static MMatrix<Tmplt> Diagonal(size_t i, Tmplt diag_value, Tmplt off_diag_value);
107 
110  ~MMatrix();
111 
114  size_t num_row() const;
115 
118  size_t num_col() const;
119 
123  Tmplt trace() const;
124 
128  Tmplt determinant() const;
129 
132  MMatrix<Tmplt> inverse() const;
133 
136  void invert();
137 
140  MMatrix<Tmplt> T() const;
141 
146  MVector<m_complex> eigenvalues() const; //return vector of eigenvalues
147  std::pair<MVector<m_complex>, MMatrix<m_complex> >
148  eigenvectors() const; //return pair of eigenvalues, eigenvectors (access using my_pair.first or my_pair.second)
149 
150  //return a submatrix, with data from min_row to max_row and min_row to max_row inclusive; indexing starts at 1
151  MMatrix<Tmplt> sub(size_t min_row, size_t max_row, size_t min_col, size_t max_col) const;
152  //create a column vector from the column^th column
153  MVector<Tmplt> get_mvector(size_t column) const;
154  //return a reference to the datum; indexing starts at 1, goes to num_row (inclusive)
155  const Tmplt& operator()(size_t row, size_t column) const;
156  Tmplt& operator()(size_t row, size_t column);
157 
159 
160  //TODO - implement iterator
161  //class iterator
162  //{
163  //}
164 
165 
174  template <class Tmplt2> friend MMatrix<Tmplt2> operator + (MMatrix<Tmplt2> m1, const MMatrix<Tmplt2> m2);
175 
176  friend class MMatrix<double>; //To do the eigenvector problem, MMatrix<double> needs to see MMatrix<complex>'s _matrix
177 
178 
179 private:
180  //build the matrix with size i,j, data not initialised
181  void build_matrix(size_t i, size_t j);
182  //build the matrix with size i,j, data initialised to data in temp
183  void build_matrix(size_t i, size_t j, Tmplt* temp);
184  //delete the matrix and set it to null
185  void delete_matrix();
186  //return the matrix overloaded to the correct type or throw
187  //if _matrix == NULL
188  static gsl_matrix* get_matrix(const MMatrix<double>& m);
189  static gsl_matrix_complex* get_matrix(const MMatrix<m_complex>& m);
190  //gsl object
191  void* _matrix;
192 };
193 
194 //Operators do what you would expect - i.e. normal mathematical functions
197 MMatrix<m_complex>& operator *=(MMatrix<m_complex>& m1, MMatrix<m_complex> m2); //M1 *= M2 returns M1 = M1*M2
198 MMatrix<double>& operator *=(MMatrix<double>& m1, MMatrix<double> m2); //M1 *= M2 returns M1 = M1*M2
203 
205 
206 template <class Tmplt> MMatrix<Tmplt> operator*(MMatrix<Tmplt>, MMatrix<Tmplt>);
207 template <class Tmplt> MMatrix<Tmplt> operator*(MMatrix<Tmplt>, MVector<Tmplt>);
208 template <class Tmplt> MMatrix<Tmplt> operator*(MMatrix<Tmplt>, Tmplt);
209 template <class Tmplt> MMatrix<Tmplt> operator/(MMatrix<Tmplt>, Tmplt);
210 
211 template <class Tmplt> MMatrix<Tmplt> operator+(MMatrix<Tmplt>, MMatrix<Tmplt>);
212 template <class Tmplt> MMatrix<Tmplt> operator-(MMatrix<Tmplt>, MMatrix<Tmplt>);
213 template <class Tmplt> MMatrix<Tmplt> operator-(MMatrix<Tmplt>); //unary minus returns matrix*(-1)
214 
215 template <class Tmplt> bool operator==(const MMatrix<Tmplt>& c1, const MMatrix<Tmplt>& c2);
216 template <class Tmplt> bool operator!=(const MMatrix<Tmplt>& c1, const MMatrix<Tmplt>& c2);
217 
218 //format is
219 // num_row num_col
220 // m11 m12 m13 ...
221 // m21 m22 m23 ...
222 // ...
223 template <class Tmplt> std::ostream& operator<<(std::ostream& out, MMatrix<Tmplt> mat);
224 template <class Tmplt> std::istream& operator>>(std::istream& in, MMatrix<Tmplt>& mat);
225 
226 //return matrix of doubles filled with either real or imaginary part of m
229 //return matrix of m_complex filled with real and imaginary parts
232 
234 
235 
236 
238 
239 template <>
240 inline size_t MMatrix<double>::num_row() const
241 {
242  if(_matrix) return ((gsl_matrix*)_matrix)->size1;
243  return 0;
244 }
245 template <>
246 inline size_t MMatrix<m_complex>::num_row() const
247 {
248  if(_matrix) return ((gsl_matrix_complex*)_matrix)->size1;
249  return 0;
250 }
251 template <>
252 inline size_t MMatrix<double>::num_col() const
253 {
254  if(_matrix) return ((gsl_matrix*)_matrix)->size2;
255  return 0;
256 }
257 
258 template <>
259 inline size_t MMatrix<m_complex>::num_col() const
260 {
261  if(_matrix) return ((gsl_matrix_complex*)_matrix)->size2;
262  return 0;
263 }
264 
266 {gsl_matrix_complex_scale((gsl_matrix_complex*)m._matrix, c); return m;}
267 
269 {gsl_matrix_scale((gsl_matrix*)m._matrix, d); return m;}
270 
272 {
273  MVector<m_complex> v0(m.num_row());
274  gsl_blas_zgemv(CblasNoTrans, m_complex_build(1.), (gsl_matrix_complex*)m._matrix, (gsl_vector_complex*)v._vector, m_complex_build(0.), (gsl_vector_complex*)v0._vector);
275  return v0;
276 }
277 
279 {
280  MVector<double> v0(m.num_row());
281  gsl_blas_dgemv(CblasNoTrans, 1., (gsl_matrix*)m._matrix, (gsl_vector*)v._vector, 0., (gsl_vector*)v0._vector);
282  return v0;
283 }
284 
286 {gsl_matrix_complex_add((gsl_matrix_complex*)m1._matrix, (gsl_matrix_complex*)m2._matrix); return m1;}
287 
289 {gsl_matrix_add((gsl_matrix*)m1._matrix, (gsl_matrix*)m2._matrix); return m1;}
290 
291 template <class Tmplt> MMatrix<Tmplt> inline operator - (MMatrix<Tmplt> m1)
292 {for(size_t i=1; i<=m1.num_row(); i++) for(size_t j=1; j<=m1.num_col(); j++) m1(i,j) = -1.*m1(i,j); return m1; }
293 
294 template <class Tmplt> MMatrix<Tmplt> inline operator*(MMatrix<Tmplt> m1, MMatrix<Tmplt> m2)
295 {return m1*=m2;}
296 
297 template <class Tmplt> MMatrix<Tmplt> inline operator*(MMatrix<Tmplt> m, MVector<Tmplt> v)
298 {return m*=v;}
299 
300 template <class Tmplt> MMatrix<Tmplt> inline operator*(MMatrix<Tmplt> m, Tmplt t)
301 {return m*=t;}
302 
303 template <class Tmplt> MMatrix<Tmplt> inline & operator/=(MMatrix<Tmplt>& m, Tmplt t)
304 {return m*=1./t;}
305 template <class Tmplt> MMatrix<Tmplt> inline operator/(MMatrix<Tmplt> m, Tmplt t)
306 {return m/=t;}
307 
308 template <class Tmplt> MMatrix<Tmplt> inline operator+(MMatrix<Tmplt> m1, MMatrix<Tmplt> m2)
309 { MMatrix<Tmplt> m3 = m1+=m2; return m3; }
310 
311 template <class Tmplt> MMatrix<Tmplt> inline & operator-=(MMatrix<Tmplt>& m1, MMatrix<Tmplt> m2)
312 {return m1 += -m2;}
313 template <class Tmplt> MMatrix<Tmplt> inline operator-(MMatrix<Tmplt> m1, MMatrix<Tmplt> m2)
314 {return m1 -= m2;}
315 
316 
317 
318 template <>
319 const double inline & MMatrix<double>::operator()(size_t i, size_t j) const
320 { return *gsl_matrix_ptr((gsl_matrix*)_matrix, i-1, j-1); }
321 template <>
322 const m_complex inline & MMatrix<m_complex>::operator()(size_t i, size_t j) const
323 { return *gsl_matrix_complex_ptr((gsl_matrix_complex*)_matrix, i-1, j-1); }
324 
325 template <>
326 double inline & MMatrix<double>::operator()(size_t i, size_t j)
327 { return *gsl_matrix_ptr((gsl_matrix*)_matrix, i-1, j-1); }
328 template <>
329 m_complex inline & MMatrix<m_complex>::operator()(size_t i, size_t j)
330 { return *gsl_matrix_complex_ptr((gsl_matrix_complex*)_matrix, i-1, j-1); }
331 
332 template <class Tmplt>
333 bool operator==(const MMatrix<Tmplt>& lhs, const MMatrix<Tmplt>& rhs)
334 {
335  if( lhs.num_row() != rhs.num_row() || lhs.num_col() != rhs.num_col() ) return false;
336  for(size_t i=1; i<=lhs.num_row(); i++)
337  for(size_t j=1; j<=lhs.num_col(); j++)
338  if(lhs(i,j) != rhs(i,j)) return false;
339  return true;
340 }
341 
342 template <class Tmplt>
343 bool inline operator!=(const MMatrix<Tmplt>& lhs, const MMatrix<Tmplt>& rhs) {return ! (lhs == rhs);}
344 
345 //iostream
346 template <class Tmplt> std::ostream& operator<<(std::ostream& out, MMatrix<Tmplt>);
347 //template <class Tmplt> std::istream& operator>>(std::istream& in, MMatrix<Tmplt>);
348 
349 template <class Tmplt>
350 gsl_matrix inline* MMatrix<Tmplt>::get_matrix(const MMatrix<double>& m)
351 {
352  if(m._matrix == NULL)
353  throw(GeneralClassicException("MMatrix::get_matrix", "Attempt to access uninitialised matrix"));
354  return (gsl_matrix*)m._matrix;
355 }
356 
357 template <class Tmplt>
358 gsl_matrix_complex inline* MMatrix<Tmplt>::get_matrix(const MMatrix<m_complex>& m)
359 {
360  if(m._matrix == NULL)
361  throw(GeneralClassicException("MMatrix::get_matrix", "Attempt to access uninitialised matrix"));
362  return (gsl_matrix_complex*)m._matrix;
363 
364 }
365 
367 }
368 
369 #endif
std::ostream & operator<<(std::ostream &os, const Attribute &attr)
Definition: Attribute.cpp:169
FLieGenerator< T, N > real(const FLieGenerator< std::complex< T >, N > &)
Take real part of a complex generator.
MMatrix< Tmplt > operator*(MMatrix< Tmplt >, Tmplt)
Definition: MMatrix.h:300
MMatrix< Tmplt > operator/(MMatrix< Tmplt >, Tmplt)
Definition: MMatrix.h:305
Mesh::Iterator & operator+=(Mesh::Iterator &lhs, const Mesh::Iterator &rhs)
m_complex m_complex_build(double r, double i)
Definition: MVector.h:62
Mesh::Iterator operator+(const Mesh::Iterator &lhs, const Mesh::Iterator &rhs)
Mesh::Iterator operator-(const Mesh::Iterator &lhs, const Mesh::Iterator &rhs)
bool operator==(const Mesh::Iterator &lhs, const Mesh::Iterator &rhs)
MMatrix< Tmplt > & operator-=(MMatrix< Tmplt > &m1, MMatrix< Tmplt > m2)
Definition: MMatrix.h:311
MMatrix< double > & operator+=(MMatrix< double > &m1, const MMatrix< double > &m2)
Definition: MMatrix.h:288
MMatrix< m_complex > complex(MMatrix< double > real)
Definition: MMatrix.cpp:396
bool operator!=(const Mesh::Iterator &lhs, const Mesh::Iterator &rhs)
Mesh::Iterator & operator-=(Mesh::Iterator &lhs, const Mesh::Iterator &rhs)
gsl_complex m_complex
Definition: MVector.h:60
MVector< m_complex > operator*(MMatrix< m_complex > m, MVector< m_complex > v)
Definition: MMatrix.h:271
MMatrix< Tmplt > & operator/=(MMatrix< Tmplt > &m, Tmplt t)
Definition: MMatrix.h:303
MMatrix< double > re(MMatrix< m_complex > mc)
Definition: MMatrix.cpp:378
MMatrix< double > & operator*=(MMatrix< double > &m, double d)
Definition: MMatrix.h:268
std::ostream & operator<<(std::ostream &out, const Mesh::Iterator &it)
Definition: Mesh.cpp:33
MMatrix< double > & operator*=(MMatrix< double > &m1, MMatrix< double > m2)
Definition: MMatrix.cpp:332
bool operator!=(const MMatrix< Tmplt > &c1, const MMatrix< Tmplt > &c2)
Definition: MMatrix.h:343
MMatrix< double > im(MMatrix< m_complex > mc)
Definition: MMatrix.cpp:387
MMatrix< Tmplt > operator+(MMatrix< Tmplt >, MMatrix< Tmplt >)
Definition: MMatrix.h:308
template std::istream & operator>>(std::istream &in, MMatrix< double > &mat)
bool operator==(const MMatrix< Tmplt > &c1, const MMatrix< Tmplt > &c2)
Definition: MMatrix.h:333
MMatrix< Tmplt > operator-(MMatrix< Tmplt >)
Definition: MMatrix.h:291
constexpr double c
The velocity of light in m/s.
Definition: Physics.h:51
void build_matrix(size_t i, size_t j, Tmplt *temp)
friend MMatrix< m_complex > & operator+=(MMatrix< m_complex > &m1, const MMatrix< m_complex > &m2)
Definition: MMatrix.h:285
friend MMatrix< Tmplt2 > operator+(MMatrix< Tmplt2 > m1, const MMatrix< Tmplt2 > m2)
~MMatrix()
destructor
Definition: MMatrix.cpp:155
MMatrix< Tmplt > & operator=(const MMatrix< Tmplt > &mm)
void invert()
turns this matrix into its inverse
static gsl_matrix * get_matrix(const MMatrix< double > &m)
Definition: MMatrix.h:350
MVector< Tmplt > get_mvector(size_t column) const
Definition: MMatrix.cpp:417
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...
Definition: MMatrix.cpp:301
void build_matrix(size_t i, size_t j)
MMatrix(const MMatrix< Tmplt > &mv)
Copy constructor makes a deep copy of mv.
MMatrix< Tmplt > T() const
returns matrix transpose T (such that M(i,j) = T(j,i))
friend MVector< m_complex > operator*(MMatrix< m_complex > m, MVector< m_complex > v)
Definition: MMatrix.h:271
size_t num_row() const
returns number of rows in the matrix
std::pair< MVector< m_complex >, MMatrix< m_complex > > eigenvectors() const
Definition: MMatrix.cpp:316
MMatrix< Tmplt > inverse() const
returns matrix inverse leaving this matrix unchanged
Definition: MMatrix.cpp:222
MMatrix< Tmplt > sub(size_t min_row, size_t max_row, size_t min_col, size_t max_col) const
Definition: MMatrix.cpp:277
static MMatrix< Tmplt > Diagonal(size_t i, Tmplt diag_value, Tmplt off_diag_value)
Construct a square matrix filling on and off diagonal values.
Definition: MMatrix.cpp:141
Tmplt & operator()(size_t row, size_t column)
Tmplt determinant() const
returns matrix determinant; throws an exception if matrix is not square
const Tmplt & operator()(size_t row, size_t column) const
static gsl_matrix_complex * get_matrix(const MMatrix< m_complex > &m)
Definition: MMatrix.h:358
MMatrix()
default constructor makes an empty MMatrix of size (0,0)
Definition: MMatrix.cpp:41
friend MMatrix< m_complex > & operator*=(MMatrix< m_complex > &m, m_complex c)
Definition: MMatrix.h:265
Tmplt trace() const
returns sum of terms with row == column, even if matrix is not square
Definition: MMatrix.cpp:290