File: MB03OD.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 (288 lines) | stat: -rw-r--r-- 10,579 bytes parent folder | download | duplicates (2)
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
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
<HTML>
<HEAD><TITLE>MB03OD - SLICOT Library Routine Documentation</TITLE>
</HEAD>
<BODY>

<H2><A Name="MB03OD">MB03OD</A></H2>
<H3>
Matrix rank determination by incremental condition estimation
</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 (optionally) a rank-revealing QR 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 QR factorization with column pivoting:
     A * P = Q * R,  where  R = [ R11 R12 ],
                                [  0  R22 ]
  with R11 defined as the largest leading submatrix whose estimated
  condition number is less than 1/RCOND.  The order of R11, RANK,
  is the effective rank of A.

  MB03OD  does not perform any scaling of the matrix A.

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

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

<B>Mode Parameters</B>
<PRE>
  JOBQR   CHARACTER*1
          = 'Q':  Perform a QR factorization with column pivoting;
          = 'N':  Do not perform the QR factorization (but assume
                  that it has been done outside).

</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 with JOBQR = 'Q', the leading M by N part of this
          array must contain the given matrix A.
          On exit with JOBQR = 'Q', the leading min(M,N) by N upper
          triangular part of A contains the triangular factor R,
          and the elements below the diagonal, with the array TAU,
          represent the orthogonal matrix Q as a product of
          min(M,N) elementary reflectors.
          On entry and on exit with JOBQR = 'N', the leading
          min(M,N) by N upper triangular part of A contains the
          triangular factor R, as determined by the QR factorization
          with pivoting.  The elements below the diagonal of A are
          not referenced.

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

  JPVT    (input/output) INTEGER array, dimension ( N )
          On entry with JOBQR = 'Q', if JPVT(i) &lt;&gt; 0, the i-th
          column of A is an initial column, otherwise it is a free
          column. Before the QR factorization of A, all initial
          columns are permuted to the leading positions; only the
          remaining free columns are moved as a result of column
          pivoting during the factorization.  For rank determination
          it is preferable that all columns be free.
          On exit with JOBQR = 'Q', if JPVT(i) = k, then the i-th
          column of A*P was the k-th column of A.
          Array JPVT is not referenced when JOBQR = 'N'.

  RCOND   (input) DOUBLE PRECISION
          RCOND is used to determine the effective rank of A, which
          is defined as the order of the largest leading triangular
          submatrix R11 in the QR factorization with pivoting of A,
          whose estimated condition number is less than 1/RCOND.
          RCOND &gt;= 0.
          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.

  TAU     (output) DOUBLE PRECISION array, dimension ( MIN( M, N ) )
          On exit with JOBQR = 'Q', the leading min(M,N) elements of
          TAU contain the scalar factors of the elementary
          reflectors.
          Array TAU is not referenced when JOBQR = 'N'.

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

  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(1:RANK,1:RANK);
          SVAL(2): smallest singular value of R(1:RANK,1:RANK);
          SVAL(3): smallest singular value of R(1:RANK+1,1:RANK+1),
                   if RANK &lt; MIN( M, N ), or of R(1:RANK,1:RANK),
                   otherwise.
          If the triangular factorization is a rank-revealing one
          (which will be the case if the leading columns 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(1:RANK,1:RANK).

</PRE>
<B>Workspace</B>
<PRE>
  DWORK   DOUBLE PRECISION array, dimension ( LDWORK )
          On exit, if  INFO = 0,  DWORK(1) returns the optimal value
          of LDWORK.

  LDWORK  INTEGER
          The length of the array DWORK.
          LDWORK &gt;= 3*N + 1,                 if JOBQR = 'Q';
          LDWORK &gt;= max( 1, 2*min( M, N ) ), if JOBQR = 'N'.
          For good performance when JOBQR = 'Q', LDWORK should be
          larger. Specifically, LDWORK &gt;= 2*N + ( N + 1 )*NB, where
          NB is the optimal block size for the LAPACK Library
          routine DGEQP3.

          If LDWORK = -1, then a workspace query is assumed;
          the routine only calculates the optimal size of the
          DWORK array, returns this value as the first entry of
          the DWORK array, and no error message related to LDWORK
          is issued by XERBLA.

</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 or uses a QR factorization with column
  pivoting of A,  A * P = Q * R,  with  R  defined above, and then
  finds the largest leading submatrix whose estimated condition
  number is less than 1/RCOND, taking the possible positive value of
  SVLMAX into account.  This is performed using the LAPACK
  incremental condition estimation scheme and a slightly modified
  rank decision test.

</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>
*     MB03OD EXAMPLE PROGRAM TEXT
*
*     .. Parameters ..
      DOUBLE PRECISION ZERO, ONE
      PARAMETER        ( ZERO = 0.0D0, ONE = 1.0D0 )
      INTEGER          NIN, NOUT
      PARAMETER        ( NIN = 5, NOUT = 6 )
      INTEGER          NMAX, MMAX
      PARAMETER        ( NMAX = 10, MMAX = 10 )
      INTEGER          LDA
      PARAMETER        ( LDA = NMAX )
      INTEGER          LDTAU
      PARAMETER        ( LDTAU = MIN(MMAX,NMAX) )
      INTEGER          LDWORK
      PARAMETER        ( LDWORK = 3*NMAX + 1 )
*     .. Local Scalars ..
      CHARACTER*1      JOBQR
      INTEGER          I, INFO, J, M, N, RANK
      DOUBLE PRECISION RCOND, SVAL(3), SVLMAX
*     ..
*     .. Local Arrays ..
      DOUBLE PRECISION A(LDA,NMAX), DWORK(LDWORK), TAU(LDTAU)
      INTEGER          JPVT(NMAX)
*     .. External Subroutines ..
      EXTERNAL         MB03OD
*     .. Intrinsic Functions ..
      INTRINSIC        MIN
*     .. Executable Statements ..
*
      WRITE ( NOUT, FMT = 99999 )
*     Skip the heading in the data file and read the data.
      READ ( NIN, FMT = '()' )
      READ ( NIN, FMT = * ) M, N, JOBQR, RCOND, SVLMAX
      IF ( N.LT.0 .OR. N.GT.NMAX ) THEN
         WRITE ( NOUT, FMT = 99972 ) N
      ELSE
         IF ( M.LT.0 .OR. M.GT.MMAX ) THEN
            WRITE ( NOUT, FMT = 99971 ) M
         ELSE
            READ ( NIN, FMT = * ) ( ( A(I,J), J = 1,N ), I = 1,M )
*           QR with column pivoting.
            DO 10 I = 1, N
               JPVT(I) = 0
   10       CONTINUE
            CALL MB03OD( JOBQR, M, N, A, LDA, JPVT, RCOND, SVLMAX, TAU,
     $                   RANK, SVAL, DWORK, LDWORK, INFO )
*
            IF ( INFO.NE.0 ) THEN
               WRITE ( NOUT, FMT = 99998 ) INFO
            ELSE
               WRITE ( NOUT, FMT = 99995 ) RANK
               WRITE ( NOUT, FMT = 99994 ) ( JPVT(I), I = 1,N )
               WRITE ( NOUT, FMT = 99993 ) ( SVAL(I), I = 1,3 )
            END IF
         END IF
      END IF
*
      STOP
*
99999 FORMAT (' MB03OD EXAMPLE PROGRAM RESULTS',/1X)
99998 FORMAT (' INFO on exit from MB03OD = ',I2)
99995 FORMAT (' The rank is ',I5)
99994 FORMAT (' Column permutations are ',/(20(I3,2X)))
99993 FORMAT (' SVAL vector is ',/(20(1X,F10.4)))
99972 FORMAT (/' N is out of range.',/' N = ',I5)
99971 FORMAT (/' M is out of range.',/' M = ',I5)
      END
</PRE>
<B>Program Data</B>
<PRE>
 MB03OD EXAMPLE PROGRAM DATA
   6     5     Q  5.D-16     0.0
   1.    2.    6.    3.    5.
  -2.   -1.   -1.    0.   -2.
   5.    5.    1.    5.    1.
  -2.   -1.   -1.    0.   -2.
   4.    8.    4.   20.    4.
  -2.   -1.   -1.    0.   -2.
</PRE>
<B>Program Results</B>
<PRE>
 MB03OD EXAMPLE PROGRAM RESULTS

 The rank is     4
 Column permutations are 
  4    3    1    5    2
 SVAL vector is 
    22.7257     1.4330     0.0000
</PRE>

<HR>
<p>
<A HREF=..\libindex.html><B>Return to index</B></A></BODY>
</HTML>