12 #ifndef EIGEN_COMPLEX_SCHUR_H 13 #define EIGEN_COMPLEX_SCHUR_H 15 #include "./HessenbergDecomposition.h" 54 typedef _MatrixType MatrixType;
56 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
57 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
58 Options = MatrixType::Options,
59 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
60 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
64 typedef typename MatrixType::Scalar
Scalar;
98 m_isInitialized(false),
99 m_matUisUptodate(false),
112 template<
typename InputType>
114 : m_matT(matrix.rows(),matrix.cols()),
115 m_matU(matrix.rows(),matrix.cols()),
116 m_hess(matrix.rows()),
117 m_isInitialized(false),
118 m_matUisUptodate(false),
121 compute(matrix.
derived(), computeU);
140 eigen_assert(m_isInitialized &&
"ComplexSchur is not initialized.");
141 eigen_assert(m_matUisUptodate &&
"The matrix U has not been computed during the ComplexSchur decomposition.");
164 eigen_assert(m_isInitialized &&
"ComplexSchur is not initialized.");
190 template<
typename InputType>
210 template<
typename HessMatrixType,
typename OrthMatrixType>
211 ComplexSchur& computeFromHessenberg(
const HessMatrixType& matrixH,
const OrthMatrixType& matrixQ,
bool computeU=
true);
219 eigen_assert(m_isInitialized &&
"ComplexSchur is not initialized.");
230 m_maxIters = maxIters;
245 static const int m_maxIterationsPerRow = 30;
248 ComplexMatrixType m_matT, m_matU;
251 bool m_isInitialized;
252 bool m_matUisUptodate;
256 bool subdiagonalEntryIsNeglegible(Index i);
257 ComplexScalar computeShift(Index iu, Index iter);
258 void reduceToTriangularForm(
bool computeU);
265 template<typename MatrixType>
266 inline bool
ComplexSchur<MatrixType>::subdiagonalEntryIsNeglegible(Index i)
268 RealScalar d = numext::norm1(m_matT.
coeff(i,i)) + numext::norm1(m_matT.
coeff(i+1,i+1));
269 RealScalar sd = numext::norm1(m_matT.
coeff(i+1,i));
272 m_matT.
coeffRef(i+1,i) = ComplexScalar(0);
280 template<
typename MatrixType>
284 if (iter == 10 || iter == 20)
287 return abs(numext::real(m_matT.
coeff(iu,iu-1))) + abs(numext::real(m_matT.
coeff(iu-1,iu-2)));
293 RealScalar normt = t.cwiseAbs().sum();
296 ComplexScalar b = t.
coeff(0,1) * t.
coeff(1,0);
297 ComplexScalar c = t.
coeff(0,0) - t.
coeff(1,1);
298 ComplexScalar disc = sqrt(c*c + RealScalar(4)*b);
299 ComplexScalar det = t.
coeff(0,0) * t.
coeff(1,1) - b;
300 ComplexScalar trace = t.
coeff(0,0) + t.
coeff(1,1);
301 ComplexScalar eival1 = (trace + disc) / RealScalar(2);
302 ComplexScalar eival2 = (trace - disc) / RealScalar(2);
304 if(numext::norm1(eival1) > numext::norm1(eival2))
305 eival2 = det / eival1;
307 eival1 = det / eival2;
310 if(numext::norm1(eival1-t.
coeff(1,1)) < numext::norm1(eival2-t.
coeff(1,1)))
311 return normt * eival1;
313 return normt * eival2;
317 template<
typename MatrixType>
318 template<
typename InputType>
321 m_matUisUptodate =
false;
322 eigen_assert(matrix.
cols() == matrix.
rows());
324 if(matrix.
cols() == 1)
326 m_matT = matrix.
derived().template cast<ComplexScalar>();
327 if(computeU) m_matU = ComplexMatrixType::Identity(1,1);
329 m_isInitialized =
true;
330 m_matUisUptodate = computeU;
335 computeFromHessenberg(m_matT, m_matU, computeU);
339 template<
typename MatrixType>
340 template<
typename HessMatrixType,
typename OrthMatrixType>
346 reduceToTriangularForm(computeU);
352 template<
typename MatrixType,
bool IsComplex>
353 struct complex_schur_reduce_to_hessenberg
359 _this.m_matT = _this.m_hess.
matrixH();
360 if(computeU) _this.m_matU = _this.m_hess.
matrixQ();
364 template<
typename MatrixType>
373 _this.m_matT = _this.m_hess.
matrixH().template cast<ComplexScalar>();
377 MatrixType Q = _this.m_hess.
matrixQ();
378 _this.m_matU = Q.template cast<ComplexScalar>();
386 template<
typename MatrixType>
389 Index maxIters = m_maxIters;
391 maxIters = m_maxIterationsPerRow * m_matT.rows();
397 Index iu = m_matT.cols() - 1;
407 if(!subdiagonalEntryIsNeglegible(iu-1))
break;
418 if(totalIter > maxIters)
break;
422 while(il > 0 && !subdiagonalEntryIsNeglegible(il-1))
431 ComplexScalar shift = computeShift(iu, iter);
434 m_matT.rightCols(m_matT.cols()-il).applyOnTheLeft(il, il+1, rot.
adjoint());
435 m_matT.topRows((std::min)(il+2,iu)+1).applyOnTheRight(il, il+1, rot);
436 if(computeU) m_matU.applyOnTheRight(il, il+1, rot);
438 for(Index i=il+1 ; i<iu ; i++)
441 m_matT.
coeffRef(i+1,i-1) = ComplexScalar(0);
442 m_matT.rightCols(m_matT.cols()-i).applyOnTheLeft(i, i+1, rot.
adjoint());
443 m_matT.topRows((std::min)(i+2,iu)+1).applyOnTheRight(i, i+1, rot);
444 if(computeU) m_matU.applyOnTheRight(i, i+1, rot);
448 if(totalIter <= maxIters)
453 m_isInitialized =
true;
454 m_matUisUptodate = computeU;
459 #endif // EIGEN_COMPLEX_SCHUR_H ComplexSchur(const EigenBase< InputType > &matrix, bool computeU=true)
Constructor; computes Schur decomposition of given matrix.
Definition: ComplexSchur.h:113
EIGEN_DEVICE_FUNC Index rows() const
Definition: EigenBase.h:58
ComplexSchur & setMaxIterations(Index maxIters)
Sets the maximum number of iterations allowed.
Definition: ComplexSchur.h:228
void makeGivens(const Scalar &p, const Scalar &q, Scalar *z=0)
Makes *this as a Givens rotation G such that applying to the left of the vector yields: ...
Definition: Jacobi.h:149
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
Definition: ForwardDeclarations.h:263
Holds information about the various numeric (i.e.
Definition: NumTraits.h:150
ComplexSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU=true)
Compute Schur decomposition from a given Hessenberg matrix.
std::complex< RealScalar > ComplexScalar
Complex scalar type for _MatrixType.
Definition: ComplexSchur.h:74
ComplexSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: ComplexSchur.h:217
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:28
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
This is an overloaded version of DenseCoeffsBase<Derived,WriteAccessors>::coeffRef(Index,Index) const provided to by-pass the creation of an evaluator of the expression, thus saving compilation efforts.
Definition: PlainObjectBase.h:177
Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime > ComplexMatrixType
Type for the matrices in the Schur decomposition.
Definition: ComplexSchur.h:81
MatrixType::Scalar Scalar
Scalar type for matrices of type _MatrixType.
Definition: ComplexSchur.h:64
HessenbergDecomposition & compute(const EigenBase< InputType > &matrix)
Computes Hessenberg decomposition of given matrix.
Definition: HessenbergDecomposition.h:152
Definition: ComplexSchur.h:20
MatrixHReturnType matrixH() const
Constructs the Hessenberg matrix H in the decomposition.
Definition: HessenbergDecomposition.h:262
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
This is an overloaded version of DenseCoeffsBase<Derived,ReadOnlyAccessors>::coeff(Index,Index) const provided to by-pass the creation of an evaluator of the expression, thus saving compilation efforts.
Definition: PlainObjectBase.h:154
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
Index getMaxIterations()
Returns the maximum number of iterations.
Definition: ComplexSchur.h:235
ComplexSchur(Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
Default constructor.
Definition: ComplexSchur.h:94
Computation was successful.
Definition: Constants.h:432
Definition: BandTriangularSolver.h:13
const ComplexMatrixType & matrixU() const
Returns the unitary matrix in the Schur decomposition.
Definition: ComplexSchur.h:138
Eigen::Index Index
Definition: ComplexSchur.h:66
EIGEN_DEVICE_FUNC Index cols() const
Definition: EigenBase.h:61
HouseholderSequenceType matrixQ() const
Reconstructs the orthogonal matrix Q in the decomposition.
Definition: HessenbergDecomposition.h:234
const int Dynamic
This value means that a positive quantity (e.g., a size) is not known at compile-time, and that instead the value is stored in some runtime variable.
Definition: Constants.h:21
Definition: ComplexSchur.h:51
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:430
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:44
Iterative procedure did not converge.
Definition: Constants.h:436
const ComplexMatrixType & matrixT() const
Returns the triangular matrix in the Schur decomposition.
Definition: ComplexSchur.h:162
JacobiRotation adjoint() const
Returns the adjoint transformation.
Definition: Jacobi.h:62