LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
dlarfx.f
Go to the documentation of this file.
1 *> \brief \b DLARFX 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 DLARFX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlarfx.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlarfx.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlarfx.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DLARFX( SIDE, M, N, V, TAU, C, LDC, WORK )
22 *
23 * .. Scalar Arguments ..
24 * CHARACTER SIDE
25 * INTEGER LDC, M, N
26 * DOUBLE PRECISION TAU
27 * ..
28 * .. Array Arguments ..
29 * DOUBLE PRECISION C( LDC, * ), V( * ), WORK( * )
30 * ..
31 *
32 *
33 *> \par Purpose:
34 * =============
35 *>
36 *> \verbatim
37 *>
38 *> DLARFX applies a real elementary reflector H to a real 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**T
43 *>
44 *> where tau is a real scalar and v is a real 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 DOUBLE PRECISION 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 DOUBLE PRECISION
83 *> The value tau in the representation of H.
84 *> \endverbatim
85 *>
86 *> \param[in,out] C
87 *> \verbatim
88 *> C is DOUBLE PRECISION 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 >= (1,M).
98 *> \endverbatim
99 *>
100 *> \param[out] WORK
101 *> \verbatim
102 *> WORK is DOUBLE PRECISION array, dimension
103 *> (N) if SIDE = 'L'
104 *> or (M) if SIDE = 'R'
105 *> WORK is not referenced if H has order < 11.
106 *> \endverbatim
107 *
108 * Authors:
109 * ========
110 *
111 *> \author Univ. of Tennessee
112 *> \author Univ. of California Berkeley
113 *> \author Univ. of Colorado Denver
114 *> \author NAG Ltd.
115 *
116 *> \date September 2012
117 *
118 *> \ingroup doubleOTHERauxiliary
119 *
120 * =====================================================================
121  SUBROUTINE dlarfx( SIDE, M, N, V, TAU, C, LDC, WORK )
122 *
123 * -- LAPACK auxiliary routine (version 3.4.2) --
124 * -- LAPACK is a software package provided by Univ. of Tennessee, --
125 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
126 * September 2012
127 *
128 * .. Scalar Arguments ..
129  CHARACTER side
130  INTEGER ldc, m, n
131  DOUBLE PRECISION tau
132 * ..
133 * .. Array Arguments ..
134  DOUBLE PRECISION c( ldc, * ), v( * ), work( * )
135 * ..
136 *
137 * =====================================================================
138 *
139 * .. Parameters ..
140  DOUBLE PRECISION zero, one
141  parameter( zero = 0.0d+0, one = 1.0d+0 )
142 * ..
143 * .. Local Scalars ..
144  INTEGER j
145  DOUBLE PRECISION 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 dlarf
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 dlarf( 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 )*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 = v( 1 )
184  t1 = tau*v1
185  v2 = v( 2 )
186  t2 = tau*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 = v( 1 )
198  t1 = tau*v1
199  v2 = v( 2 )
200  t2 = tau*v2
201  v3 = v( 3 )
202  t3 = tau*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 = v( 1 )
215  t1 = tau*v1
216  v2 = v( 2 )
217  t2 = tau*v2
218  v3 = v( 3 )
219  t3 = tau*v3
220  v4 = v( 4 )
221  t4 = tau*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 = v( 1 )
236  t1 = tau*v1
237  v2 = v( 2 )
238  t2 = tau*v2
239  v3 = v( 3 )
240  t3 = tau*v3
241  v4 = v( 4 )
242  t4 = tau*v4
243  v5 = v( 5 )
244  t5 = tau*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 = v( 1 )
260  t1 = tau*v1
261  v2 = v( 2 )
262  t2 = tau*v2
263  v3 = v( 3 )
264  t3 = tau*v3
265  v4 = v( 4 )
266  t4 = tau*v4
267  v5 = v( 5 )
268  t5 = tau*v5
269  v6 = v( 6 )
270  t6 = tau*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 = v( 1 )
287  t1 = tau*v1
288  v2 = v( 2 )
289  t2 = tau*v2
290  v3 = v( 3 )
291  t3 = tau*v3
292  v4 = v( 4 )
293  t4 = tau*v4
294  v5 = v( 5 )
295  t5 = tau*v5
296  v6 = v( 6 )
297  t6 = tau*v6
298  v7 = v( 7 )
299  t7 = tau*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 = v( 1 )
318  t1 = tau*v1
319  v2 = v( 2 )
320  t2 = tau*v2
321  v3 = v( 3 )
322  t3 = tau*v3
323  v4 = v( 4 )
324  t4 = tau*v4
325  v5 = v( 5 )
326  t5 = tau*v5
327  v6 = v( 6 )
328  t6 = tau*v6
329  v7 = v( 7 )
330  t7 = tau*v7
331  v8 = v( 8 )
332  t8 = tau*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 = v( 1 )
352  t1 = tau*v1
353  v2 = v( 2 )
354  t2 = tau*v2
355  v3 = v( 3 )
356  t3 = tau*v3
357  v4 = v( 4 )
358  t4 = tau*v4
359  v5 = v( 5 )
360  t5 = tau*v5
361  v6 = v( 6 )
362  t6 = tau*v6
363  v7 = v( 7 )
364  t7 = tau*v7
365  v8 = v( 8 )
366  t8 = tau*v8
367  v9 = v( 9 )
368  t9 = tau*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 = v( 1 )
389  t1 = tau*v1
390  v2 = v( 2 )
391  t2 = tau*v2
392  v3 = v( 3 )
393  t3 = tau*v3
394  v4 = v( 4 )
395  t4 = tau*v4
396  v5 = v( 5 )
397  t5 = tau*v5
398  v6 = v( 6 )
399  t6 = tau*v6
400  v7 = v( 7 )
401  t7 = tau*v7
402  v8 = v( 8 )
403  t8 = tau*v8
404  v9 = v( 9 )
405  t9 = tau*v9
406  v10 = v( 10 )
407  t10 = tau*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 dlarf( 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 )*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*v1
451  v2 = v( 2 )
452  t2 = tau*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*v1
465  v2 = v( 2 )
466  t2 = tau*v2
467  v3 = v( 3 )
468  t3 = tau*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*v1
482  v2 = v( 2 )
483  t2 = tau*v2
484  v3 = v( 3 )
485  t3 = tau*v3
486  v4 = v( 4 )
487  t4 = tau*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*v1
503  v2 = v( 2 )
504  t2 = tau*v2
505  v3 = v( 3 )
506  t3 = tau*v3
507  v4 = v( 4 )
508  t4 = tau*v4
509  v5 = v( 5 )
510  t5 = tau*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*v1
527  v2 = v( 2 )
528  t2 = tau*v2
529  v3 = v( 3 )
530  t3 = tau*v3
531  v4 = v( 4 )
532  t4 = tau*v4
533  v5 = v( 5 )
534  t5 = tau*v5
535  v6 = v( 6 )
536  t6 = tau*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*v1
554  v2 = v( 2 )
555  t2 = tau*v2
556  v3 = v( 3 )
557  t3 = tau*v3
558  v4 = v( 4 )
559  t4 = tau*v4
560  v5 = v( 5 )
561  t5 = tau*v5
562  v6 = v( 6 )
563  t6 = tau*v6
564  v7 = v( 7 )
565  t7 = tau*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*v1
585  v2 = v( 2 )
586  t2 = tau*v2
587  v3 = v( 3 )
588  t3 = tau*v3
589  v4 = v( 4 )
590  t4 = tau*v4
591  v5 = v( 5 )
592  t5 = tau*v5
593  v6 = v( 6 )
594  t6 = tau*v6
595  v7 = v( 7 )
596  t7 = tau*v7
597  v8 = v( 8 )
598  t8 = tau*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*v1
619  v2 = v( 2 )
620  t2 = tau*v2
621  v3 = v( 3 )
622  t3 = tau*v3
623  v4 = v( 4 )
624  t4 = tau*v4
625  v5 = v( 5 )
626  t5 = tau*v5
627  v6 = v( 6 )
628  t6 = tau*v6
629  v7 = v( 7 )
630  t7 = tau*v7
631  v8 = v( 8 )
632  t8 = tau*v8
633  v9 = v( 9 )
634  t9 = tau*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*v1
656  v2 = v( 2 )
657  t2 = tau*v2
658  v3 = v( 3 )
659  t3 = tau*v3
660  v4 = v( 4 )
661  t4 = tau*v4
662  v5 = v( 5 )
663  t5 = tau*v5
664  v6 = v( 6 )
665  t6 = tau*v6
666  v7 = v( 7 )
667  t7 = tau*v7
668  v8 = v( 8 )
669  t8 = tau*v8
670  v9 = v( 9 )
671  t9 = tau*v9
672  v10 = v( 10 )
673  t10 = tau*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 DLARFX
696 *
697  END