LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
claunhr_col_getrfnp.f
Go to the documentation of this file.
1*> \brief \b CLAUNHR_COL_GETRFNP
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CLAUNHR_COL_GETRFNP + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/claunhr_col_getrfnp.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/claunhr_col_getrfnp.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/claunhr_col_getrfnp.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CLAUNHR_COL_GETRFNP( M, N, A, LDA, D, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDA, M, N
23* ..
24* .. Array Arguments ..
25* COMPLEX A( LDA, * ), D( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> CLAUNHR_COL_GETRFNP computes the modified LU factorization without
35*> pivoting of a complex general M-by-N matrix A. The factorization has
36*> the form:
37*>
38*> A - S = L * U,
39*>
40*> where:
41*> S is a m-by-n diagonal sign matrix with the diagonal D, so that
42*> D(i) = S(i,i), 1 <= i <= min(M,N). The diagonal D is constructed
43*> as D(i)=-SIGN(A(i,i)), where A(i,i) is the value after performing
44*> i-1 steps of Gaussian elimination. This means that the diagonal
45*> element at each step of "modified" Gaussian elimination is
46*> at least one in absolute value (so that division-by-zero not
47*> not possible during the division by the diagonal element);
48*>
49*> L is a M-by-N lower triangular matrix with unit diagonal elements
50*> (lower trapezoidal if M > N);
51*>
52*> and U is a M-by-N upper triangular matrix
53*> (upper trapezoidal if M < N).
54*>
55*> This routine is an auxiliary routine used in the Householder
56*> reconstruction routine CUNHR_COL. In CUNHR_COL, this routine is
57*> applied to an M-by-N matrix A with orthonormal columns, where each
58*> element is bounded by one in absolute value. With the choice of
59*> the matrix S above, one can show that the diagonal element at each
60*> step of Gaussian elimination is the largest (in absolute value) in
61*> the column on or below the diagonal, so that no pivoting is required
62*> for numerical stability [1].
63*>
64*> For more details on the Householder reconstruction algorithm,
65*> including the modified LU factorization, see [1].
66*>
67*> This is the blocked right-looking version of the algorithm,
68*> calling Level 3 BLAS to update the submatrix. To factorize a block,
69*> this routine calls the recursive routine CLAUNHR_COL_GETRFNP2.
70*>
71*> [1] "Reconstructing Householder vectors from tall-skinny QR",
72*> G. Ballard, J. Demmel, L. Grigori, M. Jacquelin, H.D. Nguyen,
73*> E. Solomonik, J. Parallel Distrib. Comput.,
74*> vol. 85, pp. 3-31, 2015.
75*> \endverbatim
76*
77* Arguments:
78* ==========
79*
80*> \param[in] M
81*> \verbatim
82*> M is INTEGER
83*> The number of rows of the matrix A. M >= 0.
84*> \endverbatim
85*>
86*> \param[in] N
87*> \verbatim
88*> N is INTEGER
89*> The number of columns of the matrix A. N >= 0.
90*> \endverbatim
91*>
92*> \param[in,out] A
93*> \verbatim
94*> A is COMPLEX array, dimension (LDA,N)
95*> On entry, the M-by-N matrix to be factored.
96*> On exit, the factors L and U from the factorization
97*> A-S=L*U; the unit diagonal elements of L are not stored.
98*> \endverbatim
99*>
100*> \param[in] LDA
101*> \verbatim
102*> LDA is INTEGER
103*> The leading dimension of the array A. LDA >= max(1,M).
104*> \endverbatim
105*>
106*> \param[out] D
107*> \verbatim
108*> D is COMPLEX array, dimension min(M,N)
109*> The diagonal elements of the diagonal M-by-N sign matrix S,
110*> D(i) = S(i,i), where 1 <= i <= min(M,N). The elements can be
111*> only ( +1.0, 0.0 ) or (-1.0, 0.0 ).
112*> \endverbatim
113*>
114*> \param[out] INFO
115*> \verbatim
116*> INFO is INTEGER
117*> = 0: successful exit
118*> < 0: if INFO = -i, the i-th argument had an illegal value
119*> \endverbatim
120*>
121* Authors:
122* ========
123*
124*> \author Univ. of Tennessee
125*> \author Univ. of California Berkeley
126*> \author Univ. of Colorado Denver
127*> \author NAG Ltd.
128*
129*> \ingroup launhr_col_getrfnp
130*
131*> \par Contributors:
132* ==================
133*>
134*> \verbatim
135*>
136*> November 2019, Igor Kozachenko,
137*> Computer Science Division,
138*> University of California, Berkeley
139*>
140*> \endverbatim
141*
142* =====================================================================
143 SUBROUTINE claunhr_col_getrfnp( M, N, A, LDA, D, INFO )
144 IMPLICIT NONE
145*
146* -- LAPACK computational routine --
147* -- LAPACK is a software package provided by Univ. of Tennessee, --
148* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
149*
150* .. Scalar Arguments ..
151 INTEGER INFO, LDA, M, N
152* ..
153* .. Array Arguments ..
154 COMPLEX A( LDA, * ), D( * )
155* ..
156*
157* =====================================================================
158*
159* .. Parameters ..
160 COMPLEX CONE
161 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
162* ..
163* .. Local Scalars ..
164 INTEGER IINFO, J, JB, NB
165* ..
166* .. External Subroutines ..
168 $ xerbla
169* ..
170* .. External Functions ..
171 INTEGER ILAENV
172 EXTERNAL ilaenv
173* ..
174* .. Intrinsic Functions ..
175 INTRINSIC max, min
176* ..
177* .. Executable Statements ..
178*
179* Test the input parameters.
180*
181 info = 0
182 IF( m.LT.0 ) THEN
183 info = -1
184 ELSE IF( n.LT.0 ) THEN
185 info = -2
186 ELSE IF( lda.LT.max( 1, m ) ) THEN
187 info = -4
188 END IF
189 IF( info.NE.0 ) THEN
190 CALL xerbla( 'CLAUNHR_COL_GETRFNP', -info )
191 RETURN
192 END IF
193*
194* Quick return if possible
195*
196 IF( min( m, n ).EQ.0 )
197 $ RETURN
198*
199* Determine the block size for this environment.
200*
201
202 nb = ilaenv( 1, 'CLAUNHR_COL_GETRFNP', ' ', m, n, -1, -1 )
203
204 IF( nb.LE.1 .OR. nb.GE.min( m, n ) ) THEN
205*
206* Use unblocked code.
207*
208 CALL claunhr_col_getrfnp2( m, n, a, lda, d, info )
209 ELSE
210*
211* Use blocked code.
212*
213 DO j = 1, min( m, n ), nb
214 jb = min( min( m, n )-j+1, nb )
215*
216* Factor diagonal and subdiagonal blocks.
217*
218 CALL claunhr_col_getrfnp2( m-j+1, jb, a( j, j ), lda,
219 $ d( j ), iinfo )
220*
221 IF( j+jb.LE.n ) THEN
222*
223* Compute block row of U.
224*
225 CALL ctrsm( 'Left', 'Lower', 'No transpose', 'Unit',
226 $ jb,
227 $ n-j-jb+1, cone, 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 cgemm( 'No transpose', 'No transpose',
234 $ m-j-jb+1,
235 $ n-j-jb+1, jb, -cone, a( j+jb, j ), lda,
236 $ a( j, j+jb ), lda, cone, a( j+jb, j+jb ),
237 $ lda )
238 END IF
239 END IF
240 END DO
241 END IF
242 RETURN
243*
244* End of CLAUNHR_COL_GETRFNP
245*
246 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
recursive subroutine claunhr_col_getrfnp2(m, n, a, lda, d, info)
CLAUNHR_COL_GETRFNP2
subroutine claunhr_col_getrfnp(m, n, a, lda, d, info)
CLAUNHR_COL_GETRFNP
subroutine ctrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
CTRSM
Definition ctrsm.f:180