Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::ComplexSchur< _MatrixType > Class Template Reference

Performs a complex Schur decomposition of a real or complex square matrix. More...

#include <src/eigen/Eigen/src/Eigenvalues/ComplexSchur.h>

+ Inheritance diagram for Eigen::ComplexSchur< _MatrixType >:
+ Collaboration diagram for Eigen::ComplexSchur< _MatrixType >:

Public Types

enum  {
  RowsAtCompileTime = MatrixType::RowsAtCompileTime , ColsAtCompileTime = MatrixType::ColsAtCompileTime , Options = MatrixType::Options , MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime ,
  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
}
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 Scalar type for matrices of type _MatrixType.
 
typedef NumTraits< Scalar >::Real RealScalar
 
typedef Eigen::Index Index
 
typedef std::complex< RealScalarComplexScalar
 Complex scalar type for _MatrixType.
 
typedef Matrix< ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTimeComplexMatrixType
 Type for the matrices in the Schur decomposition.
 

Public Member Functions

 ComplexSchur (Index size=RowsAtCompileTime==Dynamic ? 1 :RowsAtCompileTime)
 Default constructor.
 
template<typename InputType >
 ComplexSchur (const EigenBase< InputType > &matrix, bool computeU=true)
 Constructor; computes Schur decomposition of given matrix.
 
const ComplexMatrixTypematrixU () const
 Returns the unitary matrix in the Schur decomposition.
 
const ComplexMatrixTypematrixT () const
 Returns the triangular matrix in the Schur decomposition.
 
template<typename InputType >
ComplexSchurcompute (const EigenBase< InputType > &matrix, bool computeU=true)
 Computes Schur decomposition of given matrix.
 
template<typename HessMatrixType , typename OrthMatrixType >
ComplexSchurcomputeFromHessenberg (const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU=true)
 Compute Schur decomposition from a given Hessenberg matrix.
 
ComputationInfo info () const
 Reports whether previous computation was successful.
 
ComplexSchursetMaxIterations (Index maxIters)
 Sets the maximum number of iterations allowed.
 
Index getMaxIterations ()
 Returns the maximum number of iterations.
 
template<typename InputType >
ComplexSchur< MatrixType > & compute (const EigenBase< InputType > &matrix, bool computeU)
 
template<typename HessMatrixType , typename OrthMatrixType >
ComplexSchur< MatrixType > & computeFromHessenberg (const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU)
 

Static Public Attributes

static const int m_maxIterationsPerRow = 30
 Maximum number of iterations per row.
 

Protected Attributes

ComplexMatrixType m_matT
 
ComplexMatrixType m_matU
 
HessenbergDecomposition< MatrixTypem_hess
 
ComputationInfo m_info
 
bool m_isInitialized
 
bool m_matUisUptodate
 
Index m_maxIters
 

Private Member Functions

bool subdiagonalEntryIsNeglegible (Index i)
 
ComplexScalar computeShift (Index iu, Index iter)
 
void reduceToTriangularForm (bool computeU)
 

Friends

struct internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex >
 

Detailed Description

template<typename _MatrixType>
class Eigen::ComplexSchur< _MatrixType >

Performs a complex Schur decomposition of a real or complex square matrix.

\eigenvalues_module

Template Parameters
_MatrixTypethe type of the matrix of which we are computing the Schur decomposition; this is expected to be an instantiation of the Matrix class template.

Given a real or complex square matrix A, this class computes the Schur decomposition: $ A = U T U^*$ where U is a unitary complex matrix, and T is a complex upper triangular matrix. The diagonal of the matrix T corresponds to the eigenvalues of the matrix A.

Call the function compute() to compute the Schur decomposition of a given matrix. Alternatively, you can use the ComplexSchur(const MatrixType&, bool) constructor which computes the Schur decomposition at construction time. Once the decomposition is computed, you can use the matrixU() and matrixT() functions to retrieve the matrices U and V in the decomposition.

Note
This code is inspired from Jampack
See also
class RealSchur, class EigenSolver, class ComplexEigenSolver

Member Typedef Documentation

◆ ComplexMatrixType

template<typename _MatrixType >
typedef Matrix<ComplexScalar, RowsAtCompileTime, ColsAtCompileTime, Options, MaxRowsAtCompileTime, MaxColsAtCompileTime> Eigen::ComplexSchur< _MatrixType >::ComplexMatrixType

Type for the matrices in the Schur decomposition.

This is a square matrix with entries of type ComplexScalar. The size is the same as the size of _MatrixType.

◆ ComplexScalar

template<typename _MatrixType >
typedef std::complex<RealScalar> Eigen::ComplexSchur< _MatrixType >::ComplexScalar

Complex scalar type for _MatrixType.

This is std::complex<Scalar> if Scalar is real (e.g., float or double) and just Scalar if Scalar is complex.

◆ Index

template<typename _MatrixType >
typedef Eigen::Index Eigen::ComplexSchur< _MatrixType >::Index

◆ MatrixType

template<typename _MatrixType >
typedef _MatrixType Eigen::ComplexSchur< _MatrixType >::MatrixType

◆ RealScalar

template<typename _MatrixType >
typedef NumTraits<Scalar>::Real Eigen::ComplexSchur< _MatrixType >::RealScalar

◆ Scalar

template<typename _MatrixType >
typedef MatrixType::Scalar Eigen::ComplexSchur< _MatrixType >::Scalar

Scalar type for matrices of type _MatrixType.

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType >
anonymous enum
Enumerator
RowsAtCompileTime 
ColsAtCompileTime 
Options 
MaxRowsAtCompileTime 
MaxColsAtCompileTime 
55 {
56 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
57 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
58 Options = MatrixType::Options,
59 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime,
60 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
61 };
@ MaxRowsAtCompileTime
Definition ComplexSchur.h:59
@ RowsAtCompileTime
Definition ComplexSchur.h:56
@ MaxColsAtCompileTime
Definition ComplexSchur.h:60
@ ColsAtCompileTime
Definition ComplexSchur.h:57
@ Options
Definition ComplexSchur.h:58

Constructor & Destructor Documentation

◆ ComplexSchur() [1/2]

template<typename _MatrixType >
Eigen::ComplexSchur< _MatrixType >::ComplexSchur ( Index  size = RowsAtCompileTime==Dynamic ? 1 : RowsAtCompileTime)
inlineexplicit

Default constructor.

Parameters
[in]sizePositive integer, size of the matrix whose Schur decomposition will be computed.

The default constructor is useful in cases in which the user intends to perform decompositions via compute(). The size parameter is only used as a hint. It is not an error to give a wrong size, but it may impair performance.

See also
compute() for an example.
95 : m_matT(size,size),
97 m_hess(size),
98 m_isInitialized(false),
99 m_matUisUptodate(false),
100 m_maxIters(-1)
101 {}
ComplexMatrixType m_matT
Definition ComplexSchur.h:248
Index m_maxIters
Definition ComplexSchur.h:253
bool m_isInitialized
Definition ComplexSchur.h:251
ComplexMatrixType m_matU
Definition ComplexSchur.h:248
bool m_matUisUptodate
Definition ComplexSchur.h:252
HessenbergDecomposition< MatrixType > m_hess
Definition ComplexSchur.h:249
constexpr auto size(const C &c) -> decltype(c.size())
Definition span.hpp:183

◆ ComplexSchur() [2/2]

template<typename _MatrixType >
template<typename InputType >
Eigen::ComplexSchur< _MatrixType >::ComplexSchur ( const EigenBase< InputType > &  matrix,
bool  computeU = true 
)
inlineexplicit

Constructor; computes Schur decomposition of given matrix.

Parameters
[in]matrixSquare matrix whose Schur decomposition is to be computed.
[in]computeUIf true, both T and U are computed; if false, only T is computed.

This constructor calls compute() to compute the Schur decomposition.

See also
matrixT() and matrixU() for examples.
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),
119 m_maxIters(-1)
120 {
121 compute(matrix.derived(), computeU);
122 }
ComplexSchur & compute(const EigenBase< InputType > &matrix, bool computeU=true)
Computes Schur decomposition of given matrix.

References Eigen::ComplexSchur< _MatrixType >::compute(), and Eigen::EigenBase< Derived >::derived().

+ Here is the call graph for this function:

Member Function Documentation

◆ compute() [1/2]

template<typename _MatrixType >
template<typename InputType >
ComplexSchur< MatrixType > & Eigen::ComplexSchur< _MatrixType >::compute ( const EigenBase< InputType > &  matrix,
bool  computeU 
)
320{
321 m_matUisUptodate = false;
322 eigen_assert(matrix.cols() == matrix.rows());
323
324 if(matrix.cols() == 1)
325 {
326 m_matT = matrix.derived().template cast<ComplexScalar>();
327 if(computeU) m_matU = ComplexMatrixType::Identity(1,1);
328 m_info = Success;
329 m_isInitialized = true;
330 m_matUisUptodate = computeU;
331 return *this;
332 }
333
334 internal::complex_schur_reduce_to_hessenberg<MatrixType, NumTraits<Scalar>::IsComplex>::run(*this, matrix.derived(), computeU);
336 return *this;
337}
#define eigen_assert(x)
Definition Macros.h:579
ComplexSchur & computeFromHessenberg(const HessMatrixType &matrixH, const OrthMatrixType &matrixQ, bool computeU=true)
Compute Schur decomposition from a given Hessenberg matrix.
ComputationInfo m_info
Definition ComplexSchur.h:250
@ Success
Definition Constants.h:432

References Eigen::EigenBase< Derived >::cols(), Eigen::ComplexSchur< _MatrixType >::computeFromHessenberg(), Eigen::EigenBase< Derived >::derived(), eigen_assert, Eigen::ComplexSchur< _MatrixType >::m_info, Eigen::ComplexSchur< _MatrixType >::m_isInitialized, Eigen::ComplexSchur< _MatrixType >::m_matT, Eigen::ComplexSchur< _MatrixType >::m_matU, Eigen::ComplexSchur< _MatrixType >::m_matUisUptodate, Eigen::EigenBase< Derived >::rows(), and Eigen::Success.

+ Here is the call graph for this function:

◆ compute() [2/2]

template<typename _MatrixType >
template<typename InputType >
ComplexSchur & Eigen::ComplexSchur< _MatrixType >::compute ( const EigenBase< InputType > &  matrix,
bool  computeU = true 
)

Computes Schur decomposition of given matrix.

Parameters
[in]matrixSquare matrix whose Schur decomposition is to be computed.
[in]computeUIf true, both T and U are computed; if false, only T is computed.
Returns
Reference to *this

The Schur decomposition is computed by first reducing the matrix to Hessenberg form using the class HessenbergDecomposition. The Hessenberg matrix is then reduced to triangular form by performing QR iterations with a single shift. The cost of computing the Schur decomposition depends on the number of iterations; as a rough guide, it may be taken on the number of iterations; as a rough guide, it may be taken to be $25n^3$ complex flops, or $10n^3$ complex flops if computeU is false.

Example:

Output:

See also
compute(const MatrixType&, bool, Index)

Referenced by Eigen::ComplexSchur< _MatrixType >::ComplexSchur().

+ Here is the caller graph for this function:

◆ computeFromHessenberg() [1/2]

template<typename _MatrixType >
template<typename HessMatrixType , typename OrthMatrixType >
ComplexSchur< MatrixType > & Eigen::ComplexSchur< _MatrixType >::computeFromHessenberg ( const HessMatrixType &  matrixH,
const OrthMatrixType &  matrixQ,
bool  computeU 
)
342{
343 m_matT = matrixH;
344 if(computeU)
345 m_matU = matrixQ;
346 reduceToTriangularForm(computeU);
347 return *this;
348}
void reduceToTriangularForm(bool computeU)
Definition ComplexSchur.h:387

References Eigen::ComplexSchur< _MatrixType >::m_matT, Eigen::ComplexSchur< _MatrixType >::m_matU, and Eigen::ComplexSchur< _MatrixType >::reduceToTriangularForm().

+ Here is the call graph for this function:

◆ computeFromHessenberg() [2/2]

template<typename _MatrixType >
template<typename HessMatrixType , typename OrthMatrixType >
ComplexSchur & Eigen::ComplexSchur< _MatrixType >::computeFromHessenberg ( const HessMatrixType &  matrixH,
const OrthMatrixType &  matrixQ,
bool  computeU = true 
)

Compute Schur decomposition from a given Hessenberg matrix.

Parameters
[in]matrixHMatrix in Hessenberg form H
[in]matrixQorthogonal matrix Q that transform a matrix A to H : A = Q H Q^T
computeUComputes the matriX U of the Schur vectors
Returns
Reference to *this

This routine assumes that the matrix is already reduced in Hessenberg form matrixH using either the class HessenbergDecomposition or another mean. It computes the upper quasi-triangular matrix T of the Schur decomposition of H When computeU is true, this routine computes the matrix U such that A = U T U^T = (QZ) T (QZ)^T = Q H Q^T where A is the initial matrix

NOTE Q is referenced if computeU is true; so, if the initial orthogonal matrix is not available, the user should give an identity matrix (Q.setIdentity())

See also
compute(const MatrixType&, bool)

Referenced by Eigen::ComplexSchur< _MatrixType >::compute().

+ Here is the caller graph for this function:

◆ computeShift()

template<typename MatrixType >
ComplexSchur< MatrixType >::ComplexScalar Eigen::ComplexSchur< MatrixType >::computeShift ( Index  iu,
Index  iter 
)
private

Compute the shift in the current QR iteration.

282{
283 using std::abs;
284 if (iter == 10 || iter == 20)
285 {
286 // exceptional shift, taken from http://www.netlib.org/eispack/comqr.f
287 return abs(numext::real(m_matT.coeff(iu,iu-1))) + abs(numext::real(m_matT.coeff(iu-1,iu-2)));
288 }
289
290 // compute the shift as one of the eigenvalues of t, the 2x2
291 // diagonal block on the bottom of the active submatrix
292 Matrix<ComplexScalar,2,2> t = m_matT.template block<2,2>(iu-1,iu-1);
293 RealScalar normt = t.cwiseAbs().sum();
294 t /= normt; // the normalization by sf is to avoid under/overflow
295
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);
303
304 if(numext::norm1(eival1) > numext::norm1(eival2))
305 eival2 = det / eival1;
306 else
307 eival1 = det / eival2;
308
309 // choose the eigenvalue closest to the bottom entry of the diagonal
310 if(numext::norm1(eival1-t.coeff(1,1)) < numext::norm1(eival2-t.coeff(1,1)))
311 return normt * eival1;
312 else
313 return normt * eival2;
314}
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition ArrayCwiseUnaryOps.h:152
NumTraits< Scalar >::Real RealScalar
Definition ComplexSchur.h:65
std::complex< RealScalar > ComplexScalar
Complex scalar type for _MatrixType.
Definition ComplexSchur.h:74
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
Definition PlainObjectBase.h:160
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC half abs(const half &a)
Definition Half.h:445

References Eigen::PlainObjectBase< Derived >::coeff(), Eigen::ComplexSchur< _MatrixType >::m_matT, and sqrt().

Referenced by Eigen::ComplexSchur< _MatrixType >::reduceToTriangularForm().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ getMaxIterations()

template<typename _MatrixType >
Index Eigen::ComplexSchur< _MatrixType >::getMaxIterations ( )
inline

Returns the maximum number of iterations.

236 {
237 return m_maxIters;
238 }

References Eigen::ComplexSchur< _MatrixType >::m_maxIters.

Referenced by Eigen::ComplexEigenSolver< _MatrixType >::getMaxIterations().

+ Here is the caller graph for this function:

◆ info()

template<typename _MatrixType >
ComputationInfo Eigen::ComplexSchur< _MatrixType >::info ( ) const
inline

Reports whether previous computation was successful.

Returns
Success if computation was succesful, NoConvergence otherwise.
218 {
219 eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
220 return m_info;
221 }

References eigen_assert, Eigen::ComplexSchur< _MatrixType >::m_info, and Eigen::ComplexSchur< _MatrixType >::m_isInitialized.

Referenced by Eigen::ComplexEigenSolver< _MatrixType >::info().

+ Here is the caller graph for this function:

◆ matrixT()

template<typename _MatrixType >
const ComplexMatrixType & Eigen::ComplexSchur< _MatrixType >::matrixT ( ) const
inline

Returns the triangular matrix in the Schur decomposition.

Returns
A const reference to the matrix T.

It is assumed that either the constructor ComplexSchur(const MatrixType& matrix, bool computeU) or the member function compute(const MatrixType& matrix, bool computeU) has been called before to compute the Schur decomposition of a matrix.

Note that this function returns a plain square matrix. If you want to reference only the upper triangular part, use:

schur.matrixT().triangularView<Upper>()
@ Upper
Definition Constants.h:206

Example:

Output:

 
163 {
164 eigen_assert(m_isInitialized && "ComplexSchur is not initialized.");
165 return m_matT;
166 }

References eigen_assert, Eigen::ComplexSchur< _MatrixType >::m_isInitialized, and Eigen::ComplexSchur< _MatrixType >::m_matT.

◆ matrixU()

template<typename _MatrixType >
const ComplexMatrixType & Eigen::ComplexSchur< _MatrixType >::matrixU ( ) const
inline

Returns the unitary matrix in the Schur decomposition.

Returns
A const reference to the matrix U.

It is assumed that either the constructor ComplexSchur(const MatrixType& matrix, bool computeU) or the member function compute(const MatrixType& matrix, bool computeU) has been called before to compute the Schur decomposition of a matrix, and that computeU was set to true (the default value).

Example:

Output:

 
139 {
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.");
142 return m_matU;
143 }

References eigen_assert, Eigen::ComplexSchur< _MatrixType >::m_isInitialized, Eigen::ComplexSchur< _MatrixType >::m_matU, and Eigen::ComplexSchur< _MatrixType >::m_matUisUptodate.

◆ reduceToTriangularForm()

template<typename MatrixType >
void Eigen::ComplexSchur< MatrixType >::reduceToTriangularForm ( bool  computeU)
private
388{
389 Index maxIters = m_maxIters;
390 if (maxIters == -1)
391 maxIters = m_maxIterationsPerRow * m_matT.rows();
392
393 // The matrix m_matT is divided in three parts.
394 // Rows 0,...,il-1 are decoupled from the rest because m_matT(il,il-1) is zero.
395 // Rows il,...,iu is the part we are working on (the active submatrix).
396 // Rows iu+1,...,end are already brought in triangular form.
397 Index iu = m_matT.cols() - 1;
398 Index il;
399 Index iter = 0; // number of iterations we are working on the (iu,iu) element
400 Index totalIter = 0; // number of iterations for whole matrix
401
402 while(true)
403 {
404 // find iu, the bottom row of the active submatrix
405 while(iu > 0)
406 {
407 if(!subdiagonalEntryIsNeglegible(iu-1)) break;
408 iter = 0;
409 --iu;
410 }
411
412 // if iu is zero then we are done; the whole matrix is triangularized
413 if(iu==0) break;
414
415 // if we spent too many iterations, we give up
416 iter++;
417 totalIter++;
418 if(totalIter > maxIters) break;
419
420 // find il, the top row of the active submatrix
421 il = iu-1;
422 while(il > 0 && !subdiagonalEntryIsNeglegible(il-1))
423 {
424 --il;
425 }
426
427 /* perform the QR step using Givens rotations. The first rotation
428 creates a bulge; the (il+2,il) element becomes nonzero. This
429 bulge is chased down to the bottom of the active submatrix. */
430
431 ComplexScalar shift = computeShift(iu, iter);
432 JacobiRotation<ComplexScalar> rot;
433 rot.makeGivens(m_matT.coeff(il,il) - shift, m_matT.coeff(il+1,il));
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);
437
438 for(Index i=il+1 ; i<iu ; i++)
439 {
440 rot.makeGivens(m_matT.coeffRef(i,i-1), m_matT.coeffRef(i+1,i-1), &m_matT.coeffRef(i,i-1));
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);
445 }
446 }
447
448 if(totalIter <= maxIters)
449 m_info = Success;
450 else
452
453 m_isInitialized = true;
454 m_matUisUptodate = computeU;
455}
Eigen::Index Index
Definition ComplexSchur.h:66
bool subdiagonalEntryIsNeglegible(Index i)
Definition ComplexSchur.h:266
ComplexScalar computeShift(Index iu, Index iter)
Definition ComplexSchur.h:281
static const int m_maxIterationsPerRow
Maximum number of iterations per row.
Definition ComplexSchur.h:245
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Definition PlainObjectBase.h:183
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index cols() const
Definition PlainObjectBase.h:153
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
Definition PlainObjectBase.h:151
@ NoConvergence
Definition Constants.h:436

References Eigen::JacobiRotation< Scalar >::adjoint(), Eigen::PlainObjectBase< Derived >::coeff(), Eigen::Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols >::coeffRef(), Eigen::PlainObjectBase< Derived >::cols(), Eigen::ComplexSchur< _MatrixType >::computeShift(), Eigen::ComplexSchur< _MatrixType >::m_info, Eigen::ComplexSchur< _MatrixType >::m_isInitialized, Eigen::ComplexSchur< _MatrixType >::m_matT, Eigen::ComplexSchur< _MatrixType >::m_matU, Eigen::ComplexSchur< _MatrixType >::m_matUisUptodate, Eigen::ComplexSchur< _MatrixType >::m_maxIterationsPerRow, Eigen::ComplexSchur< _MatrixType >::m_maxIters, Eigen::JacobiRotation< Scalar >::makeGivens(), Eigen::NoConvergence, Eigen::PlainObjectBase< Derived >::rows(), Eigen::ComplexSchur< _MatrixType >::subdiagonalEntryIsNeglegible(), and Eigen::Success.

Referenced by Eigen::ComplexSchur< _MatrixType >::computeFromHessenberg().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ setMaxIterations()

template<typename _MatrixType >
ComplexSchur & Eigen::ComplexSchur< _MatrixType >::setMaxIterations ( Index  maxIters)
inline

Sets the maximum number of iterations allowed.

If not specified by the user, the maximum number of iterations is m_maxIterationsPerRow times the size of the matrix.

229 {
230 m_maxIters = maxIters;
231 return *this;
232 }

References Eigen::ComplexSchur< _MatrixType >::m_maxIters.

Referenced by Eigen::ComplexEigenSolver< _MatrixType >::setMaxIterations().

+ Here is the caller graph for this function:

◆ subdiagonalEntryIsNeglegible()

template<typename MatrixType >
bool Eigen::ComplexSchur< MatrixType >::subdiagonalEntryIsNeglegible ( Index  i)
inlineprivate

If m_matT(i+1,i) is neglegible in floating point arithmetic compared to m_matT(i,i) and m_matT(j,j), then set it to zero and return true, else return false.

267{
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));
270 if (internal::isMuchSmallerThan(sd, d, NumTraits<RealScalar>::epsilon()))
271 {
272 m_matT.coeffRef(i+1,i) = ComplexScalar(0);
273 return true;
274 }
275 return false;
276}
EIGEN_DEVICE_FUNC bool isMuchSmallerThan(const Scalar &x, const OtherScalar &y, const typename NumTraits< Scalar >::Real &precision=NumTraits< Scalar >::dummy_precision())
Definition MathFunctions.h:1354

References Eigen::PlainObjectBase< Derived >::coeff(), and Eigen::Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols >::coeffRef().

Referenced by Eigen::ComplexSchur< _MatrixType >::reduceToTriangularForm().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Friends And Related Symbol Documentation

◆ internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex >

template<typename _MatrixType >
friend struct internal::complex_schur_reduce_to_hessenberg< MatrixType, NumTraits< Scalar >::IsComplex >
friend

Member Data Documentation

◆ m_hess

◆ m_info

◆ m_isInitialized

◆ m_matT

◆ m_matU

◆ m_matUisUptodate

◆ m_maxIterationsPerRow

template<typename _MatrixType >
const int Eigen::ComplexSchur< _MatrixType >::m_maxIterationsPerRow = 30
static

Maximum number of iterations per row.

If not otherwise specified, the maximum number of iterations is this number times the size of the matrix. It is currently set to 30.

Referenced by Eigen::ComplexSchur< _MatrixType >::reduceToTriangularForm().

◆ m_maxIters


The documentation for this class was generated from the following file: