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.

zlaqp3rk.f 35 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947
  1. *> \brief \b ZLAQP3RK computes a step of truncated QR factorization with column pivoting of a complex m-by-n matrix A using Level 3 BLAS and overwrites a complex m-by-nrhs matrix B with Q**H * B.
  2. *
  3. * =========== DOCUMENTATION ===========
  4. *
  5. * Online html documentation available at
  6. * http://www.netlib.org/lapack/explore-html/
  7. *
  8. *> \htmlonly
  9. *> Download ZLAQP3RK + dependencies
  10. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlaqp3rk.f">
  11. *> [TGZ]</a>
  12. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlaqp3rk.f">
  13. *> [ZIP]</a>
  14. *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlaqp3rk.f">
  15. *> [TXT]</a>
  16. *> \endhtmlonly
  17. *
  18. * Definition:
  19. * ===========
  20. *
  21. * SUBROUTINE ZLAQP3RK( M, N, NRHS, IOFFSET, NB, ABSTOL,
  22. * $ RELTOL, KP1, MAXC2NRM, A, LDA, DONE, KB,
  23. * $ MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
  24. * $ VN1, VN2, AUXV, F, LDF, IWORK, INFO )
  25. * IMPLICIT NONE
  26. * LOGICAL DONE
  27. * INTEGER INFO, IOFFSET, KB, KP1, LDA, LDF, M, N,
  28. * $ NB, NRHS
  29. * DOUBLE PRECISION ABSTOL, MAXC2NRM, MAXC2NRMK, RELMAXC2NRMK,
  30. * $ RELTOL
  31. * ..
  32. * .. Array Arguments ..
  33. * INTEGER IWORK( * ), JPIV( * )
  34. * DOUBLE PRECISION VN1( * ), VN2( * )
  35. * COMPLEX*16 A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
  36. * ..
  37. *
  38. *
  39. *> \par Purpose:
  40. * =============
  41. *>
  42. *> \verbatim
  43. *>
  44. *> ZLAQP3RK computes a step of truncated QR factorization with column
  45. *> pivoting of a complex M-by-N matrix A block A(IOFFSET+1:M,1:N)
  46. *> by using Level 3 BLAS as
  47. *>
  48. *> A * P(KB) = Q(KB) * R(KB).
  49. *>
  50. *> The routine tries to factorize NB columns from A starting from
  51. *> the row IOFFSET+1 and updates the residual matrix with BLAS 3
  52. *> xGEMM. The number of actually factorized columns is returned
  53. *> is smaller than NB.
  54. *>
  55. *> Block A(1:IOFFSET,1:N) is accordingly pivoted, but not factorized.
  56. *>
  57. *> The routine also overwrites the right-hand-sides B matrix stored
  58. *> in A(IOFFSET+1:M,1:N+1:N+NRHS) with Q(KB)**H * B.
  59. *>
  60. *> Cases when the number of factorized columns KB < NB:
  61. *>
  62. *> (1) In some cases, due to catastrophic cancellations, it cannot
  63. *> factorize all NB columns and need to update the residual matrix.
  64. *> Hence, the actual number of factorized columns in the block returned
  65. *> in KB is smaller than NB. The logical DONE is returned as FALSE.
  66. *> The factorization of the whole original matrix A_orig must proceed
  67. *> with the next block.
  68. *>
  69. *> (2) Whenever the stopping criterion ABSTOL or RELTOL is satisfied,
  70. *> the factorization of the whole original matrix A_orig is stopped,
  71. *> the logical DONE is returned as TRUE. The number of factorized
  72. *> columns which is smaller than NB is returned in KB.
  73. *>
  74. *> (3) In case both stopping criteria ABSTOL or RELTOL are not used,
  75. *> and when the residual matrix is a zero matrix in some factorization
  76. *> step KB, the factorization of the whole original matrix A_orig is
  77. *> stopped, the logical DONE is returned as TRUE. The number of
  78. *> factorized columns which is smaller than NB is returned in KB.
  79. *>
  80. *> (4) Whenever NaN is detected in the matrix A or in the array TAU,
  81. *> the factorization of the whole original matrix A_orig is stopped,
  82. *> the logical DONE is returned as TRUE. The number of factorized
  83. *> columns which is smaller than NB is returned in KB. The INFO
  84. *> parameter is set to the column index of the first NaN occurrence.
  85. *>
  86. *> \endverbatim
  87. *
  88. * Arguments:
  89. * ==========
  90. *
  91. *> \param[in] M
  92. *> \verbatim
  93. *> M is INTEGER
  94. *> The number of rows of the matrix A. M >= 0.
  95. *> \endverbatim
  96. *>
  97. *> \param[in] N
  98. *> \verbatim
  99. *> N is INTEGER
  100. *> The number of columns of the matrix A. N >= 0
  101. *> \endverbatim
  102. *>
  103. *> \param[in] NRHS
  104. *> \verbatim
  105. *> NRHS is INTEGER
  106. *> The number of right hand sides, i.e., the number of
  107. *> columns of the matrix B. NRHS >= 0.
  108. *> \endverbatim
  109. *>
  110. *> \param[in] IOFFSET
  111. *> \verbatim
  112. *> IOFFSET is INTEGER
  113. *> The number of rows of the matrix A that must be pivoted
  114. *> but not factorized. IOFFSET >= 0.
  115. *>
  116. *> IOFFSET also represents the number of columns of the whole
  117. *> original matrix A_orig that have been factorized
  118. *> in the previous steps.
  119. *> \endverbatim
  120. *>
  121. *> \param[in] NB
  122. *> \verbatim
  123. *> NB is INTEGER
  124. *> Factorization block size, i.e the number of columns
  125. *> to factorize in the matrix A. 0 <= NB
  126. *>
  127. *> If NB = 0, then the routine exits immediately.
  128. *> This means that the factorization is not performed,
  129. *> the matrices A and B and the arrays TAU, IPIV
  130. *> are not modified.
  131. *> \endverbatim
  132. *>
  133. *> \param[in] ABSTOL
  134. *> \verbatim
  135. *> ABSTOL is DOUBLE PRECISION, cannot be NaN.
  136. *>
  137. *> The absolute tolerance (stopping threshold) for
  138. *> maximum column 2-norm of the residual matrix.
  139. *> The algorithm converges (stops the factorization) when
  140. *> the maximum column 2-norm of the residual matrix
  141. *> is less than or equal to ABSTOL.
  142. *>
  143. *> a) If ABSTOL < 0.0, then this stopping criterion is not
  144. *> used, the routine factorizes columns depending
  145. *> on NB and RELTOL.
  146. *> This includes the case ABSTOL = -Inf.
  147. *>
  148. *> b) If 0.0 <= ABSTOL then the input value
  149. *> of ABSTOL is used.
  150. *> \endverbatim
  151. *>
  152. *> \param[in] RELTOL
  153. *> \verbatim
  154. *> RELTOL is DOUBLE PRECISION, cannot be NaN.
  155. *>
  156. *> The tolerance (stopping threshold) for the ratio of the
  157. *> maximum column 2-norm of the residual matrix to the maximum
  158. *> column 2-norm of the original matrix A_orig. The algorithm
  159. *> converges (stops the factorization), when this ratio is
  160. *> less than or equal to RELTOL.
  161. *>
  162. *> a) If RELTOL < 0.0, then this stopping criterion is not
  163. *> used, the routine factorizes columns depending
  164. *> on NB and ABSTOL.
  165. *> This includes the case RELTOL = -Inf.
  166. *>
  167. *> d) If 0.0 <= RELTOL then the input value of RELTOL
  168. *> is used.
  169. *> \endverbatim
  170. *>
  171. *> \param[in] KP1
  172. *> \verbatim
  173. *> KP1 is INTEGER
  174. *> The index of the column with the maximum 2-norm in
  175. *> the whole original matrix A_orig determined in the
  176. *> main routine ZGEQP3RK. 1 <= KP1 <= N_orig.
  177. *> \endverbatim
  178. *>
  179. *> \param[in] MAXC2NRM
  180. *> \verbatim
  181. *> MAXC2NRM is DOUBLE PRECISION
  182. *> The maximum column 2-norm of the whole original
  183. *> matrix A_orig computed in the main routine ZGEQP3RK.
  184. *> MAXC2NRM >= 0.
  185. *> \endverbatim
  186. *>
  187. *> \param[in,out] A
  188. *> \verbatim
  189. *> A is COMPLEX*16 array, dimension (LDA,N+NRHS)
  190. *> On entry:
  191. *> the M-by-N matrix A and M-by-NRHS matrix B, as in
  192. *>
  193. *> N NRHS
  194. *> array_A = M [ mat_A, mat_B ]
  195. *>
  196. *> On exit:
  197. *> 1. The elements in block A(IOFFSET+1:M,1:KB) below
  198. *> the diagonal together with the array TAU represent
  199. *> the orthogonal matrix Q(KB) as a product of elementary
  200. *> reflectors.
  201. *> 2. The upper triangular block of the matrix A stored
  202. *> in A(IOFFSET+1:M,1:KB) is the triangular factor obtained.
  203. *> 3. The block of the matrix A stored in A(1:IOFFSET,1:N)
  204. *> has been accordingly pivoted, but not factorized.
  205. *> 4. The rest of the array A, block A(IOFFSET+1:M,KB+1:N+NRHS).
  206. *> The left part A(IOFFSET+1:M,KB+1:N) of this block
  207. *> contains the residual of the matrix A, and,
  208. *> if NRHS > 0, the right part of the block
  209. *> A(IOFFSET+1:M,N+1:N+NRHS) contains the block of
  210. *> the right-hand-side matrix B. Both these blocks have been
  211. *> updated by multiplication from the left by Q(KB)**H.
  212. *> \endverbatim
  213. *>
  214. *> \param[in] LDA
  215. *> \verbatim
  216. *> LDA is INTEGER
  217. *> The leading dimension of the array A. LDA >= max(1,M).
  218. *> \endverbatim
  219. *>
  220. *> \param[out]
  221. *> \verbatim
  222. *> DONE is LOGICAL
  223. *> TRUE: a) if the factorization completed before processing
  224. *> all min(M-IOFFSET,NB,N) columns due to ABSTOL
  225. *> or RELTOL criterion,
  226. *> b) if the factorization completed before processing
  227. *> all min(M-IOFFSET,NB,N) columns due to the
  228. *> residual matrix being a ZERO matrix.
  229. *> c) when NaN was detected in the matrix A
  230. *> or in the array TAU.
  231. *> FALSE: otherwise.
  232. *> \endverbatim
  233. *>
  234. *> \param[out] KB
  235. *> \verbatim
  236. *> KB is INTEGER
  237. *> Factorization rank of the matrix A, i.e. the rank of
  238. *> the factor R, which is the same as the number of non-zero
  239. *> rows of the factor R. 0 <= KB <= min(M-IOFFSET,NB,N).
  240. *>
  241. *> KB also represents the number of non-zero Householder
  242. *> vectors.
  243. *> \endverbatim
  244. *>
  245. *> \param[out] MAXC2NRMK
  246. *> \verbatim
  247. *> MAXC2NRMK is DOUBLE PRECISION
  248. *> The maximum column 2-norm of the residual matrix,
  249. *> when the factorization stopped at rank KB. MAXC2NRMK >= 0.
  250. *> \endverbatim
  251. *>
  252. *> \param[out] RELMAXC2NRMK
  253. *> \verbatim
  254. *> RELMAXC2NRMK is DOUBLE PRECISION
  255. *> The ratio MAXC2NRMK / MAXC2NRM of the maximum column
  256. *> 2-norm of the residual matrix (when the factorization
  257. *> stopped at rank KB) to the maximum column 2-norm of the
  258. *> original matrix A_orig. RELMAXC2NRMK >= 0.
  259. *> \endverbatim
  260. *>
  261. *> \param[out] JPIV
  262. *> \verbatim
  263. *> JPIV is INTEGER array, dimension (N)
  264. *> Column pivot indices, for 1 <= j <= N, column j
  265. *> of the matrix A was interchanged with column JPIV(j).
  266. *> \endverbatim
  267. *>
  268. *> \param[out] TAU
  269. *> \verbatim
  270. *> TAU is COMPLEX*16 array, dimension (min(M-IOFFSET,N))
  271. *> The scalar factors of the elementary reflectors.
  272. *> \endverbatim
  273. *>
  274. *> \param[in,out] VN1
  275. *> \verbatim
  276. *> VN1 is DOUBLE PRECISION array, dimension (N)
  277. *> The vector with the partial column norms.
  278. *> \endverbatim
  279. *>
  280. *> \param[in,out] VN2
  281. *> \verbatim
  282. *> VN2 is DOUBLE PRECISION array, dimension (N)
  283. *> The vector with the exact column norms.
  284. *> \endverbatim
  285. *>
  286. *> \param[out] AUXV
  287. *> \verbatim
  288. *> AUXV is COMPLEX*16 array, dimension (NB)
  289. *> Auxiliary vector.
  290. *> \endverbatim
  291. *>
  292. *> \param[out] F
  293. *> \verbatim
  294. *> F is COMPLEX*16 array, dimension (LDF,NB)
  295. *> Matrix F**H = L*(Y**H)*A.
  296. *> \endverbatim
  297. *>
  298. *> \param[in] LDF
  299. *> \verbatim
  300. *> LDF is INTEGER
  301. *> The leading dimension of the array F. LDF >= max(1,N+NRHS).
  302. *> \endverbatim
  303. *>
  304. *> \param[out] IWORK
  305. *> \verbatim
  306. *> IWORK is INTEGER array, dimension (N-1).
  307. *> Is a work array. ( IWORK is used to store indices
  308. *> of "bad" columns for norm downdating in the residual
  309. *> matrix ).
  310. *> \endverbatim
  311. *>
  312. *> \param[out] INFO
  313. *> \verbatim
  314. *> INFO is INTEGER
  315. *> 1) INFO = 0: successful exit.
  316. *> 2) If INFO = j_1, where 1 <= j_1 <= N, then NaN was
  317. *> detected and the routine stops the computation.
  318. *> The j_1-th column of the matrix A or the j_1-th
  319. *> element of array TAU contains the first occurrence
  320. *> of NaN in the factorization step KB+1 ( when KB columns
  321. *> have been factorized ).
  322. *>
  323. *> On exit:
  324. *> KB is set to the number of
  325. *> factorized columns without
  326. *> exception.
  327. *> MAXC2NRMK is set to NaN.
  328. *> RELMAXC2NRMK is set to NaN.
  329. *> TAU(KB+1:min(M,N)) is not set and contains undefined
  330. *> elements. If j_1=KB+1, TAU(KB+1)
  331. *> may contain NaN.
  332. *> 3) If INFO = j_2, where N+1 <= j_2 <= 2*N, then no NaN
  333. *> was detected, but +Inf (or -Inf) was detected and
  334. *> the routine continues the computation until completion.
  335. *> The (j_2-N)-th column of the matrix A contains the first
  336. *> occurrence of +Inf (or -Inf) in the actorization
  337. *> step KB+1 ( when KB columns have been factorized ).
  338. *> \endverbatim
  339. *
  340. * Authors:
  341. * ========
  342. *
  343. *> \author Univ. of Tennessee
  344. *> \author Univ. of California Berkeley
  345. *> \author Univ. of Colorado Denver
  346. *> \author NAG Ltd.
  347. *
  348. *> \ingroup laqp3rk
  349. *
  350. *> \par References:
  351. * ================
  352. *> [1] A Level 3 BLAS QR factorization algorithm with column pivoting developed in 1996.
  353. *> G. Quintana-Orti, Depto. de Informatica, Universidad Jaime I, Spain.
  354. *> X. Sun, Computer Science Dept., Duke University, USA.
  355. *> C. H. Bischof, Math. and Comp. Sci. Div., Argonne National Lab, USA.
  356. *> A BLAS-3 version of the QR factorization with column pivoting.
  357. *> LAPACK Working Note 114
  358. *> \htmlonly
  359. *> <a href="https://www.netlib.org/lapack/lawnspdf/lawn114.pdf">https://www.netlib.org/lapack/lawnspdf/lawn114.pdf</a>
  360. *> \endhtmlonly
  361. *> and in
  362. *> SIAM J. Sci. Comput., 19(5):1486-1494, Sept. 1998.
  363. *> \htmlonly
  364. *> <a href="https://doi.org/10.1137/S1064827595296732">https://doi.org/10.1137/S1064827595296732</a>
  365. *> \endhtmlonly
  366. *>
  367. *> [2] A partial column norm updating strategy developed in 2006.
  368. *> Z. Drmac and Z. Bujanovic, Dept. of Math., University of Zagreb, Croatia.
  369. *> On the failure of rank revealing QR factorization software – a case study.
  370. *> LAPACK Working Note 176.
  371. *> \htmlonly
  372. *> <a href="http://www.netlib.org/lapack/lawnspdf/lawn176.pdf">http://www.netlib.org/lapack/lawnspdf/lawn176.pdf</a>
  373. *> \endhtmlonly
  374. *> and in
  375. *> ACM Trans. Math. Softw. 35, 2, Article 12 (July 2008), 28 pages.
  376. *> \htmlonly
  377. *> <a href="https://doi.org/10.1145/1377612.1377616">https://doi.org/10.1145/1377612.1377616</a>
  378. *> \endhtmlonly
  379. *
  380. *> \par Contributors:
  381. * ==================
  382. *>
  383. *> \verbatim
  384. *>
  385. *> November 2023, Igor Kozachenko, James Demmel,
  386. *> Computer Science Division,
  387. *> University of California, Berkeley
  388. *>
  389. *> \endverbatim
  390. *
  391. * =====================================================================
  392. SUBROUTINE ZLAQP3RK( M, N, NRHS, IOFFSET, NB, ABSTOL,
  393. $ RELTOL, KP1, MAXC2NRM, A, LDA, DONE, KB,
  394. $ MAXC2NRMK, RELMAXC2NRMK, JPIV, TAU,
  395. $ VN1, VN2, AUXV, F, LDF, IWORK, INFO )
  396. IMPLICIT NONE
  397. *
  398. * -- LAPACK auxiliary routine --
  399. * -- LAPACK is a software package provided by Univ. of Tennessee, --
  400. * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  401. *
  402. * .. Scalar Arguments ..
  403. LOGICAL DONE
  404. INTEGER INFO, IOFFSET, KB, KP1, LDA, LDF, M, N,
  405. $ NB, NRHS
  406. DOUBLE PRECISION ABSTOL, MAXC2NRM, MAXC2NRMK, RELMAXC2NRMK,
  407. $ RELTOL
  408. * ..
  409. * .. Array Arguments ..
  410. INTEGER IWORK( * ), JPIV( * )
  411. DOUBLE PRECISION VN1( * ), VN2( * )
  412. COMPLEX*16 A( LDA, * ), AUXV( * ), F( LDF, * ), TAU( * )
  413. * ..
  414. *
  415. * =====================================================================
  416. *
  417. * .. Parameters ..
  418. DOUBLE PRECISION ZERO, ONE
  419. PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
  420. COMPLEX*16 CZERO, CONE
  421. PARAMETER ( CZERO = ( 0.0D+0, 0.0D+0 ),
  422. $ CONE = ( 1.0D+0, 0.0D+0 ) )
  423. * ..
  424. * .. Local Scalars ..
  425. INTEGER ITEMP, J, K, MINMNFACT, MINMNUPDT,
  426. $ LSTICC, KP, I, IF
  427. DOUBLE PRECISION HUGEVAL, TAUNAN, TEMP, TEMP2, TOL3Z
  428. COMPLEX*16 AIK
  429. * ..
  430. * .. External Subroutines ..
  431. EXTERNAL ZGEMM, ZGEMV, ZLARFG, ZSWAP
  432. * ..
  433. * .. Intrinsic Functions ..
  434. INTRINSIC ABS, DBLE, DCONJG, DIMAG, MAX, MIN, SQRT
  435. * ..
  436. * .. External Functions ..
  437. LOGICAL DISNAN
  438. INTEGER IDAMAX
  439. DOUBLE PRECISION DLAMCH, DZNRM2
  440. EXTERNAL DISNAN, DLAMCH, IDAMAX, DZNRM2
  441. * ..
  442. * .. Executable Statements ..
  443. *
  444. * Initialize INFO
  445. *
  446. INFO = 0
  447. *
  448. * MINMNFACT in the smallest dimension of the submatrix
  449. * A(IOFFSET+1:M,1:N) to be factorized.
  450. *
  451. MINMNFACT = MIN( M-IOFFSET, N )
  452. MINMNUPDT = MIN( M-IOFFSET, N+NRHS )
  453. NB = MIN( NB, MINMNFACT )
  454. TOL3Z = SQRT( DLAMCH( 'Epsilon' ) )
  455. HUGEVAL = DLAMCH( 'Overflow' )
  456. *
  457. * Compute factorization in a while loop over NB columns,
  458. * K is the column index in the block A(1:M,1:N).
  459. *
  460. K = 0
  461. LSTICC = 0
  462. DONE = .FALSE.
  463. *
  464. DO WHILE ( K.LT.NB .AND. LSTICC.EQ.0 )
  465. K = K + 1
  466. I = IOFFSET + K
  467. *
  468. IF( I.EQ.1 ) THEN
  469. *
  470. * We are at the first column of the original whole matrix A_orig,
  471. * therefore we use the computed KP1 and MAXC2NRM from the
  472. * main routine.
  473. *
  474. KP = KP1
  475. *
  476. ELSE
  477. *
  478. * Determine the pivot column in K-th step, i.e. the index
  479. * of the column with the maximum 2-norm in the
  480. * submatrix A(I:M,K:N).
  481. *
  482. KP = ( K-1 ) + IDAMAX( N-K+1, VN1( K ), 1 )
  483. *
  484. * Determine the maximum column 2-norm and the relative maximum
  485. * column 2-norm of the submatrix A(I:M,K:N) in step K.
  486. *
  487. MAXC2NRMK = VN1( KP )
  488. *
  489. * ============================================================
  490. *
  491. * Check if the submatrix A(I:M,K:N) contains NaN, set
  492. * INFO parameter to the column number, where the first NaN
  493. * is found and return from the routine.
  494. * We need to check the condition only if the
  495. * column index (same as row index) of the original whole
  496. * matrix is larger than 1, since the condition for whole
  497. * original matrix is checked in the main routine.
  498. *
  499. IF( DISNAN( MAXC2NRMK ) ) THEN
  500. *
  501. DONE = .TRUE.
  502. *
  503. * Set KB, the number of factorized partial columns
  504. * that are non-zero in each step in the block,
  505. * i.e. the rank of the factor R.
  506. * Set IF, the number of processed rows in the block, which
  507. * is the same as the number of processed rows in
  508. * the original whole matrix A_orig.
  509. *
  510. KB = K - 1
  511. IF = I - 1
  512. INFO = KB + KP
  513. *
  514. * Set RELMAXC2NRMK to NaN.
  515. *
  516. RELMAXC2NRMK = MAXC2NRMK
  517. *
  518. * There is no need to apply the block reflector to the
  519. * residual of the matrix A stored in A(KB+1:M,KB+1:N),
  520. * since the submatrix contains NaN and we stop
  521. * the computation.
  522. * But, we need to apply the block reflector to the residual
  523. * right hand sides stored in A(KB+1:M,N+1:N+NRHS), if the
  524. * residual right hand sides exist. This occurs
  525. * when ( NRHS != 0 AND KB <= (M-IOFFSET) ):
  526. *
  527. * A(I+1:M,N+1:N+NRHS) := A(I+1:M,N+1:N+NRHS) -
  528. * A(I+1:M,1:KB) * F(N+1:N+NRHS,1:KB)**H.
  529. IF( NRHS.GT.0 .AND. KB.LT.(M-IOFFSET) ) THEN
  530. CALL ZGEMM( 'No transpose', 'Conjugate transpose',
  531. $ M-IF, NRHS, KB, -CONE, A( IF+1, 1 ), LDA,
  532. $ F( N+1, 1 ), LDF, CONE, A( IF+1, N+1 ), LDA )
  533. END IF
  534. *
  535. * There is no need to recompute the 2-norm of the
  536. * difficult columns, since we stop the factorization.
  537. *
  538. * Array TAU(KF+1:MINMNFACT) is not set and contains
  539. * undefined elements.
  540. *
  541. * Return from the routine.
  542. *
  543. RETURN
  544. END IF
  545. *
  546. * Quick return, if the submatrix A(I:M,K:N) is
  547. * a zero matrix. We need to check it only if the column index
  548. * (same as row index) is larger than 1, since the condition
  549. * for the whole original matrix A_orig is checked in the main
  550. * routine.
  551. *
  552. IF( MAXC2NRMK.EQ.ZERO ) THEN
  553. *
  554. DONE = .TRUE.
  555. *
  556. * Set KB, the number of factorized partial columns
  557. * that are non-zero in each step in the block,
  558. * i.e. the rank of the factor R.
  559. * Set IF, the number of processed rows in the block, which
  560. * is the same as the number of processed rows in
  561. * the original whole matrix A_orig.
  562. *
  563. KB = K - 1
  564. IF = I - 1
  565. RELMAXC2NRMK = ZERO
  566. *
  567. * There is no need to apply the block reflector to the
  568. * residual of the matrix A stored in A(KB+1:M,KB+1:N),
  569. * since the submatrix is zero and we stop the computation.
  570. * But, we need to apply the block reflector to the residual
  571. * right hand sides stored in A(KB+1:M,N+1:N+NRHS), if the
  572. * residual right hand sides exist. This occurs
  573. * when ( NRHS != 0 AND KB <= (M-IOFFSET) ):
  574. *
  575. * A(I+1:M,N+1:N+NRHS) := A(I+1:M,N+1:N+NRHS) -
  576. * A(I+1:M,1:KB) * F(N+1:N+NRHS,1:KB)**H.
  577. *
  578. IF( NRHS.GT.0 .AND. KB.LT.(M-IOFFSET) ) THEN
  579. CALL ZGEMM( 'No transpose', 'Conjugate transpose',
  580. $ M-IF, NRHS, KB, -CONE, A( IF+1, 1 ), LDA,
  581. $ F( N+1, 1 ), LDF, CONE, A( IF+1, N+1 ), LDA )
  582. END IF
  583. *
  584. * There is no need to recompute the 2-norm of the
  585. * difficult columns, since we stop the factorization.
  586. *
  587. * Set TAUs corresponding to the columns that were not
  588. * factorized to ZERO, i.e. set TAU(KB+1:MINMNFACT) = CZERO,
  589. * which is equivalent to seting TAU(K:MINMNFACT) = CZERO.
  590. *
  591. DO J = K, MINMNFACT
  592. TAU( J ) = CZERO
  593. END DO
  594. *
  595. * Return from the routine.
  596. *
  597. RETURN
  598. *
  599. END IF
  600. *
  601. * ============================================================
  602. *
  603. * Check if the submatrix A(I:M,K:N) contains Inf,
  604. * set INFO parameter to the column number, where
  605. * the first Inf is found plus N, and continue
  606. * the computation.
  607. * We need to check the condition only if the
  608. * column index (same as row index) of the original whole
  609. * matrix is larger than 1, since the condition for whole
  610. * original matrix is checked in the main routine.
  611. *
  612. IF( INFO.EQ.0 .AND. MAXC2NRMK.GT.HUGEVAL ) THEN
  613. INFO = N + K - 1 + KP
  614. END IF
  615. *
  616. * ============================================================
  617. *
  618. * Test for the second and third tolerance stopping criteria.
  619. * NOTE: There is no need to test for ABSTOL.GE.ZERO, since
  620. * MAXC2NRMK is non-negative. Similarly, there is no need
  621. * to test for RELTOL.GE.ZERO, since RELMAXC2NRMK is
  622. * non-negative.
  623. * We need to check the condition only if the
  624. * column index (same as row index) of the original whole
  625. * matrix is larger than 1, since the condition for whole
  626. * original matrix is checked in the main routine.
  627. *
  628. RELMAXC2NRMK = MAXC2NRMK / MAXC2NRM
  629. *
  630. IF( MAXC2NRMK.LE.ABSTOL .OR. RELMAXC2NRMK.LE.RELTOL ) THEN
  631. *
  632. DONE = .TRUE.
  633. *
  634. * Set KB, the number of factorized partial columns
  635. * that are non-zero in each step in the block,
  636. * i.e. the rank of the factor R.
  637. * Set IF, the number of processed rows in the block, which
  638. * is the same as the number of processed rows in
  639. * the original whole matrix A_orig;
  640. *
  641. KB = K - 1
  642. IF = I - 1
  643. *
  644. * Apply the block reflector to the residual of the
  645. * matrix A and the residual of the right hand sides B, if
  646. * the residual matrix and and/or the residual of the right
  647. * hand sides exist, i.e. if the submatrix
  648. * A(I+1:M,KB+1:N+NRHS) exists. This occurs when
  649. * KB < MINMNUPDT = min( M-IOFFSET, N+NRHS ):
  650. *
  651. * A(IF+1:M,K+1:N+NRHS) := A(IF+1:M,KB+1:N+NRHS) -
  652. * A(IF+1:M,1:KB) * F(KB+1:N+NRHS,1:KB)**H.
  653. *
  654. IF( KB.LT.MINMNUPDT ) THEN
  655. CALL ZGEMM( 'No transpose', 'Conjugate transpose',
  656. $ M-IF, N+NRHS-KB, KB,-CONE, A( IF+1, 1 ), LDA,
  657. $ F( KB+1, 1 ), LDF, CONE, A( IF+1, KB+1 ), LDA )
  658. END IF
  659. *
  660. * There is no need to recompute the 2-norm of the
  661. * difficult columns, since we stop the factorization.
  662. *
  663. * Set TAUs corresponding to the columns that were not
  664. * factorized to ZERO, i.e. set TAU(KB+1:MINMNFACT) = CZERO,
  665. * which is equivalent to seting TAU(K:MINMNFACT) = CZERO.
  666. *
  667. DO J = K, MINMNFACT
  668. TAU( J ) = CZERO
  669. END DO
  670. *
  671. * Return from the routine.
  672. *
  673. RETURN
  674. *
  675. END IF
  676. *
  677. * ============================================================
  678. *
  679. * End ELSE of IF(I.EQ.1)
  680. *
  681. END IF
  682. *
  683. * ===============================================================
  684. *
  685. * If the pivot column is not the first column of the
  686. * subblock A(1:M,K:N):
  687. * 1) swap the K-th column and the KP-th pivot column
  688. * in A(1:M,1:N);
  689. * 2) swap the K-th row and the KP-th row in F(1:N,1:K-1)
  690. * 3) copy the K-th element into the KP-th element of the partial
  691. * and exact 2-norm vectors VN1 and VN2. (Swap is not needed
  692. * for VN1 and VN2 since we use the element with the index
  693. * larger than K in the next loop step.)
  694. * 4) Save the pivot interchange with the indices relative to the
  695. * the original matrix A_orig, not the block A(1:M,1:N).
  696. *
  697. IF( KP.NE.K ) THEN
  698. CALL ZSWAP( M, A( 1, KP ), 1, A( 1, K ), 1 )
  699. CALL ZSWAP( K-1, F( KP, 1 ), LDF, F( K, 1 ), LDF )
  700. VN1( KP ) = VN1( K )
  701. VN2( KP ) = VN2( K )
  702. ITEMP = JPIV( KP )
  703. JPIV( KP ) = JPIV( K )
  704. JPIV( K ) = ITEMP
  705. END IF
  706. *
  707. * Apply previous Householder reflectors to column K:
  708. * A(I:M,K) := A(I:M,K) - A(I:M,1:K-1)*F(K,1:K-1)**H.
  709. *
  710. IF( K.GT.1 ) THEN
  711. DO J = 1, K - 1
  712. F( K, J ) = DCONJG( F( K, J ) )
  713. END DO
  714. CALL ZGEMV( 'No transpose', M-I+1, K-1, -CONE, A( I, 1 ),
  715. $ LDA, F( K, 1 ), LDF, CONE, A( I, K ), 1 )
  716. DO J = 1, K - 1
  717. F( K, J ) = DCONJG( F( K, J ) )
  718. END DO
  719. END IF
  720. *
  721. * Generate elementary reflector H(k) using the column A(I:M,K).
  722. *
  723. IF( I.LT.M ) THEN
  724. CALL ZLARFG( M-I+1, A( I, K ), A( I+1, K ), 1, TAU( K ) )
  725. ELSE
  726. TAU( K ) = CZERO
  727. END IF
  728. *
  729. * Check if TAU(K) contains NaN, set INFO parameter
  730. * to the column number where NaN is found and return from
  731. * the routine.
  732. * NOTE: There is no need to check TAU(K) for Inf,
  733. * since ZLARFG cannot produce TAU(KK) or Householder vector
  734. * below the diagonal containing Inf. Only BETA on the diagonal,
  735. * returned by ZLARFG can contain Inf, which requires
  736. * TAU(K) to contain NaN. Therefore, this case of generating Inf
  737. * by ZLARFG is covered by checking TAU(K) for NaN.
  738. *
  739. IF( DISNAN( DBLE( TAU(K) ) ) ) THEN
  740. TAUNAN = DBLE( TAU(K) )
  741. ELSE IF( DISNAN( DIMAG( TAU(K) ) ) ) THEN
  742. TAUNAN = DIMAG( TAU(K) )
  743. ELSE
  744. TAUNAN = ZERO
  745. END IF
  746. *
  747. IF( DISNAN( TAUNAN ) ) THEN
  748. *
  749. DONE = .TRUE.
  750. *
  751. * Set KB, the number of factorized partial columns
  752. * that are non-zero in each step in the block,
  753. * i.e. the rank of the factor R.
  754. * Set IF, the number of processed rows in the block, which
  755. * is the same as the number of processed rows in
  756. * the original whole matrix A_orig.
  757. *
  758. KB = K - 1
  759. IF = I - 1
  760. INFO = K
  761. *
  762. * Set MAXC2NRMK and RELMAXC2NRMK to NaN.
  763. *
  764. MAXC2NRMK = TAUNAN
  765. RELMAXC2NRMK = TAUNAN
  766. *
  767. * There is no need to apply the block reflector to the
  768. * residual of the matrix A stored in A(KB+1:M,KB+1:N),
  769. * since the submatrix contains NaN and we stop
  770. * the computation.
  771. * But, we need to apply the block reflector to the residual
  772. * right hand sides stored in A(KB+1:M,N+1:N+NRHS), if the
  773. * residual right hand sides exist. This occurs
  774. * when ( NRHS != 0 AND KB <= (M-IOFFSET) ):
  775. *
  776. * A(I+1:M,N+1:N+NRHS) := A(I+1:M,N+1:N+NRHS) -
  777. * A(I+1:M,1:KB) * F(N+1:N+NRHS,1:KB)**H.
  778. *
  779. IF( NRHS.GT.0 .AND. KB.LT.(M-IOFFSET) ) THEN
  780. CALL ZGEMM( 'No transpose', 'Conjugate transpose',
  781. $ M-IF, NRHS, KB, -CONE, A( IF+1, 1 ), LDA,
  782. $ F( N+1, 1 ), LDF, CONE, A( IF+1, N+1 ), LDA )
  783. END IF
  784. *
  785. * There is no need to recompute the 2-norm of the
  786. * difficult columns, since we stop the factorization.
  787. *
  788. * Array TAU(KF+1:MINMNFACT) is not set and contains
  789. * undefined elements.
  790. *
  791. * Return from the routine.
  792. *
  793. RETURN
  794. END IF
  795. *
  796. * ===============================================================
  797. *
  798. AIK = A( I, K )
  799. A( I, K ) = CONE
  800. *
  801. * ===============================================================
  802. *
  803. * Compute the current K-th column of F:
  804. * 1) F(K+1:N,K) := tau(K) * A(I:M,K+1:N)**H * A(I:M,K).
  805. *
  806. IF( K.LT.N+NRHS ) THEN
  807. CALL ZGEMV( 'Conjugate transpose', M-I+1, N+NRHS-K,
  808. $ TAU( K ), A( I, K+1 ), LDA, A( I, K ), 1,
  809. $ CZERO, F( K+1, K ), 1 )
  810. END IF
  811. *
  812. * 2) Zero out elements above and on the diagonal of the
  813. * column K in matrix F, i.e elements F(1:K,K).
  814. *
  815. DO J = 1, K
  816. F( J, K ) = CZERO
  817. END DO
  818. *
  819. * 3) Incremental updating of the K-th column of F:
  820. * F(1:N,K) := F(1:N,K) - tau(K) * F(1:N,1:K-1) * A(I:M,1:K-1)**H
  821. * * A(I:M,K).
  822. *
  823. IF( K.GT.1 ) THEN
  824. CALL ZGEMV( 'Conjugate Transpose', M-I+1, K-1, -TAU( K ),
  825. $ A( I, 1 ), LDA, A( I, K ), 1, CZERO,
  826. $ AUXV( 1 ), 1 )
  827. *
  828. CALL ZGEMV( 'No transpose', N+NRHS, K-1, CONE,
  829. $ F( 1, 1 ), LDF, AUXV( 1 ), 1, CONE,
  830. $ F( 1, K ), 1 )
  831. END IF
  832. *
  833. * ===============================================================
  834. *
  835. * Update the current I-th row of A:
  836. * A(I,K+1:N+NRHS) := A(I,K+1:N+NRHS)
  837. * - A(I,1:K)*F(K+1:N+NRHS,1:K)**H.
  838. *
  839. IF( K.LT.N+NRHS ) THEN
  840. CALL ZGEMM( 'No transpose', 'Conjugate transpose',
  841. $ 1, N+NRHS-K, K, -CONE, A( I, 1 ), LDA,
  842. $ F( K+1, 1 ), LDF, CONE, A( I, K+1 ), LDA )
  843. END IF
  844. *
  845. A( I, K ) = AIK
  846. *
  847. * Update the partial column 2-norms for the residual matrix,
  848. * only if the residual matrix A(I+1:M,K+1:N) exists, i.e.
  849. * when K < MINMNFACT = min( M-IOFFSET, N ).
  850. *
  851. IF( K.LT.MINMNFACT ) THEN
  852. *
  853. DO J = K + 1, N
  854. IF( VN1( J ).NE.ZERO ) THEN
  855. *
  856. * NOTE: The following lines follow from the analysis in
  857. * Lapack Working Note 176.
  858. *
  859. TEMP = ABS( A( I, J ) ) / VN1( J )
  860. TEMP = MAX( ZERO, ( ONE+TEMP )*( ONE-TEMP ) )
  861. TEMP2 = TEMP*( VN1( J ) / VN2( J ) )**2
  862. IF( TEMP2.LE.TOL3Z ) THEN
  863. *
  864. * At J-index, we have a difficult column for the
  865. * update of the 2-norm. Save the index of the previous
  866. * difficult column in IWORK(J-1).
  867. * NOTE: ILSTCC > 1, threfore we can use IWORK only
  868. * with N-1 elements, where the elements are
  869. * shifted by 1 to the left.
  870. *
  871. IWORK( J-1 ) = LSTICC
  872. *
  873. * Set the index of the last difficult column LSTICC.
  874. *
  875. LSTICC = J
  876. *
  877. ELSE
  878. VN1( J ) = VN1( J )*SQRT( TEMP )
  879. END IF
  880. END IF
  881. END DO
  882. *
  883. END IF
  884. *
  885. * End of while loop.
  886. *
  887. END DO
  888. *
  889. * Now, afler the loop:
  890. * Set KB, the number of factorized columns in the block;
  891. * Set IF, the number of processed rows in the block, which
  892. * is the same as the number of processed rows in
  893. * the original whole matrix A_orig, IF = IOFFSET + KB.
  894. *
  895. KB = K
  896. IF = I
  897. *
  898. * Apply the block reflector to the residual of the matrix A
  899. * and the residual of the right hand sides B, if the residual
  900. * matrix and and/or the residual of the right hand sides
  901. * exist, i.e. if the submatrix A(I+1:M,KB+1:N+NRHS) exists.
  902. * This occurs when KB < MINMNUPDT = min( M-IOFFSET, N+NRHS ):
  903. *
  904. * A(IF+1:M,K+1:N+NRHS) := A(IF+1:M,KB+1:N+NRHS) -
  905. * A(IF+1:M,1:KB) * F(KB+1:N+NRHS,1:KB)**H.
  906. *
  907. IF( KB.LT.MINMNUPDT ) THEN
  908. CALL ZGEMM( 'No transpose', 'Conjugate transpose',
  909. $ M-IF, N+NRHS-KB, KB, -CONE, A( IF+1, 1 ), LDA,
  910. $ F( KB+1, 1 ), LDF, CONE, A( IF+1, KB+1 ), LDA )
  911. END IF
  912. *
  913. * Recompute the 2-norm of the difficult columns.
  914. * Loop over the index of the difficult columns from the largest
  915. * to the smallest index.
  916. *
  917. DO WHILE( LSTICC.GT.0 )
  918. *
  919. * LSTICC is the index of the last difficult column is greater
  920. * than 1.
  921. * ITEMP is the index of the previous difficult column.
  922. *
  923. ITEMP = IWORK( LSTICC-1 )
  924. *
  925. * Compute the 2-norm explicilty for the last difficult column and
  926. * save it in the partial and exact 2-norm vectors VN1 and VN2.
  927. *
  928. * NOTE: The computation of VN1( LSTICC ) relies on the fact that
  929. * DZNRM2 does not fail on vectors with norm below the value of
  930. * SQRT(DLAMCH('S'))
  931. *
  932. VN1( LSTICC ) = DZNRM2( M-IF, A( IF+1, LSTICC ), 1 )
  933. VN2( LSTICC ) = VN1( LSTICC )
  934. *
  935. * Downdate the index of the last difficult column to
  936. * the index of the previous difficult column.
  937. *
  938. LSTICC = ITEMP
  939. *
  940. END DO
  941. *
  942. RETURN
  943. *
  944. * End of ZLAQP3RK
  945. *
  946. END