LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ dlaorhr_col_getrfnp2()

recursive subroutine dlaorhr_col_getrfnp2 ( integer  M,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( * )  D,
integer  INFO 
)

DLAORHR_COL_GETRFNP2

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

Purpose:
 DLAORHR_COL_GETRFNP2 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
    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 recursive version of the LU factorization algorithm.
 Denote A - S by B. The algorithm divides the matrix B into four
 submatrices:

        [  B11 | B12  ]  where B11 is n1 by n1,
    B = [ --&mdash;|--&mdash; ]        B21 is (m-n1) by n1,
        [  B21 | B22  ]        B12 is n1 by n2,
                               B22 is (m-n1) by n2,
                               with n1 = min(m,n)/2, n2 = n-n1.


 The subroutine calls itself to factor B11, solves for B21,
 solves for B12, updates B22, then calls itself to factor B22.

 For more details on the recursive LU algorithm, see [2].

 DLAORHR_COL_GETRFNP2 is called to factorize a block by the blocked
 routine DLAORHR_COL_GETRFNP, which uses blocked code calling
is self-sufficient and can be used without DLAORHR_COL_GETRFNP.

 [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.

 [2] "Recursion leads to automatic variable blocking for dense linear
     algebra algorithms", F. Gustavson, IBM J. of Res. and Dev.,
     vol. 41, no. 6, pp. 737-755, 1997.
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.
Date
November 2019
Contributors:
 November 2019, Igor Kozachenko,
                Computer Science Division,
                University of California, Berkeley

Definition at line 169 of file dlaorhr_col_getrfnp2.f.

169  IMPLICIT NONE
170 *
171 * -- LAPACK computational routine (version 3.9.0) --
172 * -- LAPACK is a software package provided by Univ. of Tennessee, --
173 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
174 * November 2019
175 *
176 * .. Scalar Arguments ..
177  INTEGER INFO, LDA, M, N
178 * ..
179 * .. Array Arguments ..
180  DOUBLE PRECISION A( LDA, * ), D( * )
181 * ..
182 *
183 * =====================================================================
184 *
185 * .. Parameters ..
186  DOUBLE PRECISION ONE
187  parameter( one = 1.0d+0 )
188 * ..
189 * .. Local Scalars ..
190  DOUBLE PRECISION SFMIN
191  INTEGER I, IINFO, N1, N2
192 * ..
193 * .. External Functions ..
194  DOUBLE PRECISION DLAMCH
195  EXTERNAL dlamch
196 * ..
197 * .. External Subroutines ..
198  EXTERNAL dgemm, dscal, dtrsm, xerbla
199 * ..
200 * .. Intrinsic Functions ..
201  INTRINSIC abs, dsign, max, min
202 * ..
203 * .. Executable Statements ..
204 *
205 * Test the input parameters
206 *
207  info = 0
208  IF( m.LT.0 ) THEN
209  info = -1
210  ELSE IF( n.LT.0 ) THEN
211  info = -2
212  ELSE IF( lda.LT.max( 1, m ) ) THEN
213  info = -4
214  END IF
215  IF( info.NE.0 ) THEN
216  CALL xerbla( 'DLAORHR_COL_GETRFNP2', -info )
217  RETURN
218  END IF
219 *
220 * Quick return if possible
221 *
222  IF( min( m, n ).EQ.0 )
223  $ RETURN
224 
225  IF ( m.EQ.1 ) THEN
226 *
227 * One row case, (also recursion termination case),
228 * use unblocked code
229 *
230 * Transfer the sign
231 *
232  d( 1 ) = -dsign( one, a( 1, 1 ) )
233 *
234 * Construct the row of U
235 *
236  a( 1, 1 ) = a( 1, 1 ) - d( 1 )
237 *
238  ELSE IF( n.EQ.1 ) THEN
239 *
240 * One column case, (also recursion termination case),
241 * use unblocked code
242 *
243 * Transfer the sign
244 *
245  d( 1 ) = -dsign( one, a( 1, 1 ) )
246 *
247 * Construct the row of U
248 *
249  a( 1, 1 ) = a( 1, 1 ) - d( 1 )
250 *
251 * Scale the elements 2:M of the column
252 *
253 * Determine machine safe minimum
254 *
255  sfmin = dlamch('S')
256 *
257 * Construct the subdiagonal elements of L
258 *
259  IF( abs( a( 1, 1 ) ) .GE. sfmin ) THEN
260  CALL dscal( m-1, one / a( 1, 1 ), a( 2, 1 ), 1 )
261  ELSE
262  DO i = 2, m
263  a( i, 1 ) = a( i, 1 ) / a( 1, 1 )
264  END DO
265  END IF
266 *
267  ELSE
268 *
269 * Divide the matrix B into four submatrices
270 *
271  n1 = min( m, n ) / 2
272  n2 = n-n1
273 
274 *
275 * Factor B11, recursive call
276 *
277  CALL dlaorhr_col_getrfnp2( n1, n1, a, lda, d, iinfo )
278 *
279 * Solve for B21
280 *
281  CALL dtrsm( 'R', 'U', 'N', 'N', m-n1, n1, one, a, lda,
282  $ a( n1+1, 1 ), lda )
283 *
284 * Solve for B12
285 *
286  CALL dtrsm( 'L', 'L', 'N', 'U', n1, n2, one, a, lda,
287  $ a( 1, n1+1 ), lda )
288 *
289 * Update B22, i.e. compute the Schur complement
290 * B22 := B22 - B21*B12
291 *
292  CALL dgemm( 'N', 'N', m-n1, n2, n1, -one, a( n1+1, 1 ), lda,
293  $ a( 1, n1+1 ), lda, one, a( n1+1, n1+1 ), lda )
294 *
295 * Factor B22, recursive call
296 *
297  CALL dlaorhr_col_getrfnp2( m-n1, n2, a( n1+1, n1+1 ), lda,
298  $ d( n1+1 ), iinfo )
299 *
300  END IF
301  RETURN
302 *
303 * End of DLAORHR_COL_GETRFNP2
304 *
Here is the call graph for this function:
Here is the caller graph for this function:
dlaorhr_col_getrfnp2
recursive subroutine dlaorhr_col_getrfnp2(M, N, A, LDA, D, INFO)
DLAORHR_COL_GETRFNP2
Definition: dlaorhr_col_getrfnp2.f:169
dtrsm
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRSM
Definition: dtrsm.f:183
dgemm
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
xerbla
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
dscal
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:81
dlamch
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:70