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