Provided by: liblapack-doc_3.12.0-3build1_all bug

NAME

       launhr_col_getrfnp - la{un,or}hr_col_getrfnp: LU factor without pivoting

SYNOPSIS

   Functions
       subroutine claunhr_col_getrfnp (m, n, a, lda, d, info)
           CLAUNHR_COL_GETRFNP
       subroutine dlaorhr_col_getrfnp (m, n, a, lda, d, info)
           DLAORHR_COL_GETRFNP
       subroutine slaorhr_col_getrfnp (m, n, a, lda, d, info)
           SLAORHR_COL_GETRFNP
       subroutine zlaunhr_col_getrfnp (m, n, a, lda, d, info)
           ZLAUNHR_COL_GETRFNP

Detailed Description

Function Documentation

   subroutine claunhr_col_getrfnp (integer m, integer n, complex, dimension( lda, * ) a, integer
       lda, complex, dimension( * ) d, integer info)
       CLAUNHR_COL_GETRFNP

       Purpose:

            CLAUNHR_COL_GETRFNP computes the modified LU factorization without
            pivoting of a complex general M-by-N matrix A. The factorization has
            the form:

                A - S = L * U,

            where:
               S is a m-by-n diagonal sign matrix with the diagonal D, so that
               D(i) = S(i,i), 1 <= i <= min(M,N). The diagonal D is constructed
               as D(i)=-SIGN(A(i,i)), where A(i,i) is the value after performing
               i-1 steps of Gaussian elimination. This means that the diagonal
               element at each step of 'modified' Gaussian elimination is
               at least one in absolute value (so that division-by-zero not
               not possible during the division by the diagonal element);

               L is a M-by-N lower triangular matrix with unit diagonal elements
               (lower trapezoidal if M > N);

               and U is a M-by-N upper triangular matrix
               (upper trapezoidal if M < N).

            This routine is an auxiliary routine used in the Householder
            reconstruction routine CUNHR_COL. In CUNHR_COL, this routine is
            applied to an M-by-N matrix A with orthonormal columns, where each
            element is bounded by one in absolute value. With the choice of
            the matrix S above, one can show that the diagonal element at each
            step of Gaussian elimination is the largest (in absolute value) in
            the column on or below the diagonal, so that no pivoting is required
            for numerical stability [1].

            For more details on the Householder reconstruction algorithm,
            including the modified LU factorization, see [1].

            This is the blocked right-looking version of the algorithm,
            calling Level 3 BLAS to update the submatrix. To factorize a block,
            this routine calls the recursive routine CLAUNHR_COL_GETRFNP2.

            [1] 'Reconstructing Householder vectors from tall-skinny QR',
                G. Ballard, J. Demmel, L. Grigori, M. Jacquelin, H.D. Nguyen,
                E. Solomonik, J. Parallel Distrib. Comput.,
                vol. 85, pp. 3-31, 2015.

       Parameters
           M

                     M is INTEGER
                     The number of rows of the matrix A.  M >= 0.

           N

                     N is INTEGER
                     The number of columns of the matrix A.  N >= 0.

           A

                     A is COMPLEX array, dimension (LDA,N)
                     On entry, the M-by-N matrix to be factored.
                     On exit, the factors L and U from the factorization
                     A-S=L*U; the unit diagonal elements of L are not stored.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,M).

           D

                     D is COMPLEX array, dimension min(M,N)
                     The diagonal elements of the diagonal M-by-N sign matrix S,
                     D(i) = S(i,i), where 1 <= i <= min(M,N). The elements can be
                     only ( +1.0, 0.0 ) or (-1.0, 0.0 ).

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

            November 2019, Igor Kozachenko,
                           Computer Science Division,
                           University of California, Berkeley

   subroutine dlaorhr_col_getrfnp (integer m, integer n, double precision, dimension( lda, * ) a,
       integer lda, double precision, dimension( * ) d, integer info)
       DLAORHR_COL_GETRFNP

       Purpose:

            DLAORHR_COL_GETRFNP computes the modified LU factorization without
            pivoting of a real general M-by-N matrix A. The factorization has
            the form:

                A - S = L * U,

            where:
               S is a m-by-n diagonal sign matrix with the diagonal D, so that
               D(i) = S(i,i), 1 <= i <= min(M,N). The diagonal D is constructed
               as D(i)=-SIGN(A(i,i)), where A(i,i) is the value after performing
               i-1 steps of Gaussian elimination. This means that the diagonal
               element at each step of 'modified' Gaussian elimination is
               at least one in absolute value (so that division-by-zero not
               not possible during the division by the diagonal element);

               L is a M-by-N lower triangular matrix with unit diagonal elements
               (lower trapezoidal if M > N);

               and U is a M-by-N upper triangular matrix
               (upper trapezoidal if M < N).

            This routine is an auxiliary routine used in the Householder
            reconstruction routine DORHR_COL. In DORHR_COL, this routine is
            applied to an M-by-N matrix A with orthonormal columns, where each
            element is bounded by one in absolute value. With the choice of
            the matrix S above, one can show that the diagonal element at each
            step of Gaussian elimination is the largest (in absolute value) in
            the column on or below the diagonal, so that no pivoting is required
            for numerical stability [1].

            For more details on the Householder reconstruction algorithm,
            including the modified LU factorization, see [1].

            This is the blocked right-looking version of the algorithm,
            calling Level 3 BLAS to update the submatrix. To factorize a block,
            this routine calls the recursive routine DLAORHR_COL_GETRFNP2.

            [1] 'Reconstructing Householder vectors from tall-skinny QR',
                G. Ballard, J. Demmel, L. Grigori, M. Jacquelin, H.D. Nguyen,
                E. Solomonik, J. Parallel Distrib. Comput.,
                vol. 85, pp. 3-31, 2015.

       Parameters
           M

                     M is INTEGER
                     The number of rows of the matrix A.  M >= 0.

           N

                     N is INTEGER
                     The number of columns of the matrix A.  N >= 0.

           A

                     A is DOUBLE PRECISION array, dimension (LDA,N)
                     On entry, the M-by-N matrix to be factored.
                     On exit, the factors L and U from the factorization
                     A-S=L*U; the unit diagonal elements of L are not stored.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,M).

           D

                     D is DOUBLE PRECISION array, dimension min(M,N)
                     The diagonal elements of the diagonal M-by-N sign matrix S,
                     D(i) = S(i,i), where 1 <= i <= min(M,N). The elements can
                     be only plus or minus one.

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

            November 2019, Igor Kozachenko,
                           Computer Science Division,
                           University of California, Berkeley

   subroutine slaorhr_col_getrfnp (integer m, integer n, real, dimension( lda, * ) a, integer
       lda, real, dimension( * ) d, integer info)
       SLAORHR_COL_GETRFNP

       Purpose:

            SLAORHR_COL_GETRFNP computes the modified LU factorization without
            pivoting of a real general M-by-N matrix A. The factorization has
            the form:

                A - S = L * U,

            where:
               S is a m-by-n diagonal sign matrix with the diagonal D, so that
               D(i) = S(i,i), 1 <= i <= min(M,N). The diagonal D is constructed
               as D(i)=-SIGN(A(i,i)), where A(i,i) is the value after performing
               i-1 steps of Gaussian elimination. This means that the diagonal
               element at each step of 'modified' Gaussian elimination is
               at least one in absolute value (so that division-by-zero not
               not possible during the division by the diagonal element);

               L is a M-by-N lower triangular matrix with unit diagonal elements
               (lower trapezoidal if M > N);

               and U is a M-by-N upper triangular matrix
               (upper trapezoidal if M < N).

            This routine is an auxiliary routine used in the Householder
            reconstruction routine SORHR_COL. In SORHR_COL, this routine is
            applied to an M-by-N matrix A with orthonormal columns, where each
            element is bounded by one in absolute value. With the choice of
            the matrix S above, one can show that the diagonal element at each
            step of Gaussian elimination is the largest (in absolute value) in
            the column on or below the diagonal, so that no pivoting is required
            for numerical stability [1].

            For more details on the Householder reconstruction algorithm,
            including the modified LU factorization, see [1].

            This is the blocked right-looking version of the algorithm,
            calling Level 3 BLAS to update the submatrix. To factorize a block,
            this routine calls the recursive routine SLAORHR_COL_GETRFNP2.

            [1] 'Reconstructing Householder vectors from tall-skinny QR',
                G. Ballard, J. Demmel, L. Grigori, M. Jacquelin, H.D. Nguyen,
                E. Solomonik, J. Parallel Distrib. Comput.,
                vol. 85, pp. 3-31, 2015.

       Parameters
           M

                     M is INTEGER
                     The number of rows of the matrix A.  M >= 0.

           N

                     N is INTEGER
                     The number of columns of the matrix A.  N >= 0.

           A

                     A is REAL array, dimension (LDA,N)
                     On entry, the M-by-N matrix to be factored.
                     On exit, the factors L and U from the factorization
                     A-S=L*U; the unit diagonal elements of L are not stored.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,M).

           D

                     D is REAL array, dimension min(M,N)
                     The diagonal elements of the diagonal M-by-N sign matrix S,
                     D(i) = S(i,i), where 1 <= i <= min(M,N). The elements can
                     be only plus or minus one.

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

            November 2019, Igor Kozachenko,
                           Computer Science Division,
                           University of California, Berkeley

   subroutine zlaunhr_col_getrfnp (integer m, integer n, complex*16, dimension( lda, * ) a,
       integer lda, complex*16, dimension( * ) d, integer info)
       ZLAUNHR_COL_GETRFNP

       Purpose:

            ZLAUNHR_COL_GETRFNP computes the modified LU factorization without
            pivoting of a complex general M-by-N matrix A. The factorization has
            the form:

                A - S = L * U,

            where:
               S is a m-by-n diagonal sign matrix with the diagonal D, so that
               D(i) = S(i,i), 1 <= i <= min(M,N). The diagonal D is constructed
               as D(i)=-SIGN(A(i,i)), where A(i,i) is the value after performing
               i-1 steps of Gaussian elimination. This means that the diagonal
               element at each step of 'modified' Gaussian elimination is
               at least one in absolute value (so that division-by-zero not
               not possible during the division by the diagonal element);

               L is a M-by-N lower triangular matrix with unit diagonal elements
               (lower trapezoidal if M > N);

               and U is a M-by-N upper triangular matrix
               (upper trapezoidal if M < N).

            This routine is an auxiliary routine used in the Householder
            reconstruction routine ZUNHR_COL. In ZUNHR_COL, this routine is
            applied to an M-by-N matrix A with orthonormal columns, where each
            element is bounded by one in absolute value. With the choice of
            the matrix S above, one can show that the diagonal element at each
            step of Gaussian elimination is the largest (in absolute value) in
            the column on or below the diagonal, so that no pivoting is required
            for numerical stability [1].

            For more details on the Householder reconstruction algorithm,
            including the modified LU factorization, see [1].

            This is the blocked right-looking version of the algorithm,
            calling Level 3 BLAS to update the submatrix. To factorize a block,
            this routine calls the recursive routine ZLAUNHR_COL_GETRFNP2.

            [1] 'Reconstructing Householder vectors from tall-skinny QR',
                G. Ballard, J. Demmel, L. Grigori, M. Jacquelin, H.D. Nguyen,
                E. Solomonik, J. Parallel Distrib. Comput.,
                vol. 85, pp. 3-31, 2015.

       Parameters
           M

                     M is INTEGER
                     The number of rows of the matrix A.  M >= 0.

           N

                     N is INTEGER
                     The number of columns of the matrix A.  N >= 0.

           A

                     A is COMPLEX*16 array, dimension (LDA,N)
                     On entry, the M-by-N matrix to be factored.
                     On exit, the factors L and U from the factorization
                     A-S=L*U; the unit diagonal elements of L are not stored.

           LDA

                     LDA is INTEGER
                     The leading dimension of the array A.  LDA >= max(1,M).

           D

                     D is COMPLEX*16 array, dimension min(M,N)
                     The diagonal elements of the diagonal M-by-N sign matrix S,
                     D(i) = S(i,i), where 1 <= i <= min(M,N). The elements can be
                     only ( +1.0, 0.0 ) or (-1.0, 0.0 ).

           INFO

                     INFO is INTEGER
                     = 0:  successful exit
                     < 0:  if INFO = -i, the i-th argument had an illegal value

       Author
           Univ. of Tennessee

           Univ. of California Berkeley

           Univ. of Colorado Denver

           NAG Ltd.

       Contributors:

            November 2019, Igor Kozachenko,
                           Computer Science Division,
                           University of California, Berkeley

Author

       Generated automatically by Doxygen for LAPACK from the source code.