Class CommonOps_CDRM

java.lang.Object
org.ejml.dense.row.CommonOps_CDRM

@Generated("org.ejml.dense.row.CommonOps_ZDRM") public class CommonOps_CDRM extends Object
Common operations on complex numbers
  • Method Details

    • identity

      public static CMatrixRMaj identity(int width)

      Creates an identity matrix of the specified size.

      aij = 0+0i if i ≠ j
      aij = 1+0i if i = j

      Parameters:
      width - The width and height of the identity matrix.
      Returns:
      A new instance of an identity matrix.
    • identity

      public static CMatrixRMaj identity(int width, int height)

      Creates a matrix with diagonal elements set to 1 and the rest 0.

      aij = 0+0i if i ≠ j
      aij = 1+0i if i = j

      Parameters:
      width - The width of the identity matrix.
      height - The height of the identity matrix.
      Returns:
      A new instance of an identity matrix.
    • diag

      public static CMatrixRMaj diag(float... data)

      Creates a new square matrix whose diagonal elements are specified by data and all the other elements are zero.

      aij = 0 if i ≤ j
      aij = diag[i] if i = j

      Parameters:
      data - Contains the values of the diagonal elements of the resulting matrix.
      Returns:
      A new complex matrix.
    • diag

      public static CMatrixRMaj diag(@Nullable @Nullable CMatrixRMaj output, int N, float... data)
    • extractDiag

      public static void extractDiag(CMatrixRMaj src, CMatrixRMaj dst)

      Extracts the diagonal elements 'src' write it to the 'dst' vector. 'dst' can either be a row or column vector.

      Parameters:
      src - Matrix whose diagonal elements are being extracted. Not modified.
      dst - A vector the results will be written into. Modified.
    • convert

      public static void convert(FMatrixD1 input, CMatrixD1 output)
      Converts the real matrix into a complex matrix.
      Parameters:
      input - Real matrix. Not modified.
      output - Complex matrix. Modified.
    • real

      public static FMatrixRMaj real(CMatrixD1 input, @Nullable @Nullable FMatrixRMaj output)
      Places the real component of the input matrix into the output matrix.
      Parameters:
      input - Complex matrix. Not modified.
      output - real matrix. Modified.
    • imaginary

      public static FMatrixRMaj imaginary(CMatrixD1 input, @Nullable @Nullable FMatrixRMaj output)
      Places the imaginary component of the input matrix into the output matrix.
      Parameters:
      input - Complex matrix. Not modified.
      output - real matrix. Modified.
    • magnitude

      public static FMatrixRMaj magnitude(CMatrixD1 input, @Nullable @Nullable FMatrixRMaj output)

      Computes the magnitude of the complex number in the input matrix and stores the results in the output matrix.

      magnitude = sqrt(real^2 + imaginary^2)
      Parameters:
      input - Complex matrix. Not modified.
      output - real matrix. Modified.
    • conjugate

      public static CMatrixD1 conjugate(CMatrixD1 input, @Nullable @Nullable CMatrixRMaj output)

      Computes the complex conjugate of the input matrix.

      reali,j = reali,j
      imaginaryi,j = -1*imaginaryi,j

      Parameters:
      input - Input matrix. Not modified.
      output - The complex conjugate of the input matrix. Modified.
    • fill

      public static void fill(CMatrixD1 a, float real, float imaginary)

      Sets every element in the matrix to the specified value.

      aij = value

      Parameters:
      a - A matrix whose elements are about to be set. Modified.
      real - The real component
      imaginary - The imaginary component
    • add

      public static void add(CMatrixD1 a, CMatrixD1 b, CMatrixD1 c)

      Performs the following operation:

      c = a + b
      cij = aij + bij

      Matrix C can be the same instance as Matrix A and/or B.

      Parameters:
      a - A Matrix. Not modified.
      b - A Matrix. Not modified.
      c - A Matrix where the results are stored. Modified.
    • subtract

      public static void subtract(CMatrixD1 a, CMatrixD1 b, CMatrixD1 c)

      Performs the following operation:

      c = a - b
      cij = aij - bij

      Matrix C can be the same instance as Matrix A and/or B.

      Parameters:
      a - A Matrix. Not modified.
      b - A Matrix. Not modified.
      c - A Matrix where the results are stored. Modified.
    • scale

      public static void scale(float alphaReal, float alphaImag, CMatrixD1 a)

      Performs an in-place element by element scalar multiplication.

      aij = α*aij

      Parameters:
      a - The matrix that is to be scaled. Modified.
      alphaReal - real component of scale factor
      alphaImag - imaginary component of scale factor
    • mult

      public static void mult(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = a * b

      cij = ∑k=1:n { * aik * bkj}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • mult

      public static void mult(float realAlpha, float imgAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = α * a * b

      cij = α ∑k=1:n { * aik * bkj}

      Parameters:
      realAlpha - real component of scaling factor.
      imgAlpha - imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAdd

      public static void multAdd(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + a * b
      cij = cij + ∑k=1:n { aik * bkj}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAdd

      public static void multAdd(float realAlpha, float imgAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + α * a * b
      cij = cij + α * ∑k=1:n { aik * bkj}

      Parameters:
      realAlpha - real component of scaling factor.
      imgAlpha - imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multTransA

      public static void multTransA(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = aH * b

      cij = ∑k=1:n { aki * bkj}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multTransA

      public static void multTransA(float realAlpha, float imagAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = α * aH * b

      cij = α ∑k=1:n { aki * bkj}

      Parameters:
      realAlpha - Real component of scaling factor.
      imagAlpha - Imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multTransB

      public static void multTransB(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = a * bH
      cij = ∑k=1:n { aik * bjk}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multTransB

      public static void multTransB(float realAlpha, float imagAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = α * a * bH
      cij = α ∑k=1:n { aik * bjk}

      Parameters:
      realAlpha - Real component of scaling factor.
      imagAlpha - Imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multTransAB

      public static void multTransAB(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = aT * bT
      cij = ∑k=1:n { aki * bjk}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multTransAB

      public static void multTransAB(float realAlpha, float imagAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = α * aH * bH
      cij = α ∑k=1:n { aki * bjk}

      Parameters:
      realAlpha - Real component of scaling factor.
      imagAlpha - Imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAddTransA

      public static void multAddTransA(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + aH * b
      cij = cij + ∑k=1:n { aki * bkj}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAddTransA

      public static void multAddTransA(float realAlpha, float imagAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + α * aH * b
      cij =cij + α * ∑k=1:n { aki * bkj}

      Parameters:
      realAlpha - Real component of scaling factor.
      imagAlpha - Imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAddTransB

      public static void multAddTransB(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + a * bH
      cij = cij + ∑k=1:n { aik * bjk}

      Parameters:
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAddTransB

      public static void multAddTransB(float realAlpha, float imagAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + α * a * bH
      cij = cij + α * ∑k=1:n { aik * bjk}

      Parameters:
      realAlpha - Real component of scaling factor.
      imagAlpha - Imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not modified.
      b - The right matrix in the multiplication operation. Not modified.
      c - Where the results of the operation are stored. Modified.
    • multAddTransAB

      public static void multAddTransAB(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + aH * bH
      cij = cij + ∑k=1:n { aki * bjk}

      Parameters:
      a - The left matrix in the multiplication operation. Not Modified.
      b - The right matrix in the multiplication operation. Not Modified.
      c - Where the results of the operation are stored. Modified.
    • multAddTransAB

      public static void multAddTransAB(float realAlpha, float imagAlpha, CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj c)

      Performs the following operation:

      c = c + α * aH * bH
      cij = cij + α * ∑k=1:n { aki * bjk}

      Parameters:
      realAlpha - Real component of scaling factor.
      imagAlpha - Imaginary component of scaling factor.
      a - The left matrix in the multiplication operation. Not Modified.
      b - The right matrix in the multiplication operation. Not Modified.
      c - Where the results of the operation are stored. Modified.
    • transpose

      public static void transpose(CMatrixRMaj mat)

      Performs an "in-place" transpose.

      For square matrices the transpose is truly in-place and does not require additional memory. For non-square matrices, internally a temporary matrix is declared and transpose(CMatrixRMaj, CMatrixRMaj) is invoked.

      Parameters:
      mat - The matrix that is to be transposed. Modified.
    • transposeConjugate

      public static void transposeConjugate(CMatrixRMaj mat)

      Performs an "in-place" conjugate transpose.

      Parameters:
      mat - The matrix that is to be transposed. Modified.
      See Also:
    • transpose

      public static CMatrixRMaj transpose(CMatrixRMaj input, @Nullable @Nullable CMatrixRMaj output)

      Transposes input matrix 'a' and stores the results in output matrix 'b':

      bij = aji
      where 'b' is the transpose of 'a'.

      Parameters:
      input - The original matrix. Not modified.
      output - Where the transpose is stored. If null a new matrix is created. Modified.
      Returns:
      The transposed matrix.
    • transposeConjugate

      public static CMatrixRMaj transposeConjugate(CMatrixRMaj input, @Nullable @Nullable CMatrixRMaj output)

      Conjugate transposes input matrix 'a' and stores the results in output matrix 'b':

      b-reali,j = a-realj,i
      b-imaginaryi,j = -1*a-imaginaryj,i
      where 'b' is the transpose of 'a'.

      Parameters:
      input - The original matrix. Not modified.
      output - Where the transpose is stored. If null a new matrix is created. Modified.
      Returns:
      The transposed matrix.
    • trace

      public static Complex_F32 trace(CMatrixRMaj input, @Nullable @Nullable Complex_F32 output)

      Computes the matrix trace:

      trace = ∑i=1:n { aii }
      where n = min(numRows,numCols)

      Parameters:
      input - (Input) A matrix
      output - (Output) Storage for the trace. Can be null.
      Returns:
      The trace's value
    • invert

      public static boolean invert(CMatrixRMaj A)

      Performs a matrix inversion operation on the specified matrix and stores the results in the same matrix.

      a = a-1

      If the algorithm could not invert the matrix then false is returned. If it returns true that just means the algorithm finished. The results could still be bad because the matrix is singular or nearly singular.

      Parameters:
      A - The matrix that is to be inverted. Results are stored here. Modified.
      Returns:
      true if it could invert the matrix false if it could not.
    • invert

      public static boolean invert(CMatrixRMaj input, CMatrixRMaj output)

      Performs a matrix inversion operation that does not modify the original and stores the results in another matrix. The two matrices must have the same dimension.

      b = a-1

      If the algorithm could not invert the matrix then false is returned. If it returns true that just means the algorithm finished. The results could still be bad because the matrix is singular or nearly singular.

      For medium to large matrices there might be a slight performance boost to using LinearSolverFactory_CDRM instead.

      Parameters:
      input - The matrix that is to be inverted. Not modified.
      output - Where the inverse matrix is stored. Modified.
      Returns:
      true if it could invert the matrix false if it could not.
    • solve

      public static boolean solve(CMatrixRMaj a, CMatrixRMaj b, CMatrixRMaj x)

      Solves for x in the following equation:

      A*x = b

      If the system could not be solved then false is returned. If it returns true that just means the algorithm finished operating, but the results could still be bad because 'A' is singular or nearly singular.

      If repeat calls to solve are being made then one should consider using LinearSolverFactory_CDRM instead.

      It is ok for 'b' and 'x' to be the same matrix.

      Parameters:
      a - A matrix that is m by n. Not modified.
      b - A matrix that is n by k. Not modified.
      x - A matrix that is m by k. Modified.
      Returns:
      true if it could invert the matrix false if it could not.
    • det

      public static Complex_F32 det(CMatrixRMaj mat)
      Returns the determinant of the matrix. If the inverse of the matrix is also needed, then using LUDecompositionAlt_CDRM directly (or any similar algorithm) can be more efficient.
      Parameters:
      mat - The matrix whose determinant is to be computed. Not modified.
      Returns:
      The determinant.
    • elementMultiply

      public static CMatrixRMaj elementMultiply(CMatrixD1 input, float real, float imaginary, @Nullable @Nullable CMatrixRMaj output)

      Performs element by element multiplication operation with a complex number

      outputij = inputij * (real + imaginary*i)

      Parameters:
      input - The left matrix in the multiplication operation. Not modified.
      real - Real component of the number it is multiplied by
      imaginary - Imaginary component of the number it is multiplied by
      output - Where the results of the operation are stored. Modified.
    • elementMultiply

      public static CMatrixRMaj elementMultiply(CMatrixD1 inputA, CMatrixD1 inputB, @Nullable @Nullable CMatrixRMaj output)

      Performs complex multiplication between two matrices with the same shape element by element.

      outputij = inputAij * inputBij

      Parameters:
      inputA - First input matrix. Not modified.
      inputB - Second input matrix. Not modified.
      output - Where the results of the operation are stored. Modified.
    • elementDivide

      public static CMatrixRMaj elementDivide(CMatrixD1 input, float real, float imaginary, @Nullable @Nullable CMatrixRMaj output)

      Performs element by element division operation with a complex number on the right

      outputij = inputij / (real + imaginary*i)

      Parameters:
      input - The left matrix in the multiplication operation. Not modified.
      real - Real component of the number it is multiplied by
      imaginary - Imaginary component of the number it is multiplied by
      output - Where the results of the operation are stored. Modified.
    • elementDivide

      public static CMatrixRMaj elementDivide(float real, float imaginary, CMatrixD1 input, @Nullable @Nullable CMatrixRMaj output)

      Performs element by element division operation with a complex number on the right

      outputij = (real + imaginary*i) / inputij

      Parameters:
      real - Real component of the number it is multiplied by
      imaginary - Imaginary component of the number it is multiplied by
      input - The right matrix in the multiplication operation. Not modified.
      output - Where the results of the operation are stored. Modified.
    • elementDivide

      public static CMatrixRMaj elementDivide(CMatrixD1 inputA, CMatrixD1 inputB, @Nullable @Nullable CMatrixRMaj output)

      Performs complex division between two matrices with the same shape element by element.

      outputij = inputAij / inputBij

      Parameters:
      inputA - First input matrix. Not modified.
      inputB - Second input matrix. Not modified.
      output - Where the results of the operation are stored. Modified.
    • elementPower

      public static CMatrixRMaj elementPower(CMatrixD1 input, float b, @Nullable @Nullable CMatrixRMaj output)

      Element by element complex power

      outputij = inputAij / inputBij

      Parameters:
      input - Input matrix. Not modified.
      b - Power
      output - Where the results of the operation are stored. Modified.
    • elementSum

      public static Complex_F32 elementSum(CMatrixD1 input, @Nullable @Nullable Complex_F32 output)

      Computes the sum of all the elements in the matrix:

      sum(i=1:m , j=1:n ; aij)

      Parameters:
      input - An m by n matrix. Not modified.
      output - (Output) Storage for the sum. Can be null.
      Returns:
      The sum of the elements.
    • elementMinReal

      public static float elementMinReal(CMatrixD1 a)

      Returns the value of the real element in the matrix that has the minimum value.

      Min{ aij } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The minimum value out of all the real values.
    • elementMinImaginary

      public static float elementMinImaginary(CMatrixD1 a)

      Returns the value of the imaginary element in the matrix that has the minimum value.

      Min{ aij } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The minimum value out of all the real values.
    • elementMaxReal

      public static float elementMaxReal(CMatrixD1 a)

      Returns the value of the real element in the matrix that has the minimum value.

      Min{ aij } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The minimum value out of all the real values.
    • elementMaxImaginary

      public static float elementMaxImaginary(CMatrixD1 a)

      Returns the value of the imaginary element in the matrix that has the minimum value.

      Min{ aij } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The minimum value out of all the real values.
    • elementMaxMagnitude2

      public static float elementMaxMagnitude2(CMatrixD1 a)

      Returns the magnitude squared of the complex element with the largest magnitude

      Max{ |aij|^2 } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The max magnitude squared
    • setIdentity

      public static void setIdentity(CMatrixRMaj mat)
      Sets all the diagonal elements equal to one and everything else equal to zero. If this is a square matrix then it will be an identity matrix.
      Parameters:
      mat - A square matrix.
    • extract

      public static CMatrixRMaj extract(CMatrixRMaj src, int srcY0, int srcY1, int srcX0, int srcX1)

      Creates a new matrix which is the specified submatrix of 'src'

      si-y0 , j-x0 = oij for all y0 ≤ i < y1 and x0 ≤ j < x1

      where 'sij' is an element in the submatrix and 'oij' is an element in the original matrix.

      Parameters:
      src - The original matrix which is to be copied. Not modified.
      srcX0 - Start column.
      srcX1 - Stop column+1.
      srcY0 - Start row.
      srcY1 - Stop row+1.
      Returns:
      Extracted submatrix.
    • extract

      public static void extract(CMatrixRMaj src, int srcY0, int srcY1, int srcX0, int srcX1, CMatrixRMaj dst, int dstY0, int dstX0)

      Extracts a submatrix from 'src' and inserts it in a submatrix in 'dst'.

      si-y0 , j-x0 = oij for all y0 ≤ i < y1 and x0 ≤ j < x1

      where 'sij' is an element in the submatrix and 'oij' is an element in the original matrix.

      Parameters:
      src - The original matrix which is to be copied. Not modified.
      srcX0 - Start column.
      srcX1 - Stop column+1.
      srcY0 - Start row.
      srcY1 - Stop row+1.
      dst - Where the submatrix are stored. Modified.
      dstY0 - Start row in dst.
      dstX0 - start column in dst.
    • columnsToVector

      public static CMatrixRMaj[] columnsToVector(CMatrixRMaj A, @Nullable @Nullable CMatrixRMaj[] v)
      Converts the columns in a matrix into a set of vectors.
      Parameters:
      A - Matrix. Not modified.
      v - Optional storage for columns.
      Returns:
      An array of vectors.
    • elementMaxAbs

      public static float elementMaxAbs(CMatrixRMaj a)

      Returns the largest absolute value of any element in the matrix.

      Max{ |aij| } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The max abs element value in the matrix.
    • elementMinAbs

      public static float elementMinAbs(CMatrixRMaj a)

      Returns the smallest absolute value of any element in the matrix.

      Min{ |aij| } for all i and j

      Parameters:
      a - A matrix. Not modified.
      Returns:
      The min abs element value in the matrix.