Mantid
Loading...
Searching...
No Matches
EigenComplexMatrix.h
Go to the documentation of this file.
1// Mantid Repository : https://github.com/mantidproject/mantid
2//
3// Copyright © 2010 ISIS Rutherford Appleton Laboratory UKRI,
4// NScD Oak Ridge National Laboratory, European Spallation Source,
5// Institut Laue - Langevin & CSNS, Institute of High Energy Physics, CAS
6// SPDX - License - Identifier: GPL - 3.0 +
7#pragma once
8
9#include "MantidCurveFitting/DllConfig.h"
13
14#include <Eigen/Dense>
15
16#include <complex>
17#include <iomanip>
18#include <stdexcept>
19#include <vector>
20
21typedef Eigen::Map<const Eigen::MatrixXcd, 0, Mantid::CurveFitting::dynamic_stride> complex_matrix_map_type;
22typedef Eigen::Map<const Eigen::VectorXcd, 0, Mantid::CurveFitting::dynamic_stride> complex_vector_map_type;
23
24namespace Mantid {
25namespace CurveFitting {
26class ComplexMatrix;
27
31class MANTID_CURVEFITTING_DLL ComplexMatrix {
32public:
36 ComplexMatrix(const size_t nx, const size_t ny);
40 ComplexMatrix(const ComplexMatrix &M, const size_t row, const size_t col, const size_t nRows, const size_t nCols);
42 explicit ComplexMatrix(ComplexMatrix &&m) noexcept;
44 explicit ComplexMatrix(Eigen::MatrixXcd &&m) noexcept;
45
47 ComplexMatrix &operator=(const ComplexMatrix &M);
49 ComplexMatrix &operator=(ComplexMatrix &&M);
51 ComplexMatrix &operator=(const Eigen::MatrixXcd M);
52
54 inline Eigen::MatrixXcd &eigen() { return m_matrix; }
56 inline const Eigen::MatrixXcd eigen() const { return m_matrix; }
57
59 bool isEmpty() const;
61 void resize(const size_t nx, const size_t ny);
63 size_t size1() const;
65 size_t size2() const;
67 void set(const size_t i, const size_t j, const ComplexType value);
69 ComplexType get(const size_t i, const size_t j) const;
71 ComplexType operator()(const size_t i, const size_t j) const;
73 ComplexType &operator()(const size_t i, const size_t j);
74
76 void identity();
78 void zero();
80 void diag(const ComplexVector &d);
90 ComplexMatrix operator*(const EigenMatrix &m) const;
93
95 ComplexVector copyRow(const size_t i) const;
97 ComplexVector copyColumn(const size_t i) const;
99 void sortColumns(const std::vector<size_t> &indices);
100
105 void solve(const ComplexVector &rhs, ComplexVector &x);
107 void invert();
109 ComplexType det();
111 void eigenSystemHermitian(EigenVector &eigenValues, ComplexMatrix &eigenVectors);
112
114 ComplexMatrix tr() const;
116 ComplexMatrix ctr() const;
119 std::vector<double> packToStdVector() const;
122 void unpackFromStdVector(const std::vector<double> &v);
123
124private:
126 Eigen::MatrixXcd m_matrix;
127};
128
130inline std::ostream &operator<<(std::ostream &ostr, const ComplexMatrix &m) {
131 std::ios::fmtflags fflags(ostr.flags());
132 ostr << std::scientific << std::setprecision(6);
133 for (size_t i = 0; i < m.size1(); ++i) {
134 for (size_t j = 0; j < m.size2(); ++j) {
135 auto value = m.get(i, j);
136 ostr << std::setw(28) << std::setprecision(13) << value.real() << "+" << value.imag() << "j ";
137 }
138 ostr << '\n';
139 }
140 ostr.flags(fflags);
141 return ostr;
142}
143
144} // namespace CurveFitting
145} // namespace Mantid
const std::vector< double > & rhs
Eigen::Map< const Eigen::VectorXcd, 0, Mantid::CurveFitting::dynamic_stride > complex_vector_map_type
Eigen::Map< const Eigen::MatrixXcd, 0, Mantid::CurveFitting::dynamic_stride > complex_matrix_map_type
double value
The value of the point.
Definition: FitMW.cpp:51
A complex-valued matrix for linear algebra computations.
Eigen::MatrixXcd m_matrix
The pointer to the complex Eigen matrix.
Eigen::MatrixXcd & eigen()
Get the reference to the eigen matrix.
const Eigen::MatrixXcd eigen() const
Get the const pointer to the GSL matrix.
A complex-valued vector for linear algebra computations.
A wrapper around Eigen::Matrix.
Definition: EigenMatrix.h:33
A wrapper around Eigen::Vector.
Definition: EigenVector.h:27
MatrixWorkspace_sptr MANTID_API_DLL operator*(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Multiply two workspaces.
MatrixWorkspace_sptr MANTID_API_DLL operator+=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Adds two workspaces.
MatrixWorkspace_sptr MANTID_API_DLL operator-=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Subtracts two workspaces.
MANTID_API_DLL std::ostream & operator<<(std::ostream &, const AlgorithmHistory &)
Prints a text representation.
MatrixWorkspace_sptr MANTID_API_DLL operator*=(const MatrixWorkspace_sptr &lhs, const MatrixWorkspace_sptr &rhs)
Multiply two workspaces.
std::complex< double > ComplexType
Helper class which provides the Collimation Length for SANS instruments.