LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dlatdf ( integer  IJOB,
integer  N,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
double precision, dimension( * )  RHS,
double precision  RDSUM,
double precision  RDSCAL,
integer, dimension( * )  IPIV,
integer, dimension( * )  JPIV 
)

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

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

Purpose:
 DLATDF uses the LU factorization of the n-by-n matrix Z computed by
 DGETC2 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 DGETC2 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 DGECON, 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 DOUBLE PRECISION array, dimension (LDZ, N)
          On entry, the LU part of the factorization of the n-by-n
          matrix Z computed by DGETC2:  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 DOUBLE PRECISION array, dimension (N)
          On entry, RHS contains contributions from other subsystems.
          On exit, RHS contains the solution of the subsystem with
          entries acoording to the value of IJOB (see above).
[in,out]RDSUM
          RDSUM is DOUBLE PRECISION
          On entry, the sum of squares of computed contributions to
          the Dif-estimate under computation by DTGSYL, 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 DTGSY2 is called by STGSYL.
[in,out]RDSCAL
          RDSCAL is DOUBLE PRECISION
          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 DTGSY2 is called by
                DTGSYL.
[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.
Date
June 2016
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 173 of file dlatdf.f.

173 *
174 * -- LAPACK auxiliary routine (version 3.6.1) --
175 * -- LAPACK is a software package provided by Univ. of Tennessee, --
176 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
177 * June 2016
178 *
179 * .. Scalar Arguments ..
180  INTEGER ijob, ldz, n
181  DOUBLE PRECISION rdscal, rdsum
182 * ..
183 * .. Array Arguments ..
184  INTEGER ipiv( * ), jpiv( * )
185  DOUBLE PRECISION rhs( * ), z( ldz, * )
186 * ..
187 *
188 * =====================================================================
189 *
190 * .. Parameters ..
191  INTEGER maxdim
192  parameter ( maxdim = 8 )
193  DOUBLE PRECISION zero, one
194  parameter ( zero = 0.0d+0, one = 1.0d+0 )
195 * ..
196 * .. Local Scalars ..
197  INTEGER i, info, j, k
198  DOUBLE PRECISION bm, bp, pmone, sminu, splus, temp
199 * ..
200 * .. Local Arrays ..
201  INTEGER iwork( maxdim )
202  DOUBLE PRECISION work( 4*maxdim ), xm( maxdim ), xp( maxdim )
203 * ..
204 * .. External Subroutines ..
205  EXTERNAL daxpy, dcopy, dgecon, dgesc2, dlassq, dlaswp,
206  $ dscal
207 * ..
208 * .. External Functions ..
209  DOUBLE PRECISION dasum, ddot
210  EXTERNAL dasum, ddot
211 * ..
212 * .. Intrinsic Functions ..
213  INTRINSIC abs, sqrt
214 * ..
215 * .. Executable Statements ..
216 *
217  IF( ijob.NE.2 ) THEN
218 *
219 * Apply permutations IPIV to RHS
220 *
221  CALL dlaswp( 1, rhs, ldz, 1, n-1, ipiv, 1 )
222 *
223 * Solve for L-part choosing RHS either to +1 or -1.
224 *
225  pmone = -one
226 *
227  DO 10 j = 1, n - 1
228  bp = rhs( j ) + one
229  bm = rhs( j ) - one
230  splus = one
231 *
232 * Look-ahead for L-part RHS(1:N-1) = + or -1, SPLUS and
233 * SMIN computed more efficiently than in BSOLVE [1].
234 *
235  splus = splus + ddot( n-j, z( j+1, j ), 1, z( j+1, j ), 1 )
236  sminu = ddot( n-j, z( j+1, j ), 1, rhs( j+1 ), 1 )
237  splus = splus*rhs( j )
238  IF( splus.GT.sminu ) THEN
239  rhs( j ) = bp
240  ELSE IF( sminu.GT.splus ) THEN
241  rhs( j ) = bm
242  ELSE
243 *
244 * In this case the updating sums are equal and we can
245 * choose RHS(J) +1 or -1. The first time this happens
246 * we choose -1, thereafter +1. This is a simple way to
247 * get good estimates of matrices like Byers well-known
248 * example (see [1]). (Not done in BSOLVE.)
249 *
250  rhs( j ) = rhs( j ) + pmone
251  pmone = one
252  END IF
253 *
254 * Compute the remaining r.h.s.
255 *
256  temp = -rhs( j )
257  CALL daxpy( n-j, temp, z( j+1, j ), 1, rhs( j+1 ), 1 )
258 *
259  10 CONTINUE
260 *
261 * Solve for U-part, look-ahead for RHS(N) = +-1. This is not done
262 * in BSOLVE and will hopefully give us a better estimate because
263 * any ill-conditioning of the original matrix is transfered to U
264 * and not to L. U(N, N) is an approximation to sigma_min(LU).
265 *
266  CALL dcopy( n-1, rhs, 1, xp, 1 )
267  xp( n ) = rhs( n ) + one
268  rhs( n ) = rhs( n ) - one
269  splus = zero
270  sminu = zero
271  DO 30 i = n, 1, -1
272  temp = one / z( i, i )
273  xp( i ) = xp( i )*temp
274  rhs( i ) = rhs( i )*temp
275  DO 20 k = i + 1, n
276  xp( i ) = xp( i ) - xp( k )*( z( i, k )*temp )
277  rhs( i ) = rhs( i ) - rhs( k )*( z( i, k )*temp )
278  20 CONTINUE
279  splus = splus + abs( xp( i ) )
280  sminu = sminu + abs( rhs( i ) )
281  30 CONTINUE
282  IF( splus.GT.sminu )
283  $ CALL dcopy( n, xp, 1, rhs, 1 )
284 *
285 * Apply the permutations JPIV to the computed solution (RHS)
286 *
287  CALL dlaswp( 1, rhs, ldz, 1, n-1, jpiv, -1 )
288 *
289 * Compute the sum of squares
290 *
291  CALL dlassq( n, rhs, 1, rdscal, rdsum )
292 *
293  ELSE
294 *
295 * IJOB = 2, Compute approximate nullvector XM of Z
296 *
297  CALL dgecon( 'I', n, z, ldz, one, temp, work, iwork, info )
298  CALL dcopy( n, work( n+1 ), 1, xm, 1 )
299 *
300 * Compute RHS
301 *
302  CALL dlaswp( 1, xm, ldz, 1, n-1, ipiv, -1 )
303  temp = one / sqrt( ddot( n, xm, 1, xm, 1 ) )
304  CALL dscal( n, temp, xm, 1 )
305  CALL dcopy( n, xm, 1, xp, 1 )
306  CALL daxpy( n, one, rhs, 1, xp, 1 )
307  CALL daxpy( n, -one, xm, 1, rhs, 1 )
308  CALL dgesc2( n, z, ldz, rhs, ipiv, jpiv, temp )
309  CALL dgesc2( n, z, ldz, xp, ipiv, jpiv, temp )
310  IF( dasum( n, xp, 1 ).GT.dasum( n, rhs, 1 ) )
311  $ CALL dcopy( n, xp, 1, rhs, 1 )
312 *
313 * Compute the sum of squares
314 *
315  CALL dlassq( n, rhs, 1, rdscal, rdsum )
316 *
317  END IF
318 *
319  RETURN
320 *
321 * End of DLATDF
322 *
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
double precision function ddot(N, DX, INCX, DY, INCY)
DDOT
Definition: ddot.f:53
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine dgesc2(N, A, LDA, RHS, IPIV, JPIV, SCALE)
DGESC2 solves a system of linear equations using the LU factorization with complete pivoting computed...
Definition: dgesc2.f:116
subroutine dlaswp(N, A, LDA, K1, K2, IPIV, INCX)
DLASWP performs a series of row interchanges on a general rectangular matrix.
Definition: dlaswp.f:116
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
double precision function dasum(N, DX, INCX)
DASUM
Definition: dasum.f:53
subroutine dlassq(N, X, INCX, SCALE, SUMSQ)
DLASSQ updates a sum of squares represented in scaled form.
Definition: dlassq.f:105
subroutine dgecon(NORM, N, A, LDA, ANORM, RCOND, WORK, IWORK, INFO)
DGECON
Definition: dgecon.f:126

Here is the call graph for this function:

Here is the caller graph for this function: