File: slarft.f

package info (click to toggle)
lapack 3.12.1-6
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 78,912 kB
  • sloc: fortran: 622,840; ansic: 217,704; f90: 6,041; makefile: 5,100; sh: 326; python: 270; xml: 182
file content (628 lines) | stat: -rw-r--r-- 19,473 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
*> \brief \b SLARFT forms the triangular factor T of a block reflector H = I - vtvH
*
*  =========== DOCUMENTATION ===========
*
* Online html documentation available at
*            http://www.netlib.org/lapack/explore-html/
*
*> \htmlonly
*> Download SLARFT + dependencies
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/slarft.f">
*> [TGZ]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/slarft.f">
*> [ZIP]</a>
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/slarft.f">
*> [TXT]</a>
*> \endhtmlonly
*
*  Definition:
*  ===========
*
*       RECURSIVE SUBROUTINE SLARFT( DIRECT, STOREV, N, K, V, LDV, TAU, T, LDT )
*
*       .. Scalar Arguments ..
*       CHARACTER          DIRECT, STOREV
*       INTEGER            K, LDT, LDV, N
*       ..
*       .. Array Arguments ..
*       REAL               T( LDT, * ), TAU( * ), V( LDV, * )
*       ..
*
*
*> \par Purpose:
*  =============
*>
*> \verbatim
*>
*> SLARFT forms the triangular factor T of a real block reflector H
*> of order n, which is defined as a product of k elementary reflectors.
*>
*> If DIRECT = 'F', H = H(1) H(2) . . . H(k) and T is upper triangular;
*>
*> If DIRECT = 'B', H = H(k) . . . H(2) H(1) and T is lower triangular.
*>
*> If STOREV = 'C', the vector which defines the elementary reflector
*> H(i) is stored in the i-th column of the array V, and
*>
*>    H  =  I - V * T * V**T
*>
*> If STOREV = 'R', the vector which defines the elementary reflector
*> H(i) is stored in the i-th row of the array V, and
*>
*>    H  =  I - V**T * T * V
*> \endverbatim
*
*  Arguments:
*  ==========
*
*> \param[in] DIRECT
*> \verbatim
*>          DIRECT is CHARACTER*1
*>          Specifies the order in which the elementary reflectors are
*>          multiplied to form the block reflector:
*>          = 'F': H = H(1) H(2) . . . H(k) (Forward)
*>          = 'B': H = H(k) . . . H(2) H(1) (Backward)
*> \endverbatim
*>
*> \param[in] STOREV
*> \verbatim
*>          STOREV is CHARACTER*1
*>          Specifies how the vectors which define the elementary
*>          reflectors are stored (see also Further Details):
*>          = 'C': columnwise
*>          = 'R': rowwise
*> \endverbatim
*>
*> \param[in] N
*> \verbatim
*>          N is INTEGER
*>          The order of the block reflector H. N >= 0.
*> \endverbatim
*>
*> \param[in] K
*> \verbatim
*>          K is INTEGER
*>          The order of the triangular factor T (= the number of
*>          elementary reflectors). K >= 1.
*> \endverbatim
*>
*> \param[in] V
*> \verbatim
*>          V is REAL array, dimension
*>                               (LDV,K) if STOREV = 'C'
*>                               (LDV,N) if STOREV = 'R'
*>          The matrix V. See further details.
*> \endverbatim
*>
*> \param[in] LDV
*> \verbatim
*>          LDV is INTEGER
*>          The leading dimension of the array V.
*>          If STOREV = 'C', LDV >= max(1,N); if STOREV = 'R', LDV >= K.
*> \endverbatim
*>
*> \param[in] TAU
*> \verbatim
*>          TAU is REAL array, dimension (K)
*>          TAU(i) must contain the scalar factor of the elementary
*>          reflector H(i).
*> \endverbatim
*>
*> \param[out] T
*> \verbatim
*>          T is REAL array, dimension (LDT,K)
*>          The k by k triangular factor T of the block reflector.
*>          If DIRECT = 'F', T is upper triangular; if DIRECT = 'B', T is
*>          lower triangular. The rest of the array is not used.
*> \endverbatim
*>
*> \param[in] LDT
*> \verbatim
*>          LDT is INTEGER
*>          The leading dimension of the array T. LDT >= K.
*> \endverbatim
*
*  Authors:
*  ========
*
*> \author Univ. of Tennessee
*> \author Univ. of California Berkeley
*> \author Johnathan Rhyne, Univ. of Colorado Denver (original author, 2024)
*> \author NAG Ltd.
*
*> \ingroup larft
*
*> \par Further Details:
*  =====================
*>
*> \verbatim
*>
*>  The shape of the matrix V and the storage of the vectors which define
*>  the H(i) is best illustrated by the following example with n = 5 and
*>  k = 3. The elements equal to 1 are not stored.
*>
*>  DIRECT = 'F' and STOREV = 'C':         DIRECT = 'F' and STOREV = 'R':
*>
*>               V = (  1       )                 V = (  1 v1 v1 v1 v1 )
*>                   ( v1  1    )                     (     1 v2 v2 v2 )
*>                   ( v1 v2  1 )                     (        1 v3 v3 )
*>                   ( v1 v2 v3 )
*>                   ( v1 v2 v3 )
*>
*>  DIRECT = 'B' and STOREV = 'C':         DIRECT = 'B' and STOREV = 'R':
*>
*>               V = ( v1 v2 v3 )                 V = ( v1 v1  1       )
*>                   ( v1 v2 v3 )                     ( v2 v2 v2  1    )
*>                   (  1 v2 v3 )                     ( v3 v3 v3 v3  1 )
*>                   (     1 v3 )
*>                   (        1 )
*> \endverbatim
*>
*  =====================================================================
      RECURSIVE SUBROUTINE SLARFT( DIRECT, STOREV, N, K, V, LDV,
     $                             TAU, T, LDT )
*
*  -- LAPACK auxiliary routine --
*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
*
*        .. Scalar Arguments
*
      CHARACTER          DIRECT, STOREV
      INTEGER            K, LDT, LDV, N
*     ..
*     .. Array Arguments ..
*
      REAL               T( LDT, * ), TAU( * ), V( LDV, * )
*     ..
*
*     .. Parameters ..
*
      REAL               ONE, NEG_ONE, ZERO
      PARAMETER(ONE=1.0E+0, ZERO = 0.0E+0, NEG_ONE=-1.0E+0)
*
*     .. Local Scalars ..
*
      INTEGER            I,J,L
      LOGICAL            QR,LQ,QL,DIRF,COLV
*
*     .. External Subroutines ..
*
      EXTERNAL           STRMM,SGEMM,SLACPY
*
*     .. External Functions..
*
      LOGICAL            LSAME
      EXTERNAL           LSAME
*     
*     The general scheme used is inspired by the approach inside DGEQRT3
*     which was (at the time of writing this code):
*     Based on the algorithm of Elmroth and Gustavson,
*     IBM J. Res. Develop. Vol 44 No. 4 July 2000.
*     ..
*     .. Executable Statements ..
*
*     Quick return if possible
*
      IF(N.EQ.0.OR.K.EQ.0) THEN
         RETURN
      END IF
*
*     Base case
*
      IF(N.EQ.1.OR.K.EQ.1) THEN
         T(1,1) = TAU(1)
         RETURN
      END IF
*
*     Beginning of executable statements
*
      L = K / 2
*
*     Determine what kind of Q we need to compute
*     We assume that if the user doesn't provide 'F' for DIRECT,
*     then they meant to provide 'B' and if they don't provide
*     'C' for STOREV, then they meant to provide 'R'
*
      DIRF = LSAME(DIRECT,'F')
      COLV = LSAME(STOREV,'C')
*
*     QR happens when we have forward direction in column storage
*
      QR = DIRF.AND.COLV
*
*     LQ happens when we have forward direction in row storage
*
      LQ = DIRF.AND.(.NOT.COLV)
*
*     QL happens when we have backward direction in column storage
*
      QL = (.NOT.DIRF).AND.COLV
*
*     The last case is RQ. Due to how we structured this, if the
*     above 3 are false, then RQ must be true, so we never store 
*     this
*     RQ happens when we have backward direction in row storage
*     RQ = (.NOT.DIRF).AND.(.NOT.COLV)
*
      IF(QR) THEN
*
*        Break V apart into 6 components
*
*        V = |---------------|
*            |V_{1,1} 0      |
*            |V_{2,1} V_{2,2}|
*            |V_{3,1} V_{3,2}|
*            |---------------|
*
*        V_{1,1}\in\R^{l,l}      unit lower triangular
*        V_{2,1}\in\R^{k-l,l}    rectangular
*        V_{3,1}\in\R^{n-k,l}    rectangular
*        
*        V_{2,2}\in\R^{k-l,k-l}  unit lower triangular
*        V_{3,2}\in\R^{n-k,k-l}  rectangular
*
*        We will construct the T matrix 
*        T = |---------------|
*            |T_{1,1} T_{1,2}|
*            |0       T_{2,2}|
*            |---------------|
*
*        T is the triangular factor obtained from block reflectors. 
*        To motivate the structure, assume we have already computed T_{1,1}
*        and T_{2,2}. Then collect the associated reflectors in V_1 and V_2
*
*        T_{1,1}\in\R^{l, l}     upper triangular
*        T_{2,2}\in\R^{k-l, k-l} upper triangular
*        T_{1,2}\in\R^{l, k-l}   rectangular
*
*        Where l = floor(k/2)
*
*        Then, consider the product:
*        
*        (I - V_1*T_{1,1}*V_1')*(I - V_2*T_{2,2}*V_2')
*        = I - V_1*T_{1,1}*V_1' - V_2*T_{2,2}*V_2' + V_1*T_{1,1}*V_1'*V_2*T_{2,2}*V_2'
*        
*        Define T_{1,2} = -T_{1,1}*V_1'*V_2*T_{2,2}
*        
*        Then, we can define the matrix V as 
*        V = |-------|
*            |V_1 V_2|
*            |-------|
*        
*        So, our product is equivalent to the matrix product
*        I - V*T*V'
*        This means, we can compute T_{1,1} and T_{2,2}, then use this information
*        to compute T_{1,2}
*
*        Compute T_{1,1} recursively
*
         CALL SLARFT(DIRECT, STOREV, N, L, V, LDV, TAU, T, LDT)
*
*        Compute T_{2,2} recursively
*
         CALL SLARFT(DIRECT, STOREV, N-L, K-L, V(L+1, L+1), LDV, 
     $               TAU(L+1), T(L+1, L+1), LDT)
*
*        Compute T_{1,2} 
*        T_{1,2} = V_{2,1}'
*
         DO J = 1, L
            DO I = 1, K-L
               T(J, L+I) = V(L+I, J)
            END DO
         END DO
*
*        T_{1,2} = T_{1,2}*V_{2,2}
*
         CALL STRMM('Right', 'Lower', 'No transpose', 'Unit', L,
     $               K-L, ONE, V(L+1, L+1), LDV, T(1, L+1), LDT)

*
*        T_{1,2} = V_{3,1}'*V_{3,2} + T_{1,2}
*        Note: We assume K <= N, and GEMM will do nothing if N=K
*
         CALL SGEMM('Transpose', 'No transpose', L, K-L, N-K, ONE, 
     $               V(K+1, 1), LDV, V(K+1, L+1), LDV, ONE, 
     $               T(1, L+1), LDT)
*
*        At this point, we have that T_{1,2} = V_1'*V_2
*        All that is left is to pre and post multiply by -T_{1,1} and T_{2,2}
*        respectively.
*
*        T_{1,2} = -T_{1,1}*T_{1,2}
*
         CALL STRMM('Left', 'Upper', 'No transpose', 'Non-unit', L,
     $               K-L, NEG_ONE, T, LDT, T(1, L+1), LDT)
*
*        T_{1,2} = T_{1,2}*T_{2,2}
*
         CALL STRMM('Right', 'Upper', 'No transpose', 'Non-unit', L, 
     $               K-L, ONE, T(L+1, L+1), LDT, T(1, L+1), LDT)

      ELSE IF(LQ) THEN
*
*        Break V apart into 6 components
*
*        V = |----------------------|
*            |V_{1,1} V_{1,2} V{1,3}|
*            |0       V_{2,2} V{2,3}|
*            |----------------------|
*
*        V_{1,1}\in\R^{l,l}      unit upper triangular
*        V_{1,2}\in\R^{l,k-l}    rectangular
*        V_{1,3}\in\R^{l,n-k}    rectangular
*        
*        V_{2,2}\in\R^{k-l,k-l}  unit upper triangular
*        V_{2,3}\in\R^{k-l,n-k}  rectangular
*
*        Where l = floor(k/2)
*
*        We will construct the T matrix 
*        T = |---------------|
*            |T_{1,1} T_{1,2}|
*            |0       T_{2,2}|
*            |---------------|
*
*        T is the triangular factor obtained from block reflectors. 
*        To motivate the structure, assume we have already computed T_{1,1}
*        and T_{2,2}. Then collect the associated reflectors in V_1 and V_2
*
*        T_{1,1}\in\R^{l, l}     upper triangular
*        T_{2,2}\in\R^{k-l, k-l} upper triangular
*        T_{1,2}\in\R^{l, k-l}   rectangular
*
*        Then, consider the product:
*        
*        (I - V_1'*T_{1,1}*V_1)*(I - V_2'*T_{2,2}*V_2)
*        = I - V_1'*T_{1,1}*V_1 - V_2'*T_{2,2}*V_2 + V_1'*T_{1,1}*V_1*V_2'*T_{2,2}*V_2
*        
*        Define T_{1,2} = -T_{1,1}*V_1*V_2'*T_{2,2}
*        
*        Then, we can define the matrix V as 
*        V = |---|
*            |V_1|
*            |V_2|
*            |---|
*        
*        So, our product is equivalent to the matrix product
*        I - V'*T*V
*        This means, we can compute T_{1,1} and T_{2,2}, then use this information
*        to compute T_{1,2}
*
*        Compute T_{1,1} recursively
*
         CALL SLARFT(DIRECT, STOREV, N, L, V, LDV, TAU, T, LDT)
*
*        Compute T_{2,2} recursively
*
         CALL SLARFT(DIRECT, STOREV, N-L, K-L, V(L+1, L+1), LDV, 
     $               TAU(L+1), T(L+1, L+1), LDT)

*
*        Compute T_{1,2}
*        T_{1,2} = V_{1,2}
*
         CALL SLACPY('All', L, K-L, V(1, L+1), LDV, T(1, L+1), LDT)
*
*        T_{1,2} = T_{1,2}*V_{2,2}'
*
         CALL STRMM('Right', 'Upper', 'Transpose', 'Unit', L, K-L,
     $               ONE, V(L+1, L+1), LDV, T(1, L+1), LDT)

*
*        T_{1,2} = V_{1,3}*V_{2,3}' + T_{1,2}
*        Note: We assume K <= N, and GEMM will do nothing if N=K
*
         CALL SGEMM('No transpose', 'Transpose', L, K-L, N-K, ONE,
     $               V(1, K+1), LDV, V(L+1, K+1), LDV, ONE, 
     $               T(1, L+1), LDT)
*
*        At this point, we have that T_{1,2} = V_1*V_2'
*        All that is left is to pre and post multiply by -T_{1,1} and T_{2,2}
*        respectively.
*
*        T_{1,2} = -T_{1,1}*T_{1,2}
*
         CALL STRMM('Left', 'Upper', 'No transpose', 'Non-unit', L,
     $               K-L, NEG_ONE, T, LDT, T(1, L+1), LDT)

*
*        T_{1,2} = T_{1,2}*T_{2,2}
*
         CALL STRMM('Right', 'Upper', 'No transpose', 'Non-unit', L,
     $               K-L, ONE, T(L+1, L+1), LDT, T(1, L+1), LDT)
      ELSE IF(QL) THEN
*
*        Break V apart into 6 components
*
*        V = |---------------|
*            |V_{1,1} V_{1,2}|
*            |V_{2,1} V_{2,2}|
*            |0       V_{3,2}|
*            |---------------|
*
*        V_{1,1}\in\R^{n-k,k-l}  rectangular
*        V_{2,1}\in\R^{k-l,k-l}  unit upper triangular
*        
*        V_{1,2}\in\R^{n-k,l}    rectangular
*        V_{2,2}\in\R^{k-l,l}    rectangular
*        V_{3,2}\in\R^{l,l}      unit upper triangular
*
*        We will construct the T matrix 
*        T = |---------------|
*            |T_{1,1} 0      |
*            |T_{2,1} T_{2,2}|
*            |---------------|
*
*        T is the triangular factor obtained from block reflectors. 
*        To motivate the structure, assume we have already computed T_{1,1}
*        and T_{2,2}. Then collect the associated reflectors in V_1 and V_2
*
*        T_{1,1}\in\R^{k-l, k-l} non-unit lower triangular
*        T_{2,2}\in\R^{l, l}     non-unit lower triangular
*        T_{2,1}\in\R^{k-l, l}   rectangular
*
*        Where l = floor(k/2)
*
*        Then, consider the product:
*        
*        (I - V_2*T_{2,2}*V_2')*(I - V_1*T_{1,1}*V_1')
*        = I - V_2*T_{2,2}*V_2' - V_1*T_{1,1}*V_1' + V_2*T_{2,2}*V_2'*V_1*T_{1,1}*V_1'
*        
*        Define T_{2,1} = -T_{2,2}*V_2'*V_1*T_{1,1}
*        
*        Then, we can define the matrix V as 
*        V = |-------|
*            |V_1 V_2|
*            |-------|
*        
*        So, our product is equivalent to the matrix product
*        I - V*T*V'
*        This means, we can compute T_{1,1} and T_{2,2}, then use this information
*        to compute T_{2,1}
*
*        Compute T_{1,1} recursively
*
         CALL SLARFT(DIRECT, STOREV, N-L, K-L, V, LDV, TAU, T, LDT)
*
*        Compute T_{2,2} recursively
*
         CALL SLARFT(DIRECT, STOREV, N, L, V(1, K-L+1), LDV,
     $               TAU(K-L+1), T(K-L+1, K-L+1), LDT)
*
*        Compute T_{2,1}
*        T_{2,1} = V_{2,2}'
*
         DO J = 1, K-L
            DO I = 1, L
               T(K-L+I, J) = V(N-K+J, K-L+I)
            END DO
         END DO
*
*        T_{2,1} = T_{2,1}*V_{2,1}
*
         CALL STRMM('Right', 'Upper', 'No transpose', 'Unit', L,
     $               K-L, ONE, V(N-K+1, 1), LDV, T(K-L+1, 1), LDT)

*
*        T_{2,1} = V_{2,2}'*V_{2,1} + T_{2,1}
*        Note: We assume K <= N, and GEMM will do nothing if N=K
*
         CALL SGEMM('Transpose', 'No transpose', L, K-L, N-K, ONE,
     $               V(1, K-L+1), LDV, V, LDV, ONE, T(K-L+1, 1), 
     $               LDT)
*
*        At this point, we have that T_{2,1} = V_2'*V_1
*        All that is left is to pre and post multiply by -T_{2,2} and T_{1,1}
*        respectively.
*
*        T_{2,1} = -T_{2,2}*T_{2,1}
*
         CALL STRMM('Left', 'Lower', 'No transpose', 'Non-unit', L,
     $               K-L, NEG_ONE, T(K-L+1, K-L+1), LDT, 
     $               T(K-L+1, 1), LDT)
*
*        T_{2,1} = T_{2,1}*T_{1,1}
*
         CALL STRMM('Right', 'Lower', 'No transpose', 'Non-unit', L,
     $               K-L, ONE, T, LDT, T(K-L+1, 1), LDT)
      ELSE
*
*        Else means RQ case
*
*        Break V apart into 6 components
*
*        V = |-----------------------|
*            |V_{1,1} V_{1,2} 0      |
*            |V_{2,1} V_{2,2} V_{2,3}|
*            |-----------------------|
*
*        V_{1,1}\in\R^{k-l,n-k}  rectangular
*        V_{1,2}\in\R^{k-l,k-l}  unit lower triangular
*
*        V_{2,1}\in\R^{l,n-k}    rectangular
*        V_{2,2}\in\R^{l,k-l}    rectangular
*        V_{2,3}\in\R^{l,l}      unit lower triangular
*
*        We will construct the T matrix 
*        T = |---------------|
*            |T_{1,1} 0      |
*            |T_{2,1} T_{2,2}|
*            |---------------|
*
*        T is the triangular factor obtained from block reflectors. 
*        To motivate the structure, assume we have already computed T_{1,1}
*        and T_{2,2}. Then collect the associated reflectors in V_1 and V_2
*
*        T_{1,1}\in\R^{k-l, k-l} non-unit lower triangular
*        T_{2,2}\in\R^{l, l}     non-unit lower triangular
*        T_{2,1}\in\R^{k-l, l}   rectangular
*
*        Where l = floor(k/2)
*
*        Then, consider the product:
*        
*        (I - V_2'*T_{2,2}*V_2)*(I - V_1'*T_{1,1}*V_1)
*        = I - V_2'*T_{2,2}*V_2 - V_1'*T_{1,1}*V_1 + V_2'*T_{2,2}*V_2*V_1'*T_{1,1}*V_1
*        
*        Define T_{2,1} = -T_{2,2}*V_2*V_1'*T_{1,1}
*        
*        Then, we can define the matrix V as 
*        V = |---|
*            |V_1|
*            |V_2|
*            |---|
*        
*        So, our product is equivalent to the matrix product
*        I - V'TV
*        This means, we can compute T_{1,1} and T_{2,2}, then use this information
*        to compute T_{2,1}
*
*        Compute T_{1,1} recursively
*
         CALL SLARFT(DIRECT, STOREV, N-L, K-L, V, LDV, TAU, T, LDT)
*
*        Compute T_{2,2} recursively
*
         CALL SLARFT(DIRECT, STOREV, N, L, V(K-L+1, 1), LDV,
     $               TAU(K-L+1), T(K-L+1, K-L+1), LDT)
*
*        Compute T_{2,1}
*        T_{2,1} = V_{2,2}
*
         CALL SLACPY('All', L, K-L, V(K-L+1, N-K+1), LDV,
     $      T(K-L+1, 1), LDT)

*
*        T_{2,1} = T_{2,1}*V_{1,2}'
*
         CALL STRMM('Right', 'Lower', 'Transpose', 'Unit', L, K-L,
     $               ONE, V(1, N-K+1), LDV, T(K-L+1, 1), LDT)

*
*        T_{2,1} = V_{2,1}*V_{1,1}' + T_{2,1} 
*        Note: We assume K <= N, and GEMM will do nothing if N=K
*
         CALL SGEMM('No transpose', 'Transpose', L, K-L, N-K, ONE, 
     $               V(K-L+1, 1), LDV, V, LDV, ONE, T(K-L+1, 1),
     $               LDT)

*
*        At this point, we have that T_{2,1} = V_2*V_1'
*        All that is left is to pre and post multiply by -T_{2,2} and T_{1,1}
*        respectively.
*
*        T_{2,1} = -T_{2,2}*T_{2,1}
*
         CALL STRMM('Left', 'Lower', 'No tranpose', 'Non-unit', L,
     $               K-L, NEG_ONE, T(K-L+1, K-L+1), LDT, 
     $               T(K-L+1, 1), LDT)

*
*        T_{2,1} = T_{2,1}*T_{1,1}
*
         CALL STRMM('Right', 'Lower', 'No tranpose', 'Non-unit', L,
     $               K-L, ONE, T, LDT, T(K-L+1, 1), LDT)
      END IF
      END SUBROUTINE