TBTK
Need a break? Support the development by playing Polarity Puzzles
Matrix.h
Go to the documentation of this file.
1 /* Copyright 2017 Kristofer Björnson
2  *
3  * Licensed under the Apache License, Version 2.0 (the "License");
4  * you may not use this file except in compliance with the License.
5  * You may obtain a copy of the License at
6  *
7  * http://www.apache.org/licenses/LICENSE-2.0
8  *
9  * Unless required by applicable law or agreed to in writing, software
10  * distributed under the License is distributed on an "AS IS" BASIS,
11  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
12  * See the License for the specific language governing permissions and
13  * limitations under the License.
14  */
15 
23 #ifndef COM_DAFER45_TBTK_MATRIX
24 #define COM_DAFER45_TBTK_MATRIX
25 
26 #include "TBTK/TBTKMacros.h"
27 
28 #include <complex>
29 
30 namespace TBTK{
31 
32 template<typename DataType, unsigned int ROWS = 0, unsigned int COLS = 0>
33 class Matrix{
34 public:
36  Matrix();
37 
39  Matrix(const Matrix<DataType, ROWS, COLS> &matrix);
40 
43 
45  ~Matrix();
46 
50  );
51 
55  );
56 
58  const DataType& at(unsigned int row, unsigned int col) const;
59 
61  DataType& at(unsigned int row, unsigned int col);
62 
64  unsigned int getNumRows() const;
65 
67  unsigned int getNumCols() const;
68 private:
70  DataType data[ROWS*COLS];
71 };
72 
73 template<typename DataType>
74 class Matrix<DataType, 0, 0>{
75 public:
77  Matrix();
78 
80  Matrix(unsigned int rows, unsigned int cols);
81 
83  Matrix(const Matrix<DataType, 0, 0> &matrix);
84 
87 
89  ~Matrix();
90 
93 
96 
98  const DataType& at(unsigned int row, unsigned int col) const;
99 
101  DataType& at(unsigned int row, unsigned int col);
102 
104  unsigned int getNumRows() const;
105 
107  unsigned int getNumCols() const;
108 
111  const Matrix<DataType, 0, 0> &rhs
112  ) const;
113 private:
115  DataType *data;
116 
118  unsigned int rows;
119 
121  unsigned int cols;
122 };
123 
124 template<>
125 class Matrix<std::complex<double>, 0, 0>{
126 public:
128  Matrix();
129 
131  Matrix(unsigned int rows, unsigned int cols);
132 
134  Matrix(const Matrix<std::complex<double>, 0, 0> &matrix);
135 
137  Matrix(Matrix<std::complex<double>, 0, 0> &&matrix);
138 
140  ~Matrix();
141 
144  const Matrix<std::complex<double>, 0, 0> &rhs
145  );
146 
149  Matrix<std::complex<double>, 0, 0> &&rhs
150  );
151 
153  const std::complex<double>& at(
154  unsigned int row,
155  unsigned int col
156  ) const;
157 
159  std::complex<double>& at(unsigned int row, unsigned int col);
160 
162  unsigned int getNumRows() const;
163 
165  unsigned int getNumCols() const;
166 
169  const Matrix<std::complex<double>, 0, 0> &rhs
170  ) const;
171 
173  void invert();
174 
176  std::complex<double> determinant();
177 private:
179  std::complex<double> *data;
180 
182  unsigned int rows;
183 
185  unsigned int cols;
186 };
187 
188 template<typename DataType, unsigned int ROWS, unsigned int COLS>
190 }
191 
192 template<typename DataType, unsigned int ROWS, unsigned int COLS>
194  for(unsigned int n = 0; n < ROWS*COLS; n++)
195  data[n] = matrix.data[n];
196 }
197 
198 template<typename DataType, unsigned int ROWS, unsigned int COLS>
200  for(unsigned int n = 0; n < ROWS*COLS; n++)
201  data[n] = matrix.data[n];
202 }
203 
204 template<typename DataType>
206  data = nullptr;
207 }
208 
209 template<typename DataType>
210 Matrix<DataType, 0, 0>::Matrix(unsigned int rows, unsigned int cols){
211  this->rows = rows;
212  this->cols = cols;
213  data = new DataType[rows*cols];
214 }
215 
216 template<typename DataType>
218  rows = matrix.rows;
219  cols = matrix.cols;
220  data = new DataType[rows*cols];
221  for(unsigned int n = 0; n < rows*cols; n++)
222  data[n] = matrix.data[n];
223 }
224 
225 template<typename DataType>
227  rows = matrix.rows;
228  cols = matrix.cols;
229  data = matrix.data;
230  matrix.data = nullptr;
231 }
232 
234  data = nullptr;
235 }
236 
237 inline Matrix<std::complex<double>, 0, 0>::Matrix(unsigned int rows, unsigned int cols){
238  this->rows = rows;
239  this->cols = cols;
240  data = new std::complex<double>[rows*cols];
241 }
242 
244  const Matrix<std::complex<double>, 0, 0> &matrix
245 ){
246  rows = matrix.rows;
247  cols = matrix.cols;
248  data = new std::complex<double>[rows*cols];
249  for(unsigned int n = 0; n < rows*cols; n++)
250  data[n] = matrix.data[n];
251 }
252 
254  Matrix<std::complex<double>, 0, 0> &&matrix
255 ){
256  rows = matrix.rows;
257  cols = matrix.cols;
258  data = matrix.data;
259  matrix.data = nullptr;
260 }
261 
262 template<typename DataType, unsigned int ROWS, unsigned int COLS>
264 }
265 
266 template<typename DataType>
268  if(data != nullptr)
269  delete [] data;
270 }
271 
272 inline Matrix<std::complex<double>, 0, 0>::~Matrix(){
273  if(data != nullptr)
274  delete [] data;
275 }
276 
277 template<typename DataType, unsigned int ROWS, unsigned int COLS>
280 ){
281  if(this != &rhs){
282  for(unsigned int n = 0; n < ROWS*COLS; n++)
283  data[n] = rhs.data[n];
284  }
285 
286  return *this;
287 }
288 
289 template<typename DataType, unsigned int ROWS, unsigned int COLS>
292 ){
293  if(this != &rhs){
294  for(unsigned int n = 0; n < ROWS*COLS; n++)
295  data[n] = rhs.data[n];
296  }
297 
298  return *this;
299 }
300 
301 template<typename DataType>
303  const Matrix<DataType, 0, 0> &rhs
304 ){
305  if(this != &rhs){
306  rows = rhs.rows;
307  cols = rhs.cols;
308 
309  if(data != nullptr)
310  delete [] data;
311 
312  data = new DataType[rows*cols];
313  for(unsigned int n = 0; n < rows*cols; n++)
314  data[n] = rhs.data[n];
315  }
316 
317  return *this;
318 }
319 
320 template<typename DataType>
323 ){
324  if(this != &rhs){
325  rows = rhs.rows;
326  cols = rhs.cols;
327 
328  if(data != nullptr)
329  delete [] data;
330 
331  data = rhs.data;
332  rhs.data = nullptr;
333  }
334 
335  return *this;
336 }
337 
339  const Matrix<std::complex<double>, 0, 0> &rhs
340 ){
341  if(this != &rhs){
342  rows = rhs.rows;
343  cols = rhs.cols;
344 
345  if(data != nullptr)
346  delete [] data;
347 
348  data = new std::complex<double>[rows*cols];
349  for(unsigned int n = 0; n < rows*cols; n++)
350  data[n] = rhs.data[n];
351  }
352 
353  return *this;
354 }
355 
357  Matrix<std::complex<double>, 0, 0> &&rhs
358 ){
359  if(this != &rhs){
360  rows = rhs.rows;
361  cols = rhs.cols;
362 
363  if(data != nullptr)
364  delete [] data;
365 
366  data = rhs.data;
367  rhs.data = nullptr;
368  }
369 
370  return *this;
371 }
372 
373 template<typename DataType, unsigned int ROWS, unsigned int COLS>
375  unsigned int row,
376  unsigned int col
377 ) const{
378  return data[row + ROWS*col];
379 }
380 
381 template<typename DataType>
383  unsigned int row,
384  unsigned int col
385 ) const{
386  return data[row + rows*col];
387 }
388 
389 inline const std::complex<double>& Matrix<std::complex<double>, 0, 0>::at(
390  unsigned int row,
391  unsigned int col
392 ) const{
393  return data[row + rows*col];
394 }
395 
396 template<typename DataType, unsigned int ROWS, unsigned int COLS>
398  unsigned int row,
399  unsigned int col
400 ){
401  return data[row + ROWS*col];
402 }
403 
404 template<typename DataType>
406  unsigned int row,
407  unsigned int col
408 ){
409  return data[row + rows*col];
410 }
411 
412 inline std::complex<double>& Matrix<std::complex<double>, 0, 0>::at(
413  unsigned int row,
414  unsigned int col
415 ){
416  return data[row + rows*col];
417 }
418 
419 template<typename DataType, unsigned int ROWS, unsigned int COLS>
421  return ROWS;
422 }
423 
424 template<typename DataType>
426  return rows;
427 }
428 
429 inline unsigned int Matrix<std::complex<double>, 0, 0>::getNumRows() const{
430  return rows;
431 }
432 
433 template<typename DataType, unsigned int ROWS, unsigned int COLS>
435  return COLS;
436 }
437 
438 template<typename DataType>
440  return cols;
441 }
442 
443 inline unsigned int Matrix<std::complex<double>, 0, 0>::getNumCols() const{
444  return cols;
445 }
446 
447 template<typename DataType>
449  const Matrix<DataType, 0, 0> &rhs
450 ) const{
451  TBTKAssert(
452  cols == rhs.rows,
453  "Matrix::operator*()",
454  "Incompatible matrix dimensions.",
455  "The matrix dimensions are " << rows << "x" << cols << " and "
456  << rhs.rows << "x" << rhs.cols << "\n"
457  );
458 
459  Matrix<DataType> result(rows, rhs.cols);
460  for(unsigned int row = 0; row < rows; row++){
461  for(unsigned int col = 0; col < rhs.cols; col++){
462  result.at(row, col) = 0.;
463 
464  for(unsigned int n = 0; n < cols; n++){
465  result.at(row, col)
466  += at(row, n)*rhs.at(n, col);
467  }
468  }
469  }
470 
471  return result;
472 }
473 
474 inline const Matrix<std::complex<double>, 0, 0> Matrix<std::complex<double>, 0, 0>::operator*(
475  const Matrix<std::complex<double>, 0, 0> &rhs
476 ) const{
477  TBTKAssert(
478  cols == rhs.rows,
479  "Matrix::operator*()",
480  "Incompatible matrix dimensions.",
481  "The matrix dimensions are " << rows << "x" << cols << " and "
482  << rhs.rows << "x" << rhs.cols << "\n"
483  );
484 
485  Matrix<std::complex<double>> result(rows, rhs.cols);
486  for(unsigned int row = 0; row < rows; row++){
487  for(unsigned int col = 0; col < rhs.cols; col++){
488  result.at(row, col) = 0.;
489 
490  for(unsigned int n = 0; n < cols; n++){
491  result.at(row, col)
492  += at(row, n)*rhs.at(n, col);
493  }
494  }
495  }
496 
497  return result;
498 }
499 
500 extern "C"{
501  void zgetrf_(
502  int *M,
503  int *N,
504  std::complex<double> *A,
505  int *lda,
506  int *ipiv,
507  int *info
508  );
509  void zgetri_(
510  int *N,
511  std::complex<double> *A,
512  int *lda,
513  int *ipiv,
514  std::complex<double> *work,
515  int *lwork,
516  int *info
517  );
518 };
519 
520 inline void Matrix<std::complex<double>, 0, 0>::invert(){
521  TBTKAssert(
522  rows == cols,
523  "Matrix::invert()",
524  "Invalid matrix dimension. Only square matrices can be"
525  << " inverted, but the matrix has size " << rows << "x" << cols
526  << "\n",
527  ""
528  );
529 
530  int *ipiv = new int[std::min(rows, cols)];
531  int lwork = rows*cols;
532  std::complex<double> *work = new std::complex<double>[lwork];
533  int info;
534 
535  zgetrf_((int*)&rows, (int*)&cols, data, (int*)&rows, ipiv, &info);
536 
537  if(info < 0){
538  TBTKExit(
539  "Matrix::invert()",
540  "Argument '" << -info << "' to zgetrf_() is invlid.",
541  "This should never happen, contact the developer."
542  );
543  }
544  else if(info > 0){
545  TBTKExit(
546  "Matrix::invert()",
547  "Unable to invert matrix since it is signular.",
548  ""
549  );
550  }
551 
552  zgetri_((int*)&rows, data, (int*)&rows, ipiv, work, &lwork, &info);
553  TBTKAssert(
554  info == 0,
555  "Matrix::invert()",
556  "Inversion failed with error code 'INFO = " << info << "'.",
557  "See the documentation for the lapack function zgetri_() for"
558  << " further information."
559  );
560 
561  delete [] ipiv;
562  delete [] work;
563 }
564 
565 inline std::complex<double> Matrix<std::complex<double>, 0, 0>::determinant(){
566  TBTKAssert(
567  rows == cols,
568  "Matrix::determinant()",
569  "Invalid matrix dimension. The determinant can only be"
570  << " calculated for square matrices, but the matrix has size "
571  << rows << "x" << cols << "\n",
572  ""
573  );
574 
575  std::complex<double> *copy = new std::complex<double>[rows*cols];
576  for(unsigned int n = 0; n < rows*cols; n++)
577  copy[n] = data[n];
578 
579  int *ipiv = new int[std::min(rows, cols)];
580 // int lwork = rows*cols;
581 // std::complex<double> *work = new std::complex<double>[lwork];
582  int info;
583 
584  zgetrf_((int*)&rows, (int*)&cols, copy, (int*)&rows, ipiv, &info);
585 
586  if(info < 0){
587  TBTKExit(
588  "Matrix::determinant()",
589  "Argument '" << -info << "' to zgetrf_() is invlid.",
590  "This should never happen, contact the developer."
591  );
592  }
593  else if(info > 0){
594  TBTKExit(
595  "Matrix::determinant()",
596  "Unable to invert matrix since it is signular.",
597  ""
598  );
599  }
600 
601  std::complex<double> det = 1.;
602  for(unsigned int n = 0; n < rows; n++){
603  Streams::out << ipiv[n] << "\n";
604  if(ipiv[n]-1 == (int)n)
605  det *= copy[rows*n + n];
606  else
607  det *= -copy[rows*n + n];
608  }
609 
610  delete [] copy;
611  delete [] ipiv;
612 
613  return det;
614 }
615 
616 }; //End namespace TBTK
617 
618 #endif
TBTK::Matrix::Matrix
Matrix()
Definition: Matrix.h:189
TBTK::Matrix< DataType, 0, 0 >
Definition: Matrix.h:74
TBTK::Matrix
Definition: Matrix.h:33
TBTK::Matrix::at
const DataType & at(unsigned int row, unsigned int col) const
Definition: Matrix.h:374
TBTK::Matrix::~Matrix
~Matrix()
Definition: Matrix.h:263
TBTK::Matrix::getNumRows
unsigned int getNumRows() const
Definition: Matrix.h:420
TBTK::Matrix< DataType, 0, 0 >::at
const DataType & at(unsigned int row, unsigned int col) const
Definition: Matrix.h:382
TBTK::Matrix::operator=
Matrix< DataType, ROWS, COLS > & operator=(const Matrix< DataType, ROWS, COLS > &rhs)
Definition: Matrix.h:278
TBTK::Matrix::getNumCols
unsigned int getNumCols() const
Definition: Matrix.h:434
TBTK::Streams::out
static std::ostream out
Definition: Streams.h:70
TBTKMacros.h
Precompiler macros.
TBTK::operator*
const Vector2d operator*(double lhs, const Vector2d &rhs)
Definition: Vector2d.h:129