File: algebraic_types

package info (click to toggle)
polymake 4.14-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 35,888 kB
  • sloc: cpp: 168,933; perl: 43,407; javascript: 31,575; ansic: 3,007; java: 2,654; python: 632; sh: 268; xml: 117; makefile: 61
file content (1497 lines) | stat: -rw-r--r-- 50,084 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
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
#  Copyright (c) 1997-2024
#  Ewgenij Gawrilow, Michael Joswig, and the polymake team
#  Technische Universität Berlin, Germany
#  https://polymake.org
#
#  This program is free software; you can redistribute it and/or modify it
#  under the terms of the GNU General Public License as published by the
#  Free Software Foundation; either version 2, or (at your option) any
#  later version: http://www.gnu.org/licenses/gpl.txt.
#
#  This program is distributed in the hope that it will be useful,
#  but WITHOUT ANY WARRANTY; without even the implied warranty of
#  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
#  GNU General Public License for more details.
#-------------------------------------------------------------------------------

# @topic category property_types/Algebraic Types
# This category contains all "algebraic" types, such as matrices, vectors, polynomials, rings, ...

# @category Algebraic Types
# A type for vectors with entries of type //Element//.
#
# You can perform algebraic operations such as addition or scalar multiplication.
#
# @example
# You can create a new Vector by entering its elements, e.g.:
# > $v = new Vector<Int>(1,2,3);
# or
# >$v = new Vector<Int>([1,2,3]);
# @tparam Element
declare property_type Vector<Element=Rational> : c++ (include => ["polymake/Vector.h"]) {

   method construct(Int) : c++;

   operator neg : c++;

   operator @arith (*:wary, *:num) : c++;

   operator @eq (*:wary, *) : c++;

   operator | |= (*&& const, *&& const) : c++;

   type_method init {
      my ($proto)=@_;
      if ($proto->params->[0] == typeof Float) {
         $proto->equal=\&compare_float_sequences;
      }
   }

   # Returns the length of the given [[Vector]].
   # @return Int
   # @example
   # > $v = new Vector<Int>(1,2,3,4);
   # > print $v->dim;
   # | 4
   user_method dim() : c++;

   # Returns a [[Vector]] containing all entries whose index is in a [[Set]] //s//.
   # @param Set s
   # @return Vector
   # @example
   # > $v = new Vector<Int>(1,2,3,4);
   # > $s = new Set(0,2);
   # > print $v->slice($s);
   # | 1 3
   user_method slice(:wary&&, *&& const) : c++ : returns(lvalue);

   # Divides every entry of [[Vector]] //v// by [[Int]] //a//. If type of //v// is [[Int]], then only the integer part of the each entry is taken into account.
   # @param Int a
   # @return Vector
   # @example
   # > $v = new Vector<Rational>(1/2,2);
   # > print $v->div_exact(2);
   # | 1/4 1
   # > $v = new Vector<Int>(1,2,3,4);
   # > print $v->div_exact(2);
   # | 0 1 1 2
   user_method div_exact(&, *) : c++;
}


# @category Arithmetic
# Compute the __greatest common divisor__ of the elements of the given vector.
# @param Vector<__Scalar__> v
# @return __Scalar__
# @example
# > $v = new Vector<Int>(3,6,9);
# > print gcd($v);
# | 3
user_function gcd(Vector) : c++ (include => "polymake/linalg.h");

# @category Arithmetic
# Compute the __least common multiple__ of the elements of the given vector.
# @param Vector<__Scalar__> v
# @return __Scalar__
# @example
# > $v = new Vector<Integer>(1,3,6,9);
# > print lcm($v);
# | 18
user_function lcm(Vector) : c++ (include => "polymake/linalg.h");

##################################################################################

# @category Algebraic Types
# A class to hold and process Plücker coordinates of a subspace.
# @tparam Scalar default: [[Rational]]
declare property_type Plucker<Scalar=Rational> : c++ (include => "polymake/Plucker.h") {

    method construct(Vector) : c++;

    method construct(Int, Int, Vector) : c++;

    method construct(Int, Int) : c++;

    operator + * : c++;

    method point() : c++;

    user_method permuted(*) : c++;

    user_method coordinates() : c++;
}

##################################################################################

# @category Algebraic Types
# A type for matrices with entries of type //Element//.
# You can create a new Matrix by entering its entries as an array of row vectors:
# @tparam Element default: [[Rational]]
# @tparam Sym not implemented - for internal type compatibility only
# @example
# > $M = new Matrix([1,3],[6,9]);
# > print $M;
# | 1 3
# | 6 9
declare property_type Matrix<Element=Rational, Sym=NonSymmetric> : c++ (name => "Matrix<%1>", include => "polymake/Matrix.h") {

   method construct(Int,Int) : c++;

   method construct(Vector+) {
      my ($proto, $vectors)=@_;
      my $M=$proto->construct->(scalar(@$vectors), $vectors->[0]->dim);
      my $i=0;
      $M->[$i++]=$_ for @$vectors;
      $M;
   }

   operator neg : c++;

   operator @arith (*:wary, *:num) : c++;

   operator @eq (*:wary, *) : c++;

   operator / | /= |= (*:wary&& const, Matrix&& const) : c++;

   operator / | /= |= (*:wary&& const, Vector&& const) : c++;

   operator / | (Vector&& const, *:wary&& const) : c++;

   type_method init {
      my ($proto)=@_;
      if ($proto->params->[0] == typeof Float) {
         $proto->equal=sub : method {
            (undef, my $m1, my $m2)=@_;
            my $vec_proto=typeof Vector<Float>;
            if ($m1->rows==$m2->rows) {
               my $i=0;
               foreach my $v1 (@$m1) {
                  $vec_proto->equal->($v1, $m2->row($i++)) or return 0;
               }
               1
            }
         };
      }
   }

   # Resizes the dimension of [[Matrix]] //M// to //r// x //c//; missing entries are filled with zero and for all i > //r//  or j > //c// the (i,j)-th entry of //M// is forgotten.
   # @param Int r new number of rows
   # @param Int c new number of columns
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > $M->resize(3,2);
   # > print $M;
   # | 1 2
   # | 4 5
   # | 0 0
   user_method resize(&, $$) : c++;

   # Change the dimensions setting all elements to 0.
   # @param Int r new number of rows
   # @param Int c new number of columns
   user_method clear(&, $$) : c++;

   # Returns the number of rows.
   # @return Int
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > print $M->rows;
   # | 2
   user_method rows() : c++;

   # Returns the number of columns.
   # @return Int
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > print $M->cols;
   # | 3
   user_method cols() : c++;

   # Returns the //i//-th row.
   # @param Int i
   # @return Vector<Element>
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > print $M->row(1);
   # | 4 5 6
   user_method row(:wary&&, $) : c++ : returns(lvalue);

   # Returns the //i//-th column.
   # @param Int i
   # @return Vector<Element>
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > print $M->col(1);
   # | 2 5
   user_method col(:wary&&, $) : c++ : returns(lvalue);

   # Returns a __minor__ of the matrix containing the rows in //r// and the columns in //c//. You can pass [[all_rows_or_cols|All]] if you want all rows or columns and ~ for the complement of a set.
   # @param Set r the rows
   # @param Set c the columns
   # @return Matrix
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > print $M->minor([0,1],[0,1]);
   # | 1 2
   # | 4 5
   # > print $M->minor(All,~[0]);
   # | 2 3
   # | 5 6
   user_method minor(:wary&&, *&& const, *&& const) : c++ : returns(lvalue);

   # Returns an element of the matrix.
   # The return value is an `lvalue', that is, it can be modified if the matrix object is mutable.
   # @param Int r the row index
   # @param Int c the column index
   # @return Element
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6]);
   # > print $M->elem(1,1);
   # | 5
   user_method elem(:wary&&, $$) : c++(name => '()') : returns(lvalue);

   # Returns the __diagonal__ of the matrix.
   # @param Int i //i//=0: the main diagonal (optional)
   #				//i//>0: the //i//-th diagonal __below__ the main diagonal
   #				//i//<0: the //i//-th diagonal __above__ the main diagonal
   # @return Vector<Element>
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6],[7,8,9]);
   # > print $M->diagonal;
   # | 1 5 9
   # > print $M->diagonal(1);
   # | 4 8
   user_method diagonal(:wary&& ; $=0) : c++ : returns(lvalue);

   # Returns the __anti-diagonal__ of the matrix.
   # @param Int i //i//=0: the main anti_diagonal (optional)
   #				//i//>0: the //i//-th anti_diagonal __below__ the main anti_diagonal
   #				//i//<0: the //i//-th anti_diagonal __above__ the main anti_diagonal
   # @return Vector<Element>
   # @example
   # > $M = new Matrix([1,2,3],[4,5,6],[7,8,9]);
   # > print $M->anti_diagonal;
   # | 3 5 7
   # > print $M->anti_diagonal(-1);
   # | 2 4
   user_method anti_diagonal(:wary&& ; $=0) : c++ : returns(lvalue);

   # Divides every entry of [[Matrix]] //M// by [[Int]] //a//. If type of //M// is [[Int]], then only the integer part of the each entry is taken into account.
   # @param Int a
   # @return Matrix
   # @example
   # > $M = new Matrix([1,2],[2,1]);
   # > print $M->div_exact(2);
   # | 1/2 1
   # | 1 1/2
   user_method div_exact(&, *) : c++;
}

# @category Linear Algebra
# Projects the row [[Matrix]] //M// of points into the [[null_space|orthogonal complement]] of a subspace given by the rows of the [[Matrix]] //N//.
# The points of //M// will be overwitten.
# @param Matrix M row matrix of the points that will be projected
# @param Matrix N row matrix of the basis elements
# @example
# > $M = new Matrix([1,0,2],[2,0,1],[3,2,3]);
# > $N = new Matrix([1,0,0],[0,1,0]);
# > project_to_orthogonal_complement($M,$N);
# > print $M;
# | 0 0 2
# | 0 0 1
# | 0 0 3
user_function project_to_orthogonal_complement(Matrix&, Matrix) : c++ (include => "polymake/linalg.h");

# find a permutation of rows between two given matrices
# matrices with [[Float]] coefficients are compared with some tolerance
# @param Matrix M1 first matrix
# @param Matrix M2 second matrix
# @param Bool expect_duplicate_rows if true, matrices may contain duplicate rows
# @return Array<Int> row indexes, or undef if there are unmatched rows
function find_matrix_row_permutation(Matrix, Matrix; $=false) : c++ (include => ["polymake/common/find_matrix_row_permutation.h"]);

# @category Algebraic Types
# Allows the use of the keyword "All" for all rows or columns, e.g. when constructing a [[minor]].
declare property_type all_rows_or_cols : c++ (name => 'pm::all_selector', include => ["polymake/Matrix.h"], builtin => enum( All ));


# @category Data Conversion
# Explicit conversion to a different element type.
# @param Vector v
# @tparam Target
# @return Vector<Target>
# @example
# > $v = new Vector<Rational>(1/2,2/3,3/4);
# > $vf = convert_to<Float>($v);
# > print $vf;
# | 0.5 0.6666666667 0.75
user_function convert_to<Element>(Vector) : c++ {
   if ($_[0]->type->params->[0]==typeof Element) {
      return $_[0];
   }
}

# @category Data Conversion
# Explicit conversion to a different element type.
# @param Matrix m
# @tparam Target
# @return Matrix<Target>
# @example
# > $M = new Matrix<Rational>([1/2,2],[3,2/3]);
# > $Mf = convert_to<Float>($M);
# > print $Mf;
# | 0.5 2
# | 3 0.6666666667
user_function convert_to<Element>(Matrix) : c++ {
   if ($_[0]->type->params->[0]==typeof Element) {
      return $_[0];
   }
}

# @category Data Conversion
# Create a [[Matrix]] by repeating the given [[Vector]] //v// as rows //i// many times.
# @param Vector v
# @param Int i
# @return Matrix
# @example
# > $v = new Vector(23,42,666);
# > $M = repeat_row($v,3);
# > print $M;
# | 23 42 666
# | 23 42 666
# | 23 42 666
user_function repeat_row(Vector:anchor,$) : c++ (include => "polymake/Matrix.h");

# @category Data Conversion
# Create a [[Matrix]] by repeating the given [[Vector]] as cols.
# @param Vector v
# @param Int i
# @return Matrix
# @example
# > $v = new Vector(23,42,666);
# > $M = repeat_col($v,3);
# > print $M;
# | 23 23 23
# | 42 42 42
# | 666 666 666
user_function repeat_col(Vector:anchor,$) : c++ (include => "polymake/Matrix.h");

# @category Data Conversion
# Convert a [[Vector]] to a [[Matrix]] with a single row.
# @param Vector v
# @return Matrix
# @example This converts a vector into a row and prints it and its type:
# > $v = new Vector([1,2,3,4]);
# > $V = vector2row($v);
# > print $V;
# | 1 2 3 4
# > print $V->type->full_name;
# | Matrix<Rational, NonSymmetric>
user_function vector2row(Vector:anchor) : c++ (include => "polymake/Matrix.h");


# @category Data Conversion
# Convert a [[Vector]] to a [[Matrix]] with a single column.
# @param Vector v
# @return Matrix
# @example This converts a vector into a column and prints it and its type:
# > $v = new Vector([1,2,3,4]);
# > $V = vector2col($v);
# > print $V;
# | 1
# | 2
# | 3
# | 4
# > print $V->type->full_name;
# | Matrix<Rational, NonSymmetric>
user_function vector2col(Vector:anchor) : c++ (include => "polymake/Matrix.h");


# @category Linear Algebra
# Produces a [[SparseMatrix]] with given [[Vector]] //d// as diagonal.
# @param Vector d the diagonal entries
# @return SparseMatrix
# @example
# > $v = new Vector(1,2,3,4);
# > $D = diag($v);
# > print $D;
# | (4) (0 1)
# | (4) (1 2)
# | (4) (2 3)
# | (4) (3 4)
# To print a more human-readable representation, use the dense() function:
# > print dense($D);
# | 1 0 0 0
# | 0 2 0 0
# | 0 0 3 0
# | 0 0 0 4
user_function diag(Vector:anchor) : c++ (include => "polymake/SparseMatrix.h");


# @category Linear Algebra
# Returns a __block diagonal matrix__ with blocks //m1// and //m2//.
# @param Matrix m1
# @param Matrix m2
# @return SparseMatrix
# @example
# > $m1 = new Matrix([1,2],[3,4]);
# > $m2 = new Matrix([1,0,2],[3,4,0]);
# > $D = diag($m1,$m2);
# > print $D;
# | (5) (0 1) (1 2)
# | (5) (0 3) (1 4)
# | 0 0 1 0 2
# | 0 0 3 4 0
# To print a more human-readable representation, use the dense() function:
# > print dense($D);
# | 1 2 0 0 0
# | 3 4 0 0 0
# | 0 0 1 0 2
# | 0 0 3 4 0
user_function diag(*:anchor,*:anchor) : c++ (include => "polymake/SparseMatrix.h");


# @category Linear Algebra
# Produces a [[SparseMatrix]] with the given [[Vector]] //d// as anti-diagonal.
# @param Vector d the anti-diagonal entries
# @return SparseMatrix
# @example
# > $M = anti_diag(new Vector([0,1,2]));
# > print $M;
# | (3) (2 2)
# | (3) (1 1)
# | (3)
# To print a more human-readable representation, use the dense() function:
# > print dense($M);
# | 0 0 2
# | 0 1 0
# | 0 0 0
user_function anti_diag(Vector:anchor) : c++ (include => "polymake/SparseMatrix.h");


# @category Linear Algebra
# Returns a __block anti-diagonal matrix__ with blocks //M1// and //M2//.
# @param Matrix M1
# @param Matrix M2
# @return SparseMatrix
# @example
# > $M = anti_diag(unit_matrix(2),unit_matrix(3));
# > print $M;
# | (5) (2 1)
# | (5) (3 1)
# | (5) (4 1)
# | (5) (0 1)
# | (5) (1 1)
# To print a more human-readable representation, use the dense() function:
# > print dense($M);
# | 0 0 1 0 0
# | 0 0 0 1 0
# | 0 0 0 0 1
# | 1 0 0 0 0
# | 0 1 0 0 0
user_function anti_diag(*:anchor,*:anchor) : c++ (include => "polymake/SparseMatrix.h");


# @category Data Conversion
# Returns a [[Container]] with the rows of the [[Matrix]] //A//.
# @param Matrix A
# @return Container<Vector>
# @example
# The following saves the rows of the vertex matrix of a square in
# the variable $w and then prints its contents using a foreach loop and concatenating
# each entry with the string "  ".
# > $w = rows(polytope::cube(2)->VERTICES);
# > foreach( @$w ){
# > print @{$_}, "  ";
# > }
# | 1-1-1  11-1  1-11  111
user_function rows(Matrix:anchor) : c++ : returns(Container);

# @category Data Conversion
# Returns a [[Container]] with the columns of the [[Matrix]] //A//.
# @param Matrix A
# @return Container<Vector>
# @example
# The following saves the columns of the vertex matrix of a square in
# the variable $w and then prints its contents using a foreach loop and concatenating
# each entry with the string "  ".
# > $w = cols(polytope::cube(2)->VERTICES);
# > foreach( @$w ){
# > print @{$_}, "  ";
# > }
# | 1111  -11-11  -1-111
user_function cols(Matrix:anchor) : c++ : returns(Container);

# @category Data Conversion
# Concatenates the rows of the [[Matrix]] //M//. If //M// is a [[SparseMatrix]], then the resulting vector is sparse as well.
# @param Matrix M
# @return Vector
# @example
# Make a vector out of the rows of the vertex matrix of a cube:
# > $v = concat_rows(polytope::cube(2)->VERTICES);
# > print $v;
# | 1 -1 -1 1 1 -1 1 -1 1 1 1 1
# @example For a sparse matrix, the resulting vector is sparse, too.
# > $vs = concat_rows(unit_matrix(3));
# > print $vs;
# | (9) (0 1) (4 1) (8 1)
user_function concat_rows(Matrix&&) : c++;


# @category Linear Algebra
# Computes the __transpose__ //M//<sup>T</sup> of a [[Matrix]] //M//, i.e., (M<sup>T</sup>)<sub>ij</sub> = M<sub>ji</sub>.
# @param Matrix M
# @return Matrix
# @example
# > $M = new Matrix([1,2,23],[23,22,21]);
# > $Mt = transpose($M);
# > print $Mt;
# | 1 23
# | 2 22
# | 23 21
user_function transpose(Matrix:anchor) : c++ (name => 'T');


# @category Linear Algebra
# Computes the __determinant__ of a [[Matrix]] using Gaussian elimination.
# If Scalar is not of field type, but element of a Euclidean ring R,
# type upgrade to element of the quotient field is performed.
# The result is recast as a Scalar, which is possible without roundoff
# since the so-computed determinant is an element of the (embedded) ring R.
# @param Matrix<__Scalar__> A
# @return __Scalar__ det(A)
# @example
# > print det(unit_matrix(3));
# | 1
# @example
# > $p = new UniPolynomial<Rational,Int>("x2+3x");
# > $M = new Matrix<UniPolynomial<Rational,Int>>([[$p, $p+1],[$p+1,$p]]);
# > print det($M);
# | -2*x^2 -6*x - 1
user_function det(Matrix:wary) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Computes the __trace__ of a [[Matrix]].
# @param Matrix<__Scalar__> M
# @return __Scalar__ trace(M)
# @example
# > $M = new Matrix([1,2,3],[23,24,25],[0,0,1]);
# > print trace($M);
# | 26
user_function trace(Matrix:wary) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Computes the __rank__ of a [[Matrix]].
# @param Matrix A
# @return Int rank($M)
# @example
# > $M = new Matrix([[1,2,3],[2,3,2],[3,4,2]]);
# > print rank($M);
# | 3
user_function rank(Matrix) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Computes the __inverse__ //M//<sup>-1</sup> of an invertible [[Matrix]] //M// using Gauss elimination.
# @param Matrix M
# @return Matrix
# @example We save the inverse of a small matrix M in the variable $iM:
# > $M = new Matrix([1,2],[3,4]);
# > $iM = inv($M);
# To print the result, type this:
# > print $iM;
# | -2 1
# | 3/2 -1/2
# As we can see, that is in fact the inverse of M.
# > print $M * $iM;
# | 1 0
# | 0 1
user_function inv(Matrix:wary) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Reduce a vector with a given [[Matrix]] using Gauss elimination.
# @param Matrix A
# @param Vector b
# @return Vector
user_function reduce(Matrix:wary, Vector:wary) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Normalize a [[Matrix]] by dividing each row by its length (l2-norm).
# @param Matrix<Float> M
# @return Matrix<Float>
# @example
# > $M = new Matrix<Float>([1.5,2],[2.5,2.5]);
# > print normalized($M);
# | 0.6 0.8
# | 0.7071067812 0.7071067812
user_function normalized(Matrix) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Calculate the average over the rows of a [[Matrix]], i.e.,  where i-th entry corrsponds to the mean of i-th column vector.
# @param Matrix<__Scalar__> M
# @return Vector<__Scalar__>
# @example
# > $M = new Matrix([2,-3,-2],[1,-1,1],[-1,2,2]);
# > print barycenter($M);
# | 2/3 -2/3 1/3
user_function barycenter(Matrix) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Return the sum of the squared entries of a [[Vector]] //v//.
# @param Vector<__Scalar__> v
# @return __Scalar__
# @example
# > $v = new Vector([1,-1,2]);
# > print(sqr($v));
# | 6
user_function sqr(Vector) : c++ (include => "polymake/Vector.h");

# @category Linear Algebra
# Computes subsets of the rows and columns of the [[Matrix]] //A// that form a basis for the linear space spanned by //A//.
# @param Matrix A
# @return Pair<Set<Int>, Set<Int>> The first set corresponds to the rows, the second to the columns.
# @example Here we have a nice matrix:
# > $M = new Matrix([[1,0,0,0],[2,0,0,0],[0,1,0,0],[0,0,1,0]]);
# Let's print bases for the row and column space:
# > ($row,$col) = basis($M);
# > print $M->minor($row,All);
# | 1 0 0 0
# | 0 1 0 0
# | 0 0 1 0
# > print $M->minor(All,$col);
# | 1 0 0
# | 2 0 0
# | 0 1 0
# | 0 0 1
user_function basis(Matrix) : returns(@) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Does the same as [[basis]] ignoring the first column of the matrix.
# @param Matrix A
# @return Pair<Set<Int>, Set<Int>> The first set corresponds to the rows, the second to the columns.
# @example Let us illustrate this using the same matrix in the example for [[basis]].
# > $M = new Matrix ([[1,0,0,0],[2,0,0,0],[0,1,0,0],[0,0,1,0]]);
# > ($row,$col) = basis_affine($M);
# > print"$row \n$col";
# | {2 3}
# | {1 2}
# > print $M->minor($row,All);
# | 0 1 0 0
# | 0 0 1 0
# > print $M->minor(All,$col);
# | 0 0
# | 0 0
# | 1 0
# | 0 1
user_function basis_affine(Matrix) : returns(@) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Computes a subset of the rows of the [[Matrix]] //M// that form a basis for the linear space spanned by the rows of //M//.
# @param Matrix M
# @return Set<Int>
# @example Here we have a nice matrix:
# > $M = new Matrix([[1,0,0,0],[2,0,0,0],[0,1,0,0],[0,0,1,0]]);
# > print(basis_rows($M));
# | {0 2 3}
# Let's print a basis of its row space:
# > print $M->minor(basis_rows($M),All);
# | 1 0 0 0
# | 0 1 0 0
# | 0 0 1 0
user_function basis_rows(Matrix) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# Computes a subset of the columns of the [[Matrix]] //M// that form a basis for the linear space spanned by //M//.
# @param Matrix M
# @return Set<Int>
# @example Here we have a nice matrix:
# > $M = new Matrix([[1,0,0,0],[2,0,0,0],[0,1,0,0],[0,0,1,0]]);
# > print(basis_cols($M));
# | {0 1 2}
# Let's print a basis of its column space:
# > print $M->minor(All,basis_cols($M));
# | 1 0 0
# | 2 0 0
# | 0 1 0
# | 0 0 1
user_function basis_cols(Matrix) : c++ (include => "polymake/linalg.h");

function numerators(Vector<Rational>) : c++ (include => "polymake/linalg.h");

function numerators(Matrix<Rational>) : c++ (include => "polymake/linalg.h");

function denominators(Vector<Rational>) : c++ (include => "polymake/linalg.h");

function denominators(Matrix<Rational>) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Remove all zero rows from a [[Matrix]].
# @param Matrix M
# @return Matrix
# @example
# > $M = new Matrix([1,2,0],[0,0,0],[1,1,0]);
# > print(remove_zero_rows($M));
# | 1 2 0
# | 1 1 0
user_function remove_zero_rows(Matrix) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Creates a unit matrix of given dimension.
# @tparam Element default: [[Rational]]
# @param Int d dimension of the matrix
# @return SparseMatrix<Element>
# @example The following stores the 3-dimensional unit matrix (ones on the diagonal and zeros otherwise) and prints it:
# > $M = unit_matrix(3);
# > print $M;
# | (3) (0 1)
# | (3) (1 1)
# | (3) (2 1)
# > print $M->type->full_name;
# | SparseMatrix<Rational, Symmetric>
# @example The following stores the 3-dimensional unit matrix (ones on the diagonal and zeros otherwise) of type Int
# in a variable and prints it:
# > $M = unit_matrix<Int>(3);
# > print $M->type->full_name;
# | SparseMatrix<Int, Symmetric>
user_function unit_matrix<Element=Rational>($) : c++ (include => "polymake/linalg.h") {
   if ($_[0] < 0) {
      croak( "unit_matrix - invalid dimension" );
   }
}

# @category Linear Algebra
# Creates a zero matrix of given dimensions.
# @tparam Element default: [[Rational]]
# @param Int i number of rows
# @param Int j number of columns
# @return Matrix<Element>
# @example The following stores a 2x3 matrix with 0 as entries (from type Rational) in a variable and prints it:
# > $M = zero_matrix(2,3);
# > print $M;
# | 0 0 0
# | 0 0 0
# > print $M->type->full_name;
# | Matrix<Rational, NonSymmetric>
# @example The following stores a 2x3 matrix with 0 as entries from type Int in a variable and prints its type:
# > $M = zero_matrix<Int>(2,3);
# > print $M->type->full_name;
# | Matrix<Int, NonSymmetric>
user_function zero_matrix<Element=Rational>($$) : c++ (include => "polymake/linalg.h") {
   if ($_[0] < 0 || $_[1] < 0) {
      croak( "zero_matrix - invalid dimension" );
   }
}

# @category Linear Algebra
# Creates a [[SparseVector]] of given length //d// with an entry 1 at position //pos// and zeroes elsewhere.
# @tparam Element default: [[Rational]]
# @param Int d the dimension of the vector
# @param Int pos the position of the 1
# @return SparseVector<Element>
# @example The following stores a vector of dimension 5 with a single 1 (as a Rational) at position 2:
# > $v = unit_vector(5,2);
# > print $v;
# | (5) (2 1)
# > print(dense($v));
# | 0 0 1 0 0
# @example The following stores a vector of dimension 5 with a single 1 (as a Int) at position 2:
# > $v = unit_vector<Int>(5,2);
# > print $v->type->full_name;
# | SparseVector<Int>
# @example The following concatenates a unit vector of dimension 3 with a 1 at position 2 and a
# unit vector of dimension 2 with a 1 at position 1:
# > $v = unit_vector(3,2) | unit_vector(2,1);
# > print $v;
# | (5) (2 1) (4 1)
user_function unit_vector<Element=Rational>($$) : c++ (include => "polymake/linalg.h") {
   if ($_[1] < 0 || $_[1] >= $_[0]) {
      croak( "unit_vector - invalid dimension or index out of range" );
   }
}

# @category Linear Algebra
# Creates a [[Vector]] with all elements equal to zero.
# @param Int d vector dimension.  If omitted, a vector of dimension 0 is created,
#               which can adjust itself when involved in a block matrix operation.
# @tparam Element default: [[Rational]]
# @return Vector<Element>
# @example The following stores a vector of dimension 5 with 0 as entries (from type Rational) in a variable and prints it:
# > $v = zero_vector(5);
# > print $v;
# | 0 0 0 0 0
# @example The following stores a vector of dimension 5 with 0 as entries from type Int in a variable and prints its type:
# > $v = zero_vector<Int>(5);
# > print $v->type->full_name;
# | Vector<Int>
# @example The following concatenates a vector of dimension 2 of ones and a vector of length 2 of zeros:
# > $v = ones_vector(2) | zero_vector(2);
# > print $v;
# | 1 1 0 0
user_function zero_vector<Element=Rational>(;$=0) : c++ (include => "polymake/linalg.h") {
   if ($_[0] < 0) {
      croak( "zero_vector - invalid dimension");
   }
}

# @category Linear Algebra
# Creates a [[Vector]] with all elements equal to 1.
# @param Int d vector dimension.  If omitted, a vector of dimension 0 is created, which can adjust itself when involved in a block matrix operation.
# @tparam Element default: [[Rational]].
# @return Vector<Element>
# @example The following stores and prints an integer vector of dimension 3 with all entries equal to 1:
# > $v = ones_vector<Int>(3);
# > print $v;
# | 1 1 1
user_function ones_vector<Element=Rational>(;$=0) : c++ (include => "polymake/linalg.h") {
   if ($_[0] < 0) {
      croak( "ones_vector - invalid dimension");
   }
}

# @category Linear Algebra
# Creates a [[Matrix]] with \\n\\ rows and \\m\\ columns such that all elements equal to 1.
# @param Int m number of rows
# @param Int n number of columns
# @tparam Element default: [[Rational]].
# @return Matrix<Element>
# @example The following creates an all-ones matrix with Rational coefficients.
# > $M = ones_matrix<Rational>(2,3);
# > print $M;
# | 1 1 1
# | 1 1 1
user_function ones_matrix<Element=Rational>($$) : c++ (include => "polymake/linalg.h") {
   if ($_[0] < 0 || $_[1] < 0) {
      croak( "ones_matrix - invalid dimension");
   }
}

# @category Linear Algebra
# Computes the __null space__ of a [[Matrix]] //M//.
# @param Matrix M
# @return Matrix
# @example
# > $M = new Matrix([1,2,0],[2,0,2]);
# > print null_space($M);
# | -1 1/2 1
user_function null_space(Matrix) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Computes the __null space__ of a [[Vector]] //v//.
# @param Vector v
# @return Matrix
# @example
# > $v = new Vector(1,2,3);
# > print null_space($v);
# | -2 1 0
# | -3 0 1
user_function null_space(Vector) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Computes the __lineality space__ of a [[Matrix]] //M//.
# @param Matrix M
# @return Matrix
# @example
# > $M = new Matrix([1,1,0,0],[1,0,1,0]);
# > print lineality_space($M);
# | 0 0 0 1
user_function lineality_space(Matrix) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Computes the [[Vector]] x that solves the system //A//x = //b// for an invertible [[Matrix]] A.
# @param Matrix A must be invertible
# @param Vector b
# @return Vector
# @example from the Wikipedia:
# > $A = new Matrix([3,2,-1],[2,-2,4],[-1,1/2,-1]);
# > $b = new Vector(1,-2,0);
# > print lin_solve($A,$b);
# | 1 -2 -2
user_function lin_solve(Matrix:wary, Vector:wary) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Computes a [[Matrix]] X that solves the system //AX// = //B//
# @param Matrix A
# @param Matrix B
# @return Matrix
# @example A non-degenerate example:
# > $A = new Matrix([[1,0,0],[1,1,0],[1,0,1],[1,1,1]]);
# > $B = new Matrix([[1,0,0],[1,0,1],[1,1,0],[1,1,1]]);
# > print solve_right($A,$B);
# | 1 0 0
# | 0 0 1
# | 0 1 0
# @example A degenerate example:
# > $A = new Matrix([[1,0,0,0,0],[0,1,0,0,0],[1,0,1,0,0],[0,1,1,0,0]]);
# > $B = new Matrix([[0,1,0,0,0],[1,0,0,0,0],[0,1,1,0,0],[1,0,1,0,0]]);
# > print solve_right($A,$B);
# | 0 1 0 0 0
# | 1 0 0 0 0
# | 0 0 1 0 0
# | 0 0 0 0 0
# | 0 0 0 0 0
user_function solve_right(Matrix:wary, Matrix:wary) : c++ (include => "polymake/matrix_linalg.h");

# @category Linear Algebra
# Computes a [[Matrix]] X that solves the system //XA// = //B//.
# @param Matrix A
# @param Matrix B
# @return Matrix
# @example This is useful, for instance, for computing the coordinates of some vectors with respect to a basis.
# The rows of the matrix solve_left(A,B) are the coordinates of the rows of //B// with respect to the rows of //A//.
# Define the matrices
# > $A = new Matrix([[-1,1,0],[0,-1,1]]);
# > $B = new Matrix([[-4,2,2],[3,-2,-1]]);
# so that the rows of //A// are a basis of the subspace of vectors with zero coordinate sum. Then the rows of
# > print solve_left($A, $B);
# | 4 2
# | -3 -1
# contain the coordinates of the rows of //B// with respect to the rows of //A//.
user_function solve_left(Matrix:wary, Matrix:wary) : c++ (include => "polymake/matrix_linalg.h");


# @category Linear Algebra
# Computes the solution of the system //A//x = //b// for a given invertible [[Matrix]] //A// and a [[Vector]] //b// by applying Cramer's rule.
# @param Matrix A must be invertible
# @param Vector b
# @return Vector
# @example from the Wikipedia:
# > $A = new Matrix([3,2,-1],[2,-2,4],[-1,1/2,-1]);
# > $b = new Vector(1,-2,0);
# > print cramer($A,$b);
# | 1 -2 -2

user_function cramer(Matrix:wary, Vector:wary) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# A naive test(exponential in the size of matrix) to check if a [[Matrix]] //M// is totally unimodular.
# The matrix //M// is totally unimodular if the determinant of each square submatrix equals 0, 1, or -1.
#
# For a better implementation try Matthias Walter's polymake extension at
#  [[https://github.com/xammy/unimodularity-test/wiki/Polymake-Extension]].
# @param Matrix M
# @return Bool
# @example
# > $M = new Matrix<Int>([-1,-1,0,0,0,1],[1,0,-1,-1,0,0],[0,1,1,0,-1,0],[0,0,0,1,1,-1]);
# > print totally_unimodular($M);
# | true

user_function totally_unimodular(Matrix) : c++ (include => "polymake/totally_unimodular.h");

# @category Linear Algebra
# Check whether both matrices are bases of the same linear subspace.
# Note: It is assumed that they are *bases* of the row space.
# @param Matrix M1
# @param Matrix M2
# @return Bool
# @example
# > $M1 = new Matrix([1,1,0],[1,0,1],[0,0,1]);
# > $M2 = new Matrix([1,0,0],[0,1,0],[0,0,1]);
# > print equal_bases($M1,$M2);
# | true
user_function equal_bases(Matrix, Matrix) {
   my ($M1, $M2)=@_;
   # we first need to check for full rank matrices as null_space gives an empty matrix in that case
   return $M1 == $M2 || $M1->rows == $M2->rows && $M1->cols == $M2->cols &&
          ($M1->rows == $M1->cols || is_zero(null_space($M1) * transpose($M2)));
}

# @category Linear Algebra
# Householder transformation of [[Vector]] //b//. Only the orthogonal matrix reflection H is returned.
# @param Vector<Float> b
# @return Matrix<Float>
user_function householder_trafo(Vector<Float>) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Returns the eigenvalues of the given [[Matrix]] //M//.
# @param Matrix<Float> M
# @return Vector<Float>
# @example
# > $M = new Matrix<Float>([[2,0,0],[0,3,4],[0,4,9]]);
# > print(eigenvalues($M));
# | 2 11 1
user_function eigenvalues<_>(Matrix<Float,_>) : c++ (include => "polymake/linalg.h");


# @category Linear Algebra
# QR decomposition of a mxn [[Matrix]] //M// with m greater than or equal to n.
# @param Matrix<Float> M
# @return Pair<Matrix,Matrix>
# @example [nocompare]
# > $M = new Matrix<Float>([23,4],[6,42]);
# > $qr = qr_decomp($M);
# > print $qr->first;
# | 0.9676172724 0.2524218971
# | 0.2524218971 -0.9676172724
# > print $qr->second;
# | 23.76972865 14.47218877
# | 0 -39.63023785
# > print $qr->first * $qr->second ;
# | 23 4
# | 6 42
user_function qr_decomp<_>(Matrix<Float,_>) : c++ (include => "polymake/linalg.h");

# @topic category property_types/Linear Algebra
# These types are needed as return types of algebraic computations.

# @category Linear Algebra
# Complete result of the __singular value decomposition__ of a [[Matrix]] //M//,
# such that left_companion * sigma * transpose(right_companion) = //M//
# Contains the following fields:
# @field Matrix<Float> sigma the diagonalized matrix
# @field Matrix<Float> left_companion matrix of left singular vectors
# @field Matrix<Float> right_companion matrix of right singular vectors
declare property_type SingularValueDecomposition : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# SVD decomposition of a Matrix. Computes the SVD of a matrix into a diagonal Marix (S), orthogonal square Matrix (U), orthogonal square Matrix (V), such that U*S*V^T=M
# The first element of the output array is S, the second U and the third V.
# @param Matrix<Float> M
# @return SingularValueDecomposition
# @example
# > $M = new Matrix<Float>([1,2],[23,24]);
# > $SVD = singular_value_decomposition($M);
# The following prints the three matrices, separated by newline characters.
# > print $SVD->left_companion ,"\n", $SVD->sigma ,"\n", $SVD->right_companion;
# | 0.06414638608 0.9979404998
# | 0.9979404998 -0.06414638608
# |
# | 33.31011547 0
# | 0 0.6604600341
# |
# | 0.6909846321 -0.7228694476
# | 0.7228694476 0.6909846321
# > print $SVD->left_companion * $SVD->sigma * transpose($SVD->right_companion);
# | 1 2
# | 23 24
user_function singular_value_decomposition<_>(Matrix<Float,_>) : c++ (include => "polymake/linalg.h");

# @category Linear Algebra
# Computes the Moore-Penrose Inverse of a [[Matrix]] //M//.
# @param Matrix<Float> M
# @return Matrix<Float>
# @example
# > $M = new Matrix<Float>([1,0],[0,1],[0,1]);
# > print(moore_penrose_inverse($M));
# | 1 0 0
# | 0 0.5 0.5
user_function moore_penrose_inverse<_>(Matrix<Float,_>) : c++ (include => "polymake/linalg.h");

# @category Data Conversion
# Return the input [[Vector]] (which is already in the dense form).
# @param Vector v
# @return Vector
user_function dense(Vector) { $_[0] }

# @category Data Conversion
# Return the input [[Matrix]] (which is already in dense form).
# @param Matrix M
# @return Matrix
user_function dense(Matrix) { $_[0] }

##################################################################################

# @category Algebraic Types
# A SparseVector is a special form of a [[Vector]] optimized for large vectors with few non-zero elements.
# Physically it is stored as a balanced binary tree with element index as a search key.
#
# The printable representation of a SparseVector looks like a sequence (l) (p<sub>1</sub> v<sub>1</sub>) ... (p<sub>k</sub> v<sub>k</sub>),
# where l is the dimension of the vector and each pair (p<sub>i</sub> v<sub>i</sub>) denotes an entry with value
# v<sub>i</sub> at position p<sub>i</sub>. All other entries are zero.
#
# Use [[dense]] to convert a sparse vector into the [[Vector|dense form]].
#
# @example construct a SparseVector from a contiguous list of scalars; zero values are filtered out automatically:
# > $v = new SparseVector<Integer>(0, 1, 2, 0, 0, 0);
# > print $v, "\n";
# | (6) (1 1) (2 2)
# @example construct a SparseVector from a hash map with non-zero values; dimension is provided explicitly:
# > $v = new SparseVector<Integer>({ 1 => 1, 2 => 2, _dim => 6 });
# > print $v, "\n";
# | (6) (1 1) (2 2)
# @tparam Element
declare property_type SparseVector<Element=Rational> : Vector<Element> : c++ (include => ["polymake/SparseVector.h"]) {

   # The number of non-zero entries of a sparse vector.
   # @return Int
   # @example
   # > $v = new SparseVector<Integer>(0, 1, 2, 0, 0, 0);
   # > print($v->size());
   # | 2
   user_method size() : c++;
}


property_type Vector {

   # accept a SparseVector of a different type as a value of Vector property but substitute the desired element type
   type_method coherent_type {
      my ($self, $value)=@_;
      instanceof SparseVector($value) || instanceof Core::PropertyType($value) && $value->name eq "SparseVector" ? typeof SparseVector($self->params->[0]) : undef;
   }
}

# @category Algebraic Types
# A SparseMatrix is a special form of a [[Matrix]] optimized for large matrices with few non-zero elements.
# Physically it is stored as a two-dimensional grid of balanced binary trees with row and column indexes as search keys.
#
# Use [[dense]] to convert a SparseMatrix into the [[Matrix|dense form]].
#
# @example construct a SparseMatrix from a list of row vectors.  Each row can be specified in a dense or sparse form.
# Number of columns is determined by the size of dense rows or by the maximal index value in sparse rows.
# > $M = new SparseMatrix<Integer>(
# >  [0, 1, 0, 0, 0],
# >  {4 => 2},
# >  {},
# >  [3, -1, 0, 0, 0]);
# > print $M;
# | (5) (1 1)
# | (5) (4 2)
# | (5)
# | (5) (0 3) (1 -1)
# @example construct a SparseMatrix with explicitly given number of columns:
# > $M = new SparseMatrix<Integer>(
# >   {1 => 10},
# >   {2 => 20},
# >   {cols => 5});
# > print $M;
# | (5) (1 10)
# | (5) (2 20)
# @tparam Element
# @tparam Sym one of [[Symmetric]] or [[NonSymmetric]], default: [[NonSymmetric]]
declare property_type SparseMatrix<Element=Rational, Sym=NonSymmetric> : Matrix<Element,Sym> : c++ (include => ["polymake/SparseMatrix.h"]) {

   # Adjusts the given [[SparseMatrix]] by removing empty rows and columns.
   # The remaining rows and columns are renumbered without gaps.
   # @example
   # > $M = new SparseMatrix<Integer>({1 =>2, _dim => 6},{5 => 3, _dim => 6}, {});
   # print(dense($M));
   # | 0 2 0 0 0 0
   # | 0 0 0 0 0 3
   # | 0 0 0 0 0 0
   # > $M->squeeze();
   # > print(dense($M));
   # | 2 0
   # | 0 3
   user_method squeeze(&) : c++;

   # Adjusts the given [[SparseMatrix]] by removing empty rows.
   # The remaining rows are renumbered without gaps.
   # @example
   # > $M = new SparseMatrix<Integer>({1 =>2, _dim => 6},{5 => 3, _dim => 6}, {});
   # print(dense($M));
   # | 0 2 0 0 0 0
   # | 0 0 0 0 0 3
   # | 0 0 0 0 0 0
   # > $M->squeeze_rows();
   # > print(dense($M));
   # | 0 2 0 0 0 0
   # | 0 0 0 0 0 3
   user_method squeeze_rows(&) : c++;

   # Adjusts the given [[SparseMatrix]] by removing empty columns.
   # The remaining columns are renumbered without gaps.
   # @example
   # > $M = new SparseMatrix<Integer>({1 =>2, _dim => 6},{5 => 3, _dim => 6}, {});
   # > print(dense($M));
   # | 0 2 0 0 0 0
   # | 0 0 0 0 0 3
   # | 0 0 0 0 0 0
   # > $M->squeeze_cols();
   # > print(dense($M));
   # | 2 0
   # | 0 3
   # | 0 0
   user_method squeeze_cols(&) : c++;

   # Resizes the [[Matrix]] according to given number of rows and columns.
   # All elements in added rows and/or columns are implicit zeros.
   # If the number of rows and/or columns of the original matrix is less than the given number or rows and/or columns, the extra entries are deleted.
   # @param Int r new number of rows
   # @param Int c new number of columns
   # @example
   # > $M = new SparseMatrix<Integer>({1 =>2, _dim => 6},{5 => 3, _dim => 6}, {});
   # print(dense($M));
   # | 0 2 0 0 0 0
   # | 0 0 0 0 0 3
   # | 0 0 0 0 0 0
   # > $M->resize(3,4);
   # > print(dense($M));
   # | 0 2 0 0
   # | 0 0 0 0
   # | 0 0 0 0
   user_method resize(&, $$) : c++;
}


property_type Matrix {

   # accept a SparseMatrix of a different type as a value of Matrix property but substitute the desired element type
   type_method coherent_type {
      my ($self, $value)=@_;
      instanceof SparseMatrix($value) || instanceof Core::PropertyType($value) && $value->name eq "SparseMatrix" ? typeof SparseMatrix(@{$self->params}) : undef;
   }
}


# @category Data Conversion
# Convert to an equivalent dense vector of the same element type.
# Returns the input [[SparseVector]] as a dense [[Vector]] of the same element type.
# @tparam Element
# @param SparseVector<Element> v
# @return Vector<Element>
# @example
# > $v = new SparseVector<Integer>([0,1,2,0,0,0]);
# > print($v);
# | (6) (1 1) (2 2)
# > print(dense($v));
# | 0 1 2 0 0 0
user_function dense<Element>(SparseVector<Element>) { new Vector<Element>(shift) }

# @category Data Conversion
# Returns the input [[SparseMatrix]] as a dense [[Matrix]] of the same element type.
# @tparam Element
# @param SparseMatrix<Element> M
# @return Matrix<Element>
# @example
# > $M = new SparseMatrix([[0,0,1],[0,0,0],[0,2,0]]);
# > print($M);
# | (3) (2 1)
# | (3)
# | (3) (1 2)
# > print(dense($M));
# | 0 0 1
# | 0 0 0
# | 0 2 0
user_function dense<Element>(SparseMatrix<Element,_>) { new Matrix<Element>(shift) }

# @category Data Conversion
# Creates a [[SparseVector]] with dimension //d// and having 1's at positions contained in the given [[Set]] //S//.
# @param Set S
# @param Int d dimension of the result
# @tparam Scalar type of apparent 1's
# @return SparseVector<Scalar>
# @example
# > $v = toVector<Float>(new Set([0,1]), 5);
# > print($v);
# | (5) (0 1) (1 1)
user_function toVector<Scalar>(Set:wary:anchor $) : c++ (name => 'same_element_sparse_vector', include => "polymake/SparseVector.h");

# @category Data Conversion
# Converts an [[IncidenceMatrix]] to a [[SparseMatrix]].
# @param IncidenceMatrix M
# @tparam Scalar
# @return SparseMatrix<Scalar>
# @example
# > $M = polytope::cube(2)->VERTICES_IN_FACETS;
# > print $M->type->full_name;
# | IncidenceMatrix<NonSymmetric>
# > print(toMatrix<Int>($M)->type->full_name);
# | SparseMatrix<Int, NonSymmetric>
user_function toMatrix<Scalar>(IncidenceMatrix:anchor) : c++ (name => 'same_element_sparse_matrix', include => "polymake/SparseMatrix.h");

# @category Data Conversion
# Converts an [[IncidenceMatrix]] to a dense 0/1 [[Matrix]].
# @param IncidenceMatrix M
# @return Matrix<Int>
# @example
# > $M = polytope::cube(2)->VERTICES_IN_FACETS;
# > print(dense($M));
# | 1 0 1 0
# | 0 1 0 1
# | 1 1 0 0
# | 0 0 1 1
user_function dense(IncidenceMatrix) { dense(toMatrix<Int>(@_)); }

# @category Data Conversion
# Converts the given [[Set]] to a dense 0/1 [[Vector]] of a given dimension.
# @param Set S
# @param Int dim
# @return Vector<Int>
# @example
# > $S = new Set([0,1]);
# > print(dense($S,3));
# | 1 1 0
user_function dense(Set $) { dense(toVector<Int>(@_)); }

# @category Data Conversion
# Get the positions of non-zero entries of a [[SparseVector]].
# @param SparseVector v
# @return Set<Int>
# @example
# > $v = new SparseVector(0,1,1,0,0,0,2,0,3);
# > print indices($v);
# | {1 2 6 8}
user_function indices(SparseVector:anchor) : c++ (include => "polymake/Set.h");

# @category Data Conversion
# Gets the positions of non-zero entries of a [[Vector]].
# @param Vector v
# @return Set<Int>
# @example
# > print support(new Vector(0,23,0,0,23,0,23,0,0,23));
# | {1 4 6 9}
user_function support(Vector:anchor) : c++ (include => "polymake/linalg.h");

# @category Data Conversion
# Gets the positions of non-zero entries of a [[SparseMatrix]].
# @param SparseMatrix M
# @return IncidenceMatrix
# @example
# > $S = new SparseMatrix([1,2,0,0,0,0],[0,0,5,0,0,32]);
# > print index_matrix($S);
# | {0 1}
# | {2 5}
user_function index_matrix(SparseMatrix:anchor) : c++ (include => "polymake/IncidenceMatrix.h");


##################################################################################

# types for tropical addition

# @category Arithmetic
# tropical addition: min
declare property_type Min : c++ (special => 'Min', include => "polymake/TropicalNumber.h") {

    user_method orientation() { return 1; }

    user_method apply {my ($o,$x,$y) = @_; return min($x,$y);}
}

# @category Arithmetic
# tropical addition: max
declare property_type Max : c++ (special => 'Max', include => "polymake/TropicalNumber.h") {

    user_method orientation() { return -1; }

    user_method apply {my ($o,$x,$y) = @_; return max($x,$y);}
}

# @topic any/tparam/Addition
# Mode of tropical addition, must be [[Min]] or [[Max]].
# There is on purpose no default value for it.

# @category Arithmetic
# @tparam Addition
# @tparam Scalar default: [[Rational]]
declare property_type TropicalNumber<Addition, Scalar=Rational> : upgrades( Scalar ) : c++ (include => "polymake/TropicalNumber.h") {

    operator ++ + += * *= / /= neg - @compare : c++;

    # The orientation of the associated addition, i.e.
    # +1 if the corresponding 0 is +inf
    # -1 if the corresponding 0 is -inf
    # @return Int
    user_method orientation() {
      return Addition->orientation();
    }

    # The zero element of the tropical semiring of this element.
    # @return Scalar
    user_method zero(:static)  : c++ (include => "polymake/TropicalNumber.h");

    type_method JSONschema { shift; (typeof Scalar)->JSONschema->(@_) }
};

# function is_tropical_addition(Min) {1}
# function is_tropical_addition(Max) {1}
# function is_tropical_addition() {0}

##################################################################################

# @topic category property_types/Linear Algebra
# These types are needed as return types of algebraic computations.

# @category Linear Algebra
# Complete result of the __Smith normal form__ computation of the input matrix //M//.
# @field SparseMatrix<Scalar> form the Smith normal form
# @field List<Pair<Scalar, Int>> torsion absolute values of the entries greater than 1 of the diagonal together with their multiplicity
# @field Int rank rank of //M//
# @field SparseMatrix<Scalar> left_companion unimodular matrix L such that //M// = L*S*R
# @field SparseMatrix<Scalar> right_companion unimodular matrix R such that //M// = L*S*R
# @tparam Scalar matrix element type
declare property_type SmithNormalForm<Scalar> : c++ (include => "polymake/Smith_normal_form.h");

# @category Linear Algebra
# Computes the __Smith normal form__ of a given [[Matrix]] //M//.
# M = L*S*R in normal case, or S = L*M*R in inverted case.
# @param Matrix M must be of integer type
# @param Bool inv if true, the companion matrices in the result will be inverted
# @return SmithNormalForm
# @example
# > $M = new Matrix<Integer>([1,2],[23,24]);
# > $SNF = smith_normal_form($M);
# The following line prints the three matrices separated by newline characters.
# > print $SNF->left_companion ,"\n", $SNF->form ,"\n", $SNF->right_companion;
# | 1 0
# | 23 1
# |
# | 1 0
# | 0 -22
# |
# | 1 2
# | 0 1
user_function smith_normal_form(Matrix; $=false) : c++ (include => "polymake/Smith_normal_form.h");

# @category Linear Algebra
# Complete result of the __Hermite normal form__ computation of the input matrix //M//.
# @field Matrix<Scalar> hnf the Hermite normal form
# @field SparseMatrix<Scalar> companion unimodular matrix R such that M*R = H
# @field Int rank rank of //M//
# @tparam Scalar matrix element type
declare property_type HermiteNormalForm<Scalar> : c++ (include => "polymake/integer_linalg.h");

# @category Linear Algebra
# Computes the (column) Hermite normal form of an integer [[Matrix]].
# Pivot entries are positive, entries to the left of a pivot are non-negative and strictly smaller than the pivot.
# @param Matrix M matrix to be transformed.
# @option Bool reduced If this is false, entries to the left of a pivot are left untouched. True by default
# @return HermiteNormalForm
# @example The following stores the result for a small matrix M in H and then prints both hnf and companion:
# > $M = new Matrix<Integer>([1,2],[2,3]);
# > $H = hermite_normal_form($M);
# > print $H->hnf;
# | 1 0
# | 0 1
# > print $H->companion;
# | -3 2
# | 2 -1
user_function hermite_normal_form(Matrix; $=true) : c++ (include => "polymake/integer_linalg.h");

# @category Linear Algebra
# Computes the __lattice null space__ of the integer [[Matrix]] //A//.
# @param Matrix A
# @return SparseMatrix Has a lattice basis of the null space as rows.
# @example
# > $M = new Matrix<Integer>([1,0,0,0],[1,0,1,0],[0,0,1,0]);
# > print null_space_integer($M);
# | (4) (1 1)
# | (4) (3 1)


user_function null_space_integer(Matrix) : c++ (include=>"polymake/integer_linalg.h");

# @category Linear Algebra
# Computes a lattice basis of the span of the rows of //A//.
# @param Matrix<Integer> A
# @return Matrix<Integer> Basis of the lattice spanned by the rows of //A//.
# @example No two of the rows of the following matrix form a basis.
# > $A = new Matrix<Integer>([[2,3],[1,3],[2,4]]);
# > print lattice_basis($A);
# | 2 3
# | -1 -1
user_function lattice_basis(Matrix) : c++ (include=>"polymake/common/lattice_tools.h");


# Local Variables:
# mode: perl
# cperl-indent-level: 3
# indent-tabs-mode:nil
# End: