-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathgq.ijs
684 lines (621 loc) · 26.3 KB
/
gq.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
NB. Generate matrix with orthonormal rows or columns from its
NB. factored form
NB.
NB. unglq Generate a matrix with orthonormal rows from
NB. output of gelqf
NB. ungql Generate a matrix with orthonormal columns
NB. from output of geqlf
NB. ungqr Generate a matrix with orthonormal columns
NB. from output of geqrf
NB. ungrq Generate a matrix with orthonormal rows from
NB. output of gerqf
NB. unglz Generate a matrix with orthonormal rows from
NB. output of tzlzf
NB. ungzl Generate a matrix with orthonormal columns
NB. from output of tzzlf
NB. ungzr Generate a matrix with orthonormal columns
NB. from output of tzzrf
NB. ungrz Generate a matrix with orthonormal rows from
NB. output of tzrzf
NB. unghrx Generate an unitary (orthogonal) matrix which
NB. is defined as the product of elementary
NB. reflectors as returned by gehrdx
NB.
NB. testungq Test ungxx by general matrix
NB. testungz Test ungxx by trapezoidal matrix
NB. testunghr Test unghrx by square matrix
NB. testgq Adv. to make verb to test ungxxx by matrix of
NB. 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
GQNB=: 32 NB. block size limit
GQNX=: 128 NB. crossover point, GQNX ≥ GQNB
NB. =========================================================
NB. Interface
NB. ---------------------------------------------------------
NB. unglq
NB.
NB. Description:
NB. Generate a matrix with orthonormal rows from output of
NB. gelqf
NB.
NB. Syntax:
NB. Q=. [k] unglq LQf
NB. where
NB. LQf - m×(n+1)-matrix, the output of gelqf
NB. Q - k×n-matrix with orthonormal rows, which is
NB. defined as the first k rows of the product of k
NB. elementary reflectors H(i) of order n:
NB. Q = Π{H(i)',i=k-1:0}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i)' * τ(i) * u(i)
NB. k ≤ n, optional, default is min(m,n)
NB.
NB. Storage layout:
NB. LQf for m=3, n=7: LQf for m=7, n=3:
NB. ( l v0 v0 v0 v0 v0 v0 τ0 ) ( l v0 v0 τ0 )
NB. ( l l v1 v1 v1 v1 v1 τ1 ) ( l l v1 τ1 )
NB. ( l l l v2 v2 v2 v2 τ2 ) ( l l l τ2 )
NB. ( l l l * )
NB. ( l l l * )
NB. ( l l l * )
NB. ( l l l * )
NB. where
NB. l - elements of m×k-matrix L
NB. vi - vector v(i)
NB. τi - scalar value conj(τ(i))
NB. (0,...,0,1,vi) - n-vector u(i)
NB. * - any value, is not used
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. L * Q = L * (Q)
NB. ((((mp unglq)~ trl ) -: (unmlqrn trlpick )) }:"1) LQf
NB. NB. I = Q * (Q^H)
NB. (( 0 idmat (<. , ])/)@(0 _1 + $) -: clean@(unmlqrc unglq)) LQf
NB.
NB. Notes:
NB. - monadic case implements LAPACK's DORGLQ, ZUNGLQ
NB. - straightforward O(k*n^3) code:
NB. Q=. k {. mp/ (idmat n) -"2 |. (+ {:"1 Qf) * (* +)~"0/~"1 }:"1 k {. Qf
unglq=: ($:~ 0 _1 <./@:+ $) :(}:"1@(([ (({."1~ -@c) larfbrcfr ]) idmat @(((] , +) #)~ c) appendr ])&:>/)@([ ( (<;.1~ (0 1:`(GQNB liso4dhs 0 , >:@(>. GQNB >.@%~ -&GQNX))`($~)} #@]))@] , <@( idmat @((- {.) , -~/@]) $)) tru1pick @(( <./ @, 0 _1 + $) {. ]))`((0 idmat , ) <:@c)@.(0 e. (, 0 _1 + $)))
NB. ---------------------------------------------------------
NB. ungql
NB.
NB. Description:
NB. Generate a matrix with orthonormal columns from output
NB. of geqlf
NB.
NB. Syntax:
NB. Q=. [k] ungql QfL
NB. where
NB. QfL - (m+1)×n-matrix, the output of geqlf
NB. Q - m×k-matrix with orthonormal columns, which is
NB. defined as the last k columns of the product of k
NB. elementary reflectors H(i) of order m:
NB. Q = Π{H(i),i=k-1:0}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i) * τ(i) * u(i)'
NB. k ≤ m, optional, default is min(m,n)
NB.
NB. Storage layout:
NB. QfL for m=3, n=7: QfL for m=7, n=3:
NB. ( * * * * τ0 τ1 τ2 ) ( τ0 τ1 τ2 )
NB. ( l l l l l v1 v2 ) ( v0 v1 v2 )
NB. ( l l l l l l v2 ) ( v0 v1 v2 )
NB. ( l l l l l l l ) ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( l v1 v2 )
NB. ( l l v2 )
NB. ( l l l )
NB. where
NB. l - elements of k×n-matrix L
NB. vi - vector v(i)
NB. τi - scalar value τ(i)
NB. (vi,1,0,...,0) - m-vector u(i)
NB. * - any value, is not used
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. Q * L = (Q) * L
NB. ((((mp~ ungql)~ (trl~ -~/@$)) -: (unmqlln (trlpick~ -~/@$))) }. ) QfL
NB. NB. I = (Q^H) * Q
NB. (((0 <. -~/) idmat ([ , <.)/)@(_1 0 + $) -: clean@(unmqllc ungql)) QfL
NB.
NB. Notes:
NB. - monadic case implements LAPACK's DORGQL, ZUNGQL
NB. - straightforward O(k*m^3) code:
NB. Q=. (-k) {."1 mp/ (idmat m) -"2 |. ( {. Qf) * (* +) "0/~"1 |: }. (-k) {."1 Qf
ungql=: ($:~ _1 0 <./@:+ $) :(}. @(([ (({. ~ #) larfblnbc ]) (idmat~ -~/)@(((] ,~ +) c)~ #) stitcht~ ])&:>/)@([ (|.@(<;.2~ '' ; (0 1:`(GQNB liso4dhs _1 , >:@(>. GQNB >.@%~ -&GQNX))`($~)} c@]))@] , <@((idmat~ -~/)@((- {:) ,~ - /@]) $)) (tru1pick~ -~/@$)@((-@(<./)@, _1 0 + $) {."1 ]))`((- idmat ,~) <:@#)@.(0 e. (, _1 0 + $)))
NB. ---------------------------------------------------------
NB. ungqr
NB.
NB. Description:
NB. Generate a matrix with orthonormal columns from output
NB. of geqrf
NB.
NB. Syntax:
NB. Q=. [k] ungqr QfR
NB. where
NB. QfR - (m+1)×n-matrix, the output of geqrf
NB. Q - m×k-matrix with orthonormal columns, which is
NB. defined as the first k columns of the product of
NB. k elementary reflectors H(i) of order m:
NB. Q = Π{H(i),i=0:k-1}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i) * τ(i) * u(i)'
NB. k ≤ m, optional, default is min(m,n)
NB.
NB. Storage layout:
NB. QfR for m=3, n=7: QfR for m=7, n=3:
NB. ( r r r r r r r ) ( r r r )
NB. ( v0 r r r r r r ) ( v0 r r )
NB. ( v0 v1 r r r r r ) ( v0 v1 r )
NB. ( τ0 τ1 τ2 * * * * ) ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( τ0 τ1 τ2 )
NB. where
NB. r - elements of k×n-matrix R
NB. vi - vector v(i)
NB. τi - scalar value τ(i)
NB. (0,...,0,1,vi) - m-vector u(i)
NB. * - any value, is not used
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. Q * R = (Q) * R
NB. ((((mp~ ungqr)~ tru ) -: (unmqrln trupick )) }: ) QfR
NB. NB. I = (Q^H) * Q
NB. (( 0 idmat ([ , <.)/)@(_1 0 + $) -: clean@(unmqrlc ungqr)) QfR
NB.
NB. Notes:
NB. - monadic case implements LAPACK's DORGQR, ZUNGQR
NB. - straightforward O(k*m^3) code:
NB. Q=. k {."1 mp/ (idmat m) -"2 ( {: Qf) * (* +) "0/~"1 |: }: k {."1 Qf
ungqr=: ($:~ _1 0 <./@:+ $) :(}: @(([ (({. ~ -@#) larfblnfc ]) idmat @(((] ,~ +) c)~ #) stitchb ])&:>/)@([ ( (<;.1~ '' ; (0 1:`(GQNB liso4dhs 0 , >:@(>. GQNB >.@%~ -&GQNX))`($~)} c@]))@] , <@( idmat @((- {:) ,~ - /@]) $)) trl1pick @(( <./ @, _1 0 + $) {."1 ]))`((0 idmat ,~) <:@#)@.(0 e. (, _1 0 + $)))
NB. ---------------------------------------------------------
NB. ungrq
NB.
NB. Description:
NB. Generate a matrix with orthonormal rows from output of
NB. gerqf
NB.
NB. Syntax:
NB. Q=. [k] ungrq RQf
NB. where
NB. RQf - m×(n+1)-matrix, the output of gerqf
NB. Q - k×n-matrix with orthonormal rows, which is
NB. defined as the last k rows of the product of k
NB. elementary reflectors H(i) of order n:
NB. Q = Π{H(i)',i=0:k-1}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i)' * τ(i) * u(i)
NB. k ≤ n, optional, default is min(m,n)
NB.
NB. Storage layout:
NB. RQf for m=3, n=7: RQf for m=7, n=3:
NB. ( τ0 v0 v0 v0 v0 r r r ) ( * r r r )
NB. ( τ1 v1 v1 v1 v1 v1 r r ) ( * r r r )
NB. ( τ2 v2 v2 v2 v2 v2 v2 r ) ( * r r r )
NB. ( * r r r )
NB. ( τ0 r r r )
NB. ( τ1 v1 r r )
NB. ( τ2 v2 v2 r )
NB. where
NB. r - elements of m×k-matrix R
NB. vi - vector v(i)
NB. τi - scalar value τ(i)
NB. (vi,1,0,...,0) - n-vector u(i)
NB. * - any value, is not used
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. R * Q = R * (Q)
NB. ((((mp ungrq)~ (tru~ -~/@$)) -: (unmrqrn (trupick~ -~/@$))) }."1) RQf
NB. NB. I = Q * (Q^H)
NB. (((0 >. -~/) idmat (<. , ])/)@(0 _1 + $) -: clean@(unmrqrc ungrq)) RQf
NB.
NB. Notes:
NB. - monadic case implements LAPACK's DORGRQ, ZUNGRQ
NB. - straightforward O(k*n^3) code:
NB. Q=. (-k) {. mp/ (idmat n) -"2 (+ {."1 Qf) * (* +)~"0/~"1 }."1 (-k) {. Qf
ungrq=: ($:~ 0 _1 <./@:+ $) :(}."1@(([ (({."1~ c) larfbrcbr ]) (idmat~ -~/)@(((] , +) #)~ c) appendl~ ])&:>/)@([ (|.@(<;.2~ (0 1:`(GQNB liso4dhs _1 , >:@(>. GQNB >.@%~ -&GQNX))`($~)} #@]))@] , <@((idmat~ -~/)@((- {.) , -~/@]) $)) (trl1pick~ -~/@$)@((-@(<./)@, 0 _1 + $) {. ]))`((-~ idmat , ) <:@c)@.(0 e. (, 0 _1 + $)))
NB. ---------------------------------------------------------
NB. unglz
NB.
NB. Description:
NB. Generate a matrix with orthonormal rows from output of
NB. tzlzf
NB.
NB. Syntax:
NB. Z=. [k] unglz LZf
NB. where
NB. LZf - m×(n+1)-matrix, the output of tzlzf
NB. Z - k×n-matrix with orthonormal rows, which is
NB. defined as the last k rows of the product of k
NB. elementary reflectors H(i) of order n:
NB. Z = Π{H(i)',i=k-1:0}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i)' * τ(i) * u(i)
NB. k ≤ n, optional, default is m
NB. m ≤ n
NB.
NB. Storage layout:
NB. LZf for m=3, n=7:
NB. ( τ0 v0 v0 v0 v0 l 0 0 )
NB. ( τ1 v1 v1 v1 v1 l l 0 )
NB. ( τ2 v2 v2 v2 v2 l l l )
NB. where
NB. l - elements of m×n-matrix L
NB. vi - (n-m)-vector v(i)
NB. τi - scalar value conj(τ(i))
NB. (vi,0,...,0,1,0,...,0) - n-vector u(i)
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. L * Z = L * (Z)
NB. (((mp unglz)~ -: [ unmlzrn ({."1~ (1 - c))~) ({."1~ -@#)) LZf
NB. NB. I = Z * Z^H
NB. (idmat@# -: clean@(mp ct))@unglz LZf
NB.
NB. Notes:
NB. - straightforward O(k*n^3) code:
NB. Z=. (-k) {. mp/ (idmat n) -"2 |. (+ {."1 Zf) * (* +)~"0/~"1 }."1 (-k) {. Zf
unglz=: ($:~ #) :(}."1@(larzbrcfr&:>/)@([ ( (<;.1~ 0 1:`( GQNB liso4dhs 0 , >:@(>. GQNB >.@%~ -&GQNX) )`($~)} #@])@] , <@(idmat~ -~/)@(, c)) ((-@(<. #) , -~/@$@]) {. ]) ,. (((0 >. -~) idmat <. , ]) #))`((-~ idmat , ) <:@c)@.(0 e. (, 0 _1 + $)))
NB. ---------------------------------------------------------
NB. ungzl
NB.
NB. Description:
NB. Generate a matrix with orthonormal columns from output
NB. of tzzlf
NB.
NB. Syntax:
NB. Z=. [k] ungzl ZfL
NB. where
NB. ZfL - (m+1)×n-matrix, the output of tzzlf
NB. Z - m×k-matrix with orthonormal columns, which is
NB. defined as the first k columns of the product of
NB. k elementary reflectors H(i) of order m:
NB. Z = Π{H(i),i=k-1:0}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i) * τ(i) * u(i)'
NB. k ≤ m, optional, default is n
NB. n ≤ m
NB.
NB. Storage layout:
NB. ZfL for m=7, n=3:
NB. ( l 0 0 )
NB. ( l l 0 )
NB. ( l l l )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( τ0 τ1 τ2 )
NB. where
NB. l - elements of m×n-matrix L
NB. vi - (n-m)-vector v(i)
NB. τi - scalar value τ(i)
NB. (0,...,0,1,0,...,0,vi) - m-vector u(i)
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. Z * L = (Z) * L
NB. (((mp~ ungzl)~ -: [ unmzlln ({. ~ (1 -~ #))~) ({. ~ c)) ZfL
NB. NB. I = Z^H * Z
NB. (idmat@c -: clean@(mp~ ct))@ungzl ZfL
NB.
NB. Notes:
NB. - straightforward O(k*m^3) code:
NB. Z=. k {."1 mp/ (idmat m) -"2 |. ( {: Zf) * (* +) "0/~"1 |: }: k {."1 Zf
ungzl=: ($:~ c) :(}: @(larzblnbc&:>/)@([ (|.@(<;.2~ '' ; 0 1:`((GQNB liso4dhs _1 , >:@(>. GQNB >.@%~ -&GQNX))`((1 0 $ 0)"_)@.(0 = ]))`($~)} c@])@] , <@ idmat @(,~ #)) (( (<. c) ,~ -~/@$@]) {. ]) , ~ (( 0 idmat <. ,~ ]) c))`((0 idmat ,~) <:@#)@.(0 e. (, _1 0 + $)))
NB. ---------------------------------------------------------
NB. ungzr
NB.
NB. Description:
NB. Generate a matrix with orthonormal columns from output
NB. of tzzrf
NB.
NB. Syntax:
NB. Z=. [k] ungzr ZfR
NB. where
NB. ZfR - (m+1)×n-matrix, the output of tzzrf
NB. Z - m×k-matrix with orthonormal columns, which is
NB. defined as the last k columns of the product of k
NB. elementary reflectors H(i) of order m:
NB. Z = Π{H(i),i=0:k-1}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i) * τ(i) * u(i)'
NB. k ≤ m, optional, default is n
NB. n ≤ m
NB.
NB. Storage layout:
NB. ZfR for m=7, n=3:
NB. ( τ0 τ1 τ2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( v0 v1 v2 )
NB. ( r r r )
NB. ( 0 r r )
NB. ( 0 0 r )
NB. where
NB. r - elements of m×n-matrix R
NB. vi - (n-m)-vector v(i)
NB. τi - scalar value τ(i)
NB. (vi,0,...,0,1,0,...,0) - m-vector u(i)
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. Z * R = (Z) * R
NB. (((mp~ ungzr)~ -: [ unmzrln ({. ~ (1 - #))~) ({. ~ -@c)) ZfR
NB. NB. I = Z^H * Z
NB. (idmat@c -: clean@(mp~ ct))@ungzr ZfR
NB.
NB. Notes:
NB. - straightforward O(k*m^3) code:
NB. Z=. (-k) {."1 mp/ (idmat m) -"2 ( {. Zf) * (* +) "0/~"1 |: }. (-k) {."1 Zf
ungzr=: ($:~ c) :(}. @(larzblnfc&:>/)@([ ( (<;.1~ '' ; 0 1:`((GQNB liso4dhs 0 , >:@(>. GQNB >.@%~ -&GQNX))`((1 0 $ 0)"_)@.(0 = ]))`($~)} c@])@] , <@(idmat~ -~/)@(,~ #)) ((-@(<. c) ,~ - /@$@]) {. ]) , (((0 <. - ) idmat <. ,~ ]) c))`((- idmat ,~) <:@#)@.(0 e. (, _1 0 + $)))
NB. ---------------------------------------------------------
NB. ungrz
NB.
NB. Description:
NB. Generate a matrix with orthonormal rows from output of
NB. tzrzf
NB.
NB. Syntax:
NB. Z=. [k] ungrz RZf
NB. where
NB. RZf - m×(n+1)-matrix, the output of tzrzf
NB. Z - k×n-matrix with orthonormal rows, which is
NB. defined as the first k rows of the product of k
NB. elementary reflectors of order n:
NB. Z = Π{H(i)',i=0:k-1}
NB. H(i) ≡ H(u(i),τ(i)) := I - u(i)' * τ(i) * u(i)
NB. k ≤ n, optional, default is m
NB. m ≤ n
NB.
NB. Storage layout:
NB. RZf for m=3, n=7:
NB. ( r r r v0 v0 v0 v0 τ0 )
NB. ( 0 r r v1 v1 v1 v1 τ1 )
NB. ( 0 0 r v2 v2 v2 v2 τ2 )
NB. where
NB. r - elements of m×n-matrix R
NB. vi - (n-m)-vector v(i)
NB. τi - scalar value conj(τ(i))
NB. (0,...,0,1,0,...,0,vi) - n-vector u(i)
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. R * Z = R * (Z)
NB. (((mp ungrz)~ -: [ unmrzrn ({."1~ (1 -~ c))~) ({."1~ #)) RZf
NB. NB. I = Z * Z^H
NB. (idmat@# -: clean@(mp ct))@ungrz RZf
NB.
NB. Notes:
NB. - straightforward O(k*n^3) code:
NB. Z=. k {. mp/ (idmat n) -"2 (+ {:"1 Zf) * (* +)~"0/~"1 }:"1 k {. Zf
ungrz=: ($:~ #) :(}:"1@(larzbrcbr&:>/)@([ (|.@(<;.2~ 0 1:`((GQNB liso4dhs _1 , >:@(>. GQNB >.@%~ -&GQNX))`((0 1 $ 0)"_)@.(0 = ]))`($~)} #@])@] , <@ idmat @(, c)) (( (<. #) , - /@$@]) {. ]) ,.~ (( 0 idmat <. , ]) #))`((0 idmat , ) <:@c)@.(0 e. (, 0 _1 + $)))
NB. ---------------------------------------------------------
NB. unghrl
NB.
NB. Description:
NB. Generate an unitary (orthogonal) matrix from output of
NB. gehrdl
NB.
NB. Syntax:
NB. Q=. unghrl HQf
NB. where
NB. HQf - n×(n+1)-matrix with packed H and Qf (see gehrdl)
NB. Q - n×n-matrix, the unitary (orthogonal), which is
NB. defined as the product of (s-1) elementary
NB. reflectors H(i) of order n:
NB. Q = Π{H(i)',i=h+s-2:h}
NB. H(i) = I - v[i]' * τ[i] * v[i]
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. I = Q * Q^H
NB. (idmat@# -: clean@(mp ct))@unghrl HQf
NB.
NB. Notes:
NB. - instead of using h and s arguments, the following
NB. product is really calculating:
NB. Q = Π{H(i)',i=n-1:0} ,
NB. where
NB. H(0:h-1) = H(h+s-1:n-1) = I .
NB. This approach delivers excessive calculations in rare
NB. case ((h>0) OR (h+s<n)).
unghrl=: unglq@(|.!.0)
NB. ---------------------------------------------------------
NB. unghru
NB.
NB. Description:
NB. Generate an unitary (orthogonal) matrix from output of
NB. gehrdu
NB.
NB. Syntax:
NB. Q=. unghru HQf
NB. where
NB. HQf - (n+1)×n-matrix with packed H and Qf (see gehrdu)
NB. Q - n×n-matrix, the unitary (orthogonal), which is
NB. defined as the product of (s-1) elementary
NB. reflectors H(i) of order n:
NB. Q = Π{H(i),i=h:h+s-2}
NB. H(i) = I - v[i] * τ[i] * v[i]'
NB.
NB. Assertions (with appropriate comparison tolerance):
NB. NB. I = Q^H * Q
NB. (idmat@# -: clean@(mp~ ct))@unghru HQf
NB.
NB. Notes:
NB. - models LAPACK's DORGHR, ZUNGHR
NB. - instead of using h and s arguments, the following
NB. product is really calculating:
NB. Q = Π{H(i),i=0:n-1} ,
NB. where
NB. H(0:h-1) = H(h+s-1:n-1) = I .
NB. This approach delivers excessive calculations in rare
NB. case ((h>0) OR (h+s<n)).
unghru=: ungqr@(0 _1&(|.!.0))
NB. =========================================================
NB. Test suite
NB. ---------------------------------------------------------
NB. testungq
NB.
NB. Description:
NB. Test:
NB. - DORGxx ZUNGxx (math/lapack2 addon)
NB. - ungxx (math/mt addon)
NB. by general matrix
NB.
NB. Syntax:
NB. log=. testungq A
NB. where
NB. A - m×n-matrix
NB. log - 6-vector of boxes, test log
testungq=: 3 : 0
_1 cocreate < 'mttmp'
load_mttmp_ 'math/mt/external/lapack2/dorglq'
load_mttmp_ 'math/mt/external/lapack2/dorgql'
load_mttmp_ 'math/mt/external/lapack2/dorgqr'
load_mttmp_ 'math/mt/external/lapack2/dorgrq'
load_mttmp_ 'math/mt/external/lapack2/zunglq'
load_mttmp_ 'math/mt/external/lapack2/zungql'
load_mttmp_ 'math/mt/external/lapack2/zungqr'
load_mttmp_ 'math/mt/external/lapack2/zungrq'
rcond=. nan`geconi@.(=/@$) y NB. meaninigful for square matrices only
ks=. /:~ ~. (, *@(>./)) 0 , (,~ <.@-:) <./ 'm n'=. $ y NB. 0,1,⌊min(m,n)/2⌋,⌊min(m,n)⌋
normw=. norm1 Awide=. |:^:(>/@$) y
normt=. norm1 Atall=. |:^:(</@$) y
LQf=. gelqf Awide
QfL=. geqlf Atall
QfR=. geqrf Atall
RQf=. gerqf Awide
argslq=. { (< Awide) ; (< normw) ; (< LQf) ; < <"0 ks
argsql=. { (< Atall) ; (< normt) ; (< QfL) ; < <"0 ks
argsqr=. { (< Atall) ; (< normt) ; (< QfR) ; < <"0 ks
argsrq=. { (< Awide) ; (< normw) ; (< RQf) ; < <"0 ks
NB. LAPACK, real datatype
log=. lcat ('dorglq_mttmp_' tmonad ( ( 3&{:: (}:"1@] ; ({. {:"1)) 2&{:: )`]`(rcond"_)`nan`lqt02))@>"0 argslq
log=. log lcat ('dorgql_mttmp_' tmonad ( (-@(3&{::) (}. @] ; ({. {. )) 2&{:: )`]`(rcond"_)`nan`qlt02))@>"0 argsql
log=. log lcat ('dorgqr_mttmp_' tmonad ( ( 3&{:: (}: @] ; ({. {: )) 2&{:: )`]`(rcond"_)`nan`qrt02))@>"0 argsqr
log=. log lcat ('dorgrq_mttmp_' tmonad ( (-@(3&{::) (}."1@] ; ({. {."1)) 2&{:: )`]`(rcond"_)`nan`rqt02))@>"0 argsrq
NB. LAPACK, complex datatype
log=. log lcat ('zunglq_mttmp_' tmonad ( ( 3&{:: (}:"1@] ; ({. {:"1)) 2&{:: )`]`(rcond"_)`nan`lqt02))@>"0 argslq
log=. log lcat ('zungql_mttmp_' tmonad ( (-@(3&{::) (}. @] ; ({. {. )) 2&{:: )`]`(rcond"_)`nan`qlt02))@>"0 argsql
log=. log lcat ('zungqr_mttmp_' tmonad ( ( 3&{:: (}: @] ; ({. {: )) 2&{:: )`]`(rcond"_)`nan`qrt02))@>"0 argsqr
log=. log lcat ('zungrq_mttmp_' tmonad ( (-@(3&{::) (}."1@] ; ({. {."1)) 2&{:: )`]`(rcond"_)`nan`rqt02))@>"0 argsrq
NB. mt, any datatype
log=. log lcat ('unglq' tmonad ( ( 2&{:: )`]`(rcond"_)`nan`lqt02)) Awide ; normw ; LQf ; m
log=. log lcat ('ungql' tmonad ( ( 2&{:: )`]`(rcond"_)`nan`qlt02)) Atall ; normt ; QfL ; n
log=. log lcat ('ungqr' tmonad ( ( 2&{:: )`]`(rcond"_)`nan`qrt02)) Atall ; normt ; QfR ; n
log=. log lcat ('ungrq' tmonad ( ( 2&{:: )`]`(rcond"_)`nan`rqt02)) Awide ; normw ; RQf ; m
log=. log lcat ('unglq' tdyad ((#@(0&{::))`( 3&{:: {. 2&{:: )`]`(rcond"_)`nan`lqt02))@>"0 argslq
log=. log lcat ('ungql' tdyad ((c@(0&{::))`(-@(3&{::) {."1 (2&{::))`]`(rcond"_)`nan`qlt02))@>"0 argsql
log=. log lcat ('ungqr' tdyad ((c@(0&{::))`( 3&{:: {."1 (2&{::))`]`(rcond"_)`nan`qrt02))@>"0 argsqr
log=. log lcat ('ungrq' tdyad ((#@(0&{::))`(-@(3&{::) {. 2&{:: )`]`(rcond"_)`nan`rqt02))@>"0 argsrq
coerase < 'mttmp'
log
)
NB. ---------------------------------------------------------
NB. testungz
NB.
NB. Description:
NB. Test ungxx by trapezoidal matrix
NB.
NB. Syntax:
NB. log=. testungz A
NB. where
NB. A - m×n-matrix
NB. log - 6-vector of boxes, test log
testungz=: 3 : 0
rcond=. nan`geconi@.(=/@$) y NB. meaninigful for square matrices only
ks=. /:~ ~. (, *@(>./)) 0 , (,~ <.@-:) <./ 'm n'=. $ y NB. 0,1,⌊min(m,n)/2⌋,⌊min(m,n)⌋
normw=. norm1 Awide=. |:^:(>/@$) y
normt=. norm1 Atall=. |:^:(</@$) y
LZf=. tzlzf Awide
ZfL=. tzzlf Atall
ZfR=. tzzrf Atall
RZf=. tzrzf Awide
log=. ('unglz' tmonad ( (2&{::)`]`(rcond"_)`nan`lzt02)) Awide ; normw ; LZf ; m
log=. log lcat ('ungzl' tmonad ( (2&{::)`]`(rcond"_)`nan`zlt02)) Atall ; normt ; ZfL ; n
log=. log lcat ('ungzr' tmonad ( (2&{::)`]`(rcond"_)`nan`zrt02)) Atall ; normt ; ZfR ; n
log=. log lcat ('ungrz' tmonad ( (2&{::)`]`(rcond"_)`nan`rzt02)) Awide ; normw ; RZf ; m
log=. log lcat ('unglz' tdyad ((3&{::)`(2&{::)`]`(rcond"_)`nan`lzt02))@>"0 { (< Awide) ; (< normw) ; (< LZf) ; < <"0 ks
log=. log lcat ('ungzl' tdyad ((3&{::)`(2&{::)`]`(rcond"_)`nan`zlt02))@>"0 { (< Atall) ; (< normt) ; (< ZfL) ; < <"0 ks
log=. log lcat ('ungzr' tdyad ((3&{::)`(2&{::)`]`(rcond"_)`nan`zrt02))@>"0 { (< Atall) ; (< normt) ; (< ZfR) ; < <"0 ks
log=. log lcat ('ungrz' tdyad ((3&{::)`(2&{::)`]`(rcond"_)`nan`rzt02))@>"0 { (< Awide) ; (< normw) ; (< RZf) ; < <"0 ks
)
NB. ---------------------------------------------------------
NB. testunghr
NB.
NB. Description:
NB. Test:
NB. - DORGHR ZUNGHR (math/lapack2 addon)
NB. - unghrx (math/mt addon)
NB. by square matrix
NB.
NB. Syntax:
NB. log=. testunghr A
NB. where
NB. A - n×n-matrix
NB. log - 6-vector of boxes, test log
testunghr=: 3 : 0
_1 cocreate < 'mttmp'
load_mttmp_ 'math/mt/external/lapack2/dorghr'
load_mttmp_ 'math/mt/external/lapack2/zunghr'
'rcondl rcondu'=. (geconi , gecon1) y
'norml normu'=. (normi , norm1) y
HlQf=. (gehrdl~ 0 , #) y
HuQf=. (gehrdu~ 0 , #) y
hst01l=: }:@[ (normi hst01) ((ct@[ mp mp~) 1 trlpick }:"1@(2&{::))~
hst01u=: }:@[ (norm1 hst01) ((ct@[ mp~ mp ) _1 trupick }: @(2&{::))~
unt01l=: (normi unt01 (mp ct))@]
unt01u=: (norm1 unt01 (mp~ ct))@]
log=. ('dorghr_mttmp_' tmonad (((1 ; c ; }: ; }:@{:)@(2&{::))`]`(rcondu"_)`nan`(hst01u >. unt01u))) y ; normu ; HuQf
log=. log lcat ('zunghr_mttmp_' tmonad (((1 ; c ; }: ; }:@{:)@(2&{::))`]`(rcondu"_)`nan`(hst01u >. unt01u))) y ; normu ; HuQf
log=. log lcat ('unghrl' tmonad (( 2&{:: )`]`(rcondl"_)`nan`(hst01l >. unt01l))) y ; norml ; HlQf
log=. log lcat ('unghru' tmonad (( 2&{:: )`]`(rcondu"_)`nan`(hst01u >. unt01u))) y ; normu ; HuQf
coerase < 'mttmp'
erase 'hst01l hst01u unt01l unt01u'
log
)
NB. ---------------------------------------------------------
NB. testgq
NB.
NB. Description:
NB. Adv. to make verb to test ungxxx by matrix of generator
NB. and shape given
NB.
NB. Syntax:
NB. log=. (mkmat testgq) (m,n)
NB. where
NB. mkmat - monad to generate a matrix; is called as:
NB. mat=. mkmat (m,n)
NB. (m,n) - 2-vector of integers, the shape of matrix mat
NB. log - 6-vector of boxes, test log
NB.
NB. Application:
NB. - test by random rectangular real matrix with elements
NB. distributed uniformly with support (0,1):
NB. log=. ?@$&0 testgq_mt_ 200 150
NB. - test by random square real matrix with elements with
NB. limited value's amplitude:
NB. log=. _1 1 0 4 _6 4&gemat_mt_ testgq_mt_ 200 200
NB. - test by random rectangular complex matrix:
NB. log=. (gemat_mt_ j. gemat_mt_) testgq_mt_ 150 200
testgq=: 1 : 'lcat_mt_@(testungq_mt_`testungz_mt_`(nolog_mt_`testunghr_mt_@.(=/@$))`:0)@u'