Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagelib
Path: blob/master/sage/geometry/integral_points.pyx
4066 views
1
r"""
2
Cython helper methods to compute integral points in polyhedra.
3
"""
4
5
#*****************************************************************************
6
# Copyright (C) 2010 Volker Braun <[email protected]>
7
#
8
# Distributed under the terms of the GNU General Public License (GPL)
9
#
10
# http://www.gnu.org/licenses/
11
#*****************************************************************************
12
13
from sage.matrix.constructor import matrix, column_matrix, vector, diagonal_matrix
14
from sage.rings.all import QQ, RR, ZZ, gcd, lcm
15
from sage.combinat.permutation import Permutation
16
from sage.combinat.cartesian_product import CartesianProduct
17
from sage.misc.all import prod
18
import copy
19
20
##############################################################################
21
# The basic idea to enumerate the lattice points in the parallelotope
22
# is to use the Smith normal form of the ray coordinate matrix to
23
# bring the lattice into a "nice" basis. Here is the straightforward
24
# implementation. Note that you do not need to reduce to the
25
# full-dimensional case, the Smith normal form takes care of that for
26
# you.
27
#
28
## def parallelotope_points(spanning_points, lattice):
29
## # compute points in the open parallelotope, see [BrunsKoch]
30
## R = matrix(spanning_points).transpose()
31
## D,U,V = R.smith_form()
32
## e = D.diagonal() # the elementary divisors
33
## d = prod(e) # the determinant
34
## u = U.inverse().columns() # generators for gp(semigroup)
35
##
36
## # "inverse" of the ray matrix as far as possible over ZZ
37
## # R*Rinv == diagonal_matrix([d]*D.ncols() + [0]*(D.nrows()-D.ncols()))
38
## # If R is full rank, this is Rinv = matrix(ZZ, R.inverse() * d)
39
## Dinv = D.transpose()
40
## for i in range(0,D.ncols()):
41
## Dinv[i,i] = d/D[i,i]
42
## Rinv = V * Dinv * U
43
##
44
## gens = []
45
## for b in CartesianProduct(*[ range(0,i) for i in e ]):
46
## # this is our generator modulo the lattice spanned by the rays
47
## gen_mod_rays = sum( b_i*u_i for b_i, u_i in zip(b,u) )
48
## q_times_d = Rinv * gen_mod_rays
49
## q_times_d = vector(ZZ,[ q_i % d for q_i in q_times_d ])
50
## gen = lattice(R*q_times_d / d)
51
## gen.set_immutable()
52
## gens.append(gen)
53
## assert(len(gens) == d)
54
## return tuple(gens)
55
#
56
# The problem with the naive implementation is that it is slow:
57
#
58
# 1. You can simplify some of the matrix multiplications
59
#
60
# 2. The inner loop keeps creating new matrices and vectors, which
61
# is slow. It needs to recycle objects. Instead of creating a new
62
# lattice point again and again, change the entries of an
63
# existing lattice point and then copy it!
64
65
66
def parallelotope_points(spanning_points, lattice):
67
r"""
68
Return integral points in the parallelotope starting at the origin
69
and spanned by the ``spanning_points``.
70
71
See :meth:`~ConvexRationalPolyhedralCone.semigroup_generators` for a description of the
72
algorithm.
73
74
INPUT:
75
76
- ``spanning_points`` -- a non-empty list of linearly independent
77
rays (`\ZZ`-vectors or :class:`toric lattice
78
<sage.geometry.toric_lattice.ToricLatticeFactory>` elements),
79
not necessarily primitive lattice points.
80
81
OUTPUT:
82
83
The tuple of all lattice points in the half-open parallelotope
84
spanned by the rays `r_i`,
85
86
.. math::
87
88
\mathop{par}(\{r_i\}) =
89
\sum_{0\leq a_i < 1} a_i r_i
90
91
By half-open parallelotope, we mean that the
92
points in the facets not meeting the origin are omitted.
93
94
EXAMPLES:
95
96
Note how the points on the outward-facing factes are omitted::
97
98
sage: from sage.geometry.integral_points import parallelotope_points
99
sage: rays = map(vector, [(2,0), (0,2)])
100
sage: parallelotope_points(rays, ZZ^2)
101
((0, 0), (1, 0), (0, 1), (1, 1))
102
103
The rays can also be toric lattice points::
104
105
sage: rays = map(ToricLattice(2), [(2,0), (0,2)])
106
sage: parallelotope_points(rays, ToricLattice(2))
107
(N(0, 0), N(1, 0), N(0, 1), N(1, 1))
108
109
A non-smooth cone::
110
111
sage: c = Cone([ (1,0), (1,2) ])
112
sage: parallelotope_points(c.rays(), c.lattice())
113
(N(0, 0), N(1, 1))
114
115
A ``ValueError`` is raised if the ``spanning_points`` are not
116
linearly independent::
117
118
sage: rays = map(ToricLattice(2), [(1,1)]*2)
119
sage: parallelotope_points(rays, ToricLattice(2))
120
Traceback (most recent call last):
121
...
122
ValueError: The spanning points are not linearly independent!
123
124
TESTS::
125
126
sage: rays = map(vector,[(-3, -2, -3, -2), (-2, -1, -8, 5), (1, 9, -7, -4), (-3, -1, -2, 2)])
127
sage: len(parallelotope_points(rays, ZZ^4))
128
967
129
"""
130
R = matrix(spanning_points).transpose()
131
e, d, VDinv = ray_matrix_normal_form(R)
132
points = loop_over_parallelotope_points(e, d, VDinv, R, lattice)
133
for p in points:
134
p.set_immutable()
135
return points
136
137
138
def ray_matrix_normal_form(R):
139
r"""
140
Compute the Smith normal form of the ray matrix for
141
:func:`parallelotope_points`.
142
143
INPUT:
144
145
- ``R`` -- `\ZZ`-matrix whose columns are the rays spanning the
146
parallelotope.
147
148
OUTPUT:
149
150
A tuple containing ``e``, ``d``, and ``VDinv``.
151
152
EXAMPLES::
153
154
sage: from sage.geometry.integral_points import ray_matrix_normal_form
155
sage: R = column_matrix(ZZ,[3,3,3])
156
sage: ray_matrix_normal_form(R)
157
([3], 3, [1])
158
"""
159
D,U,V = R.smith_form()
160
e = D.diagonal() # the elementary divisors
161
d = prod(e) # the determinant
162
if d==0:
163
raise ValueError('The spanning points are not linearly independent!')
164
Dinv = diagonal_matrix(ZZ,[ d/D[i,i] for i in range(0,D.ncols()) ])
165
VDinv = V * Dinv
166
return (e,d,VDinv)
167
168
169
170
# The optimized version avoids constructing new matrices, vectors, and lattice points
171
cpdef loop_over_parallelotope_points(e, d, VDinv, R, lattice, A=None, b=None):
172
r"""
173
The inner loop of :func:`parallelotope_points`.
174
175
INPUT:
176
177
See :meth:`parallelotope_points` for ``e``, ``d``, ``VDinv``, ``R``, ``lattice``.
178
179
- ``A``, ``b``: Either both ``None`` or a vector and number. If
180
present, only the parallelotope points satisfying `A x \leq b`
181
are returned.
182
183
OUTPUT:
184
185
The points of the half-open parallelotope as a tuple of lattice
186
points.
187
188
EXAMPLES:
189
190
sage: e = [3]
191
sage: d = prod(e)
192
sage: VDinv = matrix(ZZ, [[1]])
193
sage: R = column_matrix(ZZ,[3,3,3])
194
sage: lattice = ZZ^3
195
sage: from sage.geometry.integral_points import loop_over_parallelotope_points
196
sage: loop_over_parallelotope_points(e, d, VDinv, R, lattice)
197
((0, 0, 0), (1, 1, 1), (2, 2, 2))
198
199
sage: A = vector(ZZ, [1,0,0])
200
sage: b = 1
201
sage: loop_over_parallelotope_points(e, d, VDinv, R, lattice, A, b)
202
((0, 0, 0), (1, 1, 1))
203
"""
204
cdef int i, j
205
cdef int dim = VDinv.nrows()
206
cdef int ambient_dim = R.nrows()
207
gens = []
208
s = ZZ.zero() # summation variable
209
gen = lattice(0)
210
q_times_d = vector(ZZ, dim)
211
for base in CartesianProduct(*[ range(0,i) for i in e ]):
212
for i in range(0, dim):
213
s = ZZ.zero()
214
for j in range(0, dim):
215
s += VDinv[i,j] * base[j]
216
q_times_d[i] = s % d
217
for i in range(0, ambient_dim):
218
s = ZZ.zero()
219
for j in range(0, dim):
220
s += R[i,j] * q_times_d[j]
221
gen[i] = s / d
222
if A is not None:
223
s = ZZ.zero()
224
for i in range(0, ambient_dim):
225
s += A[i] * gen[i]
226
if s>b:
227
continue
228
gens.append(copy.copy(gen))
229
if A is None:
230
assert(len(gens) == d)
231
return tuple(gens)
232
233
234
235
##############################################################################
236
def simplex_points(vertices):
237
r"""
238
Return the integral points in a lattice simplex.
239
240
INPUT:
241
242
- ``vertices`` -- an iterable of integer coordinate vectors. The
243
indices of vertices that span the simplex under
244
consideration.
245
246
OUTPUT:
247
248
A tuple containing the integral point coordinates as `\ZZ`-vectors.
249
250
EXAMPLES::
251
252
sage: from sage.geometry.integral_points import simplex_points
253
sage: simplex_points([(1,2,3), (2,3,7), (-2,-3,-11)])
254
((-2, -3, -11), (0, 0, -2), (1, 2, 3), (2, 3, 7))
255
256
The simplex need not be full-dimensional::
257
258
sage: simplex = Polyhedron([(1,2,3,5), (2,3,7,5), (-2,-3,-11,5)])
259
sage: simplex_points(simplex.Vrepresentation())
260
((2, 3, 7, 5), (0, 0, -2, 5), (-2, -3, -11, 5), (1, 2, 3, 5))
261
262
sage: simplex_points([(2,3,7)])
263
((2, 3, 7),)
264
265
TESTS::
266
267
sage: v = [(1,0,7,-1), (-2,-2,4,-3), (-1,-1,-1,4), (2,9,0,-5), (-2,-1,5,1)]
268
sage: simplex = Polyhedron(v); simplex
269
A 4-dimensional polyhedron in QQ^4 defined as the convex hull of 5 vertices
270
sage: pts = simplex_points(simplex.Vrepresentation())
271
sage: len(pts)
272
49
273
sage: for p in pts: p.set_immutable()
274
sage: len(set(pts))
275
49
276
277
sage: all(simplex.contains(p) for p in pts)
278
True
279
280
sage: v = [(4,-1,-1,-1), (-1,4,-1,-1), (-1,-1,4,-1), (-1,-1,-1,4), (-1,-1,-1,-1)]
281
sage: P4mirror = Polyhedron(v); P4mirror
282
A 4-dimensional polyhedron in QQ^4 defined as the convex hull of 5 vertices
283
sage: len( simplex_points(P4mirror.Vrepresentation()) )
284
126
285
"""
286
rays = [ vector(ZZ, v) for v in vertices ]
287
if len(rays)==0:
288
return tuple()
289
origin = rays.pop()
290
origin.set_immutable()
291
if len(rays)==0:
292
return tuple([origin])
293
translate_points(rays, origin)
294
295
# Find equation Ax<=b that cuts out simplex from parallelotope
296
Rt = matrix(rays)
297
R = Rt.transpose()
298
if R.is_square():
299
b = abs(R.det())
300
A = R.solve_left(vector([b]*len(rays)))
301
else:
302
RtR = Rt * R
303
b = abs(RtR.det())
304
A = RtR.solve_left(vector([b]*len(rays))) * Rt
305
306
# e, d, VDinv = ray_matrix_normal_form(R)
307
# print origin
308
# print rays
309
# print parallelotope_points(rays, origin.parent())
310
# print 'A = ', A
311
# print 'b = ', b
312
313
e, d, VDinv = ray_matrix_normal_form(R)
314
lattice = origin.parent()
315
points = loop_over_parallelotope_points(e, d, VDinv, R, lattice, A, b) + tuple(rays)
316
translate_points(points, -origin)
317
return points
318
319
320
cdef translate_points(v_list, delta):
321
r"""
322
Add ``delta`` to each vector in ``v_list``.
323
"""
324
cdef int dim = delta.degree()
325
for v in v_list:
326
for i in range(0,dim):
327
v[i] -= delta[i]
328
329
330
331
##############################################################################
332
# For points with "small" coordinates (that is, fitting into a small
333
# rectangular bounding box) it is faster to naively enumerate the
334
# points. This saves the overhead of triangulating the polytope etc.
335
336
def rectangular_box_points(box_min, box_max, polyhedron=None):
337
r"""
338
Return the integral points in the lattice bounding box that are
339
also contained in the given polyhedron.
340
341
INPUT:
342
343
- ``box_min`` -- A list of integers. The minimal value for each
344
coordinate of the rectangular bounding box.
345
346
- ``box_max`` -- A list of integers. The maximal value for each
347
coordinate of the rectangular bounding box.
348
349
- ``polyhedron`` -- A :class:`~sage.geometry.polyhedra.Polyhedron`
350
or ``None`` (default).
351
352
OUTPUT:
353
354
A tuple containing the integral points of the rectangular box
355
spanned by `box_min` and `box_max` and that lie inside the
356
``polyhedron``. For sufficiently large bounding boxes, this are
357
all integral points of the polyhedron.
358
359
If no polyhedron is specified, all integral points of the
360
rectangular box are returned.
361
362
ALGORITHM:
363
364
This function implements the naive algorithm towards counting
365
integral points. Given min and max of vertex coordinates, it
366
iterates over all points in the bounding box and checks whether
367
they lie in the polyhedron. The following optimizations are
368
implemented:
369
370
* Cython: Use machine integers and optimizing C/C++ compiler
371
where possible, arbitrary precision integers where necessary.
372
Bounds checking, no compile time limits.
373
374
* Unwind inner loop (and next-to-inner loop):
375
376
.. math::
377
378
Ax\leq b
379
\quad \Leftrightarrow \quad
380
a_1 x_1 ~\leq~ b - \sum_{i=2}^d a_i x_i
381
382
so we only have to evaluate `a_1 * x_1` in the inner loop.
383
384
* Coordinates are permuted to make the longest box edge the
385
inner loop. The inner loop is optimized to run very fast, so
386
its best to do as much work as possible there.
387
388
* Continuously reorder inequalities and test the most
389
restrictive inequalities first.
390
391
* Use convexity and only find first and last allowed point in
392
the inner loop. The points in-between must be points of the
393
polyhedron, too.
394
395
EXAMPLES::
396
397
sage: from sage.geometry.integral_points import rectangular_box_points
398
sage: rectangular_box_points([0,0,0],[1,2,3])
399
((0, 0, 0), (0, 0, 1), (0, 0, 2), (0, 0, 3),
400
(0, 1, 0), (0, 1, 1), (0, 1, 2), (0, 1, 3),
401
(0, 2, 0), (0, 2, 1), (0, 2, 2), (0, 2, 3),
402
(1, 0, 0), (1, 0, 1), (1, 0, 2), (1, 0, 3),
403
(1, 1, 0), (1, 1, 1), (1, 1, 2), (1, 1, 3),
404
(1, 2, 0), (1, 2, 1), (1, 2, 2), (1, 2, 3))
405
406
sage: cell24 = polytopes.twenty_four_cell()
407
sage: rectangular_box_points([-1]*4, [1]*4, cell24)
408
((-1, 0, 0, 0), (0, -1, 0, 0), (0, 0, -1, 0), (0, 0, 0, -1),
409
(0, 0, 0, 0),
410
(0, 0, 0, 1), (0, 0, 1, 0), (0, 1, 0, 0), (1, 0, 0, 0))
411
sage: d = 3
412
sage: dilated_cell24 = d*cell24
413
sage: len( rectangular_box_points([-d]*4, [d]*4, dilated_cell24) )
414
305
415
416
sage: d = 6
417
sage: dilated_cell24 = d*cell24
418
sage: len( rectangular_box_points([-d]*4, [d]*4, dilated_cell24) )
419
3625
420
421
sage: polytope = Polyhedron([(-4,-3,-2,-1),(3,1,1,1),(1,2,1,1),(1,1,3,0),(1,3,2,4)])
422
sage: pts = rectangular_box_points([-4]*4, [4]*4, polytope); pts
423
((-4, -3, -2, -1), (-1, 0, 0, 1), (0, 1, 1, 1), (1, 1, 1, 1), (1, 1, 3, 0),
424
(1, 2, 1, 1), (1, 2, 2, 2), (1, 3, 2, 4), (2, 1, 1, 1), (3, 1, 1, 1))
425
sage: all(polytope.contains(p) for p in pts)
426
True
427
428
sage: set(map(tuple,pts)) == \
429
... set([(-4,-3,-2,-1),(3,1,1,1),(1,2,1,1),(1,1,3,0),(1,3,2,4),
430
... (0,1,1,1),(1,2,2,2),(-1,0,0,1),(1,1,1,1),(2,1,1,1)]) # computed with PALP
431
True
432
433
Long ints and non-integral polyhedra are explictly allowed::
434
435
sage: polytope = Polyhedron([[1], [10*pi.n()]], base_ring=RDF)
436
sage: len( rectangular_box_points([-100], [100], polytope) )
437
31
438
439
sage: halfplane = Polyhedron(ieqs=[(-1,1,0)])
440
sage: rectangular_box_points([0,-1+10^50], [0,1+10^50])
441
((0, 99999999999999999999999999999999999999999999999999),
442
(0, 100000000000000000000000000000000000000000000000000),
443
(0, 100000000000000000000000000000000000000000000000001))
444
sage: len( rectangular_box_points([0,-100+10^50], [1,100+10^50], halfplane) )
445
201
446
"""
447
assert len(box_min)==len(box_max)
448
cdef int d = len(box_min)
449
diameter = sorted([ (box_max[i]-box_min[i], i) for i in range(0,d) ], reverse=True)
450
diameter_value = [ x[0] for x in diameter ]
451
diameter_index = [ x[1] for x in diameter ]
452
453
sort_perm = Permutation([ i+1 for i in diameter_index])
454
orig_perm = sort_perm.inverse()
455
orig_perm_list = [ i-1 for i in orig_perm ]
456
box_min = sort_perm.action(box_min)
457
box_max = sort_perm.action(box_max)
458
inequalities = InequalityCollection(polyhedron, sort_perm, box_min, box_max)
459
460
v = vector(ZZ, d)
461
points = []
462
for p in loop_over_rectangular_box_points(box_min, box_max, inequalities, d):
463
# v = vector(ZZ, orig_perm.action(p)) # too slow
464
for i in range(0,d):
465
v[i] = p[orig_perm_list[i]]
466
points.append(copy.copy(v))
467
return tuple(points)
468
469
470
cdef loop_over_rectangular_box_points(box_min, box_max, inequalities, int d):
471
"""
472
The inner loop of :func:`rectangular_box_points`.
473
474
INPUT:
475
476
- ``box_min``, ``box_max`` -- the bounding box.
477
478
- ``inequalities`` -- a :class:`InequalityCollection` containing
479
the inequalities defining the polyhedron.
480
481
- ``d`` -- the ambient space dimension.
482
483
OUTPUT:
484
485
The integral points in the bounding box satisfying all
486
inequalities.
487
"""
488
cdef int inc
489
points = []
490
p = copy.copy(box_min)
491
inequalities.prepare_next_to_inner_loop(p)
492
while True:
493
inequalities.prepare_inner_loop(p)
494
i_min = box_min[0]
495
i_max = box_max[0]
496
# Find the lower bound for the allowed region
497
while i_min <= i_max:
498
if inequalities.are_satisfied(i_min):
499
break
500
i_min += 1
501
# Find the upper bound for the allowed region
502
while i_min <= i_max:
503
if inequalities.are_satisfied(i_max):
504
break
505
i_max -= 1
506
# The points i_min .. i_max are contained in the polyhedron
507
i = i_min
508
while i <= i_max:
509
p[0] = i
510
points.append(tuple(p))
511
i += 1
512
# finally increment the other entries in p to move on to next inner loop
513
inc = 1
514
if d==1: return points
515
while True:
516
if p[inc]==box_max[inc]:
517
p[inc] = box_min[inc]
518
inc += 1
519
if inc==d:
520
return points
521
else:
522
p[inc] += 1
523
break
524
if inc>1:
525
inequalities.prepare_next_to_inner_loop(p)
526
527
528
529
cdef class Inequality_generic:
530
"""
531
An inequality whose coefficients are arbitrary Python/Sage objects
532
533
INPUT:
534
535
- ``A`` -- list of integers.
536
537
- ``b`` -- integer
538
539
OUTPUT:
540
541
Inequality `A x + b \geq 0`.
542
543
EXAMPLES::
544
545
sage: from sage.geometry.integral_points import Inequality_generic
546
sage: Inequality_generic([2*pi,sqrt(3),7/2], -5.5)
547
generic: (2*pi, sqrt(3), 7/2) x + -5.50000000000000 >= 0
548
"""
549
550
cdef object A
551
cdef object b
552
cdef object coeff
553
cdef object cache
554
555
def __cinit__(self, A, b):
556
"""
557
The Cython constructor
558
559
INPUT:
560
561
See :class:`Inequality_generic`.
562
563
EXAMPLES::
564
565
sage: from sage.geometry.integral_points import Inequality_generic
566
sage: Inequality_generic([2*pi,sqrt(3),7/2], -5.5)
567
generic: (2*pi, sqrt(3), 7/2) x + -5.50000000000000 >= 0
568
"""
569
self.A = A
570
self.b = b
571
self.coeff = 0
572
self.cache = 0
573
574
def __repr__(self):
575
"""
576
Return a string representation.
577
578
OUTPUT:
579
580
String.
581
582
EXAMPLES::
583
584
sage: from sage.geometry.integral_points import Inequality_generic
585
sage: Inequality_generic([2,3,7], -5).__repr__()
586
'generic: (2, 3, 7) x + -5 >= 0'
587
"""
588
s = 'generic: ('
589
s += ', '.join([str(self.A[i]) for i in range(0,len(self.A))])
590
s += ') x + ' + str(self.b) + ' >= 0'
591
return s
592
593
cdef prepare_next_to_inner_loop(self, p):
594
"""
595
In :class:`Inequality_int` this method is used to peel of the
596
next-to-inner loop.
597
598
See :meth:`InequalityCollection.prepare_inner_loop` for more details.
599
"""
600
pass
601
602
cdef prepare_inner_loop(self, p):
603
"""
604
Peel off the inner loop.
605
606
See :meth:`InequalityCollection.prepare_inner_loop` for more details.
607
"""
608
cdef int j
609
self.coeff = self.A[0]
610
self.cache = self.b
611
for j in range(1,len(self.A)):
612
self.cache += self.A[j] * p[j]
613
614
cdef bint is_not_satisfied(self, inner_loop_variable):
615
r"""
616
Test the inequality, using the cached value from :meth:`prepare_inner_loop`
617
"""
618
return inner_loop_variable*self.coeff + self.cache < 0
619
620
621
622
DEF INEQ_INT_MAX_DIM = 20
623
624
cdef class Inequality_int:
625
"""
626
Fast version of inequality in the case that all coefficient fit
627
into machine ints.
628
629
INPUT:
630
631
- ``A`` -- list of integers.
632
633
- ``b`` -- integer
634
635
- ``max_abs_coordinates`` -- the maximum of the coordinates that
636
one wants to evalate the coordinates on. Used for overflow
637
checking.
638
639
OUTPUT:
640
641
Inequality `A x + b \geq 0`. A ``OverflowError`` is raised if a
642
machine integer is not long enough to hold the results. A
643
``ValueError`` is raised if some of the input is not integral.
644
645
EXAMPLES::
646
647
sage: from sage.geometry.integral_points import Inequality_int
648
sage: Inequality_int([2,3,7], -5, [10]*3)
649
integer: (2, 3, 7) x + -5 >= 0
650
651
sage: Inequality_int([1]*21, -5, [10]*21)
652
Traceback (most recent call last):
653
...
654
OverflowError: Dimension limit exceeded.
655
656
sage: Inequality_int([2,3/2,7], -5, [10]*3)
657
Traceback (most recent call last):
658
...
659
ValueError: Not integral.
660
661
sage: Inequality_int([2,3,7], -5.2, [10]*3)
662
Traceback (most recent call last):
663
...
664
ValueError: Not integral.
665
666
sage: Inequality_int([2,3,7], -5*10^50, [10]*3) # actual error message can differ between 32 and 64 bit
667
Traceback (most recent call last):
668
...
669
OverflowError: ...
670
"""
671
cdef int A[INEQ_INT_MAX_DIM]
672
cdef int b
673
cdef int dim
674
# the innermost coefficient
675
cdef int coeff
676
cdef int cache
677
# the next-to-innermost coefficient
678
cdef int coeff_next
679
cdef int cache_next
680
681
cdef int _to_int(self, x) except? -999:
682
if not x in ZZ: raise ValueError('Not integral.')
683
return int(x) # raises OverflowError in Cython if necessary
684
685
def __cinit__(self, A, b, max_abs_coordinates):
686
"""
687
The Cython constructor
688
689
See :class:`Inequality_int` for input.
690
691
EXAMPLES::
692
693
sage: from sage.geometry.integral_points import Inequality_int
694
sage: Inequality_int([2,3,7], -5, [10]*3)
695
integer: (2, 3, 7) x + -5 >= 0
696
"""
697
cdef int i
698
self.dim = self._to_int(len(A))
699
if self.dim<1 or self.dim>INEQ_INT_MAX_DIM:
700
raise OverflowError('Dimension limit exceeded.')
701
for i in range(0,self.dim):
702
self.A[i] = self._to_int(A[i])
703
self.b = self._to_int(b)
704
self.coeff = self.A[0]
705
if self.dim>0:
706
self.coeff_next = self.A[1]
707
# finally, make sure that there cannot be any overflow during the enumeration
708
self._to_int(sum( [ZZ(b)]+[ZZ(A[i])*ZZ(max_abs_coordinates[i]) for i in range(0,self.dim)] ))
709
710
def __repr__(self):
711
"""
712
Return a string representation.
713
714
OUTPUT:
715
716
String.
717
718
EXAMPLES::
719
720
sage: from sage.geometry.integral_points import Inequality_int
721
sage: Inequality_int([2,3,7], -5, [10]*3).__repr__()
722
'integer: (2, 3, 7) x + -5 >= 0'
723
"""
724
s = 'integer: ('
725
s += ', '.join([str(self.A[i]) for i in range(0,self.dim)])
726
s += ') x + ' + str(self.b) + ' >= 0'
727
return s
728
729
cdef prepare_next_to_inner_loop(Inequality_int self, p):
730
"""
731
Peel off the next-to-inner loop.
732
733
See :meth:`InequalityCollection.prepare_inner_loop` for more details.
734
"""
735
cdef int j
736
self.cache_next = self.b
737
for j in range(2,self.dim):
738
self.cache_next += self.A[j] * p[j]
739
740
cdef prepare_inner_loop(Inequality_int self, p):
741
"""
742
Peel off the inner loop.
743
744
See :meth:`InequalityCollection.prepare_inner_loop` for more details.
745
"""
746
cdef int j
747
if self.dim>1:
748
self.cache = self.cache_next + self.coeff_next*p[1]
749
else:
750
self.cache = self.cache_next
751
752
cdef bint is_not_satisfied(Inequality_int self, int inner_loop_variable):
753
return inner_loop_variable*self.coeff + self.cache < 0
754
755
756
757
cdef class InequalityCollection:
758
"""
759
A collection of inequalities.
760
761
INPUT:
762
763
- ``polyhedron`` -- a polyhedron defining the inequalities.
764
765
- ``permutation`` -- a permution of the coordinates. Will be used
766
to permute the coordinates of the inequality.
767
768
- ``box_min``, ``box_max`` -- the (not permuted) minimal and maximal
769
coordinates of the bounding box. Used for bounds checking.
770
771
EXAMPLES::
772
773
sage: from sage.geometry.integral_points import InequalityCollection
774
sage: P_QQ = Polyhedron(identity_matrix(3).columns() + [(-2, -1,-1)], base_ring=QQ)
775
sage: ieq = InequalityCollection(P_QQ, Permutation([1,2,3]), [0]*3,[1]*3); ieq
776
The collection of inequalities
777
integer: (3, -2, -2) x + 2 >= 0
778
integer: (-1, 4, -1) x + 1 >= 0
779
integer: (-1, -1, 4) x + 1 >= 0
780
integer: (-1, -1, -1) x + 1 >= 0
781
782
sage: P_RR = Polyhedron(identity_matrix(2).columns() + [(-2.7, -1)], base_ring=RDF)
783
sage: InequalityCollection(P_RR, Permutation([1,2]), [0]*2, [1]*2)
784
The collection of inequalities
785
integer: (-1, -1) x + 1 >= 0
786
generic: (-1.0, 3.7) x + 1.0 >= 0
787
generic: (1.0, -1.35) x + 1.35 >= 0
788
789
sage: line = Polyhedron(eqns=[(2,3,7)])
790
sage: InequalityCollection(line, Permutation([1,2]), [0]*2, [1]*2 )
791
The collection of inequalities
792
integer: (3, 7) x + 2 >= 0
793
integer: (-3, -7) x + -2 >= 0
794
795
TESTS::
796
797
sage: TestSuite(ieq).run(skip='_test_pickling')
798
"""
799
cdef object ineqs_int
800
cdef object ineqs_generic
801
802
def __repr__(self):
803
r"""
804
Return a string representation.
805
806
OUTPUT:
807
808
String.
809
810
EXAMPLES::
811
812
sage: from sage.geometry.integral_points import InequalityCollection
813
sage: line = Polyhedron(eqns=[(2,3,7)])
814
sage: InequalityCollection(line, Permutation([1,2]), [0]*2, [1]*2 ).__repr__()
815
'The collection of inequalities\ninteger: (3, 7) x + 2 >= 0\ninteger: (-3, -7) x + -2 >= 0'
816
"""
817
s = 'The collection of inequalities\n'
818
for ineq in self.ineqs_int:
819
s += str(<Inequality_int>ineq) + '\n'
820
for ineq in self.ineqs_generic:
821
s += str(<Inequality_generic>ineq) + '\n'
822
return s.strip()
823
824
def _make_A_b(self, Hrep_obj, permutation):
825
r"""
826
Return the coefficients and constant of the H-representation
827
object.
828
829
INPUT:
830
831
- ``Hrep_obj`` -- a H-representation object of the polyhedron.
832
833
- ``permutation`` -- the permutation of the coordinates to
834
apply to ``A``.
835
836
OUTPUT:
837
838
A pair ``(A,b)``.
839
840
EXAXMPLES::
841
842
sage: from sage.geometry.integral_points import InequalityCollection
843
sage: line = Polyhedron(eqns=[(2,3,7)])
844
sage: ieq = InequalityCollection(line, Permutation([1,2]), [0]*2, [1]*2 )
845
sage: ieq._make_A_b(line.Hrepresentation(0), Permutation([1,2]))
846
([3, 7], 2)
847
sage: ieq._make_A_b(line.Hrepresentation(0), Permutation([2,1]))
848
([7, 3], 2)
849
"""
850
v = list(Hrep_obj)
851
A = permutation.action(v[1:])
852
b = v[0]
853
try:
854
x = lcm([a.denominator() for a in A] + [b.denominator()])
855
A = [ a*x for a in A ]
856
b = b*x
857
except AttributeError:
858
pass
859
return (A,b)
860
861
def __cinit__(self, polyhedron, permutation, box_min, box_max):
862
"""
863
The Cython constructor
864
865
See the class documentation for the desrciption of the arguments.
866
867
EXAMPLES::
868
869
sage: from sage.geometry.integral_points import InequalityCollection
870
sage: line = Polyhedron(eqns=[(2,3,7)])
871
sage: InequalityCollection(line, Permutation([1,2]), [0]*2, [1]*2 )
872
The collection of inequalities
873
integer: (3, 7) x + 2 >= 0
874
integer: (-3, -7) x + -2 >= 0
875
"""
876
max_abs_coordinates = [ max(abs(c_min), abs(c_max))
877
for c_min, c_max in zip(box_min, box_max) ]
878
max_abs_coordinates = permutation.action(max_abs_coordinates)
879
self.ineqs_int = []
880
self.ineqs_generic = []
881
if not polyhedron:
882
return
883
for Hrep_obj in polyhedron.inequality_generator():
884
A, b = self._make_A_b(Hrep_obj, permutation)
885
try:
886
H = Inequality_int(A, b, max_abs_coordinates)
887
self.ineqs_int.append(H)
888
except (OverflowError, ValueError):
889
H = Inequality_generic(A, b)
890
self.ineqs_generic.append(H)
891
for Hrep_obj in polyhedron.equation_generator():
892
A, b = self._make_A_b(Hrep_obj, permutation)
893
# add inequality
894
try:
895
H = Inequality_int(A, b, max_abs_coordinates)
896
self.ineqs_int.append(H)
897
except (OverflowError, ValueError):
898
H = Inequality_generic(A, b)
899
self.ineqs_generic.append(H)
900
# add sign-reversed inequality
901
A = [ -a for a in A ]
902
b = -b
903
try:
904
H = Inequality_int(A, b, max_abs_coordinates)
905
self.ineqs_int.append(H)
906
except (OverflowError, ValueError):
907
H = Inequality_generic(A, b)
908
self.ineqs_generic.append(H)
909
910
def prepare_next_to_inner_loop(self, p):
911
r"""
912
Peel off the next-to-inner loop.
913
914
In the next-to-inner loop of :func:`rectangular_box_points`,
915
we have to repeatedly evaluate `A x-A_0 x_0+b`. To speed up
916
computation, we pre-evaluate
917
918
.. math::
919
920
c = b + sum_{i=2} A_i x_i
921
922
and only compute `A x-A_0 x_0+b = A_1 x_1 +c \geq 0` in the
923
next-to-inner loop.
924
925
INPUT:
926
927
- ``p`` -- the point coordinates. Only ``p[2:]`` coordinates
928
are potentially used by this method.
929
930
EXAMPLES::
931
932
sage: from sage.geometry.integral_points import InequalityCollection, print_cache
933
sage: P = Polyhedron(ieqs=[(2,3,7,11)])
934
sage: ieq = InequalityCollection(P, Permutation([1,2,3]), [0]*3,[1]*3); ieq
935
The collection of inequalities
936
integer: (3, 7, 11) x + 2 >= 0
937
sage: ieq.prepare_next_to_inner_loop([2,1,3])
938
sage: ieq.prepare_inner_loop([2,1,3])
939
sage: print_cache(ieq)
940
Cached inner loop: 3 * x_0 + 42 >= 0
941
Cached next-to-inner loop: 3 * x_0 + 7 * x_1 + 35 >= 0
942
"""
943
for ineq in self.ineqs_int:
944
(<Inequality_int>ineq).prepare_next_to_inner_loop(p)
945
for ineq in self.ineqs_generic:
946
(<Inequality_generic>ineq).prepare_next_to_inner_loop(p)
947
948
def prepare_inner_loop(self, p):
949
r"""
950
Peel off the inner loop.
951
952
In the inner loop of :func:`rectangular_box_points`, we have
953
to repeatedly evaluate `A x+b\geq 0`. To speed up computation, we pre-evaluate
954
955
.. math::
956
957
c = A x - A_0 x_0 +b = b + sum_{i=1} A_i x_i
958
959
and only test `A_0 x_0 +c \geq 0` in the inner loop.
960
961
You must call :meth:`prepare_next_to_inner_loop` before
962
calling this method.
963
964
INPUT:
965
966
- ``p`` -- the coordinates of the point to loop over. Only the
967
``p[1:]`` entries are used.
968
969
EXAMPLES::
970
971
sage: from sage.geometry.integral_points import InequalityCollection, print_cache
972
sage: P = Polyhedron(ieqs=[(2,3,7,11)])
973
sage: ieq = InequalityCollection(P, Permutation([1,2,3]), [0]*3,[1]*3); ieq
974
The collection of inequalities
975
integer: (3, 7, 11) x + 2 >= 0
976
sage: ieq.prepare_next_to_inner_loop([2,1,3])
977
sage: ieq.prepare_inner_loop([2,1,3])
978
sage: print_cache(ieq)
979
Cached inner loop: 3 * x_0 + 42 >= 0
980
Cached next-to-inner loop: 3 * x_0 + 7 * x_1 + 35 >= 0
981
"""
982
for ineq in self.ineqs_int:
983
(<Inequality_int>ineq).prepare_inner_loop(p)
984
for ineq in self.ineqs_generic:
985
(<Inequality_generic>ineq).prepare_inner_loop(p)
986
987
def swap_ineq_to_front(self, i):
988
r"""
989
Swap the ``i``-th entry of the list to the front of the list of inequalities.
990
991
INPUT:
992
993
- ``i`` -- Integer. The :class:`Inequality_int` to swap to the
994
beginnig of the list of integral inequalities.
995
996
EXAMPLES::
997
998
sage: from sage.geometry.integral_points import InequalityCollection
999
sage: P_QQ = Polyhedron(identity_matrix(3).columns() + [(-2, -1,-1)], base_ring=QQ)
1000
sage: iec = InequalityCollection(P_QQ, Permutation([1,2,3]), [0]*3,[1]*3)
1001
sage: iec
1002
The collection of inequalities
1003
integer: (3, -2, -2) x + 2 >= 0
1004
integer: (-1, 4, -1) x + 1 >= 0
1005
integer: (-1, -1, 4) x + 1 >= 0
1006
integer: (-1, -1, -1) x + 1 >= 0
1007
sage: iec.swap_ineq_to_front(3)
1008
sage: iec
1009
The collection of inequalities
1010
integer: (-1, -1, -1) x + 1 >= 0
1011
integer: (3, -2, -2) x + 2 >= 0
1012
integer: (-1, 4, -1) x + 1 >= 0
1013
integer: (-1, -1, 4) x + 1 >= 0
1014
"""
1015
i_th_entry = self.ineqs_int.pop(i)
1016
self.ineqs_int.insert(0, i_th_entry)
1017
1018
def are_satisfied(self, inner_loop_variable):
1019
r"""
1020
Return whether all inequalities are satisfied.
1021
1022
You must call :meth:`prepare_inner_loop` before calling this
1023
method.
1024
1025
INPUT:
1026
1027
- ``inner_loop_variable`` -- Integer. the 0-th coordinate of
1028
the lattice point.
1029
1030
OUTPUT:
1031
1032
Boolean. Whether the lattice point is in the polyhedron.
1033
1034
EXAMPLES::
1035
1036
sage: from sage.geometry.integral_points import InequalityCollection
1037
sage: line = Polyhedron(eqns=[(2,3,7)])
1038
sage: ieq = InequalityCollection(line, Permutation([1,2]), [0]*2, [1]*2 )
1039
sage: ieq.prepare_next_to_inner_loop([3,4])
1040
sage: ieq.prepare_inner_loop([3,4])
1041
sage: ieq.are_satisfied(3)
1042
False
1043
"""
1044
cdef int i
1045
for i in range(0,len(self.ineqs_int)):
1046
ineq = self.ineqs_int[i]
1047
if (<Inequality_int>ineq).is_not_satisfied(inner_loop_variable):
1048
if i>0:
1049
self.swap_ineq_to_front(i)
1050
return False
1051
for i in range(0,len(self.ineqs_generic)):
1052
ineq = self.ineqs_generic[i]
1053
if (<Inequality_generic>ineq).is_not_satisfied(inner_loop_variable):
1054
return False
1055
return True
1056
1057
1058
1059
1060
cpdef print_cache(InequalityCollection inequality_collection):
1061
r"""
1062
Print the cached values in :class:`Inequality_int` (for
1063
debugging/doctesting only).
1064
1065
EXAMPLES::
1066
1067
sage: from sage.geometry.integral_points import InequalityCollection, print_cache
1068
sage: P = Polyhedron(ieqs=[(2,3,7)])
1069
sage: ieq = InequalityCollection(P, Permutation([1,2]), [0]*2,[1]*2); ieq
1070
The collection of inequalities
1071
integer: (3, 7) x + 2 >= 0
1072
sage: ieq.prepare_next_to_inner_loop([3,5])
1073
sage: ieq.prepare_inner_loop([3,5])
1074
sage: print_cache(ieq)
1075
Cached inner loop: 3 * x_0 + 37 >= 0
1076
Cached next-to-inner loop: 3 * x_0 + 7 * x_1 + 2 >= 0
1077
"""
1078
cdef Inequality_int ieq = <Inequality_int>(inequality_collection.ineqs_int[0])
1079
print 'Cached inner loop: ' + \
1080
str(ieq.coeff) + ' * x_0 + ' + str(ieq.cache) + ' >= 0'
1081
print 'Cached next-to-inner loop: ' + \
1082
str(ieq.coeff) + ' * x_0 + ' + \
1083
str(ieq.coeff_next) + ' * x_1 + ' + str(ieq.cache_next) + ' >= 0'
1084
1085
1086
1087