LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssyevx_2stage.f
Go to the documentation of this file.
1*> \brief <b> SSYEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices</b>
2*
3* @generated from dsyevx_2stage.f, fortran d -> s, Sat Nov 5 23:55:46 2016
4*
5* =========== DOCUMENTATION ===========
6*
7* Online html documentation available at
8* http://www.netlib.org/lapack/explore-html/
9*
10*> Download SSYEVX_2STAGE + dependencies
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssyevx_2stage.f">
12*> [TGZ]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssyevx_2stage.f">
14*> [ZIP]</a>
15*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssyevx_2stage.f">
16*> [TXT]</a>
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SSYEVX_2STAGE( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
22* IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
23* LWORK, IWORK, IFAIL, INFO )
24*
25* IMPLICIT NONE
26*
27* .. Scalar Arguments ..
28* CHARACTER JOBZ, RANGE, UPLO
29* INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
30* REAL ABSTOL, VL, VU
31* ..
32* .. Array Arguments ..
33* INTEGER IFAIL( * ), IWORK( * )
34* REAL A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> SSYEVX_2STAGE computes selected eigenvalues and, optionally, eigenvectors
44*> of a real symmetric matrix A using the 2stage technique for
45*> the reduction to tridiagonal. Eigenvalues and eigenvectors can be
46*> selected by specifying either a range of values or a range of indices
47*> for the desired eigenvalues.
48*> \endverbatim
49*
50* Arguments:
51* ==========
52*
53*> \param[in] JOBZ
54*> \verbatim
55*> JOBZ is CHARACTER*1
56*> = 'N': Compute eigenvalues only;
57*> = 'V': Compute eigenvalues and eigenvectors.
58*> Not available in this release.
59*> \endverbatim
60*>
61*> \param[in] RANGE
62*> \verbatim
63*> RANGE is CHARACTER*1
64*> = 'A': all eigenvalues will be found.
65*> = 'V': all eigenvalues in the half-open interval (VL,VU]
66*> will be found.
67*> = 'I': the IL-th through IU-th eigenvalues will be found.
68*> \endverbatim
69*>
70*> \param[in] UPLO
71*> \verbatim
72*> UPLO is CHARACTER*1
73*> = 'U': Upper triangle of A is stored;
74*> = 'L': Lower triangle of A is stored.
75*> \endverbatim
76*>
77*> \param[in] N
78*> \verbatim
79*> N is INTEGER
80*> The order of the matrix A. N >= 0.
81*> \endverbatim
82*>
83*> \param[in,out] A
84*> \verbatim
85*> A is REAL array, dimension (LDA, N)
86*> On entry, the symmetric matrix A. If UPLO = 'U', the
87*> leading N-by-N upper triangular part of A contains the
88*> upper triangular part of the matrix A. If UPLO = 'L',
89*> the leading N-by-N lower triangular part of A contains
90*> the lower triangular part of the matrix A.
91*> On exit, the lower triangle (if UPLO='L') or the upper
92*> triangle (if UPLO='U') of A, including the diagonal, is
93*> destroyed.
94*> \endverbatim
95*>
96*> \param[in] LDA
97*> \verbatim
98*> LDA is INTEGER
99*> The leading dimension of the array A. LDA >= max(1,N).
100*> \endverbatim
101*>
102*> \param[in] VL
103*> \verbatim
104*> VL is REAL
105*> If RANGE='V', the lower bound of the interval to
106*> be searched for eigenvalues. VL < VU.
107*> Not referenced if RANGE = 'A' or 'I'.
108*> \endverbatim
109*>
110*> \param[in] VU
111*> \verbatim
112*> VU is REAL
113*> If RANGE='V', the upper bound of the interval to
114*> be searched for eigenvalues. VL < VU.
115*> Not referenced if RANGE = 'A' or 'I'.
116*> \endverbatim
117*>
118*> \param[in] IL
119*> \verbatim
120*> IL is INTEGER
121*> If RANGE='I', the index of the
122*> smallest eigenvalue to be returned.
123*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
124*> Not referenced if RANGE = 'A' or 'V'.
125*> \endverbatim
126*>
127*> \param[in] IU
128*> \verbatim
129*> IU is INTEGER
130*> If RANGE='I', the index of the
131*> largest eigenvalue to be returned.
132*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
133*> Not referenced if RANGE = 'A' or 'V'.
134*> \endverbatim
135*>
136*> \param[in] ABSTOL
137*> \verbatim
138*> ABSTOL is REAL
139*> The absolute error tolerance for the eigenvalues.
140*> An approximate eigenvalue is accepted as converged
141*> when it is determined to lie in an interval [a,b]
142*> of width less than or equal to
143*>
144*> ABSTOL + EPS * max( |a|,|b| ) ,
145*>
146*> where EPS is the machine precision. If ABSTOL is less than
147*> or equal to zero, then EPS*|T| will be used in its place,
148*> where |T| is the 1-norm of the tridiagonal matrix obtained
149*> by reducing A to tridiagonal form.
150*>
151*> Eigenvalues will be computed most accurately when ABSTOL is
152*> set to twice the underflow threshold 2*SLAMCH('S'), not zero.
153*> If this routine returns with INFO>0, indicating that some
154*> eigenvectors did not converge, try setting ABSTOL to
155*> 2*SLAMCH('S').
156*>
157*> See "Computing Small Singular Values of Bidiagonal Matrices
158*> with Guaranteed High Relative Accuracy," by Demmel and
159*> Kahan, LAPACK Working Note #3.
160*> \endverbatim
161*>
162*> \param[out] M
163*> \verbatim
164*> M is INTEGER
165*> The total number of eigenvalues found. 0 <= M <= N.
166*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
167*> \endverbatim
168*>
169*> \param[out] W
170*> \verbatim
171*> W is REAL array, dimension (N)
172*> On normal exit, the first M elements contain the selected
173*> eigenvalues in ascending order.
174*> \endverbatim
175*>
176*> \param[out] Z
177*> \verbatim
178*> Z is REAL array, dimension (LDZ, max(1,M))
179*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
180*> contain the orthonormal eigenvectors of the matrix A
181*> corresponding to the selected eigenvalues, with the i-th
182*> column of Z holding the eigenvector associated with W(i).
183*> If an eigenvector fails to converge, then that column of Z
184*> contains the latest approximation to the eigenvector, and the
185*> index of the eigenvector is returned in IFAIL.
186*> If JOBZ = 'N', then Z is not referenced.
187*> Note: the user must ensure that at least max(1,M) columns are
188*> supplied in the array Z; if RANGE = 'V', the exact value of M
189*> is not known in advance and an upper bound must be used.
190*> \endverbatim
191*>
192*> \param[in] LDZ
193*> \verbatim
194*> LDZ is INTEGER
195*> The leading dimension of the array Z. LDZ >= 1, and if
196*> JOBZ = 'V', LDZ >= max(1,N).
197*> \endverbatim
198*>
199*> \param[out] WORK
200*> \verbatim
201*> WORK is REAL array, dimension (MAX(1,LWORK))
202*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
203*> \endverbatim
204*>
205*> \param[in] LWORK
206*> \verbatim
207*> LWORK is INTEGER
208*> The length of the array WORK. LWORK >= 1, when N <= 1;
209*> otherwise
210*> If JOBZ = 'N' and N > 1, LWORK must be queried.
211*> LWORK = MAX(1, 8*N, dimension) where
212*> dimension = max(stage1,stage2) + (KD+1)*N + 3*N
213*> = N*KD + N*max(KD+1,FACTOPTNB)
214*> + max(2*KD*KD, KD*NTHREADS)
215*> + (KD+1)*N + 3*N
216*> where KD is the blocking size of the reduction,
217*> FACTOPTNB is the blocking used by the QR or LQ
218*> algorithm, usually FACTOPTNB=128 is a good choice
219*> NTHREADS is the number of threads used when
220*> openMP compilation is enabled, otherwise =1.
221*> If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
222*>
223*> If LWORK = -1, then a workspace query is assumed; the routine
224*> only calculates the optimal size of the WORK array, returns
225*> this value as the first entry of the WORK array, and no error
226*> message related to LWORK is issued by XERBLA.
227*> \endverbatim
228*>
229*> \param[out] IWORK
230*> \verbatim
231*> IWORK is INTEGER array, dimension (5*N)
232*> \endverbatim
233*>
234*> \param[out] IFAIL
235*> \verbatim
236*> IFAIL is INTEGER array, dimension (N)
237*> If JOBZ = 'V', then if INFO = 0, the first M elements of
238*> IFAIL are zero. If INFO > 0, then IFAIL contains the
239*> indices of the eigenvectors that failed to converge.
240*> If JOBZ = 'N', then IFAIL is not referenced.
241*> \endverbatim
242*>
243*> \param[out] INFO
244*> \verbatim
245*> INFO is INTEGER
246*> = 0: successful exit
247*> < 0: if INFO = -i, the i-th argument had an illegal value
248*> > 0: if INFO = i, then i eigenvectors failed to converge.
249*> Their indices are stored in array IFAIL.
250*> \endverbatim
251*
252* Authors:
253* ========
254*
255*> \author Univ. of Tennessee
256*> \author Univ. of California Berkeley
257*> \author Univ. of Colorado Denver
258*> \author NAG Ltd.
259*
260*> \ingroup heevx_2stage
261*
262*> \par Further Details:
263* =====================
264*>
265*> \verbatim
266*>
267*> All details about the 2stage techniques are available in:
268*>
269*> Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
270*> Parallel reduction to condensed forms for symmetric eigenvalue problems
271*> using aggregated fine-grained and memory-aware kernels. In Proceedings
272*> of 2011 International Conference for High Performance Computing,
273*> Networking, Storage and Analysis (SC '11), New York, NY, USA,
274*> Article 8 , 11 pages.
275*> http://doi.acm.org/10.1145/2063384.2063394
276*>
277*> A. Haidar, J. Kurzak, P. Luszczek, 2013.
278*> An improved parallel singular value algorithm and its implementation
279*> for multicore hardware, In Proceedings of 2013 International Conference
280*> for High Performance Computing, Networking, Storage and Analysis (SC '13).
281*> Denver, Colorado, USA, 2013.
282*> Article 90, 12 pages.
283*> http://doi.acm.org/10.1145/2503210.2503292
284*>
285*> A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
286*> A novel hybrid CPU-GPU generalized eigensolver for electronic structure
287*> calculations based on fine-grained memory aware tasks.
288*> International Journal of High Performance Computing Applications.
289*> Volume 28 Issue 2, Pages 196-209, May 2014.
290*> http://hpc.sagepub.com/content/28/2/196
291*>
292*> \endverbatim
293*
294* =====================================================================
295 SUBROUTINE ssyevx_2stage( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU,
296 $ IL, IU, ABSTOL, M, W, Z, LDZ, WORK,
297 $ LWORK, IWORK, IFAIL, INFO )
298*
299 IMPLICIT NONE
300*
301* -- LAPACK driver routine --
302* -- LAPACK is a software package provided by Univ. of Tennessee, --
303* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
304*
305* .. Scalar Arguments ..
306 CHARACTER JOBZ, RANGE, UPLO
307 INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
308 REAL ABSTOL, VL, VU
309* ..
310* .. Array Arguments ..
311 INTEGER IFAIL( * ), IWORK( * )
312 REAL A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
313* ..
314*
315* =====================================================================
316*
317* .. Parameters ..
318 REAL ZERO, ONE
319 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
320* ..
321* .. Local Scalars ..
322 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
323 $ WANTZ
324 CHARACTER ORDER
325 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
326 $ indisp, indiwo, indtau, indwkn, indwrk, iscale,
327 $ itmp1, j, jj, llwork, llwrkn,
328 $ nsplit, lwmin, lhtrd, lwtrd, kd, ib, indhous
329 REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
330 $ SIGMA, SMLNUM, TMP1, VLL, VUU
331* ..
332* .. External Functions ..
333 LOGICAL LSAME
334 INTEGER ILAENV2STAGE
335 REAL SLAMCH, SLANSY, SROUNDUP_LWORK
336 EXTERNAL lsame, slamch, slansy, ilaenv2stage,
337 $ sroundup_lwork
338* ..
339* .. External Subroutines ..
340 EXTERNAL scopy, slacpy, sorgtr, sormtr, sscal,
341 $ sstebz,
344* ..
345* .. Intrinsic Functions ..
346 INTRINSIC max, min, sqrt
347* ..
348* .. Executable Statements ..
349*
350* Test the input parameters.
351*
352 lower = lsame( uplo, 'L' )
353 wantz = lsame( jobz, 'V' )
354 alleig = lsame( range, 'A' )
355 valeig = lsame( range, 'V' )
356 indeig = lsame( range, 'I' )
357 lquery = ( lwork.EQ.-1 )
358*
359 info = 0
360 IF( .NOT.( lsame( jobz, 'N' ) ) ) THEN
361 info = -1
362 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
363 info = -2
364 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
365 info = -3
366 ELSE IF( n.LT.0 ) THEN
367 info = -4
368 ELSE IF( lda.LT.max( 1, n ) ) THEN
369 info = -6
370 ELSE
371 IF( valeig ) THEN
372 IF( n.GT.0 .AND. vu.LE.vl )
373 $ info = -8
374 ELSE IF( indeig ) THEN
375 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
376 info = -9
377 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
378 info = -10
379 END IF
380 END IF
381 END IF
382 IF( info.EQ.0 ) THEN
383 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
384 info = -15
385 END IF
386 END IF
387*
388 IF( info.EQ.0 ) THEN
389 IF( n.LE.1 ) THEN
390 lwmin = 1
391 work( 1 ) = sroundup_lwork(lwmin)
392 ELSE
393 kd = ilaenv2stage( 1, 'SSYTRD_2STAGE', jobz,
394 $ n, -1, -1, -1 )
395 ib = ilaenv2stage( 2, 'SSYTRD_2STAGE', jobz,
396 $ n, kd, -1, -1 )
397 lhtrd = ilaenv2stage( 3, 'SSYTRD_2STAGE', jobz,
398 $ n, kd, ib, -1 )
399 lwtrd = ilaenv2stage( 4, 'SSYTRD_2STAGE', jobz,
400 $ n, kd, ib, -1 )
401 lwmin = max( 8*n, 3*n + lhtrd + lwtrd )
402 work( 1 ) = real( lwmin )
403 END IF
404*
405 IF( lwork.LT.lwmin .AND. .NOT.lquery )
406 $ info = -17
407 END IF
408*
409 IF( info.NE.0 ) THEN
410 CALL xerbla( 'SSYEVX_2STAGE', -info )
411 RETURN
412 ELSE IF( lquery ) THEN
413 RETURN
414 END IF
415*
416* Quick return if possible
417*
418 m = 0
419 IF( n.EQ.0 ) THEN
420 RETURN
421 END IF
422*
423 IF( n.EQ.1 ) THEN
424 IF( alleig .OR. indeig ) THEN
425 m = 1
426 w( 1 ) = a( 1, 1 )
427 ELSE
428 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
429 m = 1
430 w( 1 ) = a( 1, 1 )
431 END IF
432 END IF
433 IF( wantz )
434 $ z( 1, 1 ) = one
435 RETURN
436 END IF
437*
438* Get machine constants.
439*
440 safmin = slamch( 'Safe minimum' )
441 eps = slamch( 'Precision' )
442 smlnum = safmin / eps
443 bignum = one / smlnum
444 rmin = sqrt( smlnum )
445 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
446*
447* Scale matrix to allowable range, if necessary.
448*
449 iscale = 0
450 abstll = abstol
451 IF( valeig ) THEN
452 vll = vl
453 vuu = vu
454 END IF
455 anrm = slansy( 'M', uplo, n, a, lda, work )
456 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
457 iscale = 1
458 sigma = rmin / anrm
459 ELSE IF( anrm.GT.rmax ) THEN
460 iscale = 1
461 sigma = rmax / anrm
462 END IF
463 IF( iscale.EQ.1 ) THEN
464 IF( lower ) THEN
465 DO 10 j = 1, n
466 CALL sscal( n-j+1, sigma, a( j, j ), 1 )
467 10 CONTINUE
468 ELSE
469 DO 20 j = 1, n
470 CALL sscal( j, sigma, a( 1, j ), 1 )
471 20 CONTINUE
472 END IF
473 IF( abstol.GT.0 )
474 $ abstll = abstol*sigma
475 IF( valeig ) THEN
476 vll = vl*sigma
477 vuu = vu*sigma
478 END IF
479 END IF
480*
481* Call SSYTRD_2STAGE to reduce symmetric matrix to tridiagonal form.
482*
483 indtau = 1
484 inde = indtau + n
485 indd = inde + n
486 indhous = indd + n
487 indwrk = indhous + lhtrd
488 llwork = lwork - indwrk + 1
489*
490 CALL ssytrd_2stage( jobz, uplo, n, a, lda, work( indd ),
491 $ work( inde ), work( indtau ), work( indhous ),
492 $ lhtrd, work( indwrk ), llwork, iinfo )
493*
494* If all eigenvalues are desired and ABSTOL is less than or equal to
495* zero, then call SSTERF or SORGTR and SSTEQR. If this fails for
496* some eigenvalue, then try SSTEBZ.
497*
498 test = .false.
499 IF( indeig ) THEN
500 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
501 test = .true.
502 END IF
503 END IF
504 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
505 CALL scopy( n, work( indd ), 1, w, 1 )
506 indee = indwrk + 2*n
507 IF( .NOT.wantz ) THEN
508 CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
509 CALL ssterf( n, w, work( indee ), info )
510 ELSE
511 CALL slacpy( 'A', n, n, a, lda, z, ldz )
512 CALL sorgtr( uplo, n, z, ldz, work( indtau ),
513 $ work( indwrk ), llwork, iinfo )
514 CALL scopy( n-1, work( inde ), 1, work( indee ), 1 )
515 CALL ssteqr( jobz, n, w, work( indee ), z, ldz,
516 $ work( indwrk ), info )
517 IF( info.EQ.0 ) THEN
518 DO 30 i = 1, n
519 ifail( i ) = 0
520 30 CONTINUE
521 END IF
522 END IF
523 IF( info.EQ.0 ) THEN
524 m = n
525 GO TO 40
526 END IF
527 info = 0
528 END IF
529*
530* Otherwise, call SSTEBZ and, if eigenvectors are desired, SSTEIN.
531*
532 IF( wantz ) THEN
533 order = 'B'
534 ELSE
535 order = 'E'
536 END IF
537 indibl = 1
538 indisp = indibl + n
539 indiwo = indisp + n
540 CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
541 $ work( indd ), work( inde ), m, nsplit, w,
542 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
543 $ iwork( indiwo ), info )
544*
545 IF( wantz ) THEN
546 CALL sstein( n, work( indd ), work( inde ), m, w,
547 $ iwork( indibl ), iwork( indisp ), z, ldz,
548 $ work( indwrk ), iwork( indiwo ), ifail, info )
549*
550* Apply orthogonal matrix used in reduction to tridiagonal
551* form to eigenvectors returned by SSTEIN.
552*
553 indwkn = inde
554 llwrkn = lwork - indwkn + 1
555 CALL sormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
556 $ z,
557 $ ldz, work( indwkn ), llwrkn, iinfo )
558 END IF
559*
560* If matrix was scaled, then rescale eigenvalues appropriately.
561*
562 40 CONTINUE
563 IF( iscale.EQ.1 ) THEN
564 IF( info.EQ.0 ) THEN
565 imax = m
566 ELSE
567 imax = info - 1
568 END IF
569 CALL sscal( imax, one / sigma, w, 1 )
570 END IF
571*
572* If eigenvalues are not in order, then sort them, along with
573* eigenvectors.
574*
575 IF( wantz ) THEN
576 DO 60 j = 1, m - 1
577 i = 0
578 tmp1 = w( j )
579 DO 50 jj = j + 1, m
580 IF( w( jj ).LT.tmp1 ) THEN
581 i = jj
582 tmp1 = w( jj )
583 END IF
584 50 CONTINUE
585*
586 IF( i.NE.0 ) THEN
587 itmp1 = iwork( indibl+i-1 )
588 w( i ) = w( j )
589 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
590 w( j ) = tmp1
591 iwork( indibl+j-1 ) = itmp1
592 CALL sswap( n, z( 1, i ), 1, z( 1, j ), 1 )
593 IF( info.NE.0 ) THEN
594 itmp1 = ifail( i )
595 ifail( i ) = ifail( j )
596 ifail( j ) = itmp1
597 END IF
598 END IF
599 60 CONTINUE
600 END IF
601*
602* Set WORK(1) to optimal workspace size.
603*
604 work( 1 ) = sroundup_lwork(lwmin)
605*
606 RETURN
607*
608* End of SSYEVX_2STAGE
609*
610 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine ssyevx_2stage(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, iwork, ifail, info)
SSYEVX_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat...
subroutine ssytrd_2stage(vect, uplo, n, a, lda, d, e, tau, hous2, lhous2, work, lwork, info)
SSYTRD_2STAGE
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
SSTEBZ
Definition sstebz.f:272
subroutine sstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
SSTEIN
Definition sstein.f:172
subroutine ssteqr(compz, n, d, e, z, ldz, work, info)
SSTEQR
Definition ssteqr.f:129
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine sorgtr(uplo, n, a, lda, tau, work, lwork, info)
SORGTR
Definition sorgtr.f:121
subroutine sormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
SORMTR
Definition sormtr.f:171