LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
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. LDA >= 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 *> \date September 2012
116 *
117 *> \ingroup complex16OTHERauxiliary
118 *
119 * =====================================================================
120  SUBROUTINE zlarfx( SIDE, M, N, V, TAU, C, LDC, WORK )
121 *
122 * -- LAPACK auxiliary routine (version 3.4.2) --
123 * -- LAPACK is a software package provided by Univ. of Tennessee, --
124 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
125 * September 2012
126 *
127 * .. Scalar Arguments ..
128  CHARACTER side
129  INTEGER ldc, m, n
130  COMPLEX*16 tau
131 * ..
132 * .. Array Arguments ..
133  COMPLEX*16 c( ldc, * ), v( * ), work( * )
134 * ..
135 *
136 * =====================================================================
137 *
138 * .. Parameters ..
139  COMPLEX*16 zero, one
140  parameter( zero = ( 0.0d+0, 0.0d+0 ),
141  $ one = ( 1.0d+0, 0.0d+0 ) )
142 * ..
143 * .. Local Scalars ..
144  INTEGER j
145  COMPLEX*16 sum, t1, t10, t2, t3, t4, t5, t6, t7, t8, t9,
146  $ v1, v10, v2, v3, v4, v5, v6, v7, v8, v9
147 * ..
148 * .. External Functions ..
149  LOGICAL lsame
150  EXTERNAL lsame
151 * ..
152 * .. External Subroutines ..
153  EXTERNAL zlarf
154 * ..
155 * .. Intrinsic Functions ..
156  INTRINSIC dconjg
157 * ..
158 * .. Executable Statements ..
159 *
160  IF( tau.EQ.zero )
161  $ return
162  IF( lsame( side, 'L' ) ) THEN
163 *
164 * Form H * C, where H has order m.
165 *
166  go to( 10, 30, 50, 70, 90, 110, 130, 150,
167  $ 170, 190 )m
168 *
169 * Code for general M
170 *
171  CALL zlarf( side, m, n, v, 1, tau, c, ldc, work )
172  go to 410
173  10 continue
174 *
175 * Special code for 1 x 1 Householder
176 *
177  t1 = one - tau*v( 1 )*dconjg( v( 1 ) )
178  DO 20 j = 1, n
179  c( 1, j ) = t1*c( 1, j )
180  20 continue
181  go to 410
182  30 continue
183 *
184 * Special code for 2 x 2 Householder
185 *
186  v1 = dconjg( v( 1 ) )
187  t1 = tau*dconjg( v1 )
188  v2 = dconjg( v( 2 ) )
189  t2 = tau*dconjg( v2 )
190  DO 40 j = 1, n
191  sum = v1*c( 1, j ) + v2*c( 2, j )
192  c( 1, j ) = c( 1, j ) - sum*t1
193  c( 2, j ) = c( 2, j ) - sum*t2
194  40 continue
195  go to 410
196  50 continue
197 *
198 * Special code for 3 x 3 Householder
199 *
200  v1 = dconjg( v( 1 ) )
201  t1 = tau*dconjg( v1 )
202  v2 = dconjg( v( 2 ) )
203  t2 = tau*dconjg( v2 )
204  v3 = dconjg( v( 3 ) )
205  t3 = tau*dconjg( v3 )
206  DO 60 j = 1, n
207  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j )
208  c( 1, j ) = c( 1, j ) - sum*t1
209  c( 2, j ) = c( 2, j ) - sum*t2
210  c( 3, j ) = c( 3, j ) - sum*t3
211  60 continue
212  go to 410
213  70 continue
214 *
215 * Special code for 4 x 4 Householder
216 *
217  v1 = dconjg( v( 1 ) )
218  t1 = tau*dconjg( v1 )
219  v2 = dconjg( v( 2 ) )
220  t2 = tau*dconjg( v2 )
221  v3 = dconjg( v( 3 ) )
222  t3 = tau*dconjg( v3 )
223  v4 = dconjg( v( 4 ) )
224  t4 = tau*dconjg( v4 )
225  DO 80 j = 1, n
226  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
227  $ v4*c( 4, j )
228  c( 1, j ) = c( 1, j ) - sum*t1
229  c( 2, j ) = c( 2, j ) - sum*t2
230  c( 3, j ) = c( 3, j ) - sum*t3
231  c( 4, j ) = c( 4, j ) - sum*t4
232  80 continue
233  go to 410
234  90 continue
235 *
236 * Special code for 5 x 5 Householder
237 *
238  v1 = dconjg( v( 1 ) )
239  t1 = tau*dconjg( v1 )
240  v2 = dconjg( v( 2 ) )
241  t2 = tau*dconjg( v2 )
242  v3 = dconjg( v( 3 ) )
243  t3 = tau*dconjg( v3 )
244  v4 = dconjg( v( 4 ) )
245  t4 = tau*dconjg( v4 )
246  v5 = dconjg( v( 5 ) )
247  t5 = tau*dconjg( v5 )
248  DO 100 j = 1, n
249  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
250  $ v4*c( 4, j ) + v5*c( 5, j )
251  c( 1, j ) = c( 1, j ) - sum*t1
252  c( 2, j ) = c( 2, j ) - sum*t2
253  c( 3, j ) = c( 3, j ) - sum*t3
254  c( 4, j ) = c( 4, j ) - sum*t4
255  c( 5, j ) = c( 5, j ) - sum*t5
256  100 continue
257  go to 410
258  110 continue
259 *
260 * Special code for 6 x 6 Householder
261 *
262  v1 = dconjg( v( 1 ) )
263  t1 = tau*dconjg( v1 )
264  v2 = dconjg( v( 2 ) )
265  t2 = tau*dconjg( v2 )
266  v3 = dconjg( v( 3 ) )
267  t3 = tau*dconjg( v3 )
268  v4 = dconjg( v( 4 ) )
269  t4 = tau*dconjg( v4 )
270  v5 = dconjg( v( 5 ) )
271  t5 = tau*dconjg( v5 )
272  v6 = dconjg( v( 6 ) )
273  t6 = tau*dconjg( v6 )
274  DO 120 j = 1, n
275  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
276  $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j )
277  c( 1, j ) = c( 1, j ) - sum*t1
278  c( 2, j ) = c( 2, j ) - sum*t2
279  c( 3, j ) = c( 3, j ) - sum*t3
280  c( 4, j ) = c( 4, j ) - sum*t4
281  c( 5, j ) = c( 5, j ) - sum*t5
282  c( 6, j ) = c( 6, j ) - sum*t6
283  120 continue
284  go to 410
285  130 continue
286 *
287 * Special code for 7 x 7 Householder
288 *
289  v1 = dconjg( v( 1 ) )
290  t1 = tau*dconjg( v1 )
291  v2 = dconjg( v( 2 ) )
292  t2 = tau*dconjg( v2 )
293  v3 = dconjg( v( 3 ) )
294  t3 = tau*dconjg( v3 )
295  v4 = dconjg( v( 4 ) )
296  t4 = tau*dconjg( v4 )
297  v5 = dconjg( v( 5 ) )
298  t5 = tau*dconjg( v5 )
299  v6 = dconjg( v( 6 ) )
300  t6 = tau*dconjg( v6 )
301  v7 = dconjg( v( 7 ) )
302  t7 = tau*dconjg( v7 )
303  DO 140 j = 1, n
304  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
305  $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
306  $ v7*c( 7, j )
307  c( 1, j ) = c( 1, j ) - sum*t1
308  c( 2, j ) = c( 2, j ) - sum*t2
309  c( 3, j ) = c( 3, j ) - sum*t3
310  c( 4, j ) = c( 4, j ) - sum*t4
311  c( 5, j ) = c( 5, j ) - sum*t5
312  c( 6, j ) = c( 6, j ) - sum*t6
313  c( 7, j ) = c( 7, j ) - sum*t7
314  140 continue
315  go to 410
316  150 continue
317 *
318 * Special code for 8 x 8 Householder
319 *
320  v1 = dconjg( v( 1 ) )
321  t1 = tau*dconjg( v1 )
322  v2 = dconjg( v( 2 ) )
323  t2 = tau*dconjg( v2 )
324  v3 = dconjg( v( 3 ) )
325  t3 = tau*dconjg( v3 )
326  v4 = dconjg( v( 4 ) )
327  t4 = tau*dconjg( v4 )
328  v5 = dconjg( v( 5 ) )
329  t5 = tau*dconjg( v5 )
330  v6 = dconjg( v( 6 ) )
331  t6 = tau*dconjg( v6 )
332  v7 = dconjg( v( 7 ) )
333  t7 = tau*dconjg( v7 )
334  v8 = dconjg( v( 8 ) )
335  t8 = tau*dconjg( v8 )
336  DO 160 j = 1, n
337  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
338  $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
339  $ v7*c( 7, j ) + v8*c( 8, j )
340  c( 1, j ) = c( 1, j ) - sum*t1
341  c( 2, j ) = c( 2, j ) - sum*t2
342  c( 3, j ) = c( 3, j ) - sum*t3
343  c( 4, j ) = c( 4, j ) - sum*t4
344  c( 5, j ) = c( 5, j ) - sum*t5
345  c( 6, j ) = c( 6, j ) - sum*t6
346  c( 7, j ) = c( 7, j ) - sum*t7
347  c( 8, j ) = c( 8, j ) - sum*t8
348  160 continue
349  go to 410
350  170 continue
351 *
352 * Special code for 9 x 9 Householder
353 *
354  v1 = dconjg( v( 1 ) )
355  t1 = tau*dconjg( v1 )
356  v2 = dconjg( v( 2 ) )
357  t2 = tau*dconjg( v2 )
358  v3 = dconjg( v( 3 ) )
359  t3 = tau*dconjg( v3 )
360  v4 = dconjg( v( 4 ) )
361  t4 = tau*dconjg( v4 )
362  v5 = dconjg( v( 5 ) )
363  t5 = tau*dconjg( v5 )
364  v6 = dconjg( v( 6 ) )
365  t6 = tau*dconjg( v6 )
366  v7 = dconjg( v( 7 ) )
367  t7 = tau*dconjg( v7 )
368  v8 = dconjg( v( 8 ) )
369  t8 = tau*dconjg( v8 )
370  v9 = dconjg( v( 9 ) )
371  t9 = tau*dconjg( v9 )
372  DO 180 j = 1, n
373  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
374  $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
375  $ v7*c( 7, j ) + v8*c( 8, j ) + v9*c( 9, j )
376  c( 1, j ) = c( 1, j ) - sum*t1
377  c( 2, j ) = c( 2, j ) - sum*t2
378  c( 3, j ) = c( 3, j ) - sum*t3
379  c( 4, j ) = c( 4, j ) - sum*t4
380  c( 5, j ) = c( 5, j ) - sum*t5
381  c( 6, j ) = c( 6, j ) - sum*t6
382  c( 7, j ) = c( 7, j ) - sum*t7
383  c( 8, j ) = c( 8, j ) - sum*t8
384  c( 9, j ) = c( 9, j ) - sum*t9
385  180 continue
386  go to 410
387  190 continue
388 *
389 * Special code for 10 x 10 Householder
390 *
391  v1 = dconjg( v( 1 ) )
392  t1 = tau*dconjg( v1 )
393  v2 = dconjg( v( 2 ) )
394  t2 = tau*dconjg( v2 )
395  v3 = dconjg( v( 3 ) )
396  t3 = tau*dconjg( v3 )
397  v4 = dconjg( v( 4 ) )
398  t4 = tau*dconjg( v4 )
399  v5 = dconjg( v( 5 ) )
400  t5 = tau*dconjg( v5 )
401  v6 = dconjg( v( 6 ) )
402  t6 = tau*dconjg( v6 )
403  v7 = dconjg( v( 7 ) )
404  t7 = tau*dconjg( v7 )
405  v8 = dconjg( v( 8 ) )
406  t8 = tau*dconjg( v8 )
407  v9 = dconjg( v( 9 ) )
408  t9 = tau*dconjg( v9 )
409  v10 = dconjg( v( 10 ) )
410  t10 = tau*dconjg( v10 )
411  DO 200 j = 1, n
412  sum = v1*c( 1, j ) + v2*c( 2, j ) + v3*c( 3, j ) +
413  $ v4*c( 4, j ) + v5*c( 5, j ) + v6*c( 6, j ) +
414  $ v7*c( 7, j ) + v8*c( 8, j ) + v9*c( 9, j ) +
415  $ v10*c( 10, j )
416  c( 1, j ) = c( 1, j ) - sum*t1
417  c( 2, j ) = c( 2, j ) - sum*t2
418  c( 3, j ) = c( 3, j ) - sum*t3
419  c( 4, j ) = c( 4, j ) - sum*t4
420  c( 5, j ) = c( 5, j ) - sum*t5
421  c( 6, j ) = c( 6, j ) - sum*t6
422  c( 7, j ) = c( 7, j ) - sum*t7
423  c( 8, j ) = c( 8, j ) - sum*t8
424  c( 9, j ) = c( 9, j ) - sum*t9
425  c( 10, j ) = c( 10, j ) - sum*t10
426  200 continue
427  go to 410
428  ELSE
429 *
430 * Form C * H, where H has order n.
431 *
432  go to( 210, 230, 250, 270, 290, 310, 330, 350,
433  $ 370, 390 )n
434 *
435 * Code for general N
436 *
437  CALL zlarf( side, m, n, v, 1, tau, c, ldc, work )
438  go to 410
439  210 continue
440 *
441 * Special code for 1 x 1 Householder
442 *
443  t1 = one - tau*v( 1 )*dconjg( v( 1 ) )
444  DO 220 j = 1, m
445  c( j, 1 ) = t1*c( j, 1 )
446  220 continue
447  go to 410
448  230 continue
449 *
450 * Special code for 2 x 2 Householder
451 *
452  v1 = v( 1 )
453  t1 = tau*dconjg( v1 )
454  v2 = v( 2 )
455  t2 = tau*dconjg( v2 )
456  DO 240 j = 1, m
457  sum = v1*c( j, 1 ) + v2*c( j, 2 )
458  c( j, 1 ) = c( j, 1 ) - sum*t1
459  c( j, 2 ) = c( j, 2 ) - sum*t2
460  240 continue
461  go to 410
462  250 continue
463 *
464 * Special code for 3 x 3 Householder
465 *
466  v1 = v( 1 )
467  t1 = tau*dconjg( v1 )
468  v2 = v( 2 )
469  t2 = tau*dconjg( v2 )
470  v3 = v( 3 )
471  t3 = tau*dconjg( v3 )
472  DO 260 j = 1, m
473  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 )
474  c( j, 1 ) = c( j, 1 ) - sum*t1
475  c( j, 2 ) = c( j, 2 ) - sum*t2
476  c( j, 3 ) = c( j, 3 ) - sum*t3
477  260 continue
478  go to 410
479  270 continue
480 *
481 * Special code for 4 x 4 Householder
482 *
483  v1 = v( 1 )
484  t1 = tau*dconjg( v1 )
485  v2 = v( 2 )
486  t2 = tau*dconjg( v2 )
487  v3 = v( 3 )
488  t3 = tau*dconjg( v3 )
489  v4 = v( 4 )
490  t4 = tau*dconjg( v4 )
491  DO 280 j = 1, m
492  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
493  $ v4*c( j, 4 )
494  c( j, 1 ) = c( j, 1 ) - sum*t1
495  c( j, 2 ) = c( j, 2 ) - sum*t2
496  c( j, 3 ) = c( j, 3 ) - sum*t3
497  c( j, 4 ) = c( j, 4 ) - sum*t4
498  280 continue
499  go to 410
500  290 continue
501 *
502 * Special code for 5 x 5 Householder
503 *
504  v1 = v( 1 )
505  t1 = tau*dconjg( v1 )
506  v2 = v( 2 )
507  t2 = tau*dconjg( v2 )
508  v3 = v( 3 )
509  t3 = tau*dconjg( v3 )
510  v4 = v( 4 )
511  t4 = tau*dconjg( v4 )
512  v5 = v( 5 )
513  t5 = tau*dconjg( v5 )
514  DO 300 j = 1, m
515  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
516  $ v4*c( j, 4 ) + v5*c( j, 5 )
517  c( j, 1 ) = c( j, 1 ) - sum*t1
518  c( j, 2 ) = c( j, 2 ) - sum*t2
519  c( j, 3 ) = c( j, 3 ) - sum*t3
520  c( j, 4 ) = c( j, 4 ) - sum*t4
521  c( j, 5 ) = c( j, 5 ) - sum*t5
522  300 continue
523  go to 410
524  310 continue
525 *
526 * Special code for 6 x 6 Householder
527 *
528  v1 = v( 1 )
529  t1 = tau*dconjg( v1 )
530  v2 = v( 2 )
531  t2 = tau*dconjg( v2 )
532  v3 = v( 3 )
533  t3 = tau*dconjg( v3 )
534  v4 = v( 4 )
535  t4 = tau*dconjg( v4 )
536  v5 = v( 5 )
537  t5 = tau*dconjg( v5 )
538  v6 = v( 6 )
539  t6 = tau*dconjg( v6 )
540  DO 320 j = 1, m
541  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
542  $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 )
543  c( j, 1 ) = c( j, 1 ) - sum*t1
544  c( j, 2 ) = c( j, 2 ) - sum*t2
545  c( j, 3 ) = c( j, 3 ) - sum*t3
546  c( j, 4 ) = c( j, 4 ) - sum*t4
547  c( j, 5 ) = c( j, 5 ) - sum*t5
548  c( j, 6 ) = c( j, 6 ) - sum*t6
549  320 continue
550  go to 410
551  330 continue
552 *
553 * Special code for 7 x 7 Householder
554 *
555  v1 = v( 1 )
556  t1 = tau*dconjg( v1 )
557  v2 = v( 2 )
558  t2 = tau*dconjg( v2 )
559  v3 = v( 3 )
560  t3 = tau*dconjg( v3 )
561  v4 = v( 4 )
562  t4 = tau*dconjg( v4 )
563  v5 = v( 5 )
564  t5 = tau*dconjg( v5 )
565  v6 = v( 6 )
566  t6 = tau*dconjg( v6 )
567  v7 = v( 7 )
568  t7 = tau*dconjg( v7 )
569  DO 340 j = 1, m
570  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
571  $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
572  $ v7*c( j, 7 )
573  c( j, 1 ) = c( j, 1 ) - sum*t1
574  c( j, 2 ) = c( j, 2 ) - sum*t2
575  c( j, 3 ) = c( j, 3 ) - sum*t3
576  c( j, 4 ) = c( j, 4 ) - sum*t4
577  c( j, 5 ) = c( j, 5 ) - sum*t5
578  c( j, 6 ) = c( j, 6 ) - sum*t6
579  c( j, 7 ) = c( j, 7 ) - sum*t7
580  340 continue
581  go to 410
582  350 continue
583 *
584 * Special code for 8 x 8 Householder
585 *
586  v1 = v( 1 )
587  t1 = tau*dconjg( v1 )
588  v2 = v( 2 )
589  t2 = tau*dconjg( v2 )
590  v3 = v( 3 )
591  t3 = tau*dconjg( v3 )
592  v4 = v( 4 )
593  t4 = tau*dconjg( v4 )
594  v5 = v( 5 )
595  t5 = tau*dconjg( v5 )
596  v6 = v( 6 )
597  t6 = tau*dconjg( v6 )
598  v7 = v( 7 )
599  t7 = tau*dconjg( v7 )
600  v8 = v( 8 )
601  t8 = tau*dconjg( v8 )
602  DO 360 j = 1, m
603  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
604  $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
605  $ v7*c( j, 7 ) + v8*c( j, 8 )
606  c( j, 1 ) = c( j, 1 ) - sum*t1
607  c( j, 2 ) = c( j, 2 ) - sum*t2
608  c( j, 3 ) = c( j, 3 ) - sum*t3
609  c( j, 4 ) = c( j, 4 ) - sum*t4
610  c( j, 5 ) = c( j, 5 ) - sum*t5
611  c( j, 6 ) = c( j, 6 ) - sum*t6
612  c( j, 7 ) = c( j, 7 ) - sum*t7
613  c( j, 8 ) = c( j, 8 ) - sum*t8
614  360 continue
615  go to 410
616  370 continue
617 *
618 * Special code for 9 x 9 Householder
619 *
620  v1 = v( 1 )
621  t1 = tau*dconjg( v1 )
622  v2 = v( 2 )
623  t2 = tau*dconjg( v2 )
624  v3 = v( 3 )
625  t3 = tau*dconjg( v3 )
626  v4 = v( 4 )
627  t4 = tau*dconjg( v4 )
628  v5 = v( 5 )
629  t5 = tau*dconjg( v5 )
630  v6 = v( 6 )
631  t6 = tau*dconjg( v6 )
632  v7 = v( 7 )
633  t7 = tau*dconjg( v7 )
634  v8 = v( 8 )
635  t8 = tau*dconjg( v8 )
636  v9 = v( 9 )
637  t9 = tau*dconjg( v9 )
638  DO 380 j = 1, m
639  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
640  $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
641  $ v7*c( j, 7 ) + v8*c( j, 8 ) + v9*c( j, 9 )
642  c( j, 1 ) = c( j, 1 ) - sum*t1
643  c( j, 2 ) = c( j, 2 ) - sum*t2
644  c( j, 3 ) = c( j, 3 ) - sum*t3
645  c( j, 4 ) = c( j, 4 ) - sum*t4
646  c( j, 5 ) = c( j, 5 ) - sum*t5
647  c( j, 6 ) = c( j, 6 ) - sum*t6
648  c( j, 7 ) = c( j, 7 ) - sum*t7
649  c( j, 8 ) = c( j, 8 ) - sum*t8
650  c( j, 9 ) = c( j, 9 ) - sum*t9
651  380 continue
652  go to 410
653  390 continue
654 *
655 * Special code for 10 x 10 Householder
656 *
657  v1 = v( 1 )
658  t1 = tau*dconjg( v1 )
659  v2 = v( 2 )
660  t2 = tau*dconjg( v2 )
661  v3 = v( 3 )
662  t3 = tau*dconjg( v3 )
663  v4 = v( 4 )
664  t4 = tau*dconjg( v4 )
665  v5 = v( 5 )
666  t5 = tau*dconjg( v5 )
667  v6 = v( 6 )
668  t6 = tau*dconjg( v6 )
669  v7 = v( 7 )
670  t7 = tau*dconjg( v7 )
671  v8 = v( 8 )
672  t8 = tau*dconjg( v8 )
673  v9 = v( 9 )
674  t9 = tau*dconjg( v9 )
675  v10 = v( 10 )
676  t10 = tau*dconjg( v10 )
677  DO 400 j = 1, m
678  sum = v1*c( j, 1 ) + v2*c( j, 2 ) + v3*c( j, 3 ) +
679  $ v4*c( j, 4 ) + v5*c( j, 5 ) + v6*c( j, 6 ) +
680  $ v7*c( j, 7 ) + v8*c( j, 8 ) + v9*c( j, 9 ) +
681  $ v10*c( j, 10 )
682  c( j, 1 ) = c( j, 1 ) - sum*t1
683  c( j, 2 ) = c( j, 2 ) - sum*t2
684  c( j, 3 ) = c( j, 3 ) - sum*t3
685  c( j, 4 ) = c( j, 4 ) - sum*t4
686  c( j, 5 ) = c( j, 5 ) - sum*t5
687  c( j, 6 ) = c( j, 6 ) - sum*t6
688  c( j, 7 ) = c( j, 7 ) - sum*t7
689  c( j, 8 ) = c( j, 8 ) - sum*t8
690  c( j, 9 ) = c( j, 9 ) - sum*t9
691  c( j, 10 ) = c( j, 10 ) - sum*t10
692  400 continue
693  go to 410
694  END IF
695  410 continue
696  return
697 *
698 * End of ZLARFX
699 *
700  END