This repository has been archived by the owner on Jan 30, 2023. It is now read-only.
-
-
Notifications
You must be signed in to change notification settings - Fork 7
/
base.pyx
2129 lines (1768 loc) · 79.6 KB
/
base.pyx
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
r"""
Combinatorial polyhedron
This module gathers algorithms for polyhedra that only depend on the
vertex-facet incidences and that are called combinatorial polyhedron.
The main class is :class:`CombinatorialPolyhedron`. Most importantly,
this class allows to iterate quickly through the faces (possibly
of given dimension) via the :class:`~sage.geometry.polyhedron.combinatorial_polyhedron.face_iterator.FaceIterator` object. The :class:`CombinatorialPolyhedron`
uses this iterator to quickly generate the f-vector, the edges,
the ridges and the face lattice.
Terminology used in this module:
- Vrepr -- ``[vertices, rays, lines]`` of the polyhedron.
- Hrepr -- inequalities and equalities of the polyhedron.
- Facets -- facets of the polyhedron.
- Vrepresentation -- represents a face by the list of Vrepr it contains.
- Hrepresentation -- represents a face by a list of Hrepr it is contained in.
- bit representation -- represents incidences as ``uint64_t``-array, where
each bit represents one incidence. There might
be trailing zeros, to fit alignment requirements.
In most instances, faces are represented by the
bit representation, where each bit corresponds to
a Vrepr or facet. Thus a bit representation can either be
a Vrepr or facet representation depending on context.
EXAMPLES:
Construction::
sage: P = polytopes.hypercube(4)
sage: C = CombinatorialPolyhedron(P); C
A 4-dimensional combinatorial polyhedron with 8 facets
Obtaining edges and ridges::
sage: C.edges()[:2]
((A vertex at (1, 1, 1, -1), A vertex at (1, 1, 1, 1)),
(A vertex at (1, 1, -1, 1), A vertex at (1, 1, 1, 1)))
sage: C.edges(names=False)[:2]
((14, 15), (13, 15))
sage: C.ridges()[:2]
((An inequality (0, 0, 1, 0) x + 1 >= 0,
An inequality (0, 0, 0, 1) x + 1 >= 0),
(An inequality (0, 1, 0, 0) x + 1 >= 0,
An inequality (0, 0, 0, 1) x + 1 >= 0))
sage: C.ridges(names=False)[:2]
((6, 7), (5, 7))
Vertex-graph and facet-graph::
sage: C.vertex_graph()
Graph on 16 vertices
sage: C.facet_graph()
Graph on 8 vertices
Face lattice::
sage: C.face_lattice()
Finite lattice containing 82 elements
Face iterator::
sage: C.face_iter()
Iterator over the proper faces of a 4-dimensional combinatorial polyhedron
sage: C.face_iter(2)
Iterator over the 2-faces of a 4-dimensional combinatorial polyhedron
AUTHOR:
- Jonathan Kliem (2019-04)
"""
# ****************************************************************************
# Copyright (C) 2019 Jonathan Kliem <jonathan.kliem@fu-berlin.de>
#
# 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 of the License, or
# (at your option) any later version.
# https://www.gnu.org/licenses/
# ****************************************************************************
from __future__ import absolute_import, division, print_function
import numbers
from sage.rings.integer import Integer
from sage.graphs.graph import Graph
from sage.graphs.digraph import DiGraph
from sage.combinat.posets.lattices import FiniteLatticePoset
from sage.geometry.polyhedron.base import Polyhedron_base
from sage.geometry.lattice_polytope import LatticePolytopeClass
from sage.structure.element import Matrix
from sage.misc.misc import is_iterator
from .conversions \
import incidence_matrix_to_bit_repr_of_facets, \
incidence_matrix_to_bit_repr_of_Vrepr, \
facets_tuple_to_bit_repr_of_facets, \
facets_tuple_to_bit_repr_of_Vrepr
from sage.rings.integer cimport smallInteger
from cysignals.signals cimport sig_check, sig_block, sig_unblock
cdef extern from "Python.h":
int unlikely(int) nogil # Defined by Cython
cdef class CombinatorialPolyhedron(SageObject):
r"""
The class of the Combinatorial Type of a Polyhedron, a Polytope.
INPUT:
- ``data`` -- an instance of
* :class:`~sage.geometry.polyhedron.parent.Polyhedron_base`
* or a :class:`~sage.geometry.lattice_polytope.LatticePolytopeClass`
* or an ``incidence_matrix`` as in
:meth:`~sage.geometry.polyhedron.base.Polyhedron_base.incidence_matrix`
In this case you should also specify the ``Vrepr`` and ``facets`` arguments
* or list of facets, each facet given as
a list of ``[vertices, rays, lines]`` if the polyhedron is unbounded,
then rays and lines and the extra argument ``nr_lines`` are required
if the polyhedron contains no lines, the rays can be thought of
as the vertices of the facets deleted from a bounded polyhedron see
:class:`~sage.geometry.polyhedron.parent.Polyhedron_base` on how to use
rays and lines.
* or an integer, representing the dimension of a polyhedron equal to its
affine hull
- ``Vrepr`` -- (optional) when ``data`` is an incidence matrix, it should
be the list of ``[vertices, rays, lines]``, if the rows in the incidence_matrix
should correspond to names
- ``facets`` -- (optional) when ``data`` is an incidence matrix or a list of facets,
it should be a list of facets that would be used instead of indices (of the columns
of the incidence matrix).
- ``unbounded`` -- value will be overwritten if ``data`` is a polyhedron;
if ``unbounded`` and ``data`` is incidence matrix or a list of facets,
need to specify ``far_face``
- ``far_face`` -- (semi-optional) when ``data` is an incidence matrix or a
list of facets and the polyhedron is unbounded this needs to be set to
the list of indices of the rays and lines
EXAMPLES:
We illustrate all possible input: a polyhedron:
sage: P = polytopes.cube()
sage: CombinatorialPolyhedron(P)
A 3-dimensional combinatorial polyhedron with 6 facets
a lattice polytope::
sage: points = [(1,0,0), (0,1,0), (0,0,1),
....: (-1,0,0), (0,-1,0), (0,0,-1)]
sage: L = LatticePolytope(points)
sage: CombinatorialPolyhedron(L)
A 3-dimensional combinatorial polyhedron with 8 facets
an incidence matrix::
sage: P = Polyhedron(rays=[[0,1]])
sage: data = P.incidence_matrix()
sage: far_face = [i for i in range(2) if not P.Vrepresentation()[i].is_vertex()]
sage: CombinatorialPolyhedron(data, unbounded=True, far_face=far_face)
A 1-dimensional combinatorial polyhedron with 1 facet
sage: C = CombinatorialPolyhedron(data, Vrepr=['myvertex'],
....: facets=['myfacet'], unbounded=True, far_face=far_face)
sage: C.Vrepresentation()
('myvertex',)
sage: C.Hrepresentation()
('myfacet',)
a list of facets::
sage: CombinatorialPolyhedron(((1,2,3),(1,2,4),(1,3,4),(2,3,4)))
A 3-dimensional combinatorial polyhedron with 4 facets
sage: facetnames = ['facet0', 'facet1', 'facet2', 'myfacet3']
sage: facetinc = ((1,2,3),(1,2,4),(1,3,4),(2,3,4))
sage: C = CombinatorialPolyhedron(facetinc, facets=facetnames)
sage: C.Vrepresentation()
(1, 2, 3, 4)
sage: C.Hrepresentation()
('facet0', 'facet1', 'facet2', 'myfacet3')
an integer::
sage: CombinatorialPolyhedron(-1).f_vector()
(1)
sage: CombinatorialPolyhedron(0).f_vector()
(1, 1)
sage: CombinatorialPolyhedron(5).f_vector()
(1, 0, 0, 0, 0, 0, 1)
Specifying that a polyhedron is unbounded is important. The following with a
polyhedron works fine::
sage: P = Polyhedron(ieqs=[[1,-1,0],[1,1,0]])
sage: C = CombinatorialPolyhedron(P) # this works fine
sage: C
A 2-dimensional combinatorial polyhedron with 2 facets
The following is incorrect, as ``unbounded`` is implicitly set to ``False``::
sage: data = P.incidence_matrix()
sage: vert = P.Vrepresentation()
sage: C = CombinatorialPolyhedron(data, Vrepr=vert)
sage: C
A 2-dimensional combinatorial polyhedron with 2 facets
sage: C.f_vector()
Traceback (most recent call last):
...
ValueError: not all vertices are intersections of facets
sage: C.vertices()
(A line in the direction (0, 1), A vertex at (1, 0), A vertex at (-1, 0))
The correct usage is::
sage: far_face = [i for i in range(3) if not P.Vrepresentation()[i].is_vertex()]
sage: C = CombinatorialPolyhedron(data, Vrepr=vert, unbounded=True, far_face=far_face)
sage: C
A 2-dimensional combinatorial polyhedron with 2 facets
sage: C.f_vector()
(1, 0, 2, 1)
sage: C.vertices()
()
TESTS:
Checking that :trac:`27987` is fixed::
sage: P1 = Polyhedron(vertices=[[0,1],[1,0]], rays=[[1,1]])
sage: P2 = Polyhedron(vertices=[[0,1],[1,0],[1,1]])
sage: P1.incidence_matrix() == P2.incidence_matrix()
True
sage: CombinatorialPolyhedron(P1).f_vector()
(1, 2, 3, 1)
sage: CombinatorialPolyhedron(P2).f_vector()
(1, 3, 3, 1)
sage: P1 = Polyhedron(vertices=[[0,1],[1,0]], rays=[[1,1]])
sage: P2 = Polyhedron(vertices=[[0,1],[1,0],[1,1]])
sage: CombinatorialPolyhedron(P1).f_vector()
(1, 2, 3, 1)
sage: CombinatorialPolyhedron(P2).f_vector()
(1, 3, 3, 1)
Some other tests regarding small polyhedra::
sage: P = Polyhedron(rays=[[1,0],[0,1]])
sage: C = CombinatorialPolyhedron(P)
sage: C
A 2-dimensional combinatorial polyhedron with 2 facets
sage: C.f_vector()
(1, 1, 2, 1)
sage: C.vertices()
(A vertex at (0, 0),)
sage: data = P.incidence_matrix()
sage: vert = P.Vrepresentation()
sage: far_face = [i for i in range(3) if not P.Vrepresentation()[i].is_vertex()]
sage: C = CombinatorialPolyhedron(data, Vrepr=vert, unbounded=True, far_face=far_face)
sage: C
A 2-dimensional combinatorial polyhedron with 2 facets
sage: C.f_vector()
(1, 1, 2, 1)
sage: C.vertices()
(A vertex at (0, 0),)
sage: CombinatorialPolyhedron(3r)
A 3-dimensional combinatorial polyhedron with 0 facets
Check that on wrong input subsequent calls of ``f_vector`` fail::
sage: data = P.incidence_matrix()
sage: vert = P.Vrepresentation()
sage: C = CombinatorialPolyhedron(data, Vrepr=vert)
sage: C.f_vector()
Traceback (most recent call last):
...
ValueError: not all vertices are intersections of facets
sage: C.f_vector()
Traceback (most recent call last):
...
ValueError: not all vertices are intersections of facets
Check that :trac:`28678` is fixed::
sage: CombinatorialPolyhedron([])
A -1-dimensional combinatorial polyhedron with 0 facets
sage: CombinatorialPolyhedron(LatticePolytope([], lattice=ToricLattice(3)))
A -1-dimensional combinatorial polyhedron with 0 facets
"""
def __init__(self, data, Vrepr=None, facets=None, unbounded=False, far_face=None):
r"""
Initialize :class:`CombinatorialPolyhedron`.
See :class:`CombinatorialPolyhedron`.
TESTS::
sage: C = CombinatorialPolyhedron([[0,1,2],[0,1,3],
....: [0,2,3],[1,2,3]]) # indirect doctest
sage: TestSuite(sage.geometry.polyhedron.combinatorial_polyhedron.base.CombinatorialPolyhedron).run()
"""
self._dimension = -2 # a "NULL" value
self._edges = NULL
self._ridges = NULL
self._face_lattice_incidences = NULL
self._equalities = ()
self._all_faces = None
self._mem_tuple = ()
# ``_length_edges_list`` should not be touched in an instance
# of :class:`CombinatorialPolyhedron`. This number can be altered,
# but should probably be a power of `2` (for memory usage).
# ``_length_edges_list`` shouldn't be too small for speed and
# shouldn't be too large, as ``ridges``, ``edges`` and ``incidences``
# each have a memory overhead of
# ``self._length_edges_list*2*sizeof(size_t *)``.
self._length_edges_list = 16348
if isinstance(data, Polyhedron_base):
# input is ``Polyhedron``
Vrepr = data.Vrepresentation()
facets = tuple(inequality for inequality in data.Hrepresentation())
self._dimension = data.dimension()
if not data.is_compact():
self._unbounded = True
far_face = tuple(i for i in range(data.n_Vrepresentation()) if not data.Vrepresentation()[i].is_vertex())
else:
self._unbounded = False
data = data.incidence_matrix()
elif isinstance(data, LatticePolytopeClass):
# input is ``LatticePolytope``
self._unbounded = False
Vrepr = data.vertices()
self._length_Vrepr = len(Vrepr)
facets = data.facets()
self._length_Hrepr = len(facets)
data = tuple(tuple(vert for vert in facet.vertices())
for facet in facets)
else:
# Input is different from ``Polyhedron`` and ``LatticePolytope``.
if not unbounded:
# bounded polyhedron
self._unbounded = False
elif not far_face:
raise ValueError("must specify far face for unbounded polyhedron")
else:
self._unbounded = True
if Vrepr:
# store vertices names
self._V = tuple(Vrepr)
Vinv = {v: i for i,v in enumerate(self._V)}
else:
self._V = None
Vinv = None
if facets:
# store facets names and compute equalities
facets = tuple(facets)
test = [1] * len(facets) # 0 if that facet is an equality
for i in range(len(facets)):
if hasattr(facets[i], "is_inequality"):
# We remove equalities.
# At the moment only equalities with this attribute ``True``
# will be detected.
if not facets[i].is_inequality():
test[i] = 0
self._H = tuple(facets[i] for i in range(len(facets)) if test[i])
self._equalities = tuple(facets[i] for i in range(len(facets)) if not test[i])
else:
self._H = None
if data == [] or data == ():
# Handling the empty polyhedron.
data = -1
if isinstance(data, Matrix):
# Input is incidence-matrix or was converted to it.
self._length_Hrepr = data.ncols()
self._length_Vrepr = data.nrows()
# Initializing the facets in their Bit-representation.
self._bitrep_facets = incidence_matrix_to_bit_repr_of_facets(data)
# Initializing the Vrepr as their Bit-representation.
self._bitrep_Vrepr = incidence_matrix_to_bit_repr_of_Vrepr(data)
self._n_facets = self.bitrep_facets().n_faces
# Initialize far_face if unbounded.
if self._unbounded:
self._far_face = facets_tuple_to_bit_repr_of_facets((tuple(far_face),), self._length_Vrepr)
else:
self._far_face = None
elif isinstance(data, numbers.Integral):
# To construct a trivial polyhedron, equal to its affine hull,
# one can give an Integer as Input.
if data < -1:
ValueError("any polyhedron must have dimension at least -1")
self._n_facets = 0
self._dimension = data
# Initializing the facets in their Bit-representation.
self._bitrep_facets = facets_tuple_to_bit_repr_of_facets((), 0)
# Initializing the Vrepr as their Bit-representation.
self._bitrep_Vrepr = facets_tuple_to_bit_repr_of_Vrepr((), 0)
self._far_face = None
else:
# Input is a "list" of facets.
# The facets given by its ``[vertices, rays, lines]``.
# Actually at least tuple, list, iterator will work.
if is_iterator(data):
data = tuple(data)
if self._V is None:
# Get the names of the Vrepr.
Vrepr = sorted(set.union(*map(set, data)))
length_Vrepr = len(Vrepr)
if Vrepr != range(len(Vrepr)):
self._V = tuple(Vrepr)
Vinv = {v: i for i,v in enumerate(self._V)}
else:
# Assuming the user gave as correct names for the vertices
# and labeled them instead by `0,...,n`.
length_Vrepr = len(self._V)
self._length_Vrepr = length_Vrepr
# Relabel the Vrepr to be `0,...,n`.
if self._V is not None:
def f(v): return Vinv[v]
else:
def f(v): return int(v)
facets = tuple(tuple(f(i) for i in j) for j in data)
self._n_facets = len(facets)
self._length_Hrepr = len(facets)
# Initializing the facets in their Bit-representation.
self._bitrep_facets = facets_tuple_to_bit_repr_of_facets(facets, length_Vrepr)
# Initializing the Vrepr as their Bit-representation.
self._bitrep_Vrepr = facets_tuple_to_bit_repr_of_Vrepr(facets, length_Vrepr)
# Initialize far_face if unbounded.
if self._unbounded:
self._far_face = facets_tuple_to_bit_repr_of_facets((tuple(far_face),), length_Vrepr)
else:
self._far_face = None
if self._unbounded:
self._far_face_tuple = tuple(far_face)
else:
self._far_face_tuple = ()
def _repr_(self):
r"""
Return a description of the combinatorial polyhedron.
EXAMPLES::
sage: P = polytopes.simplex()
sage: C = CombinatorialPolyhedron(P)
sage: C._repr_()
'A 3-dimensional combinatorial polyhedron with 4 facets'
sage: P = Polyhedron(vertices=[])
sage: C = CombinatorialPolyhedron(P)
sage: C._repr_()
'A -1-dimensional combinatorial polyhedron with 0 facets'
sage: P = Polyhedron(vertices=[[0,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C._repr_()
'A 0-dimensional combinatorial polyhedron with 1 facet'
sage: P = Polyhedron(lines=[[0,0,1],[0,1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C._repr_()
'A 2-dimensional combinatorial polyhedron with 0 facets'
sage: P = Polyhedron(rays=[[1,0,0],[0,1,0],[-1,0,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C._repr_()
'A 2-dimensional combinatorial polyhedron with 1 facet'
"""
desc = "A {}-dimensional combinatorial polyhedron with {} facet"\
.format(self.dimension(), self.n_facets())
if self.n_facets() != 1:
desc += "s"
return desc
def __reduce__(self):
r"""
Override __reduce__ to correctly pickle/unpickle.
TESTS::
sage: P = polytopes.permutahedron(4)
sage: C = CombinatorialPolyhedron(P)
sage: C1 = loads(C.dumps())
sage: it = C.face_iter()
sage: it1 = C1.face_iter()
sage: tup = tuple((face.Vrepr(), face.Hrepr()) for face in it)
sage: tup1 = tuple((face.Vrepr(), face.Hrepr()) for face in it1)
sage: tup == tup1
True
sage: P = polytopes.cyclic_polytope(4,10)
sage: C = CombinatorialPolyhedron(P)
sage: C1 = loads(C.dumps())
sage: it = C.face_iter()
sage: it1 = C1.face_iter()
sage: tup = tuple((face.Vrepr(), face.Hrepr()) for face in it)
sage: tup1 = tuple((face.Vrepr(), face.Hrepr()) for face in it1)
sage: tup == tup1
True
sage: P = Polyhedron(rays=[[1,0,0], [-1,0,0], [0,-1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C1 = loads(C.dumps())
sage: it = C.face_iter()
sage: it1 = C1.face_iter()
sage: tup = tuple((face.Vrepr(), face.Hrepr()) for face in it)
sage: tup1 = tuple((face.Vrepr(), face.Hrepr()) for face in it1)
sage: tup == tup1
True
sage: P = Polyhedron(rays=[[1,0,0], [-1,0,0],
....: [0,-1,0], [0,1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C1 = loads(C.dumps())
sage: it = C.face_iter()
sage: it1 = C1.face_iter()
sage: tup = tuple((face.Vrepr(), face.Hrepr()) for face in it)
sage: tup1 = tuple((face.Vrepr(), face.Hrepr()) for face in it1)
sage: tup == tup1
True
"""
# Give a constructor by list of facets.
if self.unbounded():
return (CombinatorialPolyhedron, (self.facets(),
self.Vrepresentation(), self.Hrepresentation(),
True, self.far_face_tuple()))
else:
return (CombinatorialPolyhedron, (self.facets(),
self.Vrepresentation(), self.Hrepresentation()))
def Vrepresentation(self):
r"""
Return a list of names of ``[vertices, rays, lines]``.
EXAMPLES::
sage: P = Polyhedron(rays=[[1,0,0], [0,1,0], \
....: [0,0,1],[0,0,-1]])
sage: C = CombinatorialPolyhedron(P)
sage: C.Vrepresentation()
(A line in the direction (0, 0, 1),
A ray in the direction (1, 0, 0),
A vertex at (0, 0, 0),
A ray in the direction (0, 1, 0))
"""
if self.V() is not None:
return self.V()
else:
return tuple(smallInteger(i) for i in range(self.length_Vrepr()))
def Hrepresentation(self):
r"""
Returns a list of names of facets and possibly some equalities.
EXAMPLES::
sage: P = polytopes.permutahedron(3)
sage: C = CombinatorialPolyhedron(P)
sage: C.Hrepresentation()
(An equation (1, 1, 1) x - 6 == 0,
An inequality (0, -1, -1) x + 5 >= 0,
An inequality (0, 0, -1) x + 3 >= 0,
An inequality (0, -1, 0) x + 3 >= 0,
An inequality (0, 1, 0) x - 1 >= 0,
An inequality (0, 1, 1) x - 3 >= 0,
An inequality (0, 0, 1) x - 1 >= 0)
"""
if self.H() is not None:
return self.equalities() + self.H()
else:
return tuple(smallInteger(i) for i in range(self.length_Hrepr()))
def dimension(self):
r"""
Return the dimension of the polyhedron.
EXAMPLES::
sage: C = CombinatorialPolyhedron([(1,2,3), (1,2,4),
....: (1,3,4), (2,3,4)])
sage: C.dimension()
3
sage: P = Polyhedron(rays=[[1,0,0],[0,1,0],[0,0,1],[0,0,-1]])
sage: CombinatorialPolyhedron(P).dimension()
3
"""
if self._dimension == -2:
# Dimension not computed yet.
if self.n_facets() == 0:
# The dimension of a trivial polyhedron is assumed to contain
# exactly one "vertex" and for each dimension one "line" as in
# :class:`~sage.geometry.polyhedron.parent.Polyhedron_base`
self._dimension = self.length_Vrepr() - 1
elif self.unbounded() or self.n_facets() <= self.length_Vrepr():
self._dimension = self.bitrep_facets().compute_dimension()
else:
# If the polyhedron has many facets,
# calculating the dimension of the dual will be faster.
# The dual exists, if the polyhedron is bounded.
self._dimension = self.bitrep_facets().compute_dimension()
return smallInteger(self._dimension)
def n_vertices(self):
r"""
Return the number of vertices.
Is equivalent to ``len(self.vertices())``.
EXAMPLES::
sage: P = polytopes.cube()
sage: C = CombinatorialPolyhedron(P)
sage: C.n_vertices()
8
sage: P = polytopes.cyclic_polytope(4,20)
sage: C = CombinatorialPolyhedron(P)
sage: C.n_vertices()
20
sage: P = Polyhedron(lines=[[0,1]], vertices=[[1,0], [-1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C.n_vertices()
0
sage: P = Polyhedron(rays=[[1,0,0], [0,1,0]], lines=[[0,0,1]])
sage: C = CombinatorialPolyhedron(P)
sage: C.n_vertices()
0
sage: C = CombinatorialPolyhedron(4)
sage: C.f_vector()
(1, 0, 0, 0, 0, 1)
sage: C.n_vertices()
0
sage: C = CombinatorialPolyhedron(0)
sage: C.f_vector()
(1, 1)
sage: C.n_vertices()
1
"""
if self.dimension() == 0:
# This specific trivial polyhedron needs special attention.
return smallInteger(1)
if self.unbounded():
# Some elements in the ``Vrepr`` might not correspond to actual combinatorial vertices.
return len(self.vertices())
else:
return smallInteger(self.length_Vrepr())
def vertices(self, names=True):
r"""
Return the elements in the Vrepresentation that are vertices.
In case of an unbounded polyhedron, there might be lines and
rays in the Vrepresentation.
If ``names`` is set to ``False``, then the vertices are given by
their indices in the Vrepresentation.
EXAMPLES::
sage: P = Polyhedron(rays=[[1,0,0],[0,1,0],[0,0,1]])
sage: C = CombinatorialPolyhedron(P)
sage: C.vertices()
(A vertex at (0, 0, 0),)
sage: C.Vrepresentation()
(A vertex at (0, 0, 0),
A ray in the direction (0, 0, 1),
A ray in the direction (0, 1, 0),
A ray in the direction (1, 0, 0))
sage: P = polytopes.cross_polytope(3)
sage: C = CombinatorialPolyhedron(P)
sage: C.vertices()
(A vertex at (-1, 0, 0),
A vertex at (0, -1, 0),
A vertex at (0, 0, -1),
A vertex at (0, 0, 1),
A vertex at (0, 1, 0),
A vertex at (1, 0, 0))
sage: C.vertices(names=False)
(0, 1, 2, 3, 4, 5)
sage: points = [(1,0,0), (0,1,0), (0,0,1),
....: (-1,0,0), (0,-1,0), (0,0,-1)]
sage: L = LatticePolytope(points)
sage: C = CombinatorialPolyhedron(L)
sage: C.vertices()
(M(1, 0, 0), M(0, 1, 0), M(0, 0, 1), M(-1, 0, 0), M(0, -1, 0), M(0, 0, -1))
sage: C.vertices(names=False)
(0, 1, 2, 3, 4, 5)
sage: P = Polyhedron(vertices=[[0,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C.vertices()
(A vertex at (0, 0),)
"""
if unlikely(self.dimension() == 0):
# Handling the case of a trivial polyhedron of dimension `0`.
if names and self.V():
return (self.V()[0],)
else:
return (smallInteger(0),)
if self.unbounded():
it = self.face_iter(0)
try:
# The Polyhedron has at least one vertex.
# In this case every element in the ``Vrepr``
# that is not contained in the far face
# is a vertex.
next(it)
except StopIteration:
# The Polyhedron has no vertex.
return ()
if names and self.V():
return tuple(self.V()[i] for i in range(self.length_Vrepr()) if not i in self.far_face_tuple())
else:
return tuple(smallInteger(i) for i in range(self.length_Vrepr()) if not i in self.far_face_tuple())
def n_facets(self):
r"""
Return the number of facets.
Is equivalent to ``len(self.facets())``.
EXAMPLES::
sage: P = polytopes.cube()
sage: C = CombinatorialPolyhedron(P)
sage: C.n_facets()
6
sage: P = polytopes.cyclic_polytope(4,20)
sage: C = CombinatorialPolyhedron(P)
sage: C.n_facets()
170
sage: P = Polyhedron(lines=[[0,1]], vertices=[[1,0], [-1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C.n_facets()
2
sage: P = Polyhedron(rays=[[1,0], [-1,0], [0,1]])
sage: C = CombinatorialPolyhedron(P)
sage: C.n_facets()
1
sage: C = CombinatorialPolyhedron(-1)
sage: C.f_vector()
(1)
sage: C.n_facets()
0
sage: C = CombinatorialPolyhedron(0)
sage: C.f_vector()
(1, 1)
sage: C.n_facets()
1
"""
if unlikely(self._dimension == 0):
# This trivial polyhedron needs special attention.
return smallInteger(1)
return smallInteger(self._n_facets)
def facets(self, names=True):
r"""
Return the facets as lists of ``[vertices, rays, lines]``.
If ``names`` is ``False``, then the Vrepresentatives in the facets
are given by their indices in the Vrepresentation.
EXAMPLES::
sage: P = polytopes.cube()
sage: C = CombinatorialPolyhedron(P)
sage: C.facets()
((A vertex at (-1, -1, 1),
A vertex at (-1, 1, 1),
A vertex at (1, -1, 1),
A vertex at (1, 1, 1)),
(A vertex at (-1, 1, -1),
A vertex at (-1, 1, 1),
A vertex at (1, 1, -1),
A vertex at (1, 1, 1)),
(A vertex at (1, -1, -1),
A vertex at (1, -1, 1),
A vertex at (1, 1, -1),
A vertex at (1, 1, 1)),
(A vertex at (-1, -1, -1),
A vertex at (-1, -1, 1),
A vertex at (-1, 1, -1),
A vertex at (-1, 1, 1)),
(A vertex at (-1, -1, -1),
A vertex at (-1, 1, -1),
A vertex at (1, -1, -1),
A vertex at (1, 1, -1)),
(A vertex at (-1, -1, -1),
A vertex at (-1, -1, 1),
A vertex at (1, -1, -1),
A vertex at (1, -1, 1)))
sage: C.facets(names=False)
((1, 3, 5, 7),
(2, 3, 6, 7),
(4, 5, 6, 7),
(0, 1, 2, 3),
(0, 2, 4, 6),
(0, 1, 4, 5))
"""
if unlikely(self.dimension() == 0):
# Special attention for this trivial case.
# There is actually one facet, but we have not initialized it.
return ((),)
# It is essential to have the facets in the exact same order as
# on input, so that pickle/unpickle by :meth:`reduce` works.
# Every facet knows its index by the facet representation.
face_iter = self.face_iter(self.dimension() - 1, dual=False)
facets = [None] * self.n_facets()
for face in face_iter:
index = face.Hrepr(names=False)[0]
verts = face.Vrepr(names=names)
facets[index] = verts
return tuple(facets)
def edges(self, names=True):
r"""
Return the edges of the polyhedron, i.e. the rank 1 faces.
If ``names`` is set to ``False``, then the Vrepresentatives in the edges
are given by their indices in the Vrepresentation.
.. NOTE::
To compute edges and f_vector, first compute the edges.
This might be faster.
EXAMPLES::
sage: P = polytopes.cyclic_polytope(3,5)
sage: C = CombinatorialPolyhedron(P)
sage: C.edges()
((A vertex at (3, 9, 27), A vertex at (4, 16, 64)),
(A vertex at (2, 4, 8), A vertex at (4, 16, 64)),
(A vertex at (1, 1, 1), A vertex at (4, 16, 64)),
(A vertex at (0, 0, 0), A vertex at (4, 16, 64)),
(A vertex at (2, 4, 8), A vertex at (3, 9, 27)),
(A vertex at (0, 0, 0), A vertex at (3, 9, 27)),
(A vertex at (1, 1, 1), A vertex at (2, 4, 8)),
(A vertex at (0, 0, 0), A vertex at (2, 4, 8)),
(A vertex at (0, 0, 0), A vertex at (1, 1, 1)))
sage: C.edges(names=False)
((3, 4), (2, 4), (1, 4), (0, 4), (2, 3), (0, 3), (1, 2), (0, 2), (0, 1))
sage: P = Polyhedron(rays=[[-1,0],[1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C.edges()
((A line in the direction (1, 0), A vertex at (0, 0)),)
sage: P = Polyhedron(vertices=[[0,0],[1,0]])
sage: C = CombinatorialPolyhedron(P)
sage: C.edges()
((A vertex at (0, 0), A vertex at (1, 0)),)
sage: from itertools import combinations
sage: N = combinations(['a','b','c','d','e'], 4)
sage: C = CombinatorialPolyhedron(N)
sage: C.edges()
(('d', 'e'),
('c', 'e'),
('b', 'e'),
('a', 'e'),
('c', 'd'),
('b', 'd'),
('a', 'd'),
('b', 'c'),
('a', 'c'),
('a', 'b'))
"""
cdef size_t len_edge_list = self._length_edges_list
if self._edges is NULL:
# compute the edges.
if self.unbounded():
self._compute_edges(dual=False)
elif self.length_Vrepr() > self.n_facets()*self.n_facets():
# This is a wild estimate
# that in this case it is better not to use the dual.
self._compute_edges(dual=False)
else:
# In most bounded cases, one should use the dual.
self._compute_ridges(dual=True)
if self._edges is NULL:
raise ValueError('could not determine edges')
# The edges are being saved in a list basically.
# The first entry represents the first vertex of the first edge.
# The second entry represents the second vertex of that edge.
# The third entry represents the first vertex of the second edge etc.
# Possibly there are many edges.
# Hence, edges are stored in an array of arrays,
# with each array containing ``len_edge_list`` of edges.
# Mapping the indices of the Vrepr to the names, if requested.
if self.V() is not None and names is True:
def f(size_t i): return self.V()[i]
else:
def f(size_t i): return smallInteger(i)
# Getting the indices of the `i`-th edge.
def vertex_one(size_t i):
return f(self._edges[i // len_edge_list][2*(i % len_edge_list)])
def vertex_two(size_t i):
return f(self._edges[i // len_edge_list][2*(i % len_edge_list)+1])
cdef size_t j
return tuple((vertex_one(j), vertex_two(j)) for j in range(self._n_edges))
def vertex_graph(self, names=True):
r"""
Return a graph in which the vertices correspond to vertices
of the polyhedron, and edges to bounded rank 1 faces.
If ``names`` is set to ``False``, the Vrepresentatives will
carry names according to the indexing of the Vrepresentation.
EXAMPLES::
sage: P = polytopes.cyclic_polytope(3,5)
sage: C = CombinatorialPolyhedron(P)
sage: C.vertex_graph()
Graph on 5 vertices
sage: G = C.vertex_graph()
sage: sorted(G.degree())
[3, 3, 4, 4, 4]
sage: P = Polyhedron(rays=[[1]])
sage: C = CombinatorialPolyhedron(P)
sage: C.graph()
Graph on 1 vertex
"""
vertices = self.vertices(names=names)
# Getting the bounded edges.
edges = tuple(edge for edge in self.edges(names=names)
if edge[0] in vertices and edge[1] in vertices)
return Graph([vertices, edges], format="vertices_and_edges")
graph = vertex_graph
def edge_graph(self, names=True):
r"""
Return the edge graph.
If ``names`` is set to ``False``, the Vrepresentatives will
carry names according to the indexing of the Vrepresentation.
EXAMPLES::
sage: P = polytopes.cyclic_polytope(3,5)
sage: C = CombinatorialPolyhedron(P)
sage: C.edge_graph()
doctest:...: DeprecationWarning: the method edge_graph of CombinatorialPolyhedron is deprecated; use vertex_graph
See https://trac.sagemath.org/28603 for details.
Graph on 5 vertices
sage: G = C.edge_graph()
sage: sorted(G.degree())
[3, 3, 4, 4, 4]
"""
from sage.misc.superseded import deprecation