LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ztfttr.f
Go to the documentation of this file.
1*> \brief \b ZTFTTR copies a triangular matrix from the rectangular full packed format (TF) to the standard full format (TR).
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZTFTTR + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ztfttr.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ztfttr.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ztfttr.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZTFTTR( TRANSR, UPLO, N, ARF, A, LDA, INFO )
20*
21* .. Scalar Arguments ..
22* CHARACTER TRANSR, UPLO
23* INTEGER INFO, N, LDA
24* ..
25* .. Array Arguments ..
26* COMPLEX*16 A( 0: LDA-1, 0: * ), ARF( 0: * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> ZTFTTR copies a triangular matrix A from rectangular full packed
36*> format (TF) to standard full format (TR).
37*> \endverbatim
38*
39* Arguments:
40* ==========
41*
42*> \param[in] TRANSR
43*> \verbatim
44*> TRANSR is CHARACTER*1
45*> = 'N': ARF is in Normal format;
46*> = 'C': ARF is in Conjugate-transpose format;
47*> \endverbatim
48*>
49*> \param[in] UPLO
50*> \verbatim
51*> UPLO is CHARACTER*1
52*> = 'U': A is upper triangular;
53*> = 'L': A is lower triangular.
54*> \endverbatim
55*>
56*> \param[in] N
57*> \verbatim
58*> N is INTEGER
59*> The order of the matrix A. N >= 0.
60*> \endverbatim
61*>
62*> \param[in] ARF
63*> \verbatim
64*> ARF is COMPLEX*16 array, dimension ( N*(N+1)/2 ),
65*> On entry, the upper or lower triangular matrix A stored in
66*> RFP format. For a further discussion see Notes below.
67*> \endverbatim
68*>
69*> \param[out] A
70*> \verbatim
71*> A is COMPLEX*16 array, dimension ( LDA, N )
72*> On exit, the triangular matrix A. If UPLO = 'U', the
73*> leading N-by-N upper triangular part of the array A contains
74*> the upper triangular matrix, and the strictly lower
75*> triangular part of A is not referenced. If UPLO = 'L', the
76*> leading N-by-N lower triangular part of the array A contains
77*> the lower triangular matrix, and the strictly upper
78*> triangular part of A is not referenced.
79*> \endverbatim
80*>
81*> \param[in] LDA
82*> \verbatim
83*> LDA is INTEGER
84*> The leading dimension of the array A. LDA >= max(1,N).
85*> \endverbatim
86*>
87*> \param[out] INFO
88*> \verbatim
89*> INFO is INTEGER
90*> = 0: successful exit
91*> < 0: if INFO = -i, the i-th argument had an illegal value
92*> \endverbatim
93*
94* Authors:
95* ========
96*
97*> \author Univ. of Tennessee
98*> \author Univ. of California Berkeley
99*> \author Univ. of Colorado Denver
100*> \author NAG Ltd.
101*
102*> \ingroup tfttr
103*
104*> \par Further Details:
105* =====================
106*>
107*> \verbatim
108*>
109*> We first consider Standard Packed Format when N is even.
110*> We give an example where N = 6.
111*>
112*> AP is Upper AP is Lower
113*>
114*> 00 01 02 03 04 05 00
115*> 11 12 13 14 15 10 11
116*> 22 23 24 25 20 21 22
117*> 33 34 35 30 31 32 33
118*> 44 45 40 41 42 43 44
119*> 55 50 51 52 53 54 55
120*>
121*>
122*> Let TRANSR = 'N'. RFP holds AP as follows:
123*> For UPLO = 'U' the upper trapezoid A(0:5,0:2) consists of the last
124*> three columns of AP upper. The lower triangle A(4:6,0:2) consists of
125*> conjugate-transpose of the first three columns of AP upper.
126*> For UPLO = 'L' the lower trapezoid A(1:6,0:2) consists of the first
127*> three columns of AP lower. The upper triangle A(0:2,0:2) consists of
128*> conjugate-transpose of the last three columns of AP lower.
129*> To denote conjugate we place -- above the element. This covers the
130*> case N even and TRANSR = 'N'.
131*>
132*> RFP A RFP A
133*>
134*> -- -- --
135*> 03 04 05 33 43 53
136*> -- --
137*> 13 14 15 00 44 54
138*> --
139*> 23 24 25 10 11 55
140*>
141*> 33 34 35 20 21 22
142*> --
143*> 00 44 45 30 31 32
144*> -- --
145*> 01 11 55 40 41 42
146*> -- -- --
147*> 02 12 22 50 51 52
148*>
149*> Now let TRANSR = 'C'. RFP A in both UPLO cases is just the conjugate-
150*> transpose of RFP A above. One therefore gets:
151*>
152*>
153*> RFP A RFP A
154*>
155*> -- -- -- -- -- -- -- -- -- --
156*> 03 13 23 33 00 01 02 33 00 10 20 30 40 50
157*> -- -- -- -- -- -- -- -- -- --
158*> 04 14 24 34 44 11 12 43 44 11 21 31 41 51
159*> -- -- -- -- -- -- -- -- -- --
160*> 05 15 25 35 45 55 22 53 54 55 22 32 42 52
161*>
162*>
163*> We next consider Standard Packed Format when N is odd.
164*> We give an example where N = 5.
165*>
166*> AP is Upper AP is Lower
167*>
168*> 00 01 02 03 04 00
169*> 11 12 13 14 10 11
170*> 22 23 24 20 21 22
171*> 33 34 30 31 32 33
172*> 44 40 41 42 43 44
173*>
174*>
175*> Let TRANSR = 'N'. RFP holds AP as follows:
176*> For UPLO = 'U' the upper trapezoid A(0:4,0:2) consists of the last
177*> three columns of AP upper. The lower triangle A(3:4,0:1) consists of
178*> conjugate-transpose of the first two columns of AP upper.
179*> For UPLO = 'L' the lower trapezoid A(0:4,0:2) consists of the first
180*> three columns of AP lower. The upper triangle A(0:1,1:2) consists of
181*> conjugate-transpose of the last two columns of AP lower.
182*> To denote conjugate we place -- above the element. This covers the
183*> case N odd and TRANSR = 'N'.
184*>
185*> RFP A RFP A
186*>
187*> -- --
188*> 02 03 04 00 33 43
189*> --
190*> 12 13 14 10 11 44
191*>
192*> 22 23 24 20 21 22
193*> --
194*> 00 33 34 30 31 32
195*> -- --
196*> 01 11 44 40 41 42
197*>
198*> Now let TRANSR = 'C'. RFP A in both UPLO cases is just the conjugate-
199*> transpose of RFP A above. One therefore gets:
200*>
201*>
202*> RFP A RFP A
203*>
204*> -- -- -- -- -- -- -- -- --
205*> 02 12 22 00 01 00 10 20 30 40 50
206*> -- -- -- -- -- -- -- -- --
207*> 03 13 23 33 11 33 11 21 31 41 51
208*> -- -- -- -- -- -- -- -- --
209*> 04 14 24 34 44 43 44 22 32 42 52
210*> \endverbatim
211*>
212* =====================================================================
213 SUBROUTINE ztfttr( TRANSR, UPLO, N, ARF, A, LDA, INFO )
214*
215* -- LAPACK computational routine --
216* -- LAPACK is a software package provided by Univ. of Tennessee, --
217* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
218*
219* .. Scalar Arguments ..
220 CHARACTER TRANSR, UPLO
221 INTEGER INFO, N, LDA
222* ..
223* .. Array Arguments ..
224 COMPLEX*16 A( 0: LDA-1, 0: * ), ARF( 0: * )
225* ..
226*
227* =====================================================================
228*
229* .. Parameters ..
230* ..
231* .. Local Scalars ..
232 LOGICAL LOWER, NISODD, NORMALTRANSR
233 INTEGER N1, N2, K, NT, NX2, NP1X2
234 INTEGER I, J, L, IJ
235* ..
236* .. External Functions ..
237 LOGICAL LSAME
238 EXTERNAL lsame
239* ..
240* .. External Subroutines ..
241 EXTERNAL xerbla
242* ..
243* .. Intrinsic Functions ..
244 INTRINSIC dconjg, max, mod
245* ..
246* .. Executable Statements ..
247*
248* Test the input parameters.
249*
250 info = 0
251 normaltransr = lsame( transr, 'N' )
252 lower = lsame( uplo, 'L' )
253 IF( .NOT.normaltransr .AND. .NOT.lsame( transr, 'C' ) ) THEN
254 info = -1
255 ELSE IF( .NOT.lower .AND. .NOT.lsame( uplo, 'U' ) ) THEN
256 info = -2
257 ELSE IF( n.LT.0 ) THEN
258 info = -3
259 ELSE IF( lda.LT.max( 1, n ) ) THEN
260 info = -6
261 END IF
262 IF( info.NE.0 ) THEN
263 CALL xerbla( 'ZTFTTR', -info )
264 RETURN
265 END IF
266*
267* Quick return if possible
268*
269 IF( n.LE.1 ) THEN
270 IF( n.EQ.1 ) THEN
271 IF( normaltransr ) THEN
272 a( 0, 0 ) = arf( 0 )
273 ELSE
274 a( 0, 0 ) = dconjg( arf( 0 ) )
275 END IF
276 END IF
277 RETURN
278 END IF
279*
280* Size of array ARF(1:2,0:nt-1)
281*
282 nt = n*( n+1 ) / 2
283*
284* set N1 and N2 depending on LOWER: for N even N1=N2=K
285*
286 IF( lower ) THEN
287 n2 = n / 2
288 n1 = n - n2
289 ELSE
290 n1 = n / 2
291 n2 = n - n1
292 END IF
293*
294* If N is odd, set NISODD = .TRUE., LDA=N+1 and A is (N+1)--by--K2.
295* If N is even, set K = N/2 and NISODD = .FALSE., LDA=N and A is
296* N--by--(N+1)/2.
297*
298 IF( mod( n, 2 ).EQ.0 ) THEN
299 k = n / 2
300 nisodd = .false.
301 IF( .NOT.lower )
302 $ np1x2 = n + n + 2
303 ELSE
304 nisodd = .true.
305 IF( .NOT.lower )
306 $ nx2 = n + n
307 END IF
308*
309 IF( nisodd ) THEN
310*
311* N is odd
312*
313 IF( normaltransr ) THEN
314*
315* N is odd and TRANSR = 'N'
316*
317 IF( lower ) THEN
318*
319* SRPA for LOWER, NORMAL and N is odd ( a(0:n-1,0:n1-1) )
320* T1 -> a(0,0), T2 -> a(0,1), S -> a(n1,0)
321* T1 -> a(0), T2 -> a(n), S -> a(n1); lda=n
322*
323 ij = 0
324 DO j = 0, n2
325 DO i = n1, n2 + j
326 a( n2+j, i ) = dconjg( arf( ij ) )
327 ij = ij + 1
328 END DO
329 DO i = j, n - 1
330 a( i, j ) = arf( ij )
331 ij = ij + 1
332 END DO
333 END DO
334*
335 ELSE
336*
337* SRPA for UPPER, NORMAL and N is odd ( a(0:n-1,0:n2-1)
338* T1 -> a(n1+1,0), T2 -> a(n1,0), S -> a(0,0)
339* T1 -> a(n2), T2 -> a(n1), S -> a(0); lda=n
340*
341 ij = nt - n
342 DO j = n - 1, n1, -1
343 DO i = 0, j
344 a( i, j ) = arf( ij )
345 ij = ij + 1
346 END DO
347 DO l = j - n1, n1 - 1
348 a( j-n1, l ) = dconjg( arf( ij ) )
349 ij = ij + 1
350 END DO
351 ij = ij - nx2
352 END DO
353*
354 END IF
355*
356 ELSE
357*
358* N is odd and TRANSR = 'C'
359*
360 IF( lower ) THEN
361*
362* SRPA for LOWER, TRANSPOSE and N is odd
363* T1 -> A(0,0) , T2 -> A(1,0) , S -> A(0,n1)
364* T1 -> A(0+0) , T2 -> A(1+0) , S -> A(0+n1*n1); lda=n1
365*
366 ij = 0
367 DO j = 0, n2 - 1
368 DO i = 0, j
369 a( j, i ) = dconjg( arf( ij ) )
370 ij = ij + 1
371 END DO
372 DO i = n1 + j, n - 1
373 a( i, n1+j ) = arf( ij )
374 ij = ij + 1
375 END DO
376 END DO
377 DO j = n2, n - 1
378 DO i = 0, n1 - 1
379 a( j, i ) = dconjg( arf( ij ) )
380 ij = ij + 1
381 END DO
382 END DO
383*
384 ELSE
385*
386* SRPA for UPPER, TRANSPOSE and N is odd
387* T1 -> A(0,n1+1), T2 -> A(0,n1), S -> A(0,0)
388* T1 -> A(n2*n2), T2 -> A(n1*n2), S -> A(0); lda = n2
389*
390 ij = 0
391 DO j = 0, n1
392 DO i = n1, n - 1
393 a( j, i ) = dconjg( arf( ij ) )
394 ij = ij + 1
395 END DO
396 END DO
397 DO j = 0, n1 - 1
398 DO i = 0, j
399 a( i, j ) = arf( ij )
400 ij = ij + 1
401 END DO
402 DO l = n2 + j, n - 1
403 a( n2+j, l ) = dconjg( arf( ij ) )
404 ij = ij + 1
405 END DO
406 END DO
407*
408 END IF
409*
410 END IF
411*
412 ELSE
413*
414* N is even
415*
416 IF( normaltransr ) THEN
417*
418* N is even and TRANSR = 'N'
419*
420 IF( lower ) THEN
421*
422* SRPA for LOWER, NORMAL, and N is even ( a(0:n,0:k-1) )
423* T1 -> a(1,0), T2 -> a(0,0), S -> a(k+1,0)
424* T1 -> a(1), T2 -> a(0), S -> a(k+1); lda=n+1
425*
426 ij = 0
427 DO j = 0, k - 1
428 DO i = k, k + j
429 a( k+j, i ) = dconjg( arf( ij ) )
430 ij = ij + 1
431 END DO
432 DO i = j, n - 1
433 a( i, j ) = arf( ij )
434 ij = ij + 1
435 END DO
436 END DO
437*
438 ELSE
439*
440* SRPA for UPPER, NORMAL, and N is even ( a(0:n,0:k-1) )
441* T1 -> a(k+1,0) , T2 -> a(k,0), S -> a(0,0)
442* T1 -> a(k+1), T2 -> a(k), S -> a(0); lda=n+1
443*
444 ij = nt - n - 1
445 DO j = n - 1, k, -1
446 DO i = 0, j
447 a( i, j ) = arf( ij )
448 ij = ij + 1
449 END DO
450 DO l = j - k, k - 1
451 a( j-k, l ) = dconjg( arf( ij ) )
452 ij = ij + 1
453 END DO
454 ij = ij - np1x2
455 END DO
456*
457 END IF
458*
459 ELSE
460*
461* N is even and TRANSR = 'C'
462*
463 IF( lower ) THEN
464*
465* SRPA for LOWER, TRANSPOSE and N is even (see paper, A=B)
466* T1 -> A(0,1) , T2 -> A(0,0) , S -> A(0,k+1) :
467* T1 -> A(0+k) , T2 -> A(0+0) , S -> A(0+k*(k+1)); lda=k
468*
469 ij = 0
470 j = k
471 DO i = k, n - 1
472 a( i, j ) = arf( ij )
473 ij = ij + 1
474 END DO
475 DO j = 0, k - 2
476 DO i = 0, j
477 a( j, i ) = dconjg( arf( ij ) )
478 ij = ij + 1
479 END DO
480 DO i = k + 1 + j, n - 1
481 a( i, k+1+j ) = arf( ij )
482 ij = ij + 1
483 END DO
484 END DO
485 DO j = k - 1, n - 1
486 DO i = 0, k - 1
487 a( j, i ) = dconjg( arf( ij ) )
488 ij = ij + 1
489 END DO
490 END DO
491*
492 ELSE
493*
494* SRPA for UPPER, TRANSPOSE and N is even (see paper, A=B)
495* T1 -> A(0,k+1) , T2 -> A(0,k) , S -> A(0,0)
496* T1 -> A(0+k*(k+1)) , T2 -> A(0+k*k) , S -> A(0+0)); lda=k
497*
498 ij = 0
499 DO j = 0, k
500 DO i = k, n - 1
501 a( j, i ) = dconjg( arf( ij ) )
502 ij = ij + 1
503 END DO
504 END DO
505 DO j = 0, k - 2
506 DO i = 0, j
507 a( i, j ) = arf( ij )
508 ij = ij + 1
509 END DO
510 DO l = k + 1 + j, n - 1
511 a( k+1+j, l ) = dconjg( arf( ij ) )
512 ij = ij + 1
513 END DO
514 END DO
515*
516* Note that here J = K-1
517*
518 DO i = 0, j
519 a( i, j ) = arf( ij )
520 ij = ij + 1
521 END DO
522*
523 END IF
524*
525 END IF
526*
527 END IF
528*
529 RETURN
530*
531* End of ZTFTTR
532*
533 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ztfttr(transr, uplo, n, arf, a, lda, info)
ZTFTTR copies a triangular matrix from the rectangular full packed format (TF) to the standard full f...
Definition ztfttr.f:214