You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

dgeqp3rk.f 39 kB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082
  1. *> \brief \b DGEQP3RK computes a truncated Householder QR factorization with column pivoting of a real m-by-n matrix A by using Level 3 BLAS and overwrites a real m-by-nrhs matrix B with Q**T * B.
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download DGEQP3RK + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgeqp3rk.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgeqp3rk.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgeqp3rk.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE DGEQP3RK( M, N, NRHS, KMAX, ABSTOL, RELTOL, A, LDA,
  22. * $ K, MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
  23. * $ WORK, LWORK, IWORK, INFO )
  24. * IMPLICIT NONE
  25. *
  26. * .. Scalar Arguments ..
  27. * INTEGER INFO, K, KMAX, LDA, LWORK, M, N, NRHS
  28. * DOUBLE PRECISION ABSTOL, MAXC2NRMK, RELMAXC2NRMK, RELTOL
  29. * ..
  30. * .. Array Arguments ..
  31. * INTEGER IWORK( * ), JPIV( * )
  32. * DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
  33. * ..
  34. *
  35. *
  36. *> \par Purpose:
  37. * =============
  38. *>
  39. *> \verbatim
  40. *>
  41. *> DGEQP3RK performs two tasks simultaneously:
  42. *>
  43. *> Task 1: The routine computes a truncated (rank K) or full rank
  44. *> Householder QR factorization with column pivoting of a real
  45. *> M-by-N matrix A using Level 3 BLAS. K is the number of columns
  46. *> that were factorized, i.e. factorization rank of the
  47. *> factor R, K <= min(M,N).
  48. *>
  49. *> A * P(K) = Q(K) * R(K) =
  50. *>
  51. *> = Q(K) * ( R11(K) R12(K) ) = Q(K) * ( R(K)_approx )
  52. *> ( 0 R22(K) ) ( 0 R(K)_residual ),
  53. *>
  54. *> where:
  55. *>
  56. *> P(K) is an N-by-N permutation matrix;
  57. *> Q(K) is an M-by-M orthogonal matrix;
  58. *> R(K)_approx = ( R11(K), R12(K) ) is a rank K approximation of the
  59. *> full rank factor R with K-by-K upper-triangular
  60. *> R11(K) and K-by-N rectangular R12(K). The diagonal
  61. *> entries of R11(K) appear in non-increasing order
  62. *> of absolute value, and absolute values of all of
  63. *> them exceed the maximum column 2-norm of R22(K)
  64. *> up to roundoff error.
  65. *> R(K)_residual = R22(K) is the residual of a rank K approximation
  66. *> of the full rank factor R. It is a
  67. *> an (M-K)-by-(N-K) rectangular matrix;
  68. *> 0 is a an (M-K)-by-K zero matrix.
  69. *>
  70. *> Task 2: At the same time, the routine overwrites a real M-by-NRHS
  71. *> matrix B with Q(K)**T * B using Level 3 BLAS.
  72. *>
  73. *> =====================================================================
  74. *>
  75. *> The matrices A and B are stored on input in the array A as
  76. *> the left and right blocks A(1:M,1:N) and A(1:M, N+1:N+NRHS)
  77. *> respectively.
  78. *>
  79. *> N NRHS
  80. *> array_A = M [ mat_A, mat_B ]
  81. *>
  82. *> The truncation criteria (i.e. when to stop the factorization)
  83. *> can be any of the following:
  84. *>
  85. *> 1) The input parameter KMAX, the maximum number of columns
  86. *> KMAX to factorize, i.e. the factorization rank is limited
  87. *> to KMAX. If KMAX >= min(M,N), the criterion is not used.
  88. *>
  89. *> 2) The input parameter ABSTOL, the absolute tolerance for
  90. *> the maximum column 2-norm of the residual matrix R22(K). This
  91. *> means that the factorization stops if this norm is less or
  92. *> equal to ABSTOL. If ABSTOL < 0.0, the criterion is not used.
  93. *>
  94. *> 3) The input parameter RELTOL, the tolerance for the maximum
  95. *> column 2-norm matrix of the residual matrix R22(K) divided
  96. *> by the maximum column 2-norm of the original matrix A, which
  97. *> is equal to abs(R(1,1)). This means that the factorization stops
  98. *> when the ratio of the maximum column 2-norm of R22(K) to
  99. *> the maximum column 2-norm of A is less than or equal to RELTOL.
  100. *> If RELTOL < 0.0, the criterion is not used.
  101. *>
  102. *> 4) In case both stopping criteria ABSTOL or RELTOL are not used,
  103. *> and when the residual matrix R22(K) is a zero matrix in some
  104. *> factorization step K. ( This stopping criterion is implicit. )
  105. *>
  106. *> The algorithm stops when any of these conditions is first
  107. *> satisfied, otherwise the whole matrix A is factorized.
  108. *>
  109. *> To factorize the whole matrix A, use the values
  110. *> KMAX >= min(M,N), ABSTOL < 0.0 and RELTOL < 0.0.
  111. *>
  112. *> The routine returns:
  113. *> a) Q(K), R(K)_approx = ( R11(K), R12(K) ),
  114. *> R(K)_residual = R22(K), P(K), i.e. the resulting matrices
  115. *> of the factorization; P(K) is represented by JPIV,
  116. *> ( if K = min(M,N), R(K)_approx is the full factor R,
  117. *> and there is no residual matrix R(K)_residual);
  118. *> b) K, the number of columns that were factorized,
  119. *> i.e. factorization rank;
  120. *> c) MAXC2NRMK, the maximum column 2-norm of the residual
  121. *> matrix R(K)_residual = R22(K),
  122. *> ( if K = min(M,N), MAXC2NRMK = 0.0 );
  123. *> d) RELMAXC2NRMK equals MAXC2NRMK divided by MAXC2NRM, the maximum
  124. *> column 2-norm of the original matrix A, which is equal
  125. *> to abs(R(1,1)), ( if K = min(M,N), RELMAXC2NRMK = 0.0 );
  126. *> e) Q(K)**T * B, the matrix B with the orthogonal
  127. *> transformation Q(K)**T applied on the left.
  128. *>
  129. *> The N-by-N permutation matrix P(K) is stored in a compact form in
  130. *> the integer array JPIV. For 1 <= j <= N, column j
  131. *> of the matrix A was interchanged with column JPIV(j).
  132. *>
  133. *> The M-by-M orthogonal matrix Q is represented as a product
  134. *> of elementary Householder reflectors
  135. *>
  136. *> Q(K) = H(1) * H(2) * . . . * H(K),
  137. *>
  138. *> where K is the number of columns that were factorized.
  139. *>
  140. *> Each H(j) has the form
  141. *>
  142. *> H(j) = I - tau * v * v**T,
  143. *>
  144. *> where 1 <= j <= K and
  145. *> I is an M-by-M identity matrix,
  146. *> tau is a real scalar,
  147. *> v is a real vector with v(1:j-1) = 0 and v(j) = 1.
  148. *>
  149. *> v(j+1:M) is stored on exit in A(j+1:M,j) and tau in TAU(j).
  150. *>
  151. *> See the Further Details section for more information.
  152. *> \endverbatim
  153. *
  154. * Arguments:
  155. * ==========
  156. *
  157. *> \param[in] M
  158. *> \verbatim
  159. *> M is INTEGER
  160. *> The number of rows of the matrix A. M >= 0.
  161. *> \endverbatim
  162. *>
  163. *> \param[in] N
  164. *> \verbatim
  165. *> N is INTEGER
  166. *> The number of columns of the matrix A. N >= 0.
  167. *> \endverbatim
  168. *>
  169. *> \param[in] NRHS
  170. *> \verbatim
  171. *> NRHS is INTEGER
  172. *> The number of right hand sides, i.e. the number of
  173. *> columns of the matrix B. NRHS >= 0.
  174. *> \endverbatim
  175. *>
  176. *> \param[in] KMAX
  177. *> \verbatim
  178. *> KMAX is INTEGER
  179. *>
  180. *> The first factorization stopping criterion. KMAX >= 0.
  181. *>
  182. *> The maximum number of columns of the matrix A to factorize,
  183. *> i.e. the maximum factorization rank.
  184. *>
  185. *> a) If KMAX >= min(M,N), then this stopping criterion
  186. *> is not used, the routine factorizes columns
  187. *> depending on ABSTOL and RELTOL.
  188. *>
  189. *> b) If KMAX = 0, then this stopping criterion is
  190. *> satisfied on input and the routine exits immediately.
  191. *> This means that the factorization is not performed,
  192. *> the matrices A and B are not modified, and
  193. *> the matrix A is itself the residual.
  194. *> \endverbatim
  195. *>
  196. *> \param[in] ABSTOL
  197. *> \verbatim
  198. *> ABSTOL is DOUBLE PRECISION
  199. *>
  200. *> The second factorization stopping criterion, cannot be NaN.
  201. *>
  202. *> The absolute tolerance (stopping threshold) for
  203. *> maximum column 2-norm of the residual matrix R22(K).
  204. *> The algorithm converges (stops the factorization) when
  205. *> the maximum column 2-norm of the residual matrix R22(K)
  206. *> is less than or equal to ABSTOL. Let SAFMIN = DLAMCH('S').
  207. *>
  208. *> a) If ABSTOL is NaN, then no computation is performed
  209. *> and an error message ( INFO = -5 ) is issued
  210. *> by XERBLA.
  211. *>
  212. *> b) If ABSTOL < 0.0, then this stopping criterion is not
  213. *> used, the routine factorizes columns depending
  214. *> on KMAX and RELTOL.
  215. *> This includes the case ABSTOL = -Inf.
  216. *>
  217. *> c) If 0.0 <= ABSTOL < 2*SAFMIN, then ABSTOL = 2*SAFMIN
  218. *> is used. This includes the case ABSTOL = -0.0.
  219. *>
  220. *> d) If 2*SAFMIN <= ABSTOL then the input value
  221. *> of ABSTOL is used.
  222. *>
  223. *> Let MAXC2NRM be the maximum column 2-norm of the
  224. *> whole original matrix A.
  225. *> If ABSTOL chosen above is >= MAXC2NRM, then this
  226. *> stopping criterion is satisfied on input and routine exits
  227. *> immediately after MAXC2NRM is computed. The routine
  228. *> returns MAXC2NRM in MAXC2NORMK,
  229. *> and 1.0 in RELMAXC2NORMK.
  230. *> This includes the case ABSTOL = +Inf. This means that the
  231. *> factorization is not performed, the matrices A and B are not
  232. *> modified, and the matrix A is itself the residual.
  233. *> \endverbatim
  234. *>
  235. *> \param[in] RELTOL
  236. *> \verbatim
  237. *> RELTOL is DOUBLE PRECISION
  238. *>
  239. *> The third factorization stopping criterion, cannot be NaN.
  240. *>
  241. *> The tolerance (stopping threshold) for the ratio
  242. *> abs(R(K+1,K+1))/abs(R(1,1)) of the maximum column 2-norm of
  243. *> the residual matrix R22(K) to the maximum column 2-norm of
  244. *> the original matrix A. The algorithm converges (stops the
  245. *> factorization), when abs(R(K+1,K+1))/abs(R(1,1)) A is less
  246. *> than or equal to RELTOL. Let EPS = DLAMCH('E').
  247. *>
  248. *> a) If RELTOL is NaN, then no computation is performed
  249. *> and an error message ( INFO = -6 ) is issued
  250. *> by XERBLA.
  251. *>
  252. *> b) If RELTOL < 0.0, then this stopping criterion is not
  253. *> used, the routine factorizes columns depending
  254. *> on KMAX and ABSTOL.
  255. *> This includes the case RELTOL = -Inf.
  256. *>
  257. *> c) If 0.0 <= RELTOL < EPS, then RELTOL = EPS is used.
  258. *> This includes the case RELTOL = -0.0.
  259. *>
  260. *> d) If EPS <= RELTOL then the input value of RELTOL
  261. *> is used.
  262. *>
  263. *> Let MAXC2NRM be the maximum column 2-norm of the
  264. *> whole original matrix A.
  265. *> If RELTOL chosen above is >= 1.0, then this stopping
  266. *> criterion is satisfied on input and routine exits
  267. *> immediately after MAXC2NRM is computed.
  268. *> The routine returns MAXC2NRM in MAXC2NORMK,
  269. *> and 1.0 in RELMAXC2NORMK.
  270. *> This includes the case RELTOL = +Inf. This means that the
  271. *> factorization is not performed, the matrices A and B are not
  272. *> modified, and the matrix A is itself the residual.
  273. *>
  274. *> NOTE: We recommend that RELTOL satisfy
  275. *> min( max(M,N)*EPS, sqrt(EPS) ) <= RELTOL
  276. *> \endverbatim
  277. *>
  278. *> \param[in,out] A
  279. *> \verbatim
  280. *> A is DOUBLE PRECISION array, dimension (LDA,N+NRHS)
  281. *>
  282. *> On entry:
  283. *>
  284. *> a) The subarray A(1:M,1:N) contains the M-by-N matrix A.
  285. *> b) The subarray A(1:M,N+1:N+NRHS) contains the M-by-NRHS
  286. *> matrix B.
  287. *>
  288. *> N NRHS
  289. *> array_A = M [ mat_A, mat_B ]
  290. *>
  291. *> On exit:
  292. *>
  293. *> a) The subarray A(1:M,1:N) contains parts of the factors
  294. *> of the matrix A:
  295. *>
  296. *> 1) If K = 0, A(1:M,1:N) contains the original matrix A.
  297. *> 2) If K > 0, A(1:M,1:N) contains parts of the
  298. *> factors:
  299. *>
  300. *> 1. The elements below the diagonal of the subarray
  301. *> A(1:M,1:K) together with TAU(1:K) represent the
  302. *> orthogonal matrix Q(K) as a product of K Householder
  303. *> elementary reflectors.
  304. *>
  305. *> 2. The elements on and above the diagonal of
  306. *> the subarray A(1:K,1:N) contain K-by-N
  307. *> upper-trapezoidal matrix
  308. *> R(K)_approx = ( R11(K), R12(K) ).
  309. *> NOTE: If K=min(M,N), i.e. full rank factorization,
  310. *> then R_approx(K) is the full factor R which
  311. *> is upper-trapezoidal. If, in addition, M>=N,
  312. *> then R is upper-triangular.
  313. *>
  314. *> 3. The subarray A(K+1:M,K+1:N) contains (M-K)-by-(N-K)
  315. *> rectangular matrix R(K)_residual = R22(K).
  316. *>
  317. *> b) If NRHS > 0, the subarray A(1:M,N+1:N+NRHS) contains
  318. *> the M-by-NRHS product Q(K)**T * B.
  319. *> \endverbatim
  320. *>
  321. *> \param[in] LDA
  322. *> \verbatim
  323. *> LDA is INTEGER
  324. *> The leading dimension of the array A. LDA >= max(1,M).
  325. *> This is the leading dimension for both matrices, A and B.
  326. *> \endverbatim
  327. *>
  328. *> \param[out] K
  329. *> \verbatim
  330. *> K is INTEGER
  331. *> Factorization rank of the matrix A, i.e. the rank of
  332. *> the factor R, which is the same as the number of non-zero
  333. *> rows of the factor R. 0 <= K <= min(M,KMAX,N).
  334. *>
  335. *> K also represents the number of non-zero Householder
  336. *> vectors.
  337. *>
  338. *> NOTE: If K = 0, a) the arrays A and B are not modified;
  339. *> b) the array TAU(1:min(M,N)) is set to ZERO,
  340. *> if the matrix A does not contain NaN,
  341. *> otherwise the elements TAU(1:min(M,N))
  342. *> are undefined;
  343. *> c) the elements of the array JPIV are set
  344. *> as follows: for j = 1:N, JPIV(j) = j.
  345. *> \endverbatim
  346. *>
  347. *> \param[out] MAXC2NRMK
  348. *> \verbatim
  349. *> MAXC2NRMK is DOUBLE PRECISION
  350. *> The maximum column 2-norm of the residual matrix R22(K),
  351. *> when the factorization stopped at rank K. MAXC2NRMK >= 0.
  352. *>
  353. *> a) If K = 0, i.e. the factorization was not performed,
  354. *> the matrix A was not modified and is itself a residual
  355. *> matrix, then MAXC2NRMK equals the maximum column 2-norm
  356. *> of the original matrix A.
  357. *>
  358. *> b) If 0 < K < min(M,N), then MAXC2NRMK is returned.
  359. *>
  360. *> c) If K = min(M,N), i.e. the whole matrix A was
  361. *> factorized and there is no residual matrix,
  362. *> then MAXC2NRMK = 0.0.
  363. *>
  364. *> NOTE: MAXC2NRMK in the factorization step K would equal
  365. *> R(K+1,K+1) in the next factorization step K+1.
  366. *> \endverbatim
  367. *>
  368. *> \param[out] RELMAXC2NRMK
  369. *> \verbatim
  370. *> RELMAXC2NRMK is DOUBLE PRECISION
  371. *> The ratio MAXC2NRMK / MAXC2NRM of the maximum column
  372. *> 2-norm of the residual matrix R22(K) (when the factorization
  373. *> stopped at rank K) to the maximum column 2-norm of the
  374. *> whole original matrix A. RELMAXC2NRMK >= 0.
  375. *>
  376. *> a) If K = 0, i.e. the factorization was not performed,
  377. *> the matrix A was not modified and is itself a residual
  378. *> matrix, then RELMAXC2NRMK = 1.0.
  379. *>
  380. *> b) If 0 < K < min(M,N), then
  381. *> RELMAXC2NRMK = MAXC2NRMK / MAXC2NRM is returned.
  382. *>
  383. *> c) If K = min(M,N), i.e. the whole matrix A was
  384. *> factorized and there is no residual matrix,
  385. *> then RELMAXC2NRMK = 0.0.
  386. *>
  387. *> NOTE: RELMAXC2NRMK in the factorization step K would equal
  388. *> abs(R(K+1,K+1))/abs(R(1,1)) in the next factorization
  389. *> step K+1.
  390. *> \endverbatim
  391. *>
  392. *> \param[out] JPIV
  393. *> \verbatim
  394. *> JPIV is INTEGER array, dimension (N)
  395. *> Column pivot indices. For 1 <= j <= N, column j
  396. *> of the matrix A was interchanged with column JPIV(j).
  397. *>
  398. *> The elements of the array JPIV(1:N) are always set
  399. *> by the routine, for example, even when no columns
  400. *> were factorized, i.e. when K = 0, the elements are
  401. *> set as JPIV(j) = j for j = 1:N.
  402. *> \endverbatim
  403. *>
  404. *> \param[out] TAU
  405. *> \verbatim
  406. *> TAU is DOUBLE PRECISION array, dimension (min(M,N))
  407. *> The scalar factors of the elementary reflectors.
  408. *>
  409. *> If 0 < K <= min(M,N), only the elements TAU(1:K) of
  410. *> the array TAU are modified by the factorization.
  411. *> After the factorization computed, if no NaN was found
  412. *> during the factorization, the remaining elements
  413. *> TAU(K+1:min(M,N)) are set to zero, otherwise the
  414. *> elements TAU(K+1:min(M,N)) are not set and therefore
  415. *> undefined.
  416. *> ( If K = 0, all elements of TAU are set to zero, if
  417. *> the matrix A does not contain NaN. )
  418. *> \endverbatim
  419. *>
  420. *> \param[out] WORK
  421. *> \verbatim
  422. *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
  423. *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
  424. *> \endverbatim
  425. *>
  426. *> \param[in] LWORK
  427. *> \verbatim
  428. *> LWORK is INTEGER
  429. *> The dimension of the array WORK.
  430. *> LWORK >= 1, if MIN(M,N) = 0, and
  431. *> LWORK >= (3*N+NRHS-1), otherwise.
  432. *> For optimal performance LWORK >= (2*N + NB*( N+NRHS+1 )),
  433. *> where NB is the optimal block size for DGEQP3RK returned
  434. *> by ILAENV. Minimal block size MINNB=2.
  435. *>
  436. *> NOTE: The decision, whether to use unblocked BLAS 2
  437. *> or blocked BLAS 3 code is based not only on the dimension
  438. *> LWORK of the availbale workspace WORK, but also also on the
  439. *> matrix A dimension N via crossover point NX returned
  440. *> by ILAENV. (For N less than NX, unblocked code should be
  441. *> used.)
  442. *>
  443. *> If LWORK = -1, then a workspace query is assumed;
  444. *> the routine only calculates the optimal size of the WORK
  445. *> array, returns this value as the first entry of the WORK
  446. *> array, and no error message related to LWORK is issued
  447. *> by XERBLA.
  448. *> \endverbatim
  449. *>
  450. *> \param[out] IWORK
  451. *> \verbatim
  452. *> IWORK is INTEGER array, dimension (N-1).
  453. *> Is a work array. ( IWORK is used to store indices
  454. *> of "bad" columns for norm downdating in the residual
  455. *> matrix in the blocked step auxiliary subroutine DLAQP3RK ).
  456. *> \endverbatim
  457. *>
  458. *> \param[out] INFO
  459. *> \verbatim
  460. *> INFO is INTEGER
  461. *> 1) INFO = 0: successful exit.
  462. *> 2) INFO < 0: if INFO = -i, the i-th argument had an
  463. *> illegal value.
  464. *> 3) If INFO = j_1, where 1 <= j_1 <= N, then NaN was
  465. *> detected and the routine stops the computation.
  466. *> The j_1-th column of the matrix A or the j_1-th
  467. *> element of array TAU contains the first occurrence
  468. *> of NaN in the factorization step K+1 ( when K columns
  469. *> have been factorized ).
  470. *>
  471. *> On exit:
  472. *> K is set to the number of
  473. *> factorized columns without
  474. *> exception.
  475. *> MAXC2NRMK is set to NaN.
  476. *> RELMAXC2NRMK is set to NaN.
  477. *> TAU(K+1:min(M,N)) is not set and contains undefined
  478. *> elements. If j_1=K+1, TAU(K+1)
  479. *> may contain NaN.
  480. *> 4) If INFO = j_2, where N+1 <= j_2 <= 2*N, then no NaN
  481. *> was detected, but +Inf (or -Inf) was detected and
  482. *> the routine continues the computation until completion.
  483. *> The (j_2-N)-th column of the matrix A contains the first
  484. *> occurrence of +Inf (or -Inf) in the factorization
  485. *> step K+1 ( when K columns have been factorized ).
  486. *> \endverbatim
  487. *
  488. * Authors:
  489. * ========
  490. *
  491. *> \author Univ. of Tennessee
  492. *> \author Univ. of California Berkeley
  493. *> \author Univ. of Colorado Denver
  494. *> \author NAG Ltd.
  495. *
  496. *> \ingroup geqp3rk
  497. *
  498. *> \par Further Details:
  499. * =====================
  500. *
  501. *> \verbatim
  502. *> DGEQP3RK is based on the same BLAS3 Householder QR factorization
  503. *> algorithm with column pivoting as in DGEQP3 routine which uses
  504. *> DLARFG routine to generate Householder reflectors
  505. *> for QR factorization.
  506. *>
  507. *> We can also write:
  508. *>
  509. *> A = A_approx(K) + A_residual(K)
  510. *>
  511. *> The low rank approximation matrix A(K)_approx from
  512. *> the truncated QR factorization of rank K of the matrix A is:
  513. *>
  514. *> A(K)_approx = Q(K) * ( R(K)_approx ) * P(K)**T
  515. *> ( 0 0 )
  516. *>
  517. *> = Q(K) * ( R11(K) R12(K) ) * P(K)**T
  518. *> ( 0 0 )
  519. *>
  520. *> The residual A_residual(K) of the matrix A is:
  521. *>
  522. *> A_residual(K) = Q(K) * ( 0 0 ) * P(K)**T =
  523. *> ( 0 R(K)_residual )
  524. *>
  525. *> = Q(K) * ( 0 0 ) * P(K)**T
  526. *> ( 0 R22(K) )
  527. *>
  528. *> The truncated (rank K) factorization guarantees that
  529. *> the maximum column 2-norm of A_residual(K) is less than
  530. *> or equal to MAXC2NRMK up to roundoff error.
  531. *>
  532. *> NOTE: An approximation of the null vectors
  533. *> of A can be easily computed from R11(K)
  534. *> and R12(K):
  535. *>
  536. *> Null( A(K) )_approx = P * ( inv(R11(K)) * R12(K) )
  537. *> ( -I )
  538. *>
  539. *> \endverbatim
  540. *
  541. *> \par References:
  542. * ================
  543. *> [1] A Level 3 BLAS QR factorization algorithm with column pivoting developed in 1996.
  544. *> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain.
  545. *> X. Sun, Computer Science Dept., Duke University, USA.
  546. *> C. H. Bischof, Math. and Comp. Sci. Div., Argonne National Lab, USA.
  547. *> A BLAS-3 version of the QR factorization with column pivoting.
  548. *> LAPACK Working Note 114
  549. *> \htmlonly
  550. *> <a href="https://www.netlib.org/lapack/lawnspdf/lawn114.pdf">https://www.netlib.org/lapack/lawnspdf/lawn114.pdf</a>
  551. *> \endhtmlonly
  552. *> and in
  553. *> SIAM J. Sci. Comput., 19(5):1486-1494, Sept. 1998.
  554. *> \htmlonly
  555. *> <a href="https://doi.org/10.1137/S1064827595296732">https://doi.org/10.1137/S1064827595296732</a>
  556. *> \endhtmlonly
  557. *>
  558. *> [2] A partial column norm updating strategy developed in 2006.
  559. *> Z. Drmac and Z. Bujanovic, Dept. of Math., University of Zagreb, Croatia.
  560. *> On the failure of rank revealing QR factorization software – a case study.
  561. *> LAPACK Working Note 176.
  562. *> \htmlonly
  563. *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">http://www.netlib.org/lapack/lawnspdf/lawn176.pdf</a>
  564. *> \endhtmlonly
  565. *> and in
  566. *> ACM Trans. Math. Softw. 35, 2, Article 12 (July 2008), 28 pages.
  567. *> \htmlonly
  568. *> <a href="https://doi.org/10.1145/1377612.1377616">https://doi.org/10.1145/1377612.1377616</a>
  569. *> \endhtmlonly
  570. *
  571. *> \par Contributors:
  572. * ==================
  573. *>
  574. *> \verbatim
  575. *>
  576. *> November 2023, Igor Kozachenko, James Demmel,
  577. *> EECS Department,
  578. *> University of California, Berkeley, USA.
  579. *>
  580. *> \endverbatim
  581. *
  582. * =====================================================================
  583. SUBROUTINE DGEQP3RK( M, N, NRHS, KMAX, ABSTOL, RELTOL, A, LDA,
  584. $ K, MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
  585. $ WORK, LWORK, IWORK, INFO )
  586. IMPLICIT NONE
  587. *
  588. * -- LAPACK computational routine --
  589. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  590. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  591. *
  592. * .. Scalar Arguments ..
  593. INTEGER INFO, K, KF, KMAX, LDA, LWORK, M, N, NRHS
  594. DOUBLE PRECISION ABSTOL, MAXC2NRMK, RELMAXC2NRMK, RELTOL
  595. * ..
  596. * .. Array Arguments ..
  597. INTEGER IWORK( * ), JPIV( * )
  598. DOUBLE PRECISION A( LDA, * ), TAU( * ), WORK( * )
  599. * ..
  600. *
  601. * =====================================================================
  602. *
  603. * .. Parameters ..
  604. INTEGER INB, INBMIN, IXOVER
  605. PARAMETER ( INB = 1, INBMIN = 2, IXOVER = 3 )
  606. DOUBLE PRECISION ZERO, ONE, TWO
  607. PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0, TWO = 2.0D+0 )
  608. * ..
  609. * .. Local Scalars ..
  610. LOGICAL LQUERY, DONE
  611. INTEGER IINFO, IOFFSET, IWS, J, JB, JBF, JMAXB, JMAX,
  612. $ JMAXC2NRM, KP1, LWKOPT, MINMN, N_SUB, NB,
  613. $ NBMIN, NX
  614. DOUBLE PRECISION EPS, HUGEVAL, MAXC2NRM, SAFMIN
  615. * ..
  616. * .. External Subroutines ..
  617. EXTERNAL DLAQP2RK, DLAQP3RK, XERBLA
  618. * ..
  619. * .. External Functions ..
  620. LOGICAL DISNAN
  621. INTEGER IDAMAX, ILAENV
  622. DOUBLE PRECISION DLAMCH, DNRM2
  623. EXTERNAL DISNAN, DLAMCH, DNRM2, IDAMAX, ILAENV
  624. * ..
  625. * .. Intrinsic Functions ..
  626. INTRINSIC DBLE, MAX, MIN
  627. * ..
  628. * .. Executable Statements ..
  629. *
  630. * Test input arguments
  631. * ====================
  632. *
  633. INFO = 0
  634. LQUERY = ( LWORK.EQ.-1 )
  635. IF( M.LT.0 ) THEN
  636. INFO = -1
  637. ELSE IF( N.LT.0 ) THEN
  638. INFO = -2
  639. ELSE IF( NRHS.LT.0 ) THEN
  640. INFO = -3
  641. ELSE IF( KMAX.LT.0 ) THEN
  642. INFO = -4
  643. ELSE IF( DISNAN( ABSTOL ) ) THEN
  644. INFO = -5
  645. ELSE IF( DISNAN( RELTOL ) ) THEN
  646. INFO = -6
  647. ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
  648. INFO = -8
  649. END IF
  650. *
  651. * If the input parameters M, N, NRHS, KMAX, LDA are valid:
  652. * a) Test the input workspace size LWORK for the minimum
  653. * size requirement IWS.
  654. * b) Determine the optimal block size NB and optimal
  655. * workspace size LWKOPT to be returned in WORK(1)
  656. * in case of (1) LWORK < IWS, (2) LQUERY = .TRUE.,
  657. * (3) when routine exits.
  658. * Here, IWS is the miminum workspace required for unblocked
  659. * code.
  660. *
  661. IF( INFO.EQ.0 ) THEN
  662. MINMN = MIN( M, N )
  663. IF( MINMN.EQ.0 ) THEN
  664. IWS = 1
  665. LWKOPT = 1
  666. ELSE
  667. *
  668. * Minimal workspace size in case of using only unblocked
  669. * BLAS 2 code in DLAQP2RK.
  670. * 1) DGEQP3RK and DLAQP2RK: 2*N to store full and partial
  671. * column 2-norms.
  672. * 2) DLAQP2RK: N+NRHS-1 to use in WORK array that is used
  673. * in DLARF subroutine inside DLAQP2RK to apply an
  674. * elementary reflector from the left.
  675. * TOTAL_WORK_SIZE = 3*N + NRHS - 1
  676. *
  677. IWS = 3*N + NRHS - 1
  678. *
  679. * Assign to NB optimal block size.
  680. *
  681. NB = ILAENV( INB, 'DGEQP3RK', ' ', M, N, -1, -1 )
  682. *
  683. * A formula for the optimal workspace size in case of using
  684. * both unblocked BLAS 2 in DLAQP2RK and blocked BLAS 3 code
  685. * in DLAQP3RK.
  686. * 1) DGEQP3RK, DLAQP2RK, DLAQP3RK: 2*N to store full and
  687. * partial column 2-norms.
  688. * 2) DLAQP2RK: N+NRHS-1 to use in WORK array that is used
  689. * in DLARF subroutine to apply an elementary reflector
  690. * from the left.
  691. * 3) DLAQP3RK: NB*(N+NRHS) to use in the work array F that
  692. * is used to apply a block reflector from
  693. * the left.
  694. * 4) DLAQP3RK: NB to use in the auxilixary array AUX.
  695. * Sizes (2) and ((3) + (4)) should intersect, therefore
  696. * TOTAL_WORK_SIZE = 2*N + NB*( N+NRHS+1 ), given NBMIN=2.
  697. *
  698. LWKOPT = 2*N + NB*( N+NRHS+1 )
  699. END IF
  700. WORK( 1 ) = DBLE( LWKOPT )
  701. *
  702. IF( ( LWORK.LT.IWS ) .AND. .NOT.LQUERY ) THEN
  703. INFO = -15
  704. END IF
  705. END IF
  706. *
  707. * NOTE: The optimal workspace size is returned in WORK(1), if
  708. * the input parameters M, N, NRHS, KMAX, LDA are valid.
  709. *
  710. IF( INFO.NE.0 ) THEN
  711. CALL XERBLA( 'DGEQP3RK', -INFO )
  712. RETURN
  713. ELSE IF( LQUERY ) THEN
  714. RETURN
  715. END IF
  716. *
  717. * Quick return if possible for M=0 or N=0.
  718. *
  719. IF( MINMN.EQ.0 ) THEN
  720. K = 0
  721. MAXC2NRMK = ZERO
  722. RELMAXC2NRMK = ZERO
  723. WORK( 1 ) = DBLE( LWKOPT )
  724. RETURN
  725. END IF
  726. *
  727. * ==================================================================
  728. *
  729. * Initialize column pivot array JPIV.
  730. *
  731. DO J = 1, N
  732. JPIV( J ) = J
  733. END DO
  734. *
  735. * ==================================================================
  736. *
  737. * Initialize storage for partial and exact column 2-norms.
  738. * a) The elements WORK(1:N) are used to store partial column
  739. * 2-norms of the matrix A, and may decrease in each computation
  740. * step; initialize to the values of complete columns 2-norms.
  741. * b) The elements WORK(N+1:2*N) are used to store complete column
  742. * 2-norms of the matrix A, they are not changed during the
  743. * computation; initialize the values of complete columns 2-norms.
  744. *
  745. DO J = 1, N
  746. WORK( J ) = DNRM2( M, A( 1, J ), 1 )
  747. WORK( N+J ) = WORK( J )
  748. END DO
  749. *
  750. * ==================================================================
  751. *
  752. * Compute the pivot column index and the maximum column 2-norm
  753. * for the whole original matrix stored in A(1:M,1:N).
  754. *
  755. KP1 = IDAMAX( N, WORK( 1 ), 1 )
  756. MAXC2NRM = WORK( KP1 )
  757. *
  758. * ==================================================================.
  759. *
  760. IF( DISNAN( MAXC2NRM ) ) THEN
  761. *
  762. * Check if the matrix A contains NaN, set INFO parameter
  763. * to the column number where the first NaN is found and return
  764. * from the routine.
  765. *
  766. K = 0
  767. INFO = KP1
  768. *
  769. * Set MAXC2NRMK and RELMAXC2NRMK to NaN.
  770. *
  771. MAXC2NRMK = MAXC2NRM
  772. RELMAXC2NRMK = MAXC2NRM
  773. *
  774. * Array TAU is not set and contains undefined elements.
  775. *
  776. WORK( 1 ) = DBLE( LWKOPT )
  777. RETURN
  778. END IF
  779. *
  780. * ===================================================================
  781. *
  782. IF( MAXC2NRM.EQ.ZERO ) THEN
  783. *
  784. * Check is the matrix A is a zero matrix, set array TAU and
  785. * return from the routine.
  786. *
  787. K = 0
  788. MAXC2NRMK = ZERO
  789. RELMAXC2NRMK = ZERO
  790. *
  791. DO J = 1, MINMN
  792. TAU( J ) = ZERO
  793. END DO
  794. *
  795. WORK( 1 ) = DBLE( LWKOPT )
  796. RETURN
  797. *
  798. END IF
  799. *
  800. * ===================================================================
  801. *
  802. HUGEVAL = DLAMCH( 'Overflow' )
  803. *
  804. IF( MAXC2NRM.GT.HUGEVAL ) THEN
  805. *
  806. * Check if the matrix A contains +Inf or -Inf, set INFO parameter
  807. * to the column number, where the first +/-Inf is found plus N,
  808. * and continue the computation.
  809. *
  810. INFO = N + KP1
  811. *
  812. END IF
  813. *
  814. * ==================================================================
  815. *
  816. * Quick return if possible for the case when the first
  817. * stopping criterion is satisfied, i.e. KMAX = 0.
  818. *
  819. IF( KMAX.EQ.0 ) THEN
  820. K = 0
  821. MAXC2NRMK = MAXC2NRM
  822. RELMAXC2NRMK = ONE
  823. DO J = 1, MINMN
  824. TAU( J ) = ZERO
  825. END DO
  826. WORK( 1 ) = DBLE( LWKOPT )
  827. RETURN
  828. END IF
  829. *
  830. * ==================================================================
  831. *
  832. EPS = DLAMCH('Epsilon')
  833. *
  834. * Adjust ABSTOL
  835. *
  836. IF( ABSTOL.GE.ZERO ) THEN
  837. SAFMIN = DLAMCH('Safe minimum')
  838. ABSTOL = MAX( ABSTOL, TWO*SAFMIN )
  839. END IF
  840. *
  841. * Adjust RELTOL
  842. *
  843. IF( RELTOL.GE.ZERO ) THEN
  844. RELTOL = MAX( RELTOL, EPS )
  845. END IF
  846. *
  847. * ===================================================================
  848. *
  849. * JMAX is the maximum index of the column to be factorized,
  850. * which is also limited by the first stopping criterion KMAX.
  851. *
  852. JMAX = MIN( KMAX, MINMN )
  853. *
  854. * ===================================================================
  855. *
  856. * Quick return if possible for the case when the second or third
  857. * stopping criterion for the whole original matrix is satified,
  858. * i.e. MAXC2NRM <= ABSTOL or RELMAXC2NRM <= RELTOL
  859. * (which is ONE <= RELTOL).
  860. *
  861. IF( MAXC2NRM.LE.ABSTOL .OR. ONE.LE.RELTOL ) THEN
  862. *
  863. K = 0
  864. MAXC2NRMK = MAXC2NRM
  865. RELMAXC2NRMK = ONE
  866. *
  867. DO J = 1, MINMN
  868. TAU( J ) = ZERO
  869. END DO
  870. *
  871. WORK( 1 ) = DBLE( LWKOPT )
  872. RETURN
  873. END IF
  874. *
  875. * ==================================================================
  876. * Factorize columns
  877. * ==================================================================
  878. *
  879. * Determine the block size.
  880. *
  881. NBMIN = 2
  882. NX = 0
  883. *
  884. IF( ( NB.GT.1 ) .AND. ( NB.LT.MINMN ) ) THEN
  885. *
  886. * Determine when to cross over from blocked to unblocked code.
  887. * (for N less than NX, unblocked code should be used).
  888. *
  889. NX = MAX( 0, ILAENV( IXOVER, 'DGEQP3RK', ' ', M, N, -1, -1 ))
  890. *
  891. IF( NX.LT.MINMN ) THEN
  892. *
  893. * Determine if workspace is large enough for blocked code.
  894. *
  895. IF( LWORK.LT.LWKOPT ) THEN
  896. *
  897. * Not enough workspace to use optimal block size that
  898. * is currently stored in NB.
  899. * Reduce NB and determine the minimum value of NB.
  900. *
  901. NB = ( LWORK-2*N ) / ( N+1 )
  902. NBMIN = MAX( 2, ILAENV( INBMIN, 'DGEQP3RK', ' ', M, N,
  903. $ -1, -1 ) )
  904. *
  905. END IF
  906. END IF
  907. END IF
  908. *
  909. * ==================================================================
  910. *
  911. * DONE is the boolean flag to rerpresent the case when the
  912. * factorization completed in the block factorization routine,
  913. * before the end of the block.
  914. *
  915. DONE = .FALSE.
  916. *
  917. * J is the column index.
  918. *
  919. J = 1
  920. *
  921. * (1) Use blocked code initially.
  922. *
  923. * JMAXB is the maximum column index of the block, when the
  924. * blocked code is used, is also limited by the first stopping
  925. * criterion KMAX.
  926. *
  927. JMAXB = MIN( KMAX, MINMN - NX )
  928. *
  929. IF( NB.GE.NBMIN .AND. NB.LT.JMAX .AND. JMAXB.GT.0 ) THEN
  930. *
  931. * Loop over the column blocks of the matrix A(1:M,1:JMAXB). Here:
  932. * J is the column index of a column block;
  933. * JB is the column block size to pass to block factorization
  934. * routine in a loop step;
  935. * JBF is the number of columns that were actually factorized
  936. * that was returned by the block factorization routine
  937. * in a loop step, JBF <= JB;
  938. * N_SUB is the number of columns in the submatrix;
  939. * IOFFSET is the number of rows that should not be factorized.
  940. *
  941. DO WHILE( J.LE.JMAXB )
  942. *
  943. JB = MIN( NB, JMAXB-J+1 )
  944. N_SUB = N-J+1
  945. IOFFSET = J-1
  946. *
  947. * Factorize JB columns among the columns A(J:N).
  948. *
  949. CALL DLAQP3RK( M, N_SUB, NRHS, IOFFSET, JB, ABSTOL,
  950. $ RELTOL, KP1, MAXC2NRM, A( 1, J ), LDA,
  951. $ DONE, JBF, MAXC2NRMK, RELMAXC2NRMK,
  952. $ JPIV( J ), TAU( J ),
  953. $ WORK( J ), WORK( N+J ),
  954. $ WORK( 2*N+1 ), WORK( 2*N+JB+1 ),
  955. $ N+NRHS-J+1, IWORK, IINFO )
  956. *
  957. * Set INFO on the first occurence of Inf.
  958. *
  959. IF( IINFO.GT.N_SUB .AND. INFO.EQ.0 ) THEN
  960. INFO = 2*IOFFSET + IINFO
  961. END IF
  962. *
  963. IF( DONE ) THEN
  964. *
  965. * Either the submatrix is zero before the end of the
  966. * column block, or ABSTOL or RELTOL criterion is
  967. * satisfied before the end of the column block, we can
  968. * return from the routine. Perform the following before
  969. * returning:
  970. * a) Set the number of factorized columns K,
  971. * K = IOFFSET + JBF from the last call of blocked
  972. * routine.
  973. * NOTE: 1) MAXC2NRMK and RELMAXC2NRMK are returned
  974. * by the block factorization routine;
  975. * 2) The remaining TAUs are set to ZERO by the
  976. * block factorization routine.
  977. *
  978. K = IOFFSET + JBF
  979. *
  980. * Set INFO on the first occurrence of NaN, NaN takes
  981. * prcedence over Inf.
  982. *
  983. IF( IINFO.LE.N_SUB .AND. IINFO.GT.0 ) THEN
  984. INFO = IOFFSET + IINFO
  985. END IF
  986. *
  987. * Return from the routine.
  988. *
  989. WORK( 1 ) = DBLE( LWKOPT )
  990. *
  991. RETURN
  992. *
  993. END IF
  994. *
  995. J = J + JBF
  996. *
  997. END DO
  998. *
  999. END IF
  1000. *
  1001. * Use unblocked code to factor the last or only block.
  1002. * J = JMAX+1 means we factorized the maximum possible number of
  1003. * columns, that is in ELSE clause we need to compute
  1004. * the MAXC2NORM and RELMAXC2NORM to return after we processed
  1005. * the blocks.
  1006. *
  1007. IF( J.LE.JMAX ) THEN
  1008. *
  1009. * N_SUB is the number of columns in the submatrix;
  1010. * IOFFSET is the number of rows that should not be factorized.
  1011. *
  1012. N_SUB = N-J+1
  1013. IOFFSET = J-1
  1014. *
  1015. CALL DLAQP2RK( M, N_SUB, NRHS, IOFFSET, JMAX-J+1,
  1016. $ ABSTOL, RELTOL, KP1, MAXC2NRM, A( 1, J ), LDA,
  1017. $ KF, MAXC2NRMK, RELMAXC2NRMK, JPIV( J ),
  1018. $ TAU( J ), WORK( J ), WORK( N+J ),
  1019. $ WORK( 2*N+1 ), IINFO )
  1020. *
  1021. * ABSTOL or RELTOL criterion is satisfied when the number of
  1022. * the factorized columns KF is smaller then the number
  1023. * of columns JMAX-J+1 supplied to be factorized by the
  1024. * unblocked routine, we can return from
  1025. * the routine. Perform the following before returning:
  1026. * a) Set the number of factorized columns K,
  1027. * b) MAXC2NRMK and RELMAXC2NRMK are returned by the
  1028. * unblocked factorization routine above.
  1029. *
  1030. K = J - 1 + KF
  1031. *
  1032. * Set INFO on the first exception occurence.
  1033. *
  1034. * Set INFO on the first exception occurence of Inf or NaN,
  1035. * (NaN takes precedence over Inf).
  1036. *
  1037. IF( IINFO.GT.N_SUB .AND. INFO.EQ.0 ) THEN
  1038. INFO = 2*IOFFSET + IINFO
  1039. ELSE IF( IINFO.LE.N_SUB .AND. IINFO.GT.0 ) THEN
  1040. INFO = IOFFSET + IINFO
  1041. END IF
  1042. *
  1043. ELSE
  1044. *
  1045. * Compute the return values for blocked code.
  1046. *
  1047. * Set the number of factorized columns if the unblocked routine
  1048. * was not called.
  1049. *
  1050. K = JMAX
  1051. *
  1052. * If there exits a residual matrix after the blocked code:
  1053. * 1) compute the values of MAXC2NRMK, RELMAXC2NRMK of the
  1054. * residual matrix, otherwise set them to ZERO;
  1055. * 2) Set TAU(K+1:MINMN) to ZERO.
  1056. *
  1057. IF( K.LT.MINMN ) THEN
  1058. JMAXC2NRM = K + IDAMAX( N-K, WORK( K+1 ), 1 )
  1059. MAXC2NRMK = WORK( JMAXC2NRM )
  1060. IF( K.EQ.0 ) THEN
  1061. RELMAXC2NRMK = ONE
  1062. ELSE
  1063. RELMAXC2NRMK = MAXC2NRMK / MAXC2NRM
  1064. END IF
  1065. *
  1066. DO J = K + 1, MINMN
  1067. TAU( J ) = ZERO
  1068. END DO
  1069. *
  1070. END IF
  1071. *
  1072. * END IF( J.LE.JMAX ) THEN
  1073. *
  1074. END IF
  1075. *
  1076. WORK( 1 ) = DBLE( LWKOPT )
  1077. *
  1078. RETURN
  1079. *
  1080. * End of DGEQP3RK
  1081. *
  1082. END