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