ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
pzgels.f
Go to the documentation of this file.
1  SUBROUTINE pzgels( TRANS, M, N, NRHS, A, IA, JA, DESCA, B, IB, JB,
2  $ DESCB, WORK, LWORK, INFO )
3 *
4 * -- ScaLAPACK routine (version 1.7) --
5 * University of Tennessee, Knoxville, Oak Ridge National Laboratory,
6 * and University of California, Berkeley.
7 * May 1, 1997
8 *
9 * .. Scalar Arguments ..
10  CHARACTER TRANS
11  INTEGER IA, IB, INFO, JA, JB, LWORK, M, N, NRHS
12 * ..
13 * .. Array Arguments ..
14  INTEGER DESCA( * ), DESCB( * )
15  COMPLEX*16 A( * ), B( * ), WORK( * )
16 * ..
17 *
18 * Purpose
19 * =======
20 *
21 * PZGELS solves overdetermined or underdetermined complex linear
22 * systems involving an M-by-N matrix sub( A ) = A(IA:IA+M-1,JA:JA+N-1),
23 * or its conjugate-transpose, using a QR or LQ factorization of
24 * sub( A ). It is assumed that sub( A ) has full rank.
25 *
26 * The following options are provided:
27 *
28 * 1. If TRANS = 'N' and m >= n: find the least squares solution of
29 * an overdetermined system, i.e., solve the least squares problem
30 * minimize || sub( B ) - sub( A )*X ||.
31 *
32 * 2. If TRANS = 'N' and m < n: find the minimum norm solution of
33 * an underdetermined system sub( A ) * X = sub( B ).
34 *
35 * 3. If TRANS = 'C' and m >= n: find the minimum norm solution of
36 * an undetermined system sub( A )**H * X = sub( B ).
37 *
38 * 4. If TRANS = 'C' and m < n: find the least squares solution of
39 * an overdetermined system, i.e., solve the least squares problem
40 * minimize || sub( B ) - sub( A )**H * X ||.
41 *
42 * where sub( B ) denotes B( IB:IB+M-1, JB:JB+NRHS-1 ) when TRANS = 'N'
43 * and B( IB:IB+N-1, JB:JB+NRHS-1 ) otherwise. Several right hand side
44 * vectors b and solution vectors x can be handled in a single call;
45 * When TRANS = 'N', the solution vectors are stored as the columns of
46 * the N-by-NRHS right hand side matrix sub( B ) and the M-by-NRHS
47 * right hand side matrix sub( B ) otherwise.
48 *
49 * Notes
50 * =====
51 *
52 * Each global data object is described by an associated description
53 * vector. This vector stores the information required to establish
54 * the mapping between an object element and its corresponding process
55 * and memory location.
56 *
57 * Let A be a generic term for any 2D block cyclicly distributed array.
58 * Such a global array has an associated description vector DESCA.
59 * In the following comments, the character _ should be read as
60 * "of the global array".
61 *
62 * NOTATION STORED IN EXPLANATION
63 * --------------- -------------- --------------------------------------
64 * DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
65 * DTYPE_A = 1.
66 * CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
67 * the BLACS process grid A is distribu-
68 * ted over. The context itself is glo-
69 * bal, but the handle (the integer
70 * value) may vary.
71 * M_A (global) DESCA( M_ ) The number of rows in the global
72 * array A.
73 * N_A (global) DESCA( N_ ) The number of columns in the global
74 * array A.
75 * MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
76 * the rows of the array.
77 * NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
78 * the columns of the array.
79 * RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
80 * row of the array A is distributed.
81 * CSRC_A (global) DESCA( CSRC_ ) The process column over which the
82 * first column of the array A is
83 * distributed.
84 * LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
85 * array. LLD_A >= MAX(1,LOCr(M_A)).
86 *
87 * Let K be the number of rows or columns of a distributed matrix,
88 * and assume that its process grid has dimension p x q.
89 * LOCr( K ) denotes the number of elements of K that a process
90 * would receive if K were distributed over the p processes of its
91 * process column.
92 * Similarly, LOCc( K ) denotes the number of elements of K that a
93 * process would receive if K were distributed over the q processes of
94 * its process row.
95 * The values of LOCr() and LOCc() may be determined via a call to the
96 * ScaLAPACK tool function, NUMROC:
97 * LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
98 * LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
99 * An upper bound for these quantities may be computed by:
100 * LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
101 * LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
102 *
103 * Arguments
104 * =========
105 *
106 * TRANS (global input) CHARACTER
107 * = 'N': the linear system involves sub( A );
108 * = 'C': the linear system involves sub( A )**H.
109 *
110 * M (global input) INTEGER
111 * The number of rows to be operated on, i.e. the number of
112 * rows of the distributed submatrix sub( A ). M >= 0.
113 *
114 * N (global input) INTEGER
115 * The number of columns to be operated on, i.e. the number of
116 * columns of the distributed submatrix sub( A ). N >= 0.
117 *
118 * NRHS (global input) INTEGER
119 * The number of right hand sides, i.e. the number of columns
120 * of the distributed submatrices sub( B ) and X. NRHS >= 0.
121 *
122 * A (local input/local output) COMPLEX*16 pointer into the
123 * local memory to an array of local dimension
124 * ( LLD_A, LOCc(JA+N-1) ). On entry, the M-by-N matrix A.
125 * if M >= N, sub( A ) is overwritten by details of its QR
126 * factorization as returned by PZGEQRF;
127 * if M < N, sub( A ) is overwritten by details of its LQ
128 * factorization as returned by PZGELQF.
129 *
130 * IA (global input) INTEGER
131 * The row index in the global array A indicating the first
132 * row of sub( A ).
133 *
134 * JA (global input) INTEGER
135 * The column index in the global array A indicating the
136 * first column of sub( A ).
137 *
138 * DESCA (global and local input) INTEGER array of dimension DLEN_.
139 * The array descriptor for the distributed matrix A.
140 *
141 * B (local input/local output) COMPLEX*16 pointer into the
142 * local memory to an array of local dimension
143 * (LLD_B, LOCc(JB+NRHS-1)). On entry, this array contains the
144 * local pieces of the distributed matrix B of right hand side
145 * vectors, stored columnwise;
146 * sub( B ) is M-by-NRHS if TRANS='N', and N-by-NRHS otherwise.
147 * On exit, sub( B ) is overwritten by the solution vectors,
148 * stored columnwise: if TRANS = 'N' and M >= N, rows 1 to N
149 * of sub( B ) contain the least squares solution vectors; the
150 * residual sum of squares for the solution in each column is
151 * given by the sum of squares of elements N+1 to M in that
152 * column; if TRANS = 'N' and M < N, rows 1 to N of sub( B )
153 * contain the minimum norm solution vectors; if TRANS = 'C'
154 * and M >= N, rows 1 to M of sub( B ) contain the minimum norm
155 * solution vectors; if TRANS = 'C' and M < N, rows 1 to M of
156 * sub( B ) contain the least squares solution vectors; the
157 * residual sum of squares for the solution in each column is
158 * given by the sum of squares of elements M+1 to N in that
159 * column.
160 *
161 * IB (global input) INTEGER
162 * The row index in the global array B indicating the first
163 * row of sub( B ).
164 *
165 * JB (global input) INTEGER
166 * The column index in the global array B indicating the
167 * first column of sub( B ).
168 *
169 * DESCB (global and local input) INTEGER array of dimension DLEN_.
170 * The array descriptor for the distributed matrix B.
171 *
172 * WORK (local workspace/local output) COMPLEX*16 array,
173 * dimension (LWORK)
174 * On exit, WORK(1) returns the minimal and optimal LWORK.
175 *
176 * LWORK (local or global input) INTEGER
177 * The dimension of the array WORK.
178 * LWORK is local input and must be at least
179 * LWORK >= LTAU + MAX( LWF, LWS ) where
180 * If M >= N, then
181 * LTAU = NUMROC( JA+MIN(M,N)-1, NB_A, MYCOL, CSRC_A, NPCOL ),
182 * LWF = NB_A * ( MpA0 + NqA0 + NB_A )
183 * LWS = MAX( (NB_A*(NB_A-1))/2, (NRHSqB0 + MpB0)*NB_A ) +
184 * NB_A * NB_A
185 * Else
186 * LTAU = NUMROC( IA+MIN(M,N)-1, MB_A, MYROW, RSRC_A, NPROW ),
187 * LWF = MB_A * ( MpA0 + NqA0 + MB_A )
188 * LWS = MAX( (MB_A*(MB_A-1))/2, ( NpB0 + MAX( NqA0 +
189 * NUMROC( NUMROC( N+IROFFB, MB_A, 0, 0, NPROW ),
190 * MB_A, 0, 0, LCMP ), NRHSqB0 ) )*MB_A ) +
191 * MB_A * MB_A
192 * End if
193 *
194 * where LCMP = LCM / NPROW with LCM = ILCM( NPROW, NPCOL ),
195 *
196 * IROFFA = MOD( IA-1, MB_A ), ICOFFA = MOD( JA-1, NB_A ),
197 * IAROW = INDXG2P( IA, MB_A, MYROW, RSRC_A, NPROW ),
198 * IACOL = INDXG2P( JA, NB_A, MYCOL, CSRC_A, NPCOL ),
199 * MpA0 = NUMROC( M+IROFFA, MB_A, MYROW, IAROW, NPROW ),
200 * NqA0 = NUMROC( N+ICOFFA, NB_A, MYCOL, IACOL, NPCOL ),
201 *
202 * IROFFB = MOD( IB-1, MB_B ), ICOFFB = MOD( JB-1, NB_B ),
203 * IBROW = INDXG2P( IB, MB_B, MYROW, RSRC_B, NPROW ),
204 * IBCOL = INDXG2P( JB, NB_B, MYCOL, CSRC_B, NPCOL ),
205 * MpB0 = NUMROC( M+IROFFB, MB_B, MYROW, IBROW, NPROW ),
206 * NpB0 = NUMROC( N+IROFFB, MB_B, MYROW, IBROW, NPROW ),
207 * NRHSqB0 = NUMROC( NRHS+ICOFFB, NB_B, MYCOL, IBCOL, NPCOL ),
208 *
209 * ILCM, INDXG2P and NUMROC are ScaLAPACK tool functions;
210 * MYROW, MYCOL, NPROW and NPCOL can be determined by calling
211 * the subroutine BLACS_GRIDINFO.
212 *
213 * If LWORK = -1, then LWORK is global input and a workspace
214 * query is assumed; the routine only calculates the minimum
215 * and optimal size for all work arrays. Each of these
216 * values is returned in the first entry of the corresponding
217 * work array, and no error message is issued by PXERBLA.
218 *
219 * INFO (global output) INTEGER
220 * = 0: successful exit
221 * < 0: If the i-th argument is an array and the j-entry had
222 * an illegal value, then INFO = -(i*100+j), if the i-th
223 * argument is a scalar and had an illegal value, then
224 * INFO = -i.
225 *
226 * =====================================================================
227 *
228 * .. Parameters ..
229  INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
230  $ lld_, mb_, m_, nb_, n_, rsrc_
231  parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
232  $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
233  $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
234  DOUBLE PRECISION ZERO, ONE
235  parameter( zero = 0.0d+0, one = 1.0d+0 )
236  COMPLEX*16 CZERO, CONE
237  parameter( czero = ( 0.0d+0, 0.0d+0 ),
238  $ cone = ( 1.0d+0, 0.0d+0 ) )
239 * ..
240 * .. Local Scalars ..
241  LOGICAL LQUERY, TPSD
242  INTEGER BROW, IACOL, IAROW, IASCL, IBCOL, IBROW, IBSCL,
243  $ icoffa, icoffb, ictxt, ipw, iroffa, iroffb,
244  $ lcm, lcmp, ltau, lwf, lwmin, lws, mpa0, mpb0,
245  $ mycol, myrow, npb0, npcol, nprow, nqa0,
246  $ nrhsqb0, scllen
247  DOUBLE PRECISION ANRM, BIGNUM, BNRM, SMLNUM
248 * ..
249 * .. Local Arrays ..
250  INTEGER IDUM1( 2 ), IDUM2( 2 )
251  DOUBLE PRECISION RWORK( 1 )
252 * ..
253 * .. External Functions ..
254  LOGICAL LSAME
255  INTEGER ILCM
256  INTEGER INDXG2P, NUMROC
257  DOUBLE PRECISION PDLAMCH, PZLANGE
258  EXTERNAL ilcm, indxg2p, lsame, numroc, pdlamch,
259  $ pzlange
260 * ..
261 * .. External Subroutines ..
262  EXTERNAL blacs_gridinfo, chk1mat, pchk2mat, pzgelqf,
264  $ pztrsm, pzunmlq, pzunmqr, pxerbla
265 * ..
266 * .. Intrinsic Functions ..
267  INTRINSIC dble, dcmplx, ichar, max, min, mod
268 * ..
269 * .. Executable Statements ..
270 *
271 * Get grid parameters
272 *
273  ictxt = desca( ctxt_ )
274  CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
275 *
276 * Test the input parameters
277 *
278  info = 0
279  IF( nprow.EQ.-1 ) THEN
280  info = -( 800 + ctxt_ )
281  ELSE
282  CALL chk1mat( m, 2, n, 3, ia, ja, desca, 8, info )
283  IF ( m .GE. n ) THEN
284  CALL chk1mat( m, 2, nrhs, 4, ib, jb, descb, 12, info )
285  ELSE
286  CALL chk1mat( n, 3, nrhs, 4, ib, jb, descb, 12, info )
287  ENDIF
288  IF( info.EQ.0 ) THEN
289  iroffa = mod( ia-1, desca( mb_ ) )
290  icoffa = mod( ja-1, desca( nb_ ) )
291  iroffb = mod( ib-1, descb( mb_ ) )
292  icoffb = mod( jb-1, descb( nb_ ) )
293  iarow = indxg2p( ia, desca( mb_ ), myrow, desca( rsrc_ ),
294  $ nprow )
295  iacol = indxg2p( ia, desca( nb_ ), mycol, desca( csrc_ ),
296  $ npcol )
297  mpa0 = numroc( m+iroffa, desca( mb_ ), myrow, iarow, nprow )
298  nqa0 = numroc( n+icoffa, desca( nb_ ), mycol, iacol, npcol )
299 *
300  ibrow = indxg2p( ib, descb( mb_ ), myrow, descb( rsrc_ ),
301  $ nprow )
302  ibcol = indxg2p( ib, descb( nb_ ), mycol, descb( csrc_ ),
303  $ npcol )
304  nrhsqb0 = numroc( nrhs+icoffb, descb( nb_ ), mycol, ibcol,
305  $ npcol )
306  IF( m.GE.n ) THEN
307  mpb0 = numroc( m+iroffb, descb( mb_ ), myrow, ibrow,
308  $ nprow )
309  ltau = numroc( ja+min(m,n)-1, desca( nb_ ), mycol,
310  $ desca( csrc_ ), npcol )
311  lwf = desca( nb_ ) * ( mpa0 + nqa0 + desca( nb_ ) )
312  lws = max( ( desca( nb_ )*( desca( nb_ ) - 1 ) ) / 2,
313  $ ( mpb0 + nrhsqb0 ) * desca( nb_ ) ) +
314  $ desca( nb_ )*desca( nb_ )
315  ELSE
316  lcm = ilcm( nprow, npcol )
317  lcmp = lcm / nprow
318  npb0 = numroc( n+iroffb, descb( mb_ ), myrow, ibrow,
319  $ nprow )
320  ltau = numroc( ia+min(m,n)-1, desca( mb_ ), myrow,
321  $ desca( rsrc_ ), nprow )
322  lwf = desca( mb_ ) * ( mpa0 + nqa0 + desca( mb_ ) )
323  lws = max( ( desca( mb_ )*( desca( mb_ ) - 1 ) ) / 2,
324  $ ( npb0 + max( nqa0 + numroc( numroc( n+iroffb,
325  $ desca( mb_ ), 0, 0, nprow ), desca( mb_ ), 0, 0,
326  $ lcmp ), nrhsqb0 ) )*desca( mb_ ) ) +
327  $ desca( mb_ ) * desca( mb_ )
328  END IF
329  lwmin = ltau + max( lwf, lws )
330  work( 1 ) = dcmplx( dble( lwmin ) )
331  lquery = ( lwork.EQ.-1 )
332 *
333  tpsd = .true.
334  IF( lsame( trans, 'N' ) )
335  $ tpsd = .false.
336 *
337  IF( .NOT.( lsame( trans, 'N' ) .OR.
338  $ lsame( trans, 'C' ) ) ) THEN
339  info = -1
340  ELSE IF( m.LT.0 ) THEN
341  info = -2
342  ELSE IF( n.LT.0 ) THEN
343  info = -3
344  ELSE IF( nrhs.LT.0 ) THEN
345  info = -4
346  ELSE IF( m.GE.n .AND. iroffa.NE.iroffb ) THEN
347  info = -10
348  ELSE IF( m.GE.n .AND. iarow.NE.ibrow ) THEN
349  info = -10
350  ELSE IF( m.LT.n .AND. icoffa.NE.iroffb ) THEN
351  info = -10
352  ELSE IF( m.GE.n .AND. desca( mb_ ).NE.descb( mb_ ) ) THEN
353  info = -( 1200 + mb_ )
354  ELSE IF( m.LT.n .AND. desca( nb_ ).NE.descb( mb_ ) ) THEN
355  info = -( 1200 + mb_ )
356  ELSE IF( ictxt.NE.descb( ctxt_ ) ) THEN
357  info = -( 1200 + ctxt_ )
358  ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
359  info = -14
360  END IF
361  END IF
362 *
363  IF( .NOT.tpsd ) THEN
364  idum1( 1 ) = ichar( 'N' )
365  ELSE
366  idum1( 1 ) = ichar( 'C' )
367  END IF
368  idum2( 1 ) = 1
369  IF( lwork.EQ.-1 ) THEN
370  idum1( 2 ) = -1
371  ELSE
372  idum1( 2 ) = 1
373  END IF
374  idum2( 2 ) = 14
375  CALL pchk2mat( m, 2, n, 3, ia, ja, desca, 8, n, 3, nrhs, 4,
376  $ ib, jb, descb, 12, 2, idum1, idum2, info )
377  END IF
378 *
379  IF( info.NE.0 ) THEN
380  CALL pxerbla( ictxt, 'PZGELS', -info )
381  RETURN
382  ELSE IF( lquery ) THEN
383  RETURN
384  END IF
385 *
386 * Quick return if possible
387 *
388  IF( min( m, n, nrhs ).EQ.0 ) THEN
389  CALL pzlaset( 'Full', max( m, n ), nrhs, czero, czero, b,
390  $ ib, jb, descb )
391  RETURN
392  END IF
393 *
394 * Get machine parameters
395 *
396  smlnum = pdlamch( ictxt, 'S' )
397  smlnum = smlnum / pdlamch( ictxt, 'P' )
398  bignum = one / smlnum
399  CALL pdlabad( ictxt, smlnum, bignum )
400 *
401 * Scale A, B if max entry outside range [SMLNUM,BIGNUM]
402 *
403  anrm = pzlange( 'M', m, n, a, ia, ja, desca, rwork )
404  iascl = 0
405  IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
406 *
407 * Scale matrix norm up to SMLNUM
408 *
409  CALL pzlascl( 'G', anrm, smlnum, m, n, a, ia, ja, desca,
410  $ info )
411  iascl = 1
412  ELSE IF( anrm.GT.bignum ) THEN
413 *
414 * Scale matrix norm down to BIGNUM
415 *
416  CALL pzlascl( 'G', anrm, bignum, m, n, a, ia, ja, desca,
417  $ info )
418  iascl = 2
419  ELSE IF( anrm.EQ.zero ) THEN
420 *
421 * Matrix all zero. Return zero solution.
422 *
423  CALL pzlaset( 'F', max( m, n ), nrhs, czero, czero, b, ib,
424  $ jb, descb )
425  GO TO 10
426  END IF
427 *
428  brow = m
429  IF( tpsd )
430  $ brow = n
431 *
432  bnrm = pzlange( 'M', brow, nrhs, b, ib, jb, descb, rwork )
433 *
434  ibscl = 0
435  IF( bnrm.GT.zero .AND. bnrm.LT.smlnum ) THEN
436 *
437 * Scale matrix norm up to SMLNUM
438 *
439  CALL pzlascl( 'G', bnrm, smlnum, brow, nrhs, b, ib, jb,
440  $ descb, info )
441  ibscl = 1
442  ELSE IF( bnrm.GT.bignum ) THEN
443 *
444 * Scale matrix norm down to BIGNUM
445 *
446  CALL pzlascl( 'G', bnrm, bignum, brow, nrhs, b, ib, jb,
447  $ descb, info )
448  ibscl = 2
449  END IF
450 *
451  ipw = ltau + 1
452 *
453  IF( m.GE.n ) THEN
454 *
455 * compute QR factorization of A
456 *
457  CALL pzgeqrf( m, n, a, ia, ja, desca, work, work( ipw ),
458  $ lwork-ltau, info )
459 *
460 * workspace at least N, optimally N*NB
461 *
462  IF( .NOT.tpsd ) THEN
463 *
464 * Least-Squares Problem min || A * X - B ||
465 *
466 * B(IB:IB+M-1,JB:JB+NRHS-1) := Q' * B(IB:IB+M-1,JB:JB+NRHS-1)
467 *
468  CALL pzunmqr( 'Left', 'Conjugate transpose', m, nrhs, n, a,
469  $ ia, ja, desca, work, b, ib, jb, descb,
470  $ work( ipw ), lwork-ltau, info )
471 *
472 * workspace at least NRHS, optimally NRHS*NB
473 *
474 * B(IB:IB+N-1,JB:JB+NRHS-1) := inv(R) *
475 * B(IB:IB+N-1,JB:JB+NRHS-1)
476 *
477  CALL pztrsm( 'Left', 'Upper', 'No transpose', 'Non-unit', n,
478  $ nrhs, cone, a, ia, ja, desca, b, ib, jb,
479  $ descb )
480 *
481  scllen = n
482 *
483  ELSE
484 *
485 * Overdetermined system of equations sub( A )' * X = sub( B )
486 *
487 * sub( B ) := inv(R') * sub( B )
488 *
489  CALL pztrsm( 'Left', 'Upper', 'Conjugate transpose',
490  $ 'Non-unit', n, nrhs, cone, a, ia, ja, desca,
491  $ b, ib, jb, descb )
492 *
493 * B(IB+N:IB+M-1,JB:JB+NRHS-1) = ZERO
494 *
495  CALL pzlaset( 'All', m-n, nrhs, czero, czero, b, ib+n, jb,
496  $ descb )
497 *
498 * B(IB:IB+M-1,JB:JB+NRHS-1) := Q(1:N,:) *
499 * B(IB:IB+N-1,JB:JB+NRHS-1)
500 *
501  CALL pzunmqr( 'Left', 'No transpose', m, nrhs, n, a, ia, ja,
502  $ desca, work, b, ib, jb, descb, work( ipw ),
503  $ lwork-ltau, info )
504 *
505 * workspace at least NRHS, optimally NRHS*NB
506 *
507  scllen = m
508 *
509  END IF
510 *
511  ELSE
512 *
513 * Compute LQ factorization of sub( A )
514 *
515  CALL pzgelqf( m, n, a, ia, ja, desca, work, work( ipw ),
516  $ lwork-ltau, info )
517 *
518 * workspace at least M, optimally M*NB.
519 *
520  IF( .NOT.tpsd ) THEN
521 *
522 * underdetermined system of equations sub( A ) * X = sub( B )
523 *
524 * B(IB:IB+M-1,JB:JB+NRHS-1) := inv(L) *
525 * B(IB:IB+M-1,JB:JB+NRHS-1)
526 *
527  CALL pztrsm( 'Left', 'Lower', 'No transpose', 'Non-unit', m,
528  $ nrhs, cone, a, ia, ja, desca, b, ib, jb,
529  $ descb )
530 *
531 * B(IB+M:IB+N-1,JB:JB+NRHS-1) = 0
532 *
533  CALL pzlaset( 'All', n-m, nrhs, czero, czero, b, ib+m, jb,
534  $ descb )
535 *
536 * B(IB:IB+N-1,JB:JB+NRHS-1) := Q(1:N,:)' *
537 * B(IB:IB+M-1,JB:JB+NRHS-1)
538 *
539  CALL pzunmlq( 'Left', 'Conjugate transpose', n, nrhs, m, a,
540  $ ia, ja, desca, work, b, ib, jb, descb,
541  $ work( ipw ), lwork-ltau, info )
542 *
543 * workspace at least NRHS, optimally NRHS*NB
544 *
545  scllen = n
546 *
547  ELSE
548 *
549 * overdetermined system min || A' * X - B ||
550 *
551 * B(IB:IB+N-1,JB:JB+NRHS-1) := Q * B(IB:IB+N-1,JB:JB+NRHS-1)
552 *
553  CALL pzunmlq( 'Left', 'No transpose', n, nrhs, m, a, ia, ja,
554  $ desca, work, b, ib, jb, descb, work( ipw ),
555  $ lwork-ltau, info )
556 *
557 * workspace at least NRHS, optimally NRHS*NB
558 *
559 * B(IB:IB+M-1,JB:JB+NRHS-1) := inv(L') *
560 * B(IB:IB+M-1,JB:JB+NRHS-1)
561 *
562  CALL pztrsm( 'Left', 'Lower', 'Conjugate transpose',
563  $ 'Non-unit', m, nrhs, cone, a, ia, ja, desca,
564  $ b, ib, jb, descb )
565 *
566  scllen = m
567 *
568  END IF
569 *
570  END IF
571 *
572 * Undo scaling
573 *
574  IF( iascl.EQ.1 ) THEN
575  CALL pzlascl( 'G', anrm, smlnum, scllen, nrhs, b, ib, jb,
576  $ descb, info )
577  ELSE IF( iascl.EQ.2 ) THEN
578  CALL pzlascl( 'G', anrm, bignum, scllen, nrhs, b, ib, jb,
579  $ descb, info )
580  END IF
581  IF( ibscl.EQ.1 ) THEN
582  CALL pzlascl( 'G', smlnum, bnrm, scllen, nrhs, b, ib, jb,
583  $ descb, info )
584  ELSE IF( ibscl.EQ.2 ) THEN
585  CALL pzlascl( 'G', bignum, bnrm, scllen, nrhs, b, ib, jb,
586  $ descb, info )
587  END IF
588 *
589  10 CONTINUE
590 *
591  work( 1 ) = dcmplx( dble( lwmin ) )
592 *
593  RETURN
594 *
595 * End of PZGELS
596 *
597  END
pzgeqrf
subroutine pzgeqrf(M, N, A, IA, JA, DESCA, TAU, WORK, LWORK, INFO)
Definition: pzgeqrf.f:3
pzunmlq
subroutine pzunmlq(SIDE, TRANS, M, N, K, A, IA, JA, DESCA, TAU, C, IC, JC, DESCC, WORK, LWORK, INFO)
Definition: pzunmlq.f:3
max
#define max(A, B)
Definition: pcgemr.c:180
pzunmqr
subroutine pzunmqr(SIDE, TRANS, M, N, K, A, IA, JA, DESCA, TAU, C, IC, JC, DESCC, WORK, LWORK, INFO)
Definition: pzunmqr.f:3
pzgelqf
subroutine pzgelqf(M, N, A, IA, JA, DESCA, TAU, WORK, LWORK, INFO)
Definition: pzgelqf.f:3
pdlabad
subroutine pdlabad(ICTXT, SMALL, LARGE)
Definition: pdlabad.f:2
pzlascl
subroutine pzlascl(TYPE, CFROM, CTO, M, N, A, IA, JA, DESCA, INFO)
Definition: pzlascl.f:3
pzlaset
subroutine pzlaset(UPLO, M, N, ALPHA, BETA, A, IA, JA, DESCA)
Definition: pzblastst.f:7509
pchk2mat
subroutine pchk2mat(MA, MAPOS0, NA, NAPOS0, IA, JA, DESCA, DESCAPOS0, MB, MBPOS0, NB, NBPOS0, IB, JB, DESCB, DESCBPOS0, NEXTRA, EX, EXPOS, INFO)
Definition: pchkxmat.f:175
pzgels
subroutine pzgels(TRANS, M, N, NRHS, A, IA, JA, DESCA, B, IB, JB, DESCB, WORK, LWORK, INFO)
Definition: pzgels.f:3
chk1mat
subroutine chk1mat(MA, MAPOS0, NA, NAPOS0, IA, JA, DESCA, DESCAPOS0, INFO)
Definition: chk1mat.f:3
pxerbla
subroutine pxerbla(ICTXT, SRNAME, INFO)
Definition: pxerbla.f:2
min
#define min(A, B)
Definition: pcgemr.c:181