-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtrf.ijs
1860 lines (1787 loc) · 68.5 KB
/
trf.ijs
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
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
NB. Triangular factorization
NB.
NB. getrfxxxx Triangular factorization with partial pivoting
NB. of a general matrix
NB. hetrfpx Triangular factorization with full pivoting of
NB. a Hermitian (symmetric) matrix
NB. potrfx Cholesky factorization of a Hermitian
NB. (symmetric) positive definite matrix
NB. pttrfx Triangular factorization of a Hermitian
NB. (symmetric) positive definite tridiagonal
NB. matrix
NB.
NB. testgetrf Test getrfxxxx by general matrix
NB. testhetrf Test hetrfpx by Hermitian (symmetric) matrix
NB. testpotrf Test potrfx by Hermitian (symmetric) positive
NB. definite matrix
NB. testpttrf Test pttrfx by Hermitian (symmetric) positive
NB. definite tridiagonal matrix
NB. testtrf Adv. to make verb to test xxtrfxxxx by matrix
NB. of generator and shape given
NB.
NB. Copyright 2010,2011,2013,2017,2018,2020,2021,2023,2024,
NB. 2025 Igor Zhuravlov
NB.
NB. This file is part of mt
NB.
NB. mt is free software: you can redistribute it and/or
NB. modify it under the terms of the GNU Lesser General
NB. Public License as published by the Free Software
NB. Foundation, either version 3 of the License, or (at your
NB. option) any later version.
NB.
NB. mt is distributed in the hope that it will be useful, but
NB. WITHOUT ANY WARRANTY; without even the implied warranty
NB. of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
NB. See the GNU Lesser General Public License for more
NB. details.
NB.
NB. You should have received a copy of the GNU Lesser General
NB. Public License along with mt. If not, see
NB. <http://www.gnu.org/licenses/>.
NB. =========================================================
NB. Configuration
coclass 'mt'
NB. =========================================================
NB. Local definitions
NB. ---------------------------------------------------------
NB. Blocked code constants
TRFNB=: 64 NB. block size limit
NB. ---------------------------------------------------------
NB. lahefpl
NB.
NB. Description:
NB. Partial triangular factorization of a Hermitian
NB. (symmetric) matrix:
NB. subP * subL1 * subT * subL1^H * subP^H = subA
NB. using the combination of Parlett and Reid, and Aasen
NB. methods [1].
NB.
NB. Syntax:
NB. 'ipo B lto t0o t1o'=. lahefpl ipi;subA;lti;t0i;t1i
NB. where
NB. ipi -:i. (n-i)
NB. subA - (n-i)×(n-i)-matrix to factorize, the
NB. Hermitian (symmetric), the bottom right
NB. part of A
NB. lti - (n-i)-vector, which is defined as:
NB. lti=. ti 0} li
NB. t0i - i-vector, real, leading elements of main
NB. diagonal of T
NB. t1i - max(0,min(i,n-1))-vector, leading elements
NB. of subdiagonal of T
NB. ipo - (n-i)-vector of integers, the inversed full
NB. permutation, corresponding to subP
NB. B - (n-i)×(n-i)-matrix, contains subL110,
NB. subH10^H, strict lower triangle of subL100,
NB. strict upper triangle of subH00^H, and
NB. subA11upd (see layout 1 below)
NB. lto - max(0,n-i-TRFNB)-vector, aimed to be lti at
NB. the next call from hetrfpl
NB. t0o - min(n,i+TRFNB)-vector, real, leading
NB. elements of main diagonal of T
NB. t1o - min(n-1,i+TRFNB)-vector, leading elements
NB. of subdiagonal of T
NB. li - (n-i)-vector, the 1st scaled column of
NB. subL1, having 1 in the 1st element
NB. ti - scalar, the max from li before it was
NB. scaled, can be any value when i=0
NB. subL100 - min(TRFNB,n-i)×min(TRFNB,n-i)-matrix, the
NB. unit lower triangular, top left part of
NB. subL1
NB. subL110 - max(0,n-i-TRFNB)×min(TRFNB,n-i)-matrix,
NB. the bottom left part of subL1
NB. subH00 - min(TRFNB,n-i)×min(TRFNB,n-i)-matrix, the
NB. lower Hessenberg, top left part of subH
NB. subH10 - max(0,n-i-TRFNB)×min(TRFNB,n-i)-matrix,
NB. the bottom left part of subH
NB. subA11upd - max(0,n-i-TRFNB)×max(0,n-i-TRFNB), not yet
NB. factored bottom right part of subA
NB. subP - (n-i)×(n-i)-matrix, boolean, the full
NB. permutation of subA
NB. subL1 - (n-i)×(n-i)-matrix, the unit lower
NB. triangular, being columns i:min(i+TRFNB,n)
NB. of L1
NB. subH - (n-i)×(n-i)-matrix, the lower Hessenberg,
NB. which is defined as:
NB. subH=. subL1 mp subT
NB. subT - (n-i)×(n-i)-matrix, the Hermitian
NB. (symmetric) tridiagonal, being bottom right
NB. part of T
NB. i ∈ {(0:⌈n/TRFNB⌉)*TRFNB}, lIO subA's 1st row
NB. and column within A
NB. A - n×n-matrix, the Hermitian (symmetric)
NB. L1 - n×n-matrix, the unit lower triangular (unit
NB. diagonal is not stored)
NB. T - n×n-matrix, the Hermitian (symmetric)
NB. tridiagonal
NB.
NB. Storage layout 1 (assume TRFNB<n):
NB. input subA: output B:
NB. ( subA00 subA10^H ) TRFNB ( subL100H00h subH10^H ) TRFNB
NB. ( subA10 subA11 ) n-i-TRFNB ( subL110 subA11upd ) n-i-TRFNB
NB. TRFNB n-i-TRFNB TRFNB n-i-TRFNB
NB. where
NB. subL100H00h - TRFNB×TRFNB-matrix, a strict lower
NB. triangle of L100 combined with a strict
NB. upper triangle of subH00^H, the diagonal
NB. is undefined
NB.
NB. Example for output B when TRFNB=2, n=5:
NB. ( * h00h h10h h10h h10h )
NB. ( l100 * h10h h10h h10h )
NB. ( l110 l110 a11u a11u a11u )
NB. ( l110 l110 a11u a11u a11u )
NB. ( l110 l110 a11u a11u a11u )
NB. where
NB. * - value is undefined
NB. l100,l110 - elements of 1st TRFNB columns of subL1
NB. l100 - elements of subL100
NB. l110 - elements of subL110
NB. h00h,h10h - elements of 1st TRFNB rows of subH^H
NB. h00h - elements of subH00^H
NB. h10h - elements of subH10^H
NB. a11u - elements of subA11upd
NB.
NB. Storage layout 2 (assume TRFNB<n):
NB. ISO from head: ISO from tail:
NB. ipi subA lti ipi subA lti
NB. 0 0 -(n-i) -(n-i)
NB. ... ... ... ...
NB. * j j 0 -(n-i-j) -(n-i-j) -(n-i-j)
NB. ... ... ... ... ... ...
NB. * j+p j+p p p-(n-i-j) p-(n-i-j) p-(n-i-j)
NB. ... ... ... ... ... ...
NB. n-i-1 n-i-1 n-i-j-1 -1 -1 -1
NB. where
NB. ipi - (n-i)-vector
NB. subA - (n-i)×(n-i)-matrix
NB. lti - (n-i-j)-vector
NB. j - integer in range 0:TRFNB-1
NB. * - marked elements to be permuted
NB.
NB. Algorithm (assume 0<i<n-TRFNB):
NB. In:
NB. ipi -:i. (n-i)
NB. subA - A[i:n-1,i:n-1]
NB. lti -:ti 0} li
NB. t0i - (T[0,0],T[1,1],...,T[i-1,i-1])
NB. t1i - (T[1,0],T[2,1],...,T[i,i-1])
NB. li - vector to write into subL1[0:n-i-1,0] i.e.
NB. L1[i:n-1,i]
NB. ti - any scalar
NB. Out:
NB. ipo - ip[i:n-1]
NB. B - see layout 1
NB. lto -:to 0} lo
NB. t0o - (T[0,0],T[1,1],...,T[i+TRFNB-1,i+TRFNB-1]) =
NB. (t0i[0],t0i[1],...,t0i[i-1],subT[0,0],subT[1,1],...,subT[TRFNB-1,TRFNB-1])
NB. t1o - (T[1,0],T[2,1],...,T[i+TRFNB,i+TRFNB-1]) =
NB. (t1i[0],t1i[1],...,t1i[i-1],subT[1,0],subT[2,1],...,subT[TRFNB,TRFNB-1])
NB. lo - vector to write into subL1[TRFNB:n-i-1,TRFNB]
NB. i.e. L1[i+TRFNB:n-1,i+TRFNB]
NB. to -:T[i+TRFNB+1,i+TRFNB]
NB. 1) count amount of iterations to do:
NB. J := min(TRFNB,n-i)
NB. 2) do iterations j=0:J-1 :
NB. 2.1) reconstruct (subH[j:n-i-1,j])^H, the tail part
NB. of conjugated j-th column of subH:
NB. w := conj(subH[j:n-i-1,j]) :=
NB. := subA[j,j:n-i-1] - subL1[j,0:j-1]*(subH[j:n-i-1,0:j-1])^H ,
NB. note1: following omitted steps will reconstruct
NB. full subH[0:n-i-1,j] (for j>0 only):
NB. subH[0:j-2,j] := 0
NB. subH[j-1,j] := subT[j,j-1]
NB. note2: both vectors used are the only two parts
NB. of subA j-th row, matrix used is located
NB. above tail part
NB. 2.2) batch write w and lti into subA:
NB. subA[j,j:n-i-1] := w
NB. subA[j:n-i-1,j] := lti
NB. note: heads of both w and lti will be writed in
NB. the same place in the subA main diagonal,
NB. so resulting value will be undefined, btw
NB. it won't be used
NB. 2.3) try to update w:
NB. w := w - subT[j,j-1] * conj(subL1[j:n-i-1,j-1])
NB. w[0] := Re(w[0])
NB. note: now w[0] contains T[j,j], which must be
NB. real
NB. 2.3.1) if failed (i.e. if j=0) then leave w
NB. unchanged:
NB. w := w - 0
NB. 2.4) build new non-transposed non-scaled lti:
NB. lti := conj(w[1:]) - subL1[j+1:n-i-1,j]*T[j,j]
NB. 2.5) prepare non-standard transposition dip:
NB. 2.5.1) find lIO 1st element with max value:
NB. p := liofmax lti
NB. note: to force L1 to be truly diagonally
NB. dominant replace sorim by soris in
NB. liofmax definition
NB. 2.5.2) remap lISO to measure from tail, making
NB. possible to apply the same dip to arrays
NB. of different lengths (see layout 2):
NB. dip := (-(n-i-j),p-(n-i-j))
NB. note: don't waste time for standardizing the
NB. transposition, use Adverse (::) later
NB. instead
NB. 2.6) try to apply the transposition dip to ipi:
NB. ipi[dip[0]] ↔ ipi[dip[1]]
NB. 2.6.1) if failed (i.e. if p=0), then leave ipi
NB. unchanged
NB. 2.7) try to apply the transposition dip to rows and
NB. columns of subA:
NB. subA[dip[0],:] ↔ subA[dip[1],:]
NB. subA[:,dip[0]] ↔ subA[:,dip[1]]
NB. 2.7.1) if failed (i.e. if p=0), then leave subA
NB. unchanged
NB. 2.8) try to apply the transposition dip to lti:
NB. lti[dip[0]] ↔ lti[dip[1]]
NB. note: now lti[0] contains T[j+1,j]
NB. 2.8.1) if failed (i.e. if p=0), then leave lti
NB. unchanged
NB. 2.9) try to scale lti by head, leaving head itself
NB. unchanged:
NB. lti[1:] := lti / lti[0]
NB. note: now lti[1:] contains subL1[j+1:n-i-1,j+1]
NB. i.e. L1[i+j+1:n-1,i+j+1]
NB. 2.9.1) if failed (i.e. if i+j=n-1 , the last
NB. iteration in the last partition takes
NB. place, and lti is empty vector), then
NB. leave lti unchanged
NB. 2.10) update t0i by appending subT[j,j]:
NB. t0i := t0i , T[j,j]
NB. 2.11) try to update t1i by appending subT[j+1,j]):
NB. t1i := t1i , T[j+1,j]
NB. 2.11.1) if failed (i.e. if i+j=n-1 , the last
NB. iteration in the last partition takes
NB. place, and lti is empty vector), then
NB. leave t1i unchanged
NB. 2.12) link modified arrays to form output of current
NB. iteration and input for next one:
NB. ipi ; subA ; lti ; t0i ; t1i
NB. 3) use output from the last iteration as algorithm's
NB. output
NB.
NB. Notes:
NB. - for n>TRFNB models LAPACK's DLASYF_AA('L') and
NB. ZLAHEF_AA('L')
NB. - for n<:TRFNB is similar to LAPACK's DSYTF2('L') and
NB. ZHETF2('L'), but uses another factorization
NB. - diagonals 0 and 1 of subH (main diagonal and
NB. subdiagonal of subH^H) aren't reconstructed since
NB. aren't used
NB.
NB. References:
NB. [1] Miroslav Rozloznik, Gil Shklarski, Sivan Toledo.
NB. Partitioned triangular tridiagonalization.
NB. 26 September 2007.
NB. http://www.cs.cas.cz/miro/rst08.pdf
lahefpl=: (3 : 0)^:(TRFNB<.#@(0&{::))
'ip A lt t0 t1'=. y
w=. lt ((-~&# { ]) ((}.~ #) - ({.~ #) mp ]) ((-~,-@[)&# {. ])) A
A=. (w,:lt) ((0 lisoE),:(0 lisoS))&c} A
w=. (9&o.&.(0&{)) w - lt ((({.@[) * +@((_1 lisoS)&# ({,) ])) :: 0) A
lt=. lt (+@}.@] - ((* }.)~ {.)) w
dip=. (liofmax (-@] <@, -) #) lt
ip=. dip (C. :: ]) ip
A=. dip fp :: ] A
lt=. dip (C. :: ]) lt
lt=. 0 (({`[`(] % {))} :: ]) lt
t0=. t0 , {. w
t1=. t1 , 0 ({ :: ]) lt
ip ; A ; lt ; t0 ; t1
)
NB. ---------------------------------------------------------
NB. lahefpu
NB.
NB. Description:
NB. Partial triangular factorization of a Hermitian
NB. (symmetric) matrix:
NB. subP * subU1 * subT * subU1^H * subP^H = subA
NB. using the combination of Parlett and Reid, and Aasen
NB. methods, inspired by [1].
NB.
NB. Syntax:
NB. 'ipo B uto t0o t1o'=. lahefpu ipi;subA;uti;t0i;t1i
NB. where
NB. ipi -:i. (n+i+1)
NB. subA - (n+i+1)×(n+i+1)-matrix to factorize,
NB. Hermitian (symmetric), the top left part of
NB. A
NB. uti - (n+i+1)-vector, which is defined as:
NB. uti=. ti _1} ui
NB. t0i - (i+1)-vector, real, tail elements of main
NB. diagonal of T
NB. t1i - max(0,min(i+1,n-1))-vector, tail elements
NB. of superdiagonal of T
NB. ipo - (n+i+1)-vector of integers, the inversed
NB. full permutation, corresponding to subP
NB. B - (n+i+1)×(n+i+1)-matrix, contains subU101,
NB. subH01^H, strict upper triangle of subU111,
NB. strict lower triangle of subH11^H, and
NB. subA00upd (see layout 1 below)
NB. uto - max(0,n+i+1-TRFNB)-vector, aimed to be uti
NB. at the next call from hetrfpu
NB. t0o - min(n,i+1+TRFNB)-vector, real, tail
NB. elements of main diagonal of T
NB. t1o - min(n-1,i+1+TRFNB)-vector, tail elements of
NB. superdiagonal of T
NB. ui - (n+i+1)-vector, the last scaled column of
NB. subU1, having 1 in the last element
NB. ti - scalar, the max from ui before it was
NB. scaled, can be any value when i=_1
NB. subU111 - min(TRFNB,n+i+1)×min(TRFNB,n+i+1)-matrix,
NB. the unit upper triangular, bottom right
NB. part of subU1
NB. subU101 - max(0,n+i+1-TRFNB)×min(TRFNB,n+i+1)-matrix,
NB. the top right part of subU1
NB. subH11 - min(TRFNB,n+i+1)×min(TRFNB,n+i+1)-matrix,
NB. the upper Hessenberg, bottom right part of
NB. subH
NB. subH01 - max(0,n+i+1-TRFNB)×min(TRFNB,n+i+1)-matrix,
NB. the top right part of subH
NB. subA00upd - max(0,n+i+1-TRFNB)×max(0,n+i+1-TRFNB), not
NB. yet factored top left part of subA
NB. subP - (n+i+1)×(n+i+1)-matrix, boolean, the full
NB. permutation of subA
NB. subU1 - (n+i+1)×(n+i+1)-matrix, the unit upper
NB. triangular, being columns 0:n+i of U1
NB. subH - (n+i+1)×(n+i+1)-matrix, the upper
NB. Hessenberg, is defined as:
NB. subH=. subU1 mp subT
NB. subT - (n+i+1)×(n+i+1)-matrix, the Hermitian
NB. (symmetric) tridiagonal, being top left
NB. part of T
NB. i ∈ {_1-(0:⌈n/TRFNB⌉)*TRFNB}, lIO subA's last
NB. row and column within A
NB. A - n×n-matrix, the Hermitian (symmetric)
NB. U1 - n×n-matrix, the unit upper triangular (unit
NB. diagonal is not stored)
NB. T - n×n-matrix, the Hermitian (symmetric)
NB. tridiagonal
NB.
NB. Storage layout 1 (assume TRFNB<n):
NB. input subA: output B:
NB. ( subA00 subA10^H ) n+i-TRFNB ( subA00upd subU101 ) n+i-TRFNB
NB. ( subA10 subA11 ) TRFNB ( subH01^H subU111H11h ) TRFNB
NB. n+i-TRFNB TRFNB n+i-TRFNB TRFNB
NB. where
NB. subU111H11h - TRFNB×TRFNB-matrix, a strict upper
NB. triangle of U111 combined with a strict
NB. lower triangle of subH11^H, the diagonal
NB. is undefined
NB.
NB. Example for output B when TRFNB=2, n=5:
NB. ( a00u a00u a00u u101 u101 )
NB. ( a00u a00u a00u u101 u101 )
NB. ( a00u a00u a00u u101 u101 )
NB. ( h01h h01h h01h * u111 )
NB. ( h01h h01h h01h h11h * )
NB. where
NB. * - value is undefined
NB. u101,u111 - elements of last TRFNB columns of subU1
NB. u101 - elements of subU101
NB. u111 - elements of subU111
NB. h01h,h11h - elements of last TRFNB rows of subH^H
NB. h01h - elements of subH01^H
NB. h11h - elements of subH11^H
NB. a00u - elements of subA00upd
NB.
NB. Storage layout 2 (assume TRFNB<n):
NB. ISO from head:
NB. ipi subA uti
NB. 0 0 0
NB. ... ... ...
NB. * p p p
NB. ... ... ...
NB. * n-i-j-1 n-i-j-1 n-i-j-1
NB. ... ...
NB. n-i-1 n-i-1
NB. where
NB. ipi - (n-i)-vector
NB. subA - (n-i)×(n-i)-matrix
NB. uti - (n-i-j)-vector
NB. j - integer in range 0:TRFNB-1
NB. * - marked elements to be permuted
NB.
NB. Algorithm (assume TRFNB-n-1<i<-1):
NB. In:
NB. ipi -:i. (n+i+1)
NB. subA - A[0:n+i,0:n+i]
NB. uti -:ti _1} ui
NB. t0i - (T[n+i+1,n+i+1],T[n+i+2,n+i+2],...,T[n-1,n-1])
NB. t1i - (T[n+i,n+i+1],T[n+i+1,n+i+2],...,T[n-2,n-1])
NB. ui - vector to write into subU1[0:n+i,n+i] i.e.
NB. U1[0:n+i,n+i]
NB. ti - any scalar
NB. Out:
NB. ipo - ip[0:n+i]
NB. B - see layout
NB. uto -:to _1} uo
NB. t0o - (T[n+i-TRFNB+1,n+i-TRFNB+1],T[n+i-TRFNB+2,n+i-TRFNB+2],...,T[n-1,n-1]) =
NB. (subT[0,0],subT[1,1],...,subT[TRFNB-1,TRFNB-1],t0i[0],t0i[1],...,t0i[-i-2])
NB. t1o - (T[n+i-TRFNB,n+i-TRFNB+1],T[n+i-TRFNB+1,n+i-TRFNB+2],...,T[n-2,n-1]) =
NB. (subT[0,1],subT[1,2],...,subT[TRFNB-1,TRFNB],t1i[0],t1i[1],...,t1i[-i-2])
NB. uo - vector to write into subU1[0:n+i-TRFNB,n+i-TRFNB]
NB. i.e. U1[0:n+i-TRFNB,n+i-TRFNB]
NB. to -:T[n+i-TRFNB-1,n+i-TRFNB]
NB. 1) count amount of iterations to do:
NB. J := min(TRFNB,n+i+1)
NB. 2) do iterations j={n+i,n+i-1,...,n+i-J+1}:
NB. 2.1) reconstruct (subH[0:j,j])^H, the head part
NB. of conjugated j-th column of subH:
NB. w := conj(subH[0:j,j]) :=
NB. := subA[j,0:j] - subU1[j,j+1:n-i]*(subH[0:j,j+1:n-i])^H ,
NB. note1: following omitted steps will reconstruct
NB. full subH[0:n+i+1,j] (for j<n+i only):
NB. subH[j+1,j] := subT[j-1,j]
NB. subH[j+2:n+i+1,j] := 0
NB. note2: both vectors used are the only two parts
NB. of subA j-th row, matrix used is located
NB. below head part
NB. 2.2) batch write uti and w into subA:
NB. subA[0:j,j] := uti
NB. subA[j,0:j] := w
NB. note: tails of both uti and w will be writed in
NB. the same place in the subA main diagonal,
NB. so resulting value will be undefined, btw
NB. it won't be used
NB. 2.3) try to update w:
NB. w := w - subT[j-1,j] * conj(subU1[0:j,j+1])
NB. w[j] := Re(w[j])
NB. note: now w[j] contains T[j,j], which must be
NB. real
NB. 2.3.1) if failed (i.e. if j=n+i-J) then leave w
NB. unchanged:
NB. w := w - 0
NB. 2.4) build new non-transposed non-scaled uti:
NB. uti := conj(w[0:j-1]) - subU1[0:j-1,j]*T[j,j]
NB. 2.5) prepare non-standard transposition dip:
NB. 2.5.1) find lIO last element with max value:
NB. p := liolmax uti
NB. note: to force U1 to be truly diagonally
NB. dominant replace sorim by soris in
NB. liolmax definition
NB. 2.5.2) form lISO to measure from head, making
NB. possible to apply the same dip to arrays
NB. of different lengths (see layout 2):
NB. dip := (n-i-j-1,p)
NB. note: don't waste time for standardizing the
NB. transposition, use Adverse (::) later
NB. instead
NB. 2.6) try to apply the transposition dip to ipi:
NB. ipi[dip[0]] ↔ ipi[dip[1]]
NB. 2.6.1) if failed (i.e. if p=n-i-j-1), then
NB. leave ipi unchanged
NB. 2.7) try to apply the transposition dip to rows and
NB. columns of subA:
NB. subA[dip[0],:] ↔ subA[dip[1],:]
NB. subA[:,dip[0]] ↔ subA[:,dip[1]]
NB. 2.7.1) if failed (i.e. if p=n-i-j-1), then
NB. leave subA unchanged
NB. 2.8) try to apply the transposition dip to uti:
NB. uti[dip[0]] ↔ uti[dip[1]]
NB. note: now uti[j] contains T[j,j+1]
NB. 2.8.1) if failed (i.e. if p=n-i-j-1), then
NB. leave uti unchanged
NB. 2.9) try to scale uti by tail, leaving tail itself
NB. unchanged:
NB. uti[0:j-2] := uti / uti[j-1]
NB. note: now uti[0:j-2] contains subU1[0:j-2,j-1]
NB. i.e. U1[0:j-2,j-1]
NB. 2.9.1) if failed (i.e. if j=0 , the last
NB. iteration in the last partition takes
NB. place, and uti is empty vector), then
NB. leave uti unchanged
NB. 2.10) update t0i by prepending subT[j,j]:
NB. t0i := T[j,j] , t0i
NB. 2.11) try to update t1i by prepending subT[j,j+1]):
NB. t1i := T[j,j+1] , t1i
NB. 2.11.1) if failed (i.e. if j=0 , the last
NB. iteration in the last partition takes
NB. place, and uti is empty vector), then
NB. leave t1i unchanged
NB. 2.12) link modified arrays to form output of current
NB. iteration and input for next one:
NB. ipi ; subA ; uti ; t0i ; t1i
NB. 3) use output from the last iteration as algorithm's
NB. output
NB.
NB. Notes:
NB. - diagonals 0 and _1 of subH (main diagonal and
NB. superdiagonal of subH^H) aren't reconstructed since
NB. aren't used
NB.
NB. References:
NB. [1] Miroslav Rozloznik, Gil Shklarski, Sivan Toledo.
NB. Partitioned triangular tridiagonalization.
NB. 26 September 2007.
NB. http://www.cs.cas.cz/miro/rst08.pdf
lahefpu=: (3 : 0)^:(TRFNB<.#@(0&{::))
'ip A ut t0 t1'=. y
w=. ut ((<:@-&# { ]) (({.~ c) - (}.~ c) mp ]) ((-,[)&# {. ])) A
A=. (ut,:w) ((0 lisoN),:(0 lisoW))&c} A
w=. (9&o.&.(_1&{)) w - ut ((({:@[) * +@((1 lisoN)&# ({,) ])) :: 0) A
ut=. ut (+@}:@] - ((* }:)~ {:)) w
dip=. ((<:@(1>.#)) (<@,) liolmax) ut
ip=. dip (C. :: ]) ip
A=. dip fp :: ] A
ut=. dip (C. :: ]) ut
ut=. _1 (({`[`(] % {))} :: ]) ut
t0=. ({: w) , t0
t1=. (_1 ({ :: ]) ut) , t1
ip ; A ; ut ; t0 ; t1
)
NB. =========================================================
NB. Interface
NB. ---------------------------------------------------------
NB. getrflu1p
NB.
NB. Description:
NB. Triangular factorization with partial pivoting of a
NB. general matrix:
NB. L * U1 * P = A
NB.
NB. Syntax:
NB. 'ip LU1'=. getrflu1p A
NB. where
NB. A - m×n-matrix to factorize
NB. ip - n-vector, columns inversed permutation of A
NB. LU1 - m×n-matrix, the lower triangle contains L, and
NB. the strict upper triangle contains U1 without
NB. unit diagonal
NB. P - n×n-matrix, columns permutation of A
NB. L - m×min(m,n)-matrix, the lower trapezoidal
NB. U1 - min(m,n)×n-matrix, the unit upper trapezoidal
NB. (unit diagonal is not stored)
NB.
NB. Storage layout:
NB. A's partitioning: L's partitioning:
NB. k ( Aaa Aab ) := ( Aa ) := A k ( Laa ) := L
NB. m-k ( Aba Abb ) ( Ab ) m-k ( Lba Lbb )
NB. k n-k n k n-k
NB. U1's partitioning: LU1's partitioning:
NB. k ( U1aa U1ab ) := ( U1a ) := U1 k ( LaaU1aa U1ab ) := LU1
NB. m-k ( U1bb ) ( U1b ) m-k ( Lba LbbU1bb )
NB. k n-k n k n-k
NB. where
NB. LaaU1aa - Laa and U1aa combined, U1aa's unit diagonal
NB. is not stored
NB. LbbU1bb - Lbb and U1bb combined, U1bb's unit diagonal
NB. is not stored
NB.
NB. Algorithm:
NB. In: A
NB. Out: ip LU1
NB. 1) acquire geometry of A:
NB. sh := shape
NB. m := quantity of rows
NB. n := quantity of columns
NB. 2) if m=0 or n=0 then:
NB. 2.1) set output:
NB. ip := i. n
NB. LU1 := A
NB. 3) elseif m=1 then:
NB. 3.1) prepare non-standard transposition dip:
NB. 3.1.1) find lIO 1st element with max value in
NB. 1st row of A:
NB. p := liofmax A[0,:]
NB. note: to force U1 to be truly diagonally
NB. dominant replace sorim by soris in
NB. liofmax definition
NB. 3.1.2) compose non-standard transposition:
NB. dip := (0,p)
NB. note: don't waste time for standardizing the
NB. transposition, use Adverse (::) later
NB. instead
NB. 3.2) prepare ip:
NB. 3.2.1) init ip:
NB. ip=. i. n
NB. 3.2.2) try to apply the transposition dip to
NB. ip:
NB. ip[0] ↔ ip[p]
NB. 3.2.3) if failed (i.e. if p=0), then leave ip
NB. unchanged
NB. 3.3) prepare LU1:
NB. 3.3.1) try to apply the transposition dip to
NB. columns of A:
NB. A[:,0] ↔ A[:,p]
NB. 3.3.2) if failed (i.e. if p=0), then leave A
NB. unchanged
NB. 3.3.3) factorize single row:
NB. L := A[0,0]
NB. U1 := L^_1 * A
NB. 4) else:
NB. 4.1) find split point:
NB. k := min(n,⌈m/2⌉)
NB. 4.2) factorize Aa recursively:
NB. Laa * U1a * P = Aa
NB. note1: P is represented by vector of inversed
NB. permutation ip
NB. note2: Laa and U1a are stored in the same
NB. matrix LaaU1a, U1a's unit diagonal isn't
NB. stored
NB. 4.3) permute columns of Ab according to P^H :
NB. Ab := Ab * P^H
NB. note: purge original A, reuse name 'y' to
NB. store Ab
NB. 4.4) compute Lba:
NB. 4.4.1) extract LaaU1aa:
NB. LaaU1aa := LaaU1a[:,0:k-1]
NB. 4.4.2) extract Aba:
NB. Aba := Ab[:,0:k-1]
NB. 4.4.3) solve:
NB. Lba * U1aa = Aba
NB. 4.5) update and factorize Abb recursively:
NB. 4.5.1) extract Abb:
NB. Abb := Ab[k:n-1,:]
NB. 4.5.2) extract U1ab:
NB. U1ab := LaaU1a[k:n-1,:]
NB. 4.5.3) update Abb:
NB. Abb := Abb - Lba * U1ab
NB. 4.5.4) factorize Abb recursively:
NB. Lbb * U1bb * Pb = Abb
NB. note: Pb is represented by vector of
NB. inversed permutation ipb
NB. 4.6) prepare delta of inversed permutation dipb,
NB. which defines inversed permutation of tail part
NB. 4.7) assemble output
NB. 4.7.1) permute tail part of ip according to
NB. Pb^H :
NB. ip[k:n-1] := Pb^H * ip[k:n-1]
NB. 4.7.2) permute columns of U1ab according to
NB. Pb^H :
NB. U1ab := U1ab * Pb^H
NB. 4.7.3) assemble trapezoidal matrices L and U1
NB.
NB. Assertions:
NB. P -: %. iP
NB. P -: |: iP
NB. P -: P4ip ip
NB. A -: clean p {"1 L mp U1
NB. A -: clean p C."1 L mp U1
NB. A -: clean ip C.^:_1"1 L mp U1
NB. A -: clean L mp U1 mp iP NB. apply p to columns
NB. where
NB. 'ip LU1'=. getrflu1p A
NB. p=. /: ip
NB. iP=. P4p ip
NB. P=. P4p p
NB. L=. trl LU1
NB. U1=. tru1 LU1
getrflu1p=: 3 : 0
'm n'=. sh=. $ y
if. 0 e. sh do.
(i. n) ; y
elseif. 1 = m do.
dip=. < 0 , liofmax {. y
ip=. dip (C. :: ]) i. n
y=. ((] 0:} %) 0&({,)) dip (C."1 :: ]) y
ip ; y
else.
k=. n (<. >.@-:) m
'ip LaaU1a'=. getrflu1p k {. y
y=. ip (C."1) k }. y
Lba=. LaaU1a trsmrunu&(k&({."1)) y
'ipb LbbU1bb'=. getrflu1p y (- Lba&mp)&(k&(}."1)) LaaU1a
dipb=. (i. k) , (k + ipb)
(dipb C. ip) ; (dipb (C."1) LaaU1a) , (Lba ,. LbbU1bb)
end.
)
NB. ---------------------------------------------------------
NB. getrfpl1u
NB.
NB. Description:
NB. Triangular factorization with partial pivoting of a
NB. general matrix:
NB. P * L1 * U = A
NB.
NB. Syntax:
NB. 'ip L1U'=. getrfpl1u A
NB. where
NB. A - m×n-matrix to factorize
NB. ip - m-vector, rows inversed permutation of A
NB. L1U - m×n-matrix, the upper triangle contains U, and
NB. the strict lower triangle contains L1 without
NB. unit diagonal
NB. P - n×n-matrix, rows permutation of A
NB. L1 - m×min(m,n)-matrix, the unit lower trapezoidal
NB. (unit diagonal is not stored)
NB. U - min(m,n)×n-matrix, the upper trapezoidal
NB.
NB. Storage layout:
NB. A's partitioning: U's partitioning:
NB. k ( Aaa Aba ) := ( Aa Ab ) := A k ( Uaa Uba ) := U
NB. m-k ( Aab Abb ) m-k ( Ubb )
NB. k n-k k n-k k n-k
NB. L1's partitioning: L1U's partitioning:
NB. k ( L1aa ) := ( L1a L1b ) := L1 k ( L1aaUaa Uba ) := L1U
NB. m-k ( L1ab L1bb ) m-k ( L1ab L1bbUbb )
NB. k n-k k n-k k n-k
NB. where
NB. L1aaUaa - L1aa and Uaa combined, L1aa's unit diagonal
NB. is not stored
NB. L1bbUbb - L1bb and Ubb combined, L1bb's unit diagonal
NB. is not stored
NB.
NB. Algorithm:
NB. In: A
NB. Out: ip L1U
NB. 1) acquire geometry of A:
NB. sh := shape
NB. m := quantity of rows
NB. n := quantity of columns
NB. 2) if m=0 or n=0 then:
NB. 2.1) set output:
NB. ip := i. m
NB. L1U := A
NB. 3) elseif n=1 then:
NB. 3.1) prepare non-standard transposition dip:
NB. 3.1.1) find lIO 1st element with max value in
NB. 1st column of A:
NB. p := liofmax A[:,0]
NB. note: to force L1 to be truly diagonally
NB. dominant replace sorim by soris in
NB. liofmax definition
NB. 3.1.2) compose non-standard transposition:
NB. dip := (0,p)
NB. note: don't waste time for standardizing the
NB. transposition, use Adverse (::) later
NB. instead
NB. 3.2) prepare ip:
NB. 3.2.1) init ip:
NB. ip=. i. m
NB. 3.2.2) try to apply the transposition dip to
NB. ip:
NB. ip[0] ↔ ip[p]
NB. 3.2.3) if failed (i.e. if p=0), then leave ip
NB. unchanged
NB. 3.3) prepare L1U:
NB. 3.3.1) try to apply the transposition dip to
NB. rows of A:
NB. A[0,:] ↔ A[p,:]
NB. 3.3.2) if failed (i.e. if p=0), then leave A
NB. unchanged
NB. 3.3.3) factorize single column:
NB. U := A[0,0]
NB. L1 := U^_1 * A
NB. 4) else:
NB. 4.1) find split point:
NB. k := min(m,⌈n/2⌉)
NB. 4.2) factorize Aa recursively:
NB. P * L1a * Uaa = Aa
NB. note1: P is represented by vector of inversed
NB. permutation ip
NB. note2: L1a and Uaa are stored in the same
NB. matrix L1aUaa, L1a's unit diagonal isn't
NB. stored
NB. 4.3) permute rows of Ab according to P^H :
NB. Ab := P^H * Ab
NB. note: purge original A, reuse name 'y' to
NB. store Ab
NB. 4.4) compute Uab:
NB. 4.4.1) extract L1aaUaa:
NB. L1aaUaa := L1aUaa[0:k-1,:]
NB. 4.4.2) extract Aba:
NB. Aba := Ab[0:k-1,:]
NB. 4.4.3) solve:
NB. L1aa * Uba = Aba
NB. 4.5) update and factorize Abb recursively:
NB. 4.5.1) extract Abb:
NB. Abb := Ab[k:m-1,:]
NB. 4.5.2) extract L1ab:
NB. L1ab := L1aUaa[k:m-1,:]
NB. 4.5.3) update Abb:
NB. Abb := Abb - L1ab * Uba
NB. 4.5.4) factorize Abb recursively:
NB. Pb * L1bb * Ubb = Abb
NB. note: Pb is represented by vector of
NB. inversed permutation ipb
NB. 4.6) prepare delta of inversed permutation dipb,
NB. which defines inversed permutation of tail part
NB. 4.7) assemble output
NB. 4.7.1) permute tail part of ip according to
NB. Pb^H :
NB. ip[k:m-1] := Pb^H * ip[k:m-1]
NB. 4.7.2) permute rows of L1ab according to
NB. Pb^H :
NB. L1ab := Pb^H * L1ab
NB. 4.7.3) assemble trapezoidal matrices L1 and U
NB.
NB. Assertions:
NB. P -: %. iP
NB. P -: |: iP
NB. P -: P4ip ip
NB. A -: clean p { L1 mp U
NB. A -: clean p C. L1 mp U
NB. A -: clean ip C.^:_1 L1 mp U
NB. A -: clean P mp L1 mp U
NB. where
NB. 'ip L1U'=. getrfpl1u A
NB. p=. /: ip
NB. iP=. P4p ip
NB. P=. P4p p
NB. L1=. trl1 L1U
NB. U=. tru L1U
NB.
NB. Notes:
NB. - implements LAPACK's xGETRF
getrfpl1u=: 3 : 0
'm n'=. sh=. $ y
if. 0 e. sh do.
(i. m) ; y
elseif. 1 = n do.
dip=. < 0 , liofmax y
ip=. dip C. :: ] i. m
y=. ((] 0:} %) 0&({,)) dip (C. :: ]) y
ip ; y
else.
k=. m (<. >.@-:) n
'pi L1aUaa'=. getrfpl1u k ({."1) y
y=. pi C. k (}."1) y
Uab=. L1aUaa trsmllnu&(k&{.) y
'ipb L1bbUbb'=. getrfpl1u y (- mp&Uab)&(k&}.) L1aUaa
dipb=. (i. k) , (k + ipb)
(dipb C. pi) ; (dipb C. L1aUaa) ,. (Uab , L1bbUbb)
end.
)
NB. ---------------------------------------------------------
NB. getrfpu1l
NB.
NB. Description:
NB. Triangular factorization with partial pivoting of a
NB. general matrix:
NB. P * U1 * L = A
NB.
NB. Syntax:
NB. 'ip U1L'=. getrfpu1l A
NB. where
NB. ip - m-vector, rows inversed permutation of A
NB. A - m×n-matrix to factorize
NB. U1L - m×n-matrix, the lower triangle contains L, and
NB. the strict upper triangle contains U1 without
NB. unit diagonal
NB. P - n×n-matrix, rows permutation of A
NB. L - m×min(m,n)-matrix, the lower trapezoidal
NB. U1 - min(m,n)×n-matrix, the unit upper trapezoidal
NB. (unit diagonal is not stored)
NB.
NB. Storage layout:
NB. A's partitioning: L's partitioning:
NB. m-k ( Aaa Aba ) := ( Aa Ab ) := A m-k ( Laa ) := L
NB. k ( Aab Abb ) k ( Lab Lbb )
NB. n-k k n-k k n-k k
NB. U1's partitioning: U1L's partitioning:
NB. m-k ( U1aa U1ba ) := ( U1a U1b ) := U1 m-k ( U1aaLaa U1ba ) := U1L
NB. k ( U1bb ) k ( Lab U1bbLbb )
NB. n-k k n-k k n-k k
NB. where
NB. U1aaLaa - U1aa and Laa combined, U1aa's unit diagonal
NB. is not stored
NB. U1bbLbb - U1bb and Lbb combined, U1bb's unit diagonal
NB. is not stored
NB.
NB. Algorithm:
NB. In: A
NB. Out: ip U1L
NB. 1) acquire geometry of A:
NB. sh := shape
NB. m := quantity of rows
NB. n := quantity of columns
NB. 2) if m=0 or n=0 then:
NB. 2.1) set output:
NB. ip := i. m
NB. U1L := A
NB. 3) elseif n=1 then:
NB. 3.1) prepare non-standard transposition dip:
NB. 3.1.1) find lIO last element with max value in
NB. last column of A:
NB. p := liolmax A[:,_1]
NB. note: to force U1 to be truly diagonally
NB. dominant replace sorim by soris in
NB. liolmax definition
NB. 3.1.2) compose non-standard transposition:
NB. dip := (_1,p)
NB. note: don't waste time for standardizing the
NB. transposition, use Adverse (::) later
NB. instead
NB. 3.2) prepare ip:
NB. 3.2.1) init ip:
NB. ip=. i. m
NB. 3.2.2) try to apply the transposition dip to
NB. ip:
NB. ip[_1] ↔ ip[p]
NB. 3.2.3) if failed (i.e. if p=m-1), then leave ip
NB. unchanged
NB. 3.3) prepare U1L:
NB. 3.3.1) try to apply the transposition dip to
NB. rows of A:
NB. A[_1,:] ↔ A[p,:]
NB. 3.3.2) if failed (i.e. if p=m-1), then leave A
NB. unchanged
NB. 3.3.3) factorize single column:
NB. L := A[_1,_1]
NB. U1 := A * L^_1
NB. 4) else:
NB. 4.1) find split point:
NB. k := min(m,⌈n/2⌉)
NB. 4.2) factorize Ab recursively:
NB. P * U1b * Lbb = Ab
NB. note1: P is represented by vector of inversed
NB. permutation ip
NB. note2: Lbb and U1b are stored in the same
NB. matrix U1bLbb, U1b's unit diagonal isn't
NB. stored
NB. 4.3) permute rows of Aa according to P^H :
NB. Aa := P^H * Aa
NB. note: purge original A, reuse name 'y' to
NB. store Aa
NB. 4.4) compute Lab:
NB. 4.4.1) extract U1bbLbb:
NB. U1bbLbb := U1bLbb[m-k:m-1,:]
NB. 4.4.2) extract Aab:
NB. Aab := Aa[m-k:m-1,:]
NB. 4.4.3) solve:
NB. U1bb * Lab = Aab
NB. 4.5) update and factorize Aaa recursively:
NB. 4.5.1) extract Aaa:
NB. Aaa := Aa[0:m-k-1,:]
NB. 4.5.2) extract U1ba:
NB. U1ba := U1bLbb[0:m-k-1,:]
NB. 4.5.3) update Aaa:
NB. Aaa := Aaa - U1ba * Lab
NB. 4.5.4) factorize Aaa recursively:
NB. Pa * U1aa * Laa = Aaa
NB. note: Pa is represented by vector of
NB. inversed permutation ipa
NB. 4.6) prepare delta of inversed permutation dipa,
NB. which defines inversed permutation of head part
NB. 4.7) assemble output
NB. 4.7.1) permute head part of ip according to
NB. Pa^H :
NB. ip[0:m-k-1] := Pa^H * ip[0:m-k-1]
NB. 4.7.2) permute rows of U1ba according to
NB. Pa^H :
NB. U1ba := Pa^H * U1ba
NB. 4.7.3) assemble trapezoidal matrices U1 and L
NB.
NB. Assertions:
NB. P -: %. iP
NB. P -: |: iP
NB. P -: P4ip ip
NB. A -: clean p { U1 mp L
NB. A -: clean p C. U1 mp L
NB. A -: clean ip C.^:_1 U1 mp L
NB. A -: clean P mp U1 mp L
NB. where
NB. 'ip U1L'=. getrfpu1l A
NB. p=. /: ip
NB. iP=. P4p ip
NB. P=. P4p p
NB. U1=. (tru1~ -~/@$) U1L
NB. L=. (trl~ -~/@$) U1L
getrfpu1l=: 3 : 0