LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zlarfx.f
Go to the documentation of this file.
1*> \brief \b ZLARFX applies an elementary reflector to a general rectangular matrix, with loop unrolling when the reflector has order ≤ 10.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZLARFX + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zlarfx.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zlarfx.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zlarfx.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZLARFX( SIDE, M, N, V, TAU, C, LDC, WORK )
22*
23* .. Scalar Arguments ..
24* CHARACTER SIDE
25* INTEGER LDC, M, N
26* COMPLEX*16 TAU
27* ..
28* .. Array Arguments ..
29* COMPLEX*16 C( LDC, * ), V( * ), WORK( * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> ZLARFX applies a complex elementary reflector H to a complex m by n
39*> matrix C, from either the left or the right. H is represented in the
40*> form
41*>
42*> H = I - tau * v * v**H
43*>
44*> where tau is a complex scalar and v is a complex vector.
45*>
46*> If tau = 0, then H is taken to be the unit matrix
47*>
48*> This version uses inline code if H has order < 11.
49*> \endverbatim
50*
51* Arguments:
52* ==========
53*
54*> \param[in] SIDE
55*> \verbatim
56*> SIDE is CHARACTER*1
57*> = 'L': form H * C
58*> = 'R': form C * H
59*> \endverbatim
60*>
61*> \param[in] M
62*> \verbatim
63*> M is INTEGER
64*> The number of rows of the matrix C.
65*> \endverbatim
66*>
67*> \param[in] N
68*> \verbatim
69*> N is INTEGER
70*> The number of columns of the matrix C.
71*> \endverbatim
72*>
73*> \param[in] V
74*> \verbatim
75*> V is COMPLEX*16 array, dimension (M) if SIDE = 'L'
76*> or (N) if SIDE = 'R'
77*> The vector v in the representation of H.
78*> \endverbatim
79*>
80*> \param[in] TAU
81*> \verbatim
82*> TAU is COMPLEX*16
83*> The value tau in the representation of H.
84*> \endverbatim
85*>
86*> \param[in,out] C
87*> \verbatim
88*> C is COMPLEX*16 array, dimension (LDC,N)
89*> On entry, the m by n matrix C.
90*> On exit, C is overwritten by the matrix H * C if SIDE = 'L',
91*> or C * H if SIDE = 'R'.
92*> \endverbatim
93*>
94*> \param[in] LDC
95*> \verbatim
96*> LDC is INTEGER
97*> The leading dimension of the array C. LDC >= max(1,M).
98*> \endverbatim
99*>
100*> \param[out] WORK
101*> \verbatim
102*> WORK is COMPLEX*16 array, dimension (N) if SIDE = 'L'
103*> or (M) if SIDE = 'R'
104*> WORK is not referenced if H has order < 11.
105*> \endverbatim
106*
107* Authors:
108* ========
109*
110*> \author Univ. of Tennessee
111*> \author Univ. of California Berkeley
112*> \author Univ. of Colorado Denver
113*> \author NAG Ltd.
114*
115*> \ingroup larfx
116*
117* =====================================================================
118 SUBROUTINE zlarfx( SIDE, M, N, V, TAU, C, LDC, WORK )
119*
120* -- LAPACK auxiliary routine --
121* -- LAPACK is a software package provided by Univ. of Tennessee, --
122* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
123*
124* .. Scalar Arguments ..
125 CHARACTER SIDE
126 INTEGER LDC, M, N
127 COMPLEX*16 TAU
128* ..
129* .. Array Arguments ..
130 COMPLEX*16 C( LDC, * ), V( * ), WORK( * )
131* ..
132*
133* =====================================================================
134*
135* .. Parameters ..
136 COMPLEX*16 ZERO, ONE
137 parameter( zero = ( 0.0d+0, 0.0d+0 ),
138 $ one = ( 1.0d+0, 0.0d+0 ) )
139* ..
140* .. Local Scalars ..
141 INTEGER J
142 COMPLEX*16 SUM, T1, T10, T2, T3, T4, T5, T6, T7, T8, T9,
143 $ V1, V10, V2, V3, V4, V5, V6, V7, V8, V9
144* ..
145* .. External Functions ..
146 LOGICAL LSAME
147 EXTERNAL lsame
148* ..
149* .. External Subroutines ..
150 EXTERNAL zlarf
151* ..
152* .. Intrinsic Functions ..
153 INTRINSIC dconjg
154* ..
155* .. Executable Statements ..
156*
157 IF( tau.EQ.zero )
158 $ RETURN
159 IF( lsame( side, 'L' ) ) THEN
160*
161* Form H * C, where H has order m.
162*
163 GO TO ( 10, 30, 50, 70, 90, 110, 130, 150,
164 $ 170, 190 )m
165*
166* Code for general M
167*
168 CALL zlarf( side, m, n, v, 1, tau, c, ldc, work )
169 GO TO 410
170 10 CONTINUE
171*
172* Special code for 1 x 1 Householder
173*
174 t1 = one - tau*v( 1 )*dconjg( v( 1 ) )
175 DO 20 j = 1, n
176 c( 1, j ) = t1*c( 1, j )
177 20 CONTINUE
178 GO TO 410
179 30 CONTINUE
180*
181* Special code for 2 x 2 Householder
182*
183 v1 = dconjg( v( 1 ) )
184 t1 = tau*dconjg( v1 )
185 v2 = dconjg( v( 2 ) )
186 t2 = tau*dconjg( v2 )
187 DO 40 j = 1, n
188 sum = v1*c( 1, j ) + v2*c( 2, j )
189 c( 1, j ) = c( 1, j ) - sum*t1
190 c( 2, j ) = c( 2, j ) - sum*t2
191 40 CONTINUE
192 GO TO 410
193 50 CONTINUE
194*
195* Special code for 3 x 3 Householder
196*
197 v1 = dconjg( v( 1 ) )
198 t1 = tau*dconjg( v1 )
199 v2 = dconjg( v( 2 ) )
200 t2 = tau*dconjg( v2 )
201 v3 = dconjg( v( 3 ) )
202 t3 = tau*dconjg( v3 )
203 DO 60 j = 1, n
204 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j )
205 c( 1, j ) = c( 1, j ) - sum*t1
206 c( 2, j ) = c( 2, j ) - sum*t2
207 c( 3, j ) = c( 3, j ) - sum*t3
208 60 CONTINUE
209 GO TO 410
210 70 CONTINUE
211*
212* Special code for 4 x 4 Householder
213*
214 v1 = dconjg( v( 1 ) )
215 t1 = tau*dconjg( v1 )
216 v2 = dconjg( v( 2 ) )
217 t2 = tau*dconjg( v2 )
218 v3 = dconjg( v( 3 ) )
219 t3 = tau*dconjg( v3 )
220 v4 = dconjg( v( 4 ) )
221 t4 = tau*dconjg( v4 )
222 DO 80 j = 1, n
223 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
224 $ v4*c( 4, j )
225 c( 1, j ) = c( 1, j ) - sum*t1
226 c( 2, j ) = c( 2, j ) - sum*t2
227 c( 3, j ) = c( 3, j ) - sum*t3
228 c( 4, j ) = c( 4, j ) - sum*t4
229 80 CONTINUE
230 GO TO 410
231 90 CONTINUE
232*
233* Special code for 5 x 5 Householder
234*
235 v1 = dconjg( v( 1 ) )
236 t1 = tau*dconjg( v1 )
237 v2 = dconjg( v( 2 ) )
238 t2 = tau*dconjg( v2 )
239 v3 = dconjg( v( 3 ) )
240 t3 = tau*dconjg( v3 )
241 v4 = dconjg( v( 4 ) )
242 t4 = tau*dconjg( v4 )
243 v5 = dconjg( v( 5 ) )
244 t5 = tau*dconjg( v5 )
245 DO 100 j = 1, n
246 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
247 $ v4*c( 4, j ) + v5*c( 5, j )
248 c( 1, j ) = c( 1, j ) - sum*t1
249 c( 2, j ) = c( 2, j ) - sum*t2
250 c( 3, j ) = c( 3, j ) - sum*t3
251 c( 4, j ) = c( 4, j ) - sum*t4
252 c( 5, j ) = c( 5, j ) - sum*t5
253 100 CONTINUE
254 GO TO 410
255 110 CONTINUE
256*
257* Special code for 6 x 6 Householder
258*
259 v1 = dconjg( v( 1 ) )
260 t1 = tau*dconjg( v1 )
261 v2 = dconjg( v( 2 ) )
262 t2 = tau*dconjg( v2 )
263 v3 = dconjg( v( 3 ) )
264 t3 = tau*dconjg( v3 )
265 v4 = dconjg( v( 4 ) )
266 t4 = tau*dconjg( v4 )
267 v5 = dconjg( v( 5 ) )
268 t5 = tau*dconjg( v5 )
269 v6 = dconjg( v( 6 ) )
270 t6 = tau*dconjg( v6 )
271 DO 120 j = 1, n
272 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
273 $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j )
274 c( 1, j ) = c( 1, j ) - sum*t1
275 c( 2, j ) = c( 2, j ) - sum*t2
276 c( 3, j ) = c( 3, j ) - sum*t3
277 c( 4, j ) = c( 4, j ) - sum*t4
278 c( 5, j ) = c( 5, j ) - sum*t5
279 c( 6, j ) = c( 6, j ) - sum*t6
280 120 CONTINUE
281 GO TO 410
282 130 CONTINUE
283*
284* Special code for 7 x 7 Householder
285*
286 v1 = dconjg( v( 1 ) )
287 t1 = tau*dconjg( v1 )
288 v2 = dconjg( v( 2 ) )
289 t2 = tau*dconjg( v2 )
290 v3 = dconjg( v( 3 ) )
291 t3 = tau*dconjg( v3 )
292 v4 = dconjg( v( 4 ) )
293 t4 = tau*dconjg( v4 )
294 v5 = dconjg( v( 5 ) )
295 t5 = tau*dconjg( v5 )
296 v6 = dconjg( v( 6 ) )
297 t6 = tau*dconjg( v6 )
298 v7 = dconjg( v( 7 ) )
299 t7 = tau*dconjg( v7 )
300 DO 140 j = 1, n
301 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
302 $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
303 $ v7*c( 7, j )
304 c( 1, j ) = c( 1, j ) - sum*t1
305 c( 2, j ) = c( 2, j ) - sum*t2
306 c( 3, j ) = c( 3, j ) - sum*t3
307 c( 4, j ) = c( 4, j ) - sum*t4
308 c( 5, j ) = c( 5, j ) - sum*t5
309 c( 6, j ) = c( 6, j ) - sum*t6
310 c( 7, j ) = c( 7, j ) - sum*t7
311 140 CONTINUE
312 GO TO 410
313 150 CONTINUE
314*
315* Special code for 8 x 8 Householder
316*
317 v1 = dconjg( v( 1 ) )
318 t1 = tau*dconjg( v1 )
319 v2 = dconjg( v( 2 ) )
320 t2 = tau*dconjg( v2 )
321 v3 = dconjg( v( 3 ) )
322 t3 = tau*dconjg( v3 )
323 v4 = dconjg( v( 4 ) )
324 t4 = tau*dconjg( v4 )
325 v5 = dconjg( v( 5 ) )
326 t5 = tau*dconjg( v5 )
327 v6 = dconjg( v( 6 ) )
328 t6 = tau*dconjg( v6 )
329 v7 = dconjg( v( 7 ) )
330 t7 = tau*dconjg( v7 )
331 v8 = dconjg( v( 8 ) )
332 t8 = tau*dconjg( v8 )
333 DO 160 j = 1, n
334 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
335 $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
336 $ v7*c( 7, j ) + v8*c( 8, j )
337 c( 1, j ) = c( 1, j ) - sum*t1
338 c( 2, j ) = c( 2, j ) - sum*t2
339 c( 3, j ) = c( 3, j ) - sum*t3
340 c( 4, j ) = c( 4, j ) - sum*t4
341 c( 5, j ) = c( 5, j ) - sum*t5
342 c( 6, j ) = c( 6, j ) - sum*t6
343 c( 7, j ) = c( 7, j ) - sum*t7
344 c( 8, j ) = c( 8, j ) - sum*t8
345 160 CONTINUE
346 GO TO 410
347 170 CONTINUE
348*
349* Special code for 9 x 9 Householder
350*
351 v1 = dconjg( v( 1 ) )
352 t1 = tau*dconjg( v1 )
353 v2 = dconjg( v( 2 ) )
354 t2 = tau*dconjg( v2 )
355 v3 = dconjg( v( 3 ) )
356 t3 = tau*dconjg( v3 )
357 v4 = dconjg( v( 4 ) )
358 t4 = tau*dconjg( v4 )
359 v5 = dconjg( v( 5 ) )
360 t5 = tau*dconjg( v5 )
361 v6 = dconjg( v( 6 ) )
362 t6 = tau*dconjg( v6 )
363 v7 = dconjg( v( 7 ) )
364 t7 = tau*dconjg( v7 )
365 v8 = dconjg( v( 8 ) )
366 t8 = tau*dconjg( v8 )
367 v9 = dconjg( v( 9 ) )
368 t9 = tau*dconjg( v9 )
369 DO 180 j = 1, n
370 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
371 $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
372 $ v7*c( 7, j ) + v8*c( 8, j ) + v9*c( 9, j )
373 c( 1, j ) = c( 1, j ) - sum*t1
374 c( 2, j ) = c( 2, j ) - sum*t2
375 c( 3, j ) = c( 3, j ) - sum*t3
376 c( 4, j ) = c( 4, j ) - sum*t4
377 c( 5, j ) = c( 5, j ) - sum*t5
378 c( 6, j ) = c( 6, j ) - sum*t6
379 c( 7, j ) = c( 7, j ) - sum*t7
380 c( 8, j ) = c( 8, j ) - sum*t8
381 c( 9, j ) = c( 9, j ) - sum*t9
382 180 CONTINUE
383 GO TO 410
384 190 CONTINUE
385*
386* Special code for 10 x 10 Householder
387*
388 v1 = dconjg( v( 1 ) )
389 t1 = tau*dconjg( v1 )
390 v2 = dconjg( v( 2 ) )
391 t2 = tau*dconjg( v2 )
392 v3 = dconjg( v( 3 ) )
393 t3 = tau*dconjg( v3 )
394 v4 = dconjg( v( 4 ) )
395 t4 = tau*dconjg( v4 )
396 v5 = dconjg( v( 5 ) )
397 t5 = tau*dconjg( v5 )
398 v6 = dconjg( v( 6 ) )
399 t6 = tau*dconjg( v6 )
400 v7 = dconjg( v( 7 ) )
401 t7 = tau*dconjg( v7 )
402 v8 = dconjg( v( 8 ) )
403 t8 = tau*dconjg( v8 )
404 v9 = dconjg( v( 9 ) )
405 t9 = tau*dconjg( v9 )
406 v10 = dconjg( v( 10 ) )
407 t10 = tau*dconjg( v10 )
408 DO 200 j = 1, n
409 sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
410 $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
411 $ v7*c( 7, j ) + v8*c( 8, j ) + v9*c( 9, j ) +
412 $ v10*c( 10, j )
413 c( 1, j ) = c( 1, j ) - sum*t1
414 c( 2, j ) = c( 2, j ) - sum*t2
415 c( 3, j ) = c( 3, j ) - sum*t3
416 c( 4, j ) = c( 4, j ) - sum*t4
417 c( 5, j ) = c( 5, j ) - sum*t5
418 c( 6, j ) = c( 6, j ) - sum*t6
419 c( 7, j ) = c( 7, j ) - sum*t7
420 c( 8, j ) = c( 8, j ) - sum*t8
421 c( 9, j ) = c( 9, j ) - sum*t9
422 c( 10, j ) = c( 10, j ) - sum*t10
423 200 CONTINUE
424 GO TO 410
425 ELSE
426*
427* Form C * H, where H has order n.
428*
429 GO TO ( 210, 230, 250, 270, 290, 310, 330, 350,
430 $ 370, 390 )n
431*
432* Code for general N
433*
434 CALL zlarf( side, m, n, v, 1, tau, c, ldc, work )
435 GO TO 410
436 210 CONTINUE
437*
438* Special code for 1 x 1 Householder
439*
440 t1 = one - tau*v( 1 )*dconjg( v( 1 ) )
441 DO 220 j = 1, m
442 c( j, 1 ) = t1*c( j, 1 )
443 220 CONTINUE
444 GO TO 410
445 230 CONTINUE
446*
447* Special code for 2 x 2 Householder
448*
449 v1 = v( 1 )
450 t1 = tau*dconjg( v1 )
451 v2 = v( 2 )
452 t2 = tau*dconjg( v2 )
453 DO 240 j = 1, m
454 sum = v1*c( j, 1 ) + v2*c( j, 2 )
455 c( j, 1 ) = c( j, 1 ) - sum*t1
456 c( j, 2 ) = c( j, 2 ) - sum*t2
457 240 CONTINUE
458 GO TO 410
459 250 CONTINUE
460*
461* Special code for 3 x 3 Householder
462*
463 v1 = v( 1 )
464 t1 = tau*dconjg( v1 )
465 v2 = v( 2 )
466 t2 = tau*dconjg( v2 )
467 v3 = v( 3 )
468 t3 = tau*dconjg( v3 )
469 DO 260 j = 1, m
470 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 )
471 c( j, 1 ) = c( j, 1 ) - sum*t1
472 c( j, 2 ) = c( j, 2 ) - sum*t2
473 c( j, 3 ) = c( j, 3 ) - sum*t3
474 260 CONTINUE
475 GO TO 410
476 270 CONTINUE
477*
478* Special code for 4 x 4 Householder
479*
480 v1 = v( 1 )
481 t1 = tau*dconjg( v1 )
482 v2 = v( 2 )
483 t2 = tau*dconjg( v2 )
484 v3 = v( 3 )
485 t3 = tau*dconjg( v3 )
486 v4 = v( 4 )
487 t4 = tau*dconjg( v4 )
488 DO 280 j = 1, m
489 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
490 $ v4*c( j, 4 )
491 c( j, 1 ) = c( j, 1 ) - sum*t1
492 c( j, 2 ) = c( j, 2 ) - sum*t2
493 c( j, 3 ) = c( j, 3 ) - sum*t3
494 c( j, 4 ) = c( j, 4 ) - sum*t4
495 280 CONTINUE
496 GO TO 410
497 290 CONTINUE
498*
499* Special code for 5 x 5 Householder
500*
501 v1 = v( 1 )
502 t1 = tau*dconjg( v1 )
503 v2 = v( 2 )
504 t2 = tau*dconjg( v2 )
505 v3 = v( 3 )
506 t3 = tau*dconjg( v3 )
507 v4 = v( 4 )
508 t4 = tau*dconjg( v4 )
509 v5 = v( 5 )
510 t5 = tau*dconjg( v5 )
511 DO 300 j = 1, m
512 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
513 $ v4*c( j, 4 ) + v5*c( j, 5 )
514 c( j, 1 ) = c( j, 1 ) - sum*t1
515 c( j, 2 ) = c( j, 2 ) - sum*t2
516 c( j, 3 ) = c( j, 3 ) - sum*t3
517 c( j, 4 ) = c( j, 4 ) - sum*t4
518 c( j, 5 ) = c( j, 5 ) - sum*t5
519 300 CONTINUE
520 GO TO 410
521 310 CONTINUE
522*
523* Special code for 6 x 6 Householder
524*
525 v1 = v( 1 )
526 t1 = tau*dconjg( v1 )
527 v2 = v( 2 )
528 t2 = tau*dconjg( v2 )
529 v3 = v( 3 )
530 t3 = tau*dconjg( v3 )
531 v4 = v( 4 )
532 t4 = tau*dconjg( v4 )
533 v5 = v( 5 )
534 t5 = tau*dconjg( v5 )
535 v6 = v( 6 )
536 t6 = tau*dconjg( v6 )
537 DO 320 j = 1, m
538 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
539 $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 )
540 c( j, 1 ) = c( j, 1 ) - sum*t1
541 c( j, 2 ) = c( j, 2 ) - sum*t2
542 c( j, 3 ) = c( j, 3 ) - sum*t3
543 c( j, 4 ) = c( j, 4 ) - sum*t4
544 c( j, 5 ) = c( j, 5 ) - sum*t5
545 c( j, 6 ) = c( j, 6 ) - sum*t6
546 320 CONTINUE
547 GO TO 410
548 330 CONTINUE
549*
550* Special code for 7 x 7 Householder
551*
552 v1 = v( 1 )
553 t1 = tau*dconjg( v1 )
554 v2 = v( 2 )
555 t2 = tau*dconjg( v2 )
556 v3 = v( 3 )
557 t3 = tau*dconjg( v3 )
558 v4 = v( 4 )
559 t4 = tau*dconjg( v4 )
560 v5 = v( 5 )
561 t5 = tau*dconjg( v5 )
562 v6 = v( 6 )
563 t6 = tau*dconjg( v6 )
564 v7 = v( 7 )
565 t7 = tau*dconjg( v7 )
566 DO 340 j = 1, m
567 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
568 $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
569 $ v7*c( j, 7 )
570 c( j, 1 ) = c( j, 1 ) - sum*t1
571 c( j, 2 ) = c( j, 2 ) - sum*t2
572 c( j, 3 ) = c( j, 3 ) - sum*t3
573 c( j, 4 ) = c( j, 4 ) - sum*t4
574 c( j, 5 ) = c( j, 5 ) - sum*t5
575 c( j, 6 ) = c( j, 6 ) - sum*t6
576 c( j, 7 ) = c( j, 7 ) - sum*t7
577 340 CONTINUE
578 GO TO 410
579 350 CONTINUE
580*
581* Special code for 8 x 8 Householder
582*
583 v1 = v( 1 )
584 t1 = tau*dconjg( v1 )
585 v2 = v( 2 )
586 t2 = tau*dconjg( v2 )
587 v3 = v( 3 )
588 t3 = tau*dconjg( v3 )
589 v4 = v( 4 )
590 t4 = tau*dconjg( v4 )
591 v5 = v( 5 )
592 t5 = tau*dconjg( v5 )
593 v6 = v( 6 )
594 t6 = tau*dconjg( v6 )
595 v7 = v( 7 )
596 t7 = tau*dconjg( v7 )
597 v8 = v( 8 )
598 t8 = tau*dconjg( v8 )
599 DO 360 j = 1, m
600 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
601 $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
602 $ v7*c( j, 7 ) + v8*c( j, 8 )
603 c( j, 1 ) = c( j, 1 ) - sum*t1
604 c( j, 2 ) = c( j, 2 ) - sum*t2
605 c( j, 3 ) = c( j, 3 ) - sum*t3
606 c( j, 4 ) = c( j, 4 ) - sum*t4
607 c( j, 5 ) = c( j, 5 ) - sum*t5
608 c( j, 6 ) = c( j, 6 ) - sum*t6
609 c( j, 7 ) = c( j, 7 ) - sum*t7
610 c( j, 8 ) = c( j, 8 ) - sum*t8
611 360 CONTINUE
612 GO TO 410
613 370 CONTINUE
614*
615* Special code for 9 x 9 Householder
616*
617 v1 = v( 1 )
618 t1 = tau*dconjg( v1 )
619 v2 = v( 2 )
620 t2 = tau*dconjg( v2 )
621 v3 = v( 3 )
622 t3 = tau*dconjg( v3 )
623 v4 = v( 4 )
624 t4 = tau*dconjg( v4 )
625 v5 = v( 5 )
626 t5 = tau*dconjg( v5 )
627 v6 = v( 6 )
628 t6 = tau*dconjg( v6 )
629 v7 = v( 7 )
630 t7 = tau*dconjg( v7 )
631 v8 = v( 8 )
632 t8 = tau*dconjg( v8 )
633 v9 = v( 9 )
634 t9 = tau*dconjg( v9 )
635 DO 380 j = 1, m
636 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
637 $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
638 $ v7*c( j, 7 ) + v8*c( j, 8 ) + v9*c( j, 9 )
639 c( j, 1 ) = c( j, 1 ) - sum*t1
640 c( j, 2 ) = c( j, 2 ) - sum*t2
641 c( j, 3 ) = c( j, 3 ) - sum*t3
642 c( j, 4 ) = c( j, 4 ) - sum*t4
643 c( j, 5 ) = c( j, 5 ) - sum*t5
644 c( j, 6 ) = c( j, 6 ) - sum*t6
645 c( j, 7 ) = c( j, 7 ) - sum*t7
646 c( j, 8 ) = c( j, 8 ) - sum*t8
647 c( j, 9 ) = c( j, 9 ) - sum*t9
648 380 CONTINUE
649 GO TO 410
650 390 CONTINUE
651*
652* Special code for 10 x 10 Householder
653*
654 v1 = v( 1 )
655 t1 = tau*dconjg( v1 )
656 v2 = v( 2 )
657 t2 = tau*dconjg( v2 )
658 v3 = v( 3 )
659 t3 = tau*dconjg( v3 )
660 v4 = v( 4 )
661 t4 = tau*dconjg( v4 )
662 v5 = v( 5 )
663 t5 = tau*dconjg( v5 )
664 v6 = v( 6 )
665 t6 = tau*dconjg( v6 )
666 v7 = v( 7 )
667 t7 = tau*dconjg( v7 )
668 v8 = v( 8 )
669 t8 = tau*dconjg( v8 )
670 v9 = v( 9 )
671 t9 = tau*dconjg( v9 )
672 v10 = v( 10 )
673 t10 = tau*dconjg( v10 )
674 DO 400 j = 1, m
675 sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
676 $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
677 $ v7*c( j, 7 ) + v8*c( j, 8 ) + v9*c( j, 9 ) +
678 $ v10*c( j, 10 )
679 c( j, 1 ) = c( j, 1 ) - sum*t1
680 c( j, 2 ) = c( j, 2 ) - sum*t2
681 c( j, 3 ) = c( j, 3 ) - sum*t3
682 c( j, 4 ) = c( j, 4 ) - sum*t4
683 c( j, 5 ) = c( j, 5 ) - sum*t5
684 c( j, 6 ) = c( j, 6 ) - sum*t6
685 c( j, 7 ) = c( j, 7 ) - sum*t7
686 c( j, 8 ) = c( j, 8 ) - sum*t8
687 c( j, 9 ) = c( j, 9 ) - sum*t9
688 c( j, 10 ) = c( j, 10 ) - sum*t10
689 400 CONTINUE
690 GO TO 410
691 END IF
692 410 CONTINUE
693 RETURN
694*
695* End of ZLARFX
696*
697 END
subroutine zlarf(side, m, n, v, incv, tau, c, ldc, work)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition zlarf.f:128
subroutine zlarfx(side, m, n, v, tau, c, ldc, work)
ZLARFX applies an elementary reflector to a general rectangular matrix, with loop unrolling when the ...
Definition zlarfx.f:119