Class ComplexHess

Direct Known Subclasses:
HermHess

public class ComplexHess extends ComplexUnitaryDecomposition

Computes the Hessenberg decomposition of a complex dense square matrix.

The Hessenberg decomposition decomposes a given square matrix \( A \) into the product: \[ A = QHQ^{H} \] where \( Q \) is a unitary matrix and \( H \) is an upper Hessenberg matrix, which (i.e., it has the same eigenvalues).

A matrix \( H \) is in upper Hessenburg form if it is nearly upper triangular. Specifically, if \( H \) has all zeros below the first sub-diagonal.

For example, the following matrix is in upper Hessenburg form where each '\( \times \)' may hold a different value: \[ \begin{bmatrix} \times & \times & \times & \times & \times \\ \times & \times & \times & \times & \times \\ 0 & \times & \times & \times & \times \\ 0 & 0 & \times & \times & \times \\ 0 & 0 & 0 & \times & \times \end{bmatrix} \]

Efficiency Considerations:

  • If the unitary matrix \( Q \) is not required, setting computeQ = false in the constructor may improve performance.
  • Support for in-place decomposition to reduce memory usage.
  • Support for decomposition of matrix sub-blocks, enabling efficient eigenvalue computations.

Usage:

The decomposition workflow typically follows these steps:
  1. Instantiate an instance of ComplexHess.
  2. Call decompose(CMatrix) to perform the factorization.
  3. Retrieve the resulting matrices using getH() and ComplexUnitaryDecomposition.getQ().
See Also:
  • Constructor Details

    • ComplexHess

      public ComplexHess()

      Creates a complex Hessenburg decomposer. This decomposer will compute the Hessenburg decomposition for complex dense matrices.

      By default, the unitary matrix will be computed. To specify if the unitary matrix should be computed, use ComplexHess(boolean).

      See Also:
    • ComplexHess

      public ComplexHess(boolean computeQ)

      Creates a complex Hessenburg decomposer. This decomposer will compute the Hessenburg decomposition for complex dense matrices.

      Parameters:
      computeQ - Flag indicating if the unitary matrix in the Hessenburg decomposition should be computed. If it is not needed, setting this to false may yield an increase in performance.
      See Also:
    • ComplexHess

      public ComplexHess(boolean computeQ, boolean inPlace)

      Creates a complex Hessenburg decomposer. This decomposer will compute the Hessenburg decomposition for complex dense matrices.

      Parameters:
      computeQ - Flag indicating if the unitary matrix in the Hessenburg decomposition should be computed. If it is not needed, setting this to false may yield an increase in performance.
      inPlace - Flag indicating if the decomposition should be done in-place.
      • If true, then the decomposition will be done in place.
      • If false, then the decomposition will be done out-of-place.
      See Also:
  • Method Details

    • decompose

      public ComplexHess decompose(CMatrix src)

      Computes the Hessenberg decomposition of the specified matrix.

      Note, the computation of the orthogonal matrix \( Q \) in the decomposition is deferred until ComplexUnitaryDecomposition.getQ() is explicitly called. This allows for efficient decompositions when \( Q \) is not needed.

      Overrides:
      decompose in class UnitaryDecomposition<CMatrix,Complex128[]>
      Parameters:
      src - The source matrix to decompose.
      Returns:
      A reference to this decomposer.
      Throws:
      LinearAlgebraException - If src is not a square matrix.
    • decompose

      public ComplexHess decompose(CMatrix src, int iLow, int iHigh)

      Applies decomposition to the source matrix. Note, the computation of the unitary matrix in the decomposition is deferred until ComplexUnitaryDecomposition.getQ() is explicitly called. This allows for efficient decompositions when \( Q \) is not needed.

      This method can be used specify that only a sub-block within the full matrix needs to be reduced. This is useful when you know that an upper and lower diagonal block of the matrix is already in the correct form, and you only need to reduce an inner sub-block of the full matrix. Most commonly this would be useful after balancing a matrix using ComplexBalancer, which results in the form \[ \begin{bmatrix} T_1 & X & Y \\ \mathbf{0} & B & Z \\ \mathbf{0} & \mathbf{0} & T_1 \end{bmatrix} \] where \( T_{1} \) and \( T_{2} \) are in upper-triangular form. As such, only the \( B \) block needs to be reduced. The staring row/column index of \( B \) (inclusive) is specified by iLow and the ending row/column index (exclusive) is specified by iHigh. It should be noted that the blocks \( X \) and \( Z \) will also be updated during the reduction of \( B \) so the full matrix must still be passed.

      Overrides:
      decompose in class UnitaryDecomposition<CMatrix,Complex128[]>
      Parameters:
      src - The source matrix to decompose.
      iLow - Lower bound (inclusive) of the sub-matrix to reduce to upper Hessenburg form.
      iHigh - Upper bound (exclusive) of the sub-matrix to reduce to upper Hessenburg form.
      Returns:
      A reference to this decomposer.
      Throws:
      LinearAlgebraException - If src is not a square matrix.
    • initQ

      protected CMatrix initQ()
      Creates and initializes \( Q \) to the appropriately sized identity matrix.
      Specified by:
      initQ in class UnitaryDecomposition<CMatrix,Complex128[]>
      Returns:
      An identity matrix with the appropriate size.
    • getUpper

      public CMatrix getUpper()
      Gets the upper Hessenburg matrix, \( H \), from the last decomposition. Same as getH()
      Specified by:
      getUpper in class UnitaryDecomposition<CMatrix,Complex128[]>
      Returns:
      The upper Hessenburg matrix from the last decomposition.
    • getH

      public CMatrix getH()
      Gets the upper Hessenburg matrix \( H \) from the Hessenburg decomposition.
      Returns:
      The upper Hessenburg matrix \( H \) from the Hessenburg decomposition.