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

◆ slatdf()

subroutine slatdf ( integer ijob,
integer n,
real, dimension( ldz, * ) z,
integer ldz,
real, dimension( * ) rhs,
real rdsum,
real rdscal,
integer, dimension( * ) ipiv,
integer, dimension( * ) jpiv )

SLATDF uses the LU factorization of the n-by-n matrix computed by sgetc2 and computes a contribution to the reciprocal Dif-estimate.

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

Purpose:
!>
!> SLATDF uses the LU factorization of the n-by-n matrix Z computed by
!> SGETC2 and computes a contribution to the reciprocal Dif-estimate
!> by solving Z * x = b for x, and choosing the r.h.s. b such that
!> the norm of x is as large as possible. On entry RHS = b holds the
!> contribution from earlier solved sub-systems, and on return RHS = x.
!>
!> The factorization of Z returned by SGETC2 has the form Z = P*L*U*Q,
!> where P and Q are permutation matrices. L is lower triangular with
!> unit diagonal elements and U is upper triangular.
!> 
Parameters
[in]IJOB
!>          IJOB is INTEGER
!>          IJOB = 2: First compute an approximative null-vector e
!>              of Z using SGECON, e is normalized and solve for
!>              Zx = +-e - f with the sign giving the greater value
!>              of 2-norm(x). About 5 times as expensive as Default.
!>          IJOB .ne. 2: Local look ahead strategy where all entries of
!>              the r.h.s. b is chosen as either +1 or -1 (Default).
!> 
[in]N
!>          N is INTEGER
!>          The number of columns of the matrix Z.
!> 
[in]Z
!>          Z is REAL array, dimension (LDZ, N)
!>          On entry, the LU part of the factorization of the n-by-n
!>          matrix Z computed by SGETC2:  Z = P * L * U * Q
!> 
[in]LDZ
!>          LDZ is INTEGER
!>          The leading dimension of the array Z.  LDA >= max(1, N).
!> 
[in,out]RHS
!>          RHS is REAL array, dimension N.
!>          On entry, RHS contains contributions from other subsystems.
!>          On exit, RHS contains the solution of the subsystem with
!>          entries according to the value of IJOB (see above).
!> 
[in,out]RDSUM
!>          RDSUM is REAL
!>          On entry, the sum of squares of computed contributions to
!>          the Dif-estimate under computation by STGSYL, where the
!>          scaling factor RDSCAL (see below) has been factored out.
!>          On exit, the corresponding sum of squares updated with the
!>          contributions from the current sub-system.
!>          If TRANS = 'T' RDSUM is not touched.
!>          NOTE: RDSUM only makes sense when STGSY2 is called by STGSYL.
!> 
[in,out]RDSCAL
!>          RDSCAL is REAL
!>          On entry, scaling factor used to prevent overflow in RDSUM.
!>          On exit, RDSCAL is updated w.r.t. the current contributions
!>          in RDSUM.
!>          If TRANS = 'T', RDSCAL is not touched.
!>          NOTE: RDSCAL only makes sense when STGSY2 is called by
!>                STGSYL.
!> 
[in]IPIV
!>          IPIV is INTEGER array, dimension (N).
!>          The pivot indices; for 1 <= i <= N, row i of the
!>          matrix has been interchanged with row IPIV(i).
!> 
[in]JPIV
!>          JPIV is INTEGER array, dimension (N).
!>          The pivot indices; for 1 <= j <= N, column j of the
!>          matrix has been interchanged with column JPIV(j).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
This routine is a further developed implementation of algorithm BSOLVE in [1] using complete pivoting in the LU factorization.
Contributors:
Bo Kagstrom and Peter Poromaa, Department of Computing Science, Umea University, S-901 87 Umea, Sweden.
References:
!>
!>
!>  [1] Bo Kagstrom and Lars Westin,
!>      Generalized Schur Methods with Condition Estimators for
!>      Solving the Generalized Sylvester Equation, IEEE Transactions
!>      on Automatic Control, Vol. 34, No. 7, July 1989, pp 745-751.
!>
!>  [2] Peter Poromaa,
!>      On Efficient and Robust Estimators for the Separation
!>      between two Regular Matrix Pairs with Applications in
!>      Condition Estimation. Report IMINF-95.05, Departement of
!>      Computing Science, Umea University, S-901 87 Umea, Sweden, 1995.
!> 

Definition at line 167 of file slatdf.f.

169*
170* -- LAPACK auxiliary routine --
171* -- LAPACK is a software package provided by Univ. of Tennessee, --
172* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
173*
174* .. Scalar Arguments ..
175 INTEGER IJOB, LDZ, N
176 REAL RDSCAL, RDSUM
177* ..
178* .. Array Arguments ..
179 INTEGER IPIV( * ), JPIV( * )
180 REAL RHS( * ), Z( LDZ, * )
181* ..
182*
183* =====================================================================
184*
185* .. Parameters ..
186 INTEGER MAXDIM
187 parameter( maxdim = 8 )
188 REAL ZERO, ONE
189 parameter( zero = 0.0e+0, one = 1.0e+0 )
190* ..
191* .. Local Scalars ..
192 INTEGER I, INFO, J, K
193 REAL BM, BP, PMONE, SMINU, SPLUS, TEMP
194* ..
195* .. Local Arrays ..
196 INTEGER IWORK( MAXDIM )
197 REAL WORK( 4*MAXDIM ), XM( MAXDIM ), XP( MAXDIM )
198* ..
199* .. External Subroutines ..
200 EXTERNAL saxpy, scopy, sgecon, sgesc2, slassq,
201 $ slaswp,
202 $ sscal
203* ..
204* .. External Functions ..
205 REAL SASUM, SDOT
206 EXTERNAL sasum, sdot
207* ..
208* .. Intrinsic Functions ..
209 INTRINSIC abs, sqrt
210* ..
211* .. Executable Statements ..
212*
213 IF( ijob.NE.2 ) THEN
214*
215* Apply permutations IPIV to RHS
216*
217 CALL slaswp( 1, rhs, ldz, 1, n-1, ipiv, 1 )
218*
219* Solve for L-part choosing RHS either to +1 or -1.
220*
221 pmone = -one
222*
223 DO 10 j = 1, n - 1
224 bp = rhs( j ) + one
225 bm = rhs( j ) - one
226 splus = one
227*
228* Look-ahead for L-part RHS(1:N-1) = + or -1, SPLUS and
229* SMIN computed more efficiently than in BSOLVE [1].
230*
231 splus = splus + sdot( n-j, z( j+1, j ), 1, z( j+1, j ),
232 $ 1 )
233 sminu = sdot( n-j, z( j+1, j ), 1, rhs( j+1 ), 1 )
234 splus = splus*rhs( j )
235 IF( splus.GT.sminu ) THEN
236 rhs( j ) = bp
237 ELSE IF( sminu.GT.splus ) THEN
238 rhs( j ) = bm
239 ELSE
240*
241* In this case the updating sums are equal and we can
242* choose RHS(J) +1 or -1. The first time this happens
243* we choose -1, thereafter +1. This is a simple way to
244* get good estimates of matrices like Byers well-known
245* example (see [1]). (Not done in BSOLVE.)
246*
247 rhs( j ) = rhs( j ) + pmone
248 pmone = one
249 END IF
250*
251* Compute the remaining r.h.s.
252*
253 temp = -rhs( j )
254 CALL saxpy( n-j, temp, z( j+1, j ), 1, rhs( j+1 ), 1 )
255*
256 10 CONTINUE
257*
258* Solve for U-part, look-ahead for RHS(N) = +-1. This is not done
259* in BSOLVE and will hopefully give us a better estimate because
260* any ill-conditioning of the original matrix is transferred to U
261* and not to L. U(N, N) is an approximation to sigma_min(LU).
262*
263 CALL scopy( n-1, rhs, 1, xp, 1 )
264 xp( n ) = rhs( n ) + one
265 rhs( n ) = rhs( n ) - one
266 splus = zero
267 sminu = zero
268 DO 30 i = n, 1, -1
269 temp = one / z( i, i )
270 xp( i ) = xp( i )*temp
271 rhs( i ) = rhs( i )*temp
272 DO 20 k = i + 1, n
273 xp( i ) = xp( i ) - xp( k )*( z( i, k )*temp )
274 rhs( i ) = rhs( i ) - rhs( k )*( z( i, k )*temp )
275 20 CONTINUE
276 splus = splus + abs( xp( i ) )
277 sminu = sminu + abs( rhs( i ) )
278 30 CONTINUE
279 IF( splus.GT.sminu )
280 $ CALL scopy( n, xp, 1, rhs, 1 )
281*
282* Apply the permutations JPIV to the computed solution (RHS)
283*
284 CALL slaswp( 1, rhs, ldz, 1, n-1, jpiv, -1 )
285*
286* Compute the sum of squares
287*
288 CALL slassq( n, rhs, 1, rdscal, rdsum )
289*
290 ELSE
291*
292* IJOB = 2, Compute approximate nullvector XM of Z
293*
294 CALL sgecon( 'I', n, z, ldz, one, temp, work, iwork, info )
295 CALL scopy( n, work( n+1 ), 1, xm, 1 )
296*
297* Compute RHS
298*
299 CALL slaswp( 1, xm, ldz, 1, n-1, ipiv, -1 )
300 temp = one / sqrt( sdot( n, xm, 1, xm, 1 ) )
301 CALL sscal( n, temp, xm, 1 )
302 CALL scopy( n, xm, 1, xp, 1 )
303 CALL saxpy( n, one, rhs, 1, xp, 1 )
304 CALL saxpy( n, -one, xm, 1, rhs, 1 )
305 CALL sgesc2( n, z, ldz, rhs, ipiv, jpiv, temp )
306 CALL sgesc2( n, z, ldz, xp, ipiv, jpiv, temp )
307 IF( sasum( n, xp, 1 ).GT.sasum( n, rhs, 1 ) )
308 $ CALL scopy( n, xp, 1, rhs, 1 )
309*
310* Compute the sum of squares
311*
312 CALL slassq( n, rhs, 1, rdscal, rdsum )
313*
314 END IF
315*
316 RETURN
317*
318* End of SLATDF
319*
real function sasum(n, sx, incx)
SASUM
Definition sasum.f:72
subroutine saxpy(n, sa, sx, incx, sy, incy)
SAXPY
Definition saxpy.f:89
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
real function sdot(n, sx, incx, sy, incy)
SDOT
Definition sdot.f:82
subroutine sgecon(norm, n, a, lda, anorm, rcond, work, iwork, info)
SGECON
Definition sgecon.f:130
subroutine sgesc2(n, a, lda, rhs, ipiv, jpiv, scale)
SGESC2 solves a system of linear equations using the LU factorization with complete pivoting computed...
Definition sgesc2.f:112
subroutine slassq(n, x, incx, scale, sumsq)
SLASSQ updates a sum of squares represented in scaled form.
Definition slassq.f90:122
subroutine slaswp(n, a, lda, k1, k2, ipiv, incx)
SLASWP performs a series of row interchanges on a general rectangular matrix.
Definition slaswp.f:113
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
Here is the call graph for this function:
Here is the caller graph for this function: