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

◆ clatdf()

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

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

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

Purpose:
 CLATDF computes the contribution to the reciprocal Dif-estimate
 by solving for x in Z * x = b, where b is chosen such that the norm
 of x is as large as possible. It is assumed that LU decomposition
 of Z has been computed by CGETC2. On entry RHS = f holds the
 contribution from earlier solved sub-systems, and on return RHS = x.

 The factorization of Z returned by CGETC2 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 CGECON, 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 COMPLEX array, dimension (LDZ, N)
          On entry, the LU part of the factorization of the n-by-n
          matrix Z computed by CGETC2:  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 COMPLEX 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 CTGSYL, 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 CTGSY2 is called by CTGSYL.
[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 CTGSY2 is called by
          CTGSYL.
[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 UMINF-95.05, Department of Computing Science, Umea University, S-901 87 Umea, Sweden, 1995.

Definition at line 167 of file clatdf.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 COMPLEX RHS( * ), Z( LDZ, * )
181* ..
182*
183* =====================================================================
184*
185* .. Parameters ..
186 INTEGER MAXDIM
187 parameter( maxdim = 2 )
188 REAL ZERO, ONE
189 parameter( zero = 0.0e+0, one = 1.0e+0 )
190 COMPLEX CONE
191 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
192* ..
193* .. Local Scalars ..
194 INTEGER I, INFO, J, K
195 REAL RTEMP, SCALE, SMINU, SPLUS
196 COMPLEX BM, BP, PMONE, TEMP
197* ..
198* .. Local Arrays ..
199 REAL RWORK( MAXDIM )
200 COMPLEX WORK( 4*MAXDIM ), XM( MAXDIM ), XP( MAXDIM )
201* ..
202* .. External Subroutines ..
203 EXTERNAL caxpy, ccopy, cgecon, cgesc2, classq, claswp,
204 $ cscal
205* ..
206* .. External Functions ..
207 REAL SCASUM
208 COMPLEX CDOTC
209 EXTERNAL scasum, cdotc
210* ..
211* .. Intrinsic Functions ..
212 INTRINSIC abs, real, sqrt
213* ..
214* .. Executable Statements ..
215*
216 IF( ijob.NE.2 ) THEN
217*
218* Apply permutations IPIV to RHS
219*
220 CALL claswp( 1, rhs, ldz, 1, n-1, ipiv, 1 )
221*
222* Solve for L-part choosing RHS either to +1 or -1.
223*
224 pmone = -cone
225 DO 10 j = 1, n - 1
226 bp = rhs( j ) + cone
227 bm = rhs( j ) - cone
228 splus = one
229*
230* Look-ahead for L- part RHS(1:N-1) = +-1
231* SPLUS and SMIN computed more efficiently than in BSOLVE[1].
232*
233 splus = splus + real( cdotc( n-j, z( j+1, j ), 1, z( j+1,
234 $ j ), 1 ) )
235 sminu = real( cdotc( n-j, z( j+1, j ), 1, rhs( j+1 ), 1 ) )
236 splus = splus*real( rhs( j ) )
237 IF( splus.GT.sminu ) THEN
238 rhs( j ) = bp
239 ELSE IF( sminu.GT.splus ) THEN
240 rhs( j ) = bm
241 ELSE
242*
243* In this case the updating sums are equal and we can
244* choose RHS(J) +1 or -1. The first time this happens we
245* choose -1, thereafter +1. This is a simple way to get
246* good estimates of matrices like Byers well-known example
247* (see [1]). (Not done in BSOLVE.)
248*
249 rhs( j ) = rhs( j ) + pmone
250 pmone = cone
251 END IF
252*
253* Compute the remaining r.h.s.
254*
255 temp = -rhs( j )
256 CALL caxpy( n-j, temp, z( j+1, j ), 1, rhs( j+1 ), 1 )
257 10 CONTINUE
258*
259* Solve for U- part, lockahead for RHS(N) = +-1. This is not done
260* In BSOLVE and will hopefully give us a better estimate because
261* any ill-conditioning of the original matrix is transferred to U
262* and not to L. U(N, N) is an approximation to sigma_min(LU).
263*
264 CALL ccopy( n-1, rhs, 1, work, 1 )
265 work( n ) = rhs( n ) + cone
266 rhs( n ) = rhs( n ) - cone
267 splus = zero
268 sminu = zero
269 DO 30 i = n, 1, -1
270 temp = cone / z( i, i )
271 work( i ) = work( i )*temp
272 rhs( i ) = rhs( i )*temp
273 DO 20 k = i + 1, n
274 work( i ) = work( i ) - work( k )*( z( i, k )*temp )
275 rhs( i ) = rhs( i ) - rhs( k )*( z( i, k )*temp )
276 20 CONTINUE
277 splus = splus + abs( work( i ) )
278 sminu = sminu + abs( rhs( i ) )
279 30 CONTINUE
280 IF( splus.GT.sminu )
281 $ CALL ccopy( n, work, 1, rhs, 1 )
282*
283* Apply the permutations JPIV to the computed solution (RHS)
284*
285 CALL claswp( 1, rhs, ldz, 1, n-1, jpiv, -1 )
286*
287* Compute the sum of squares
288*
289 CALL classq( n, rhs, 1, rdscal, rdsum )
290 RETURN
291 END IF
292*
293* ENTRY IJOB = 2
294*
295* Compute approximate nullvector XM of Z
296*
297 CALL cgecon( 'I', n, z, ldz, one, rtemp, work, rwork, info )
298 CALL ccopy( n, work( n+1 ), 1, xm, 1 )
299*
300* Compute RHS
301*
302 CALL claswp( 1, xm, ldz, 1, n-1, ipiv, -1 )
303 temp = cone / sqrt( cdotc( n, xm, 1, xm, 1 ) )
304 CALL cscal( n, temp, xm, 1 )
305 CALL ccopy( n, xm, 1, xp, 1 )
306 CALL caxpy( n, cone, rhs, 1, xp, 1 )
307 CALL caxpy( n, -cone, xm, 1, rhs, 1 )
308 CALL cgesc2( n, z, ldz, rhs, ipiv, jpiv, scale )
309 CALL cgesc2( n, z, ldz, xp, ipiv, jpiv, scale )
310 IF( scasum( n, xp, 1 ).GT.scasum( n, rhs, 1 ) )
311 $ CALL ccopy( n, xp, 1, rhs, 1 )
312*
313* Compute the sum of squares
314*
315 CALL classq( n, rhs, 1, rdscal, rdsum )
316 RETURN
317*
318* End of CLATDF
319*
real function scasum(n, cx, incx)
SCASUM
Definition scasum.f:72
subroutine caxpy(n, ca, cx, incx, cy, incy)
CAXPY
Definition caxpy.f:88
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
complex function cdotc(n, cx, incx, cy, incy)
CDOTC
Definition cdotc.f:83
subroutine cgecon(norm, n, a, lda, anorm, rcond, work, rwork, info)
CGECON
Definition cgecon.f:132
subroutine cgesc2(n, a, lda, rhs, ipiv, jpiv, scale)
CGESC2 solves a system of linear equations using the LU factorization with complete pivoting computed...
Definition cgesc2.f:115
subroutine classq(n, x, incx, scale, sumsq)
CLASSQ updates a sum of squares represented in scaled form.
Definition classq.f90:124
subroutine claswp(n, a, lda, k1, k2, ipiv, incx)
CLASWP performs a series of row interchanges on a general rectangular matrix.
Definition claswp.f:115
subroutine cscal(n, ca, cx, incx)
CSCAL
Definition cscal.f:78
Here is the call graph for this function:
Here is the caller graph for this function: