LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ claunhr_col_getrfnp()

subroutine claunhr_col_getrfnp ( integer  M,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( * )  D,
integer  INFO 
)

CLAUNHR_COL_GETRFNP

Download CLAUNHR_COL_GETRFNP + dependencies [TGZ] [ZIP] [TXT]

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
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in,out]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.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[out]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 ).
[out]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.
Date
November 2019
Contributors:
 November 2019, Igor Kozachenko,
                Computer Science Division,
                University of California, Berkeley

Definition at line 148 of file claunhr_col_getrfnp.f.

148  IMPLICIT NONE
149 *
150 * -- LAPACK computational routine (version 3.9.0) --
151 * -- LAPACK is a software package provided by Univ. of Tennessee, --
152 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
153 * November 2019
154 *
155 * .. Scalar Arguments ..
156  INTEGER INFO, LDA, M, N
157 * ..
158 * .. Array Arguments ..
159  COMPLEX A( LDA, * ), D( * )
160 * ..
161 *
162 * =====================================================================
163 *
164 * .. Parameters ..
165  COMPLEX CONE
166  parameter( cone = ( 1.0e+0, 0.0e+0 ) )
167 * ..
168 * .. Local Scalars ..
169  INTEGER IINFO, J, JB, NB
170 * ..
171 * .. External Subroutines ..
173 * ..
174 * .. External Functions ..
175  INTEGER ILAENV
176  EXTERNAL ilaenv
177 * ..
178 * .. Intrinsic Functions ..
179  INTRINSIC max, min
180 * ..
181 * .. Executable Statements ..
182 *
183 * Test the input parameters.
184 *
185  info = 0
186  IF( m.LT.0 ) THEN
187  info = -1
188  ELSE IF( n.LT.0 ) THEN
189  info = -2
190  ELSE IF( lda.LT.max( 1, m ) ) THEN
191  info = -4
192  END IF
193  IF( info.NE.0 ) THEN
194  CALL xerbla( 'CLAUNHR_COL_GETRFNP', -info )
195  RETURN
196  END IF
197 *
198 * Quick return if possible
199 *
200  IF( min( m, n ).EQ.0 )
201  $ RETURN
202 *
203 * Determine the block size for this environment.
204 *
205 
206  nb = ilaenv( 1, 'CLAUNHR_COL_GETRFNP', ' ', m, n, -1, -1 )
207 
208  IF( nb.LE.1 .OR. nb.GE.min( m, n ) ) THEN
209 *
210 * Use unblocked code.
211 *
212  CALL claunhr_col_getrfnp2( m, n, a, lda, d, info )
213  ELSE
214 *
215 * Use blocked code.
216 *
217  DO j = 1, min( m, n ), nb
218  jb = min( min( m, n )-j+1, nb )
219 *
220 * Factor diagonal and subdiagonal blocks.
221 *
222  CALL claunhr_col_getrfnp2( m-j+1, jb, a( j, j ), lda,
223  $ d( j ), iinfo )
224 *
225  IF( j+jb.LE.n ) THEN
226 *
227 * Compute block row of U.
228 *
229  CALL ctrsm( 'Left', 'Lower', 'No transpose', 'Unit', jb,
230  $ n-j-jb+1, cone, a( j, j ), lda, a( j, j+jb ),
231  $ lda )
232  IF( j+jb.LE.m ) THEN
233 *
234 * Update trailing submatrix.
235 *
236  CALL cgemm( 'No transpose', 'No transpose', m-j-jb+1,
237  $ n-j-jb+1, jb, -cone, a( j+jb, j ), lda,
238  $ a( j, j+jb ), lda, cone, a( j+jb, j+jb ),
239  $ lda )
240  END IF
241  END IF
242  END DO
243  END IF
244  RETURN
245 *
246 * End of CLAUNHR_COL_GETRFNP
247 *
Here is the call graph for this function:
Here is the caller graph for this function:
cgemm
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189
ctrsm
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:182
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
ilaenv
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: tstiee.f:83
claunhr_col_getrfnp2
recursive subroutine claunhr_col_getrfnp2(M, N, A, LDA, D, INFO)
CLAUNHR_COL_GETRFNP2
Definition: claunhr_col_getrfnp2.f:169