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