LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dlaorhr_col_getrfnp()

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

DLAORHR_COL_GETRFNP

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

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

Definition at line 145 of file dlaorhr_col_getrfnp.f.

146 IMPLICIT NONE
147*
148* -- LAPACK computational routine --
149* -- LAPACK is a software package provided by Univ. of Tennessee, --
150* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
151*
152* .. Scalar Arguments ..
153 INTEGER INFO, LDA, M, N
154* ..
155* .. Array Arguments ..
156 DOUBLE PRECISION A( LDA, * ), D( * )
157* ..
158*
159* =====================================================================
160*
161* .. Parameters ..
162 DOUBLE PRECISION ONE
163 parameter( one = 1.0d+0 )
164* ..
165* .. Local Scalars ..
166 INTEGER IINFO, J, JB, NB
167* ..
168* .. External Subroutines ..
170* ..
171* .. External Functions ..
172 INTEGER ILAENV
173 EXTERNAL ilaenv
174* ..
175* .. Intrinsic Functions ..
176 INTRINSIC max, min
177* ..
178* .. Executable Statements ..
179*
180* Test the input parameters.
181*
182 info = 0
183 IF( m.LT.0 ) THEN
184 info = -1
185 ELSE IF( n.LT.0 ) THEN
186 info = -2
187 ELSE IF( lda.LT.max( 1, m ) ) THEN
188 info = -4
189 END IF
190 IF( info.NE.0 ) THEN
191 CALL xerbla( 'DLAORHR_COL_GETRFNP', -info )
192 RETURN
193 END IF
194*
195* Quick return if possible
196*
197 IF( min( m, n ).EQ.0 )
198 $ RETURN
199*
200* Determine the block size for this environment.
201*
202
203 nb = ilaenv( 1, 'DLAORHR_COL_GETRFNP', ' ', m, n, -1, -1 )
204
205 IF( nb.LE.1 .OR. nb.GE.min( m, n ) ) THEN
206*
207* Use unblocked code.
208*
209 CALL dlaorhr_col_getrfnp2( m, n, a, lda, d, info )
210 ELSE
211*
212* Use blocked code.
213*
214 DO j = 1, min( m, n ), nb
215 jb = min( min( m, n )-j+1, nb )
216*
217* Factor diagonal and subdiagonal blocks.
218*
219 CALL dlaorhr_col_getrfnp2( m-j+1, jb, a( j, j ), lda,
220 $ d( j ), iinfo )
221*
222 IF( j+jb.LE.n ) THEN
223*
224* Compute block row of U.
225*
226 CALL dtrsm( 'Left', 'Lower', 'No transpose', 'Unit', jb,
227 $ n-j-jb+1, one, a( j, j ), lda, a( j, j+jb ),
228 $ lda )
229 IF( j+jb.LE.m ) THEN
230*
231* Update trailing submatrix.
232*
233 CALL dgemm( 'No transpose', 'No transpose', m-j-jb+1,
234 $ n-j-jb+1, jb, -one, a( j+jb, j ), lda,
235 $ a( j, j+jb ), lda, one, a( j+jb, j+jb ),
236 $ lda )
237 END IF
238 END IF
239 END DO
240 END IF
241 RETURN
242*
243* End of DLAORHR_COL_GETRFNP
244*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
recursive subroutine dlaorhr_col_getrfnp2(m, n, a, lda, d, info)
DLAORHR_COL_GETRFNP2
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: