File: MB03PY.html

package info (click to toggle)
slicot 5.9.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 23,528 kB
  • sloc: fortran: 148,076; makefile: 964; sh: 57
file content (211 lines) | stat: -rw-r--r-- 7,563 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
<HTML>
<HEAD><TITLE>MB03PY - SLICOT Library Routine Documentation</TITLE>
</HEAD>
<BODY>

<H2><A Name="MB03PY">MB03PY</A></H2>
<H3>
Matrix rank determination by incremental condition estimation, during the pivoted RQ factorization process (row pivoting)
</H3>
<A HREF ="#Specification"><B>[Specification]</B></A>
<A HREF ="#Arguments"><B>[Arguments]</B></A>
<A HREF ="#Method"><B>[Method]</B></A>
<A HREF ="#References"><B>[References]</B></A>
<A HREF ="#Comments"><B>[Comments]</B></A>
<A HREF ="#Example"><B>[Example]</B></A>

<P>
<B><FONT SIZE="+1">Purpose</FONT></B>
<PRE>
  To compute a rank-revealing RQ factorization of a real general
  M-by-N matrix  A,  which may be rank-deficient, and estimate its
  effective rank using incremental condition estimation.

  The routine uses a truncated RQ factorization with row pivoting:
                                [ R11 R12 ]
     P * A = R * Q,  where  R = [         ],
                                [  0  R22 ]
  with R22 defined as the largest trailing upper triangular
  submatrix whose estimated condition number is less than 1/RCOND.
  The order of R22, RANK, is the effective rank of A.  Condition
  estimation is performed during the RQ factorization process.
  Matrix R11 is full (but of small norm), or empty.

  MB03PY  does not perform any scaling of the matrix A.

</PRE>
<A name="Specification"><B><FONT SIZE="+1">Specification</FONT></B></A>
<PRE>
      SUBROUTINE MB03PY( M, N, A, LDA, RCOND, SVLMAX, RANK, SVAL, JPVT,
     $                   TAU, DWORK, INFO )
C     .. Scalar Arguments ..
      INTEGER            INFO, LDA, M, N, RANK
      DOUBLE PRECISION   RCOND, SVLMAX
C     .. Array Arguments ..
      INTEGER            JPVT( * )
      DOUBLE PRECISION   A( LDA, * ), DWORK( * ), SVAL( 3 ), TAU( * )

</PRE>
<A name="Arguments"><B><FONT SIZE="+1">Arguments</FONT></B></A>
<P>

</PRE>
<B>Input/Output Parameters</B>
<PRE>
  M       (input) INTEGER
          The number of rows of the matrix A.  M &gt;= 0.

  N       (input) INTEGER
          The number of columns of the matrix A.  N &gt;= 0.

  A       (input/output) DOUBLE PRECISION array, dimension
          ( LDA, N )
          On entry, the leading M-by-N part of this array must
          contain the given matrix A.
          On exit, the upper triangle of the subarray
          A(M-RANK+1:M,N-RANK+1:N) contains the RANK-by-RANK upper
          triangular matrix R22;  the remaining elements in the last
          RANK  rows, with the array TAU, represent the orthogonal
          matrix Q as a product of  RANK  elementary reflectors
          (see METHOD).  The first  M-RANK  rows contain the result
          of the RQ factorization process used.

  LDA     INTEGER
          The leading dimension of the array A.  LDA &gt;= max(1,M).

  RCOND   (input) DOUBLE PRECISION
          RCOND is used to determine the effective rank of A, which
          is defined as the order of the largest trailing triangular
          submatrix R22 in the RQ factorization with pivoting of A,
          whose estimated condition number is less than 1/RCOND.
          0 &lt;= RCOND &lt;= 1.
          NOTE that when SVLMAX &gt; 0, the estimated rank could be
          less than that defined above (see SVLMAX).

  SVLMAX  (input) DOUBLE PRECISION
          If A is a submatrix of another matrix B, and the rank
          decision should be related to that matrix, then SVLMAX
          should be an estimate of the largest singular value of B
          (for instance, the Frobenius norm of B).  If this is not
          the case, the input value SVLMAX = 0 should work.
          SVLMAX &gt;= 0.

  RANK    (output) INTEGER
          The effective (estimated) rank of A, i.e., the order of
          the submatrix R22.

  SVAL    (output) DOUBLE PRECISION array, dimension ( 3 )
          The estimates of some of the singular values of the
          triangular factor R:
          SVAL(1): largest singular value of
                   R(M-RANK+1:M,N-RANK+1:N);
          SVAL(2): smallest singular value of
                   R(M-RANK+1:M,N-RANK+1:N);
          SVAL(3): smallest singular value of R(M-RANK:M,N-RANK:N),
                   if RANK &lt; MIN( M, N ), or of
                   R(M-RANK+1:M,N-RANK+1:N), otherwise.
          If the triangular factorization is a rank-revealing one
          (which will be the case if the trailing rows were well-
          conditioned), then SVAL(1) will also be an estimate for
          the largest singular value of A, and SVAL(2) and SVAL(3)
          will be estimates for the RANK-th and (RANK+1)-st singular
          values of A, respectively.
          By examining these values, one can confirm that the rank
          is well defined with respect to the chosen value of RCOND.
          The ratio SVAL(1)/SVAL(2) is an estimate of the condition
          number of R(M-RANK+1:M,N-RANK+1:N).

  JPVT    (output) INTEGER array, dimension ( M )
          If JPVT(i) = k, then the i-th row of P*A was the k-th row
          of A.

  TAU     (output) DOUBLE PRECISION array, dimension ( MIN( M, N ) )
          The trailing  RANK  elements of TAU contain the scalar
          factors of the elementary reflectors.

</PRE>
<B>Workspace</B>
<PRE>
  DWORK   DOUBLE PRECISION array, dimension ( 3*M-1 )

</PRE>
<B>Error Indicator</B>
<PRE>
  INFO    INTEGER
          = 0:  successful exit;
          &lt; 0:  if INFO = -i, the i-th argument had an illegal
                value.

</PRE>
<A name="Method"><B><FONT SIZE="+1">Method</FONT></B></A>
<PRE>
  The routine computes a truncated RQ factorization with row
  pivoting of A,  P * A = R * Q,  with  R  defined above, and,
  during this process, finds the largest trailing submatrix whose
  estimated condition number is less than 1/RCOND, taking the
  possible positive value of SVLMAX into account.  This is performed
  using an adaptation of the LAPACK incremental condition estimation
  scheme and a slightly modified rank decision test.  The
  factorization process stops when  RANK  has been determined.

  The matrix Q is represented as a product of elementary reflectors

     Q = H(k-rank+1) H(k-rank+2) . . . H(k), where k = min(m,n).

  Each H(i) has the form

     H = I - tau * v * v'

  where tau is a real scalar, and v is a real vector with
  v(n-k+i+1:n) = 0 and v(n-k+i) = 1; v(1:n-k+i-1) is stored on exit
  in A(m-k+i,1:n-k+i-1), and tau in TAU(i).

  The matrix P is represented in jpvt as follows: If
     jpvt(j) = i
  then the jth row of P is the ith canonical unit vector.

</PRE>
<A name="References"><B><FONT SIZE="+1">References</FONT></B></A>
<PRE>
  [1] Bischof, C.H. and P. Tang.
      Generalizing Incremental Condition Estimation.
      LAPACK Working Notes 32, Mathematics and Computer Science
      Division, Argonne National Laboratory, UT, CS-91-132,
      May 1991.

  [2] Bischof, C.H. and P. Tang.
      Robust Incremental Condition Estimation.
      LAPACK Working Notes 33, Mathematics and Computer Science
      Division, Argonne National Laboratory, UT, CS-91-133,
      May 1991.

</PRE>
<A name="Numerical Aspects"><B><FONT SIZE="+1">Numerical Aspects</FONT></B></A>
<PRE>
  The algorithm is backward stable.

</PRE>

<A name="Comments"><B><FONT SIZE="+1">Further Comments</FONT></B></A>
<PRE>
  None
</PRE>

<A name="Example"><B><FONT SIZE="+1">Example</FONT></B></A>
<P>
<B>Program Text</B>
<PRE>
  None
</PRE>
<B>Program Data</B>
<PRE>
  None
</PRE>
<B>Program Results</B>
<PRE>
  None
</PRE>

<HR>
<A HREF=support.html><B>Return to Supporting Routines index</B></A></BODY>
</HTML>