Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

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

#cython: wraparound=False, boundscheck=False 

r""" 

Cython helper methods to compute integral points in polyhedra. 

""" 

  

#***************************************************************************** 

# Copyright (C) 2010 Volker Braun <vbraun.name@gmail.com> 

# 

# 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. 

# http://www.gnu.org/licenses/ 

#***************************************************************************** 

from __future__ import print_function, absolute_import 

  

from cysignals.signals cimport sig_check 

import copy 

import itertools 

  

from sage.matrix.constructor import matrix, column_matrix, vector, diagonal_matrix 

from sage.rings.all import QQ, RR, ZZ 

from sage.rings.integer cimport Integer 

from sage.arith.all import gcd, lcm 

from sage.combinat.permutation import Permutation 

from sage.misc.all import prod, uniq 

from sage.modules.free_module import FreeModule 

from sage.modules.vector_integer_dense cimport Vector_integer_dense 

from sage.matrix.matrix_integer_dense cimport Matrix_integer_dense 

  

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

# The basic idea to enumerate the lattice points in the parallelotope 

# is to use the Smith normal form of the ray coordinate matrix to 

# bring the lattice into a "nice" basis. Here is the straightforward 

# implementation. Note that you do not need to reduce to the 

# full-dimensional case, the Smith normal form takes care of that for 

# you. 

# 

## def parallelotope_points(spanning_points, lattice): 

## # compute points in the open parallelotope, see [BK2001] 

## R = matrix(spanning_points).transpose() 

## D,U,V = R.smith_form() 

## e = D.diagonal() # the elementary divisors 

## d = prod(e) # the determinant 

## u = U.inverse().columns() # generators for gp(semigroup) 

## 

## # "inverse" of the ray matrix as far as possible over ZZ 

## # R*Rinv == diagonal_matrix([d]*D.ncols() + [0]*(D.nrows()-D.ncols())) 

## # If R is full rank, this is Rinv = matrix(ZZ, R.inverse() * d) 

## Dinv = D.transpose() 

## for i in range(0,D.ncols()): 

## Dinv[i,i] = d/D[i,i] 

## Rinv = V * Dinv * U 

## 

## gens = [] 

## for b in CartesianProduct(*[ range(0,i) for i in e ]): 

## # this is our generator modulo the lattice spanned by the rays 

## gen_mod_rays = sum( b_i*u_i for b_i, u_i in zip(b,u) ) 

## q_times_d = Rinv * gen_mod_rays 

## q_times_d = vector(ZZ,[ q_i % d for q_i in q_times_d ]) 

## gen = lattice(R*q_times_d / d) 

## gen.set_immutable() 

## gens.append(gen) 

## assert(len(gens) == d) 

## return tuple(gens) 

# 

# The problem with the naive implementation is that it is slow: 

# 

# 1. You can simplify some of the matrix multiplications 

# 

# 2. The inner loop keeps creating new matrices and vectors, which 

# is slow. It needs to recycle objects. Instead of creating a new 

# lattice point again and again, change the entries of an 

# existing lattice point and then copy it! 

  

  

cpdef tuple parallelotope_points(spanning_points, lattice): 

r""" 

Return integral points in the parallelotope starting at the origin 

and spanned by the ``spanning_points``. 

  

See :meth:`~ConvexRationalPolyhedralCone.semigroup_generators` for a description of the 

algorithm. 

  

INPUT: 

  

- ``spanning_points`` -- a non-empty list of linearly independent 

rays (`\ZZ`-vectors or :class:`toric lattice 

<sage.geometry.toric_lattice.ToricLatticeFactory>` elements), 

not necessarily primitive lattice points. 

  

OUTPUT: 

  

The tuple of all lattice points in the half-open parallelotope 

spanned by the rays `r_i`, 

  

.. MATH:: 

  

\mathop{par}(\{r_i\}) = 

\sum_{0\leq a_i < 1} a_i r_i 

  

By half-open parallelotope, we mean that the 

points in the facets not meeting the origin are omitted. 

  

EXAMPLES: 

  

Note how the points on the outward-facing factes are omitted:: 

  

sage: from sage.geometry.integral_points import parallelotope_points 

sage: rays = list(map(vector, [(2,0), (0,2)])) 

sage: parallelotope_points(rays, ZZ^2) 

((0, 0), (1, 0), (0, 1), (1, 1)) 

  

The rays can also be toric lattice points:: 

  

sage: rays = list(map(ToricLattice(2), [(2,0), (0,2)])) 

sage: parallelotope_points(rays, ToricLattice(2)) 

(N(0, 0), N(1, 0), N(0, 1), N(1, 1)) 

  

A non-smooth cone:: 

  

sage: c = Cone([ (1,0), (1,2) ]) 

sage: parallelotope_points(c.rays(), c.lattice()) 

(N(0, 0), N(1, 1)) 

  

A ``ValueError`` is raised if the ``spanning_points`` are not 

linearly independent:: 

  

sage: rays = list(map(ToricLattice(2), [(1,1)]*2)) 

sage: parallelotope_points(rays, ToricLattice(2)) 

Traceback (most recent call last): 

... 

ValueError: The spanning points are not linearly independent! 

  

TESTS:: 

  

sage: rays = list(map(vector,[(-3, -2, -3, -2), (-2, -1, -8, 5), (1, 9, -7, -4), (-3, -1, -2, 2)])) 

sage: len(parallelotope_points(rays, ZZ^4)) 

967 

""" 

cdef Matrix_integer_dense VDinv 

cdef Matrix_integer_dense R = matrix(spanning_points).transpose() 

e, d, VDinv = ray_matrix_normal_form(R) 

cdef tuple points = loop_over_parallelotope_points(e, d, VDinv, R, lattice) 

for p in points: 

p.set_immutable() 

return points 

  

  

cpdef tuple ray_matrix_normal_form(R): 

r""" 

Compute the Smith normal form of the ray matrix for 

:func:`parallelotope_points`. 

  

INPUT: 

  

- ``R`` -- `\ZZ`-matrix whose columns are the rays spanning the 

parallelotope. 

  

OUTPUT: 

  

A tuple containing ``e``, ``d``, and ``VDinv``. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import ray_matrix_normal_form 

sage: R = column_matrix(ZZ,[3,3,3]) 

sage: ray_matrix_normal_form(R) 

([3], 3, [1]) 

""" 

D,U,V = R.smith_form() 

e = D.diagonal() # the elementary divisors 

cdef Integer d = prod(e) # the determinant 

if d == ZZ.zero(): 

raise ValueError('The spanning points are not linearly independent!') 

cdef int i 

Dinv = diagonal_matrix(ZZ, [ d // e[i] for i in range(0,D.ncols()) ]) 

VDinv = V * Dinv 

return (e, d, VDinv) 

  

  

  

# The optimized version avoids constructing new matrices, vectors, and lattice points 

cpdef tuple loop_over_parallelotope_points(e, d, Matrix_integer_dense VDinv, 

Matrix_integer_dense R, lattice, 

A=None, b=None): 

r""" 

The inner loop of :func:`parallelotope_points`. 

  

INPUT: 

  

See :meth:`parallelotope_points` for ``e``, ``d``, ``VDinv``, ``R``, ``lattice``. 

  

- ``A``, ``b``: Either both ``None`` or a vector and number. If 

present, only the parallelotope points satisfying `A x \leq b` 

are returned. 

  

OUTPUT: 

  

The points of the half-open parallelotope as a tuple of lattice 

points. 

  

EXAMPLES:: 

  

sage: e = [3] 

sage: d = prod(e) 

sage: VDinv = matrix(ZZ, [[1]]) 

sage: R = column_matrix(ZZ,[3,3,3]) 

sage: lattice = ZZ^3 

sage: from sage.geometry.integral_points import loop_over_parallelotope_points 

sage: loop_over_parallelotope_points(e, d, VDinv, R, lattice) 

((0, 0, 0), (1, 1, 1), (2, 2, 2)) 

  

sage: A = vector(ZZ, [1,0,0]) 

sage: b = 1 

sage: loop_over_parallelotope_points(e, d, VDinv, R, lattice, A, b) 

((0, 0, 0), (1, 1, 1)) 

""" 

cdef int i, j 

cdef int dim = VDinv.nrows() 

cdef int ambient_dim = R.nrows() 

s = ZZ.zero() # summation variable 

cdef list gens = [] 

gen = lattice(ZZ.zero()) 

cdef Vector_integer_dense q_times_d = vector(ZZ, dim) 

for base in itertools.product(*[ range(0,i) for i in e ]): 

for i in range(0, dim): 

s = ZZ.zero() 

for j in range(0, dim): 

s += VDinv.get_unsafe(i,j) * base[j] 

q_times_d.set_unsafe(i, s % d) 

for i in range(0, ambient_dim): 

s = ZZ.zero() 

for j in range(0, dim): 

s += R.get_unsafe(i,j) * q_times_d.get_unsafe(j) 

gen[i] = s / d 

if A is not None: 

s = ZZ.zero() 

for i in range(0, ambient_dim): 

s += A[i] * gen[i] 

if s > b: 

continue 

gens.append(copy.copy(gen)) 

if A is None: 

assert(len(gens) == d) 

return tuple(gens) 

  

  

  

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

cpdef tuple simplex_points(vertices): 

r""" 

Return the integral points in a lattice simplex. 

  

INPUT: 

  

- ``vertices`` -- an iterable of integer coordinate vectors. The 

indices of vertices that span the simplex under 

consideration. 

  

OUTPUT: 

  

A tuple containing the integral point coordinates as `\ZZ`-vectors. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import simplex_points 

sage: simplex_points([(1,2,3), (2,3,7), (-2,-3,-11)]) 

((-2, -3, -11), (0, 0, -2), (1, 2, 3), (2, 3, 7)) 

  

The simplex need not be full-dimensional:: 

  

sage: simplex = Polyhedron([(1,2,3,5), (2,3,7,5), (-2,-3,-11,5)]) 

sage: simplex_points(simplex.Vrepresentation()) 

((2, 3, 7, 5), (0, 0, -2, 5), (-2, -3, -11, 5), (1, 2, 3, 5)) 

  

sage: simplex_points([(2,3,7)]) 

((2, 3, 7),) 

  

TESTS:: 

  

sage: v = [(1,0,7,-1), (-2,-2,4,-3), (-1,-1,-1,4), (2,9,0,-5), (-2,-1,5,1)] 

sage: simplex = Polyhedron(v); simplex 

A 4-dimensional polyhedron in ZZ^4 defined as the convex hull of 5 vertices 

sage: pts = simplex_points(simplex.Vrepresentation()) 

sage: len(pts) 

49 

sage: for p in pts: p.set_immutable() 

sage: len(set(pts)) 

49 

  

sage: all(simplex.contains(p) for p in pts) 

True 

  

sage: v = [(4,-1,-1,-1), (-1,4,-1,-1), (-1,-1,4,-1), (-1,-1,-1,4), (-1,-1,-1,-1)] 

sage: P4mirror = Polyhedron(v); P4mirror 

A 4-dimensional polyhedron in ZZ^4 defined as the convex hull of 5 vertices 

sage: len(simplex_points(P4mirror.Vrepresentation())) 

126 

  

sage: vertices = list(map(vector, [(1,2,3), (2,3,7), (-2,-3,-11)])) 

sage: for v in vertices: v.set_immutable() 

sage: simplex_points(vertices) 

((-2, -3, -11), (0, 0, -2), (1, 2, 3), (2, 3, 7)) 

""" 

cdef list rays = [vector(ZZ, list(v)) for v in vertices] 

if not rays: 

return () 

cdef Vector_integer_dense origin = rays.pop() 

origin.set_immutable() 

if not rays: 

return (origin,) 

translate_points(rays, origin) 

  

# Find equation Ax<=b that cuts out simplex from parallelotope 

cdef Matrix_integer_dense Rt = matrix(ZZ, rays) 

cdef Matrix_integer_dense R = Rt.transpose() 

cdef int nrays = len(rays) 

cdef Integer b 

M = FreeModule(ZZ, nrays) 

if R.is_square(): 

b = abs(R.det()) 

A = R.solve_left(M([b]*nrays)) 

else: 

RtR = Rt * R 

b = abs(RtR.det()) 

A = RtR.solve_left(M([b]*nrays)) * Rt 

  

e, d, VDinv = ray_matrix_normal_form(R) 

lattice = origin.parent() 

points = loop_over_parallelotope_points(e, d, VDinv, R, lattice, A, b) + tuple(rays) 

translate_points(points, -origin) 

return points 

  

  

cdef translate_points(v_list, Vector_integer_dense delta): 

r""" 

Add ``delta`` to each vector in ``v_list``. 

""" 

cdef int dim = delta.degree() 

cdef int i 

for v in v_list: 

for i in range(0,dim): 

v[i] -= delta.get_unsafe(i) 

  

  

  

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

# For points with "small" coordinates (that is, fitting into a small 

# rectangular bounding box) it is faster to naively enumerate the 

# points. This saves the overhead of triangulating the polytope etc. 

  

cpdef rectangular_box_points(list box_min, list box_max, 

polyhedron=None, count_only=False, 

return_saturated=False): 

r""" 

Return the integral points in the lattice bounding box that are 

also contained in the given polyhedron. 

  

INPUT: 

  

- ``box_min`` -- A list of integers. The minimal value for each 

coordinate of the rectangular bounding box. 

  

- ``box_max`` -- A list of integers. The maximal value for each 

coordinate of the rectangular bounding box. 

  

- ``polyhedron`` -- A 

:class:`~sage.geometry.polyhedron.base.Polyhedron_base`, a PPL 

:class:`~sage.libs.ppl.C_Polyhedron`, or ``None`` (default). 

  

- ``count_only`` -- Boolean (default: ``False``). Whether to 

return only the total number of vertices, and not their 

coordinates. Enabling this option speeds up the 

enumeration. Cannot be combined with the ``return_saturated`` 

option. 

  

- ``return_saturated`` -- Boolean (default: ``False``. Whether to 

also return which inequalities are saturated for each point of 

the polyhedron. Enabling this slows down the enumeration. Cannot 

be combined with the ``count_only`` option. 

  

OUTPUT: 

  

By default, this function returns a tuple containing the integral 

points of the rectangular box spanned by ``box_min`` and ``box_max`` 

and that lie inside the ``polyhedron``. For sufficiently large 

bounding boxes, this are all integral points of the polyhedron. 

  

If no polyhedron is specified, all integral points of the 

rectangular box are returned. 

  

If ``count_only`` is specified, only the total number (an integer) 

of found lattice points is returned. 

  

If ``return_saturated`` is enabled, then for each integral point a 

pair ``(point, Hrep)`` is returned where ``point`` is the point 

and ``Hrep`` is the set of indices of the H-representation objects 

that are saturated at the point. 

  

ALGORITHM: 

  

This function implements the naive algorithm towards counting 

integral points. Given min and max of vertex coordinates, it 

iterates over all points in the bounding box and checks whether 

they lie in the polyhedron. The following optimizations are 

implemented: 

  

* Cython: Use machine integers and optimizing C/C++ compiler 

where possible, arbitrary precision integers where necessary. 

Bounds checking, no compile time limits. 

  

* Unwind inner loop (and next-to-inner loop): 

  

.. MATH:: 

  

Ax \leq b 

\quad \Leftrightarrow \quad 

a_1 x_1 ~\leq~ b - \sum_{i=2}^d a_i x_i 

  

so we only have to evaluate `a_1 * x_1` in the inner loop. 

  

* Coordinates are permuted to make the longest box edge the 

inner loop. The inner loop is optimized to run very fast, so 

its best to do as much work as possible there. 

  

* Continuously reorder inequalities and test the most 

restrictive inequalities first. 

  

* Use convexity and only find first and last allowed point in 

the inner loop. The points in-between must be points of the 

polyhedron, too. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import rectangular_box_points 

sage: rectangular_box_points([0,0,0],[1,2,3]) 

((0, 0, 0), (0, 0, 1), (0, 0, 2), (0, 0, 3), 

(0, 1, 0), (0, 1, 1), (0, 1, 2), (0, 1, 3), 

(0, 2, 0), (0, 2, 1), (0, 2, 2), (0, 2, 3), 

(1, 0, 0), (1, 0, 1), (1, 0, 2), (1, 0, 3), 

(1, 1, 0), (1, 1, 1), (1, 1, 2), (1, 1, 3), 

(1, 2, 0), (1, 2, 1), (1, 2, 2), (1, 2, 3)) 

  

sage: from sage.geometry.integral_points import rectangular_box_points 

sage: rectangular_box_points([0,0,0],[1,2,3], count_only=True) 

24 

  

sage: cell24 = polytopes.twenty_four_cell() 

sage: rectangular_box_points([-1]*4, [1]*4, cell24) 

((-1, 0, 0, 0), (0, -1, 0, 0), (0, 0, -1, 0), (0, 0, 0, -1), 

(0, 0, 0, 0), 

(0, 0, 0, 1), (0, 0, 1, 0), (0, 1, 0, 0), (1, 0, 0, 0)) 

sage: d = 3 

sage: dilated_cell24 = d*cell24 

sage: len( rectangular_box_points([-d]*4, [d]*4, dilated_cell24) ) 

305 

  

sage: d = 6 

sage: dilated_cell24 = d*cell24 

sage: len( rectangular_box_points([-d]*4, [d]*4, dilated_cell24) ) 

3625 

  

sage: rectangular_box_points([-d]*4, [d]*4, dilated_cell24, count_only=True) 

3625 

  

sage: polytope = Polyhedron([(-4,-3,-2,-1),(3,1,1,1),(1,2,1,1),(1,1,3,0),(1,3,2,4)]) 

sage: pts = rectangular_box_points([-4]*4, [4]*4, polytope); pts 

((-4, -3, -2, -1), (-1, 0, 0, 1), (0, 1, 1, 1), (1, 1, 1, 1), (1, 1, 3, 0), 

(1, 2, 1, 1), (1, 2, 2, 2), (1, 3, 2, 4), (2, 1, 1, 1), (3, 1, 1, 1)) 

sage: all(polytope.contains(p) for p in pts) 

True 

  

sage: set(map(tuple,pts)) == \ 

....: set([(-4,-3,-2,-1),(3,1,1,1),(1,2,1,1),(1,1,3,0),(1,3,2,4), 

....: (0,1,1,1),(1,2,2,2),(-1,0,0,1),(1,1,1,1),(2,1,1,1)]) # computed with PALP 

True 

  

Long ints and non-integral polyhedra are explicitly allowed:: 

  

sage: polytope = Polyhedron([[1], [10*pi.n()]], base_ring=RDF) 

sage: len( rectangular_box_points([-100], [100], polytope) ) 

31 

  

sage: halfplane = Polyhedron(ieqs=[(-1,1,0)]) 

sage: rectangular_box_points([0,-1+10^50], [0,1+10^50]) 

((0, 99999999999999999999999999999999999999999999999999), 

(0, 100000000000000000000000000000000000000000000000000), 

(0, 100000000000000000000000000000000000000000000000001)) 

sage: len( rectangular_box_points([0,-100+10^50], [1,100+10^50], halfplane) ) 

201 

  

Using a PPL polyhedron:: 

  

sage: from sage.libs.ppl import Variable, Generator_System, C_Polyhedron, point 

sage: gs = Generator_System() 

sage: x = Variable(0); y = Variable(1); z = Variable(2) 

sage: gs.insert(point(0*x + 1*y + 0*z)) 

sage: gs.insert(point(0*x + 1*y + 3*z)) 

sage: gs.insert(point(3*x + 1*y + 0*z)) 

sage: gs.insert(point(3*x + 1*y + 3*z)) 

sage: poly = C_Polyhedron(gs) 

sage: rectangular_box_points([0]*3, [3]*3, poly) 

((0, 1, 0), (0, 1, 1), (0, 1, 2), (0, 1, 3), (1, 1, 0), (1, 1, 1), (1, 1, 2), (1, 1, 3), 

(2, 1, 0), (2, 1, 1), (2, 1, 2), (2, 1, 3), (3, 1, 0), (3, 1, 1), (3, 1, 2), (3, 1, 3)) 

  

Optionally, return the information about the saturated inequalities as well:: 

  

sage: cube = polytopes.cube() 

sage: cube.Hrepresentation(0) 

An inequality (0, 0, -1) x + 1 >= 0 

sage: cube.Hrepresentation(1) 

An inequality (0, -1, 0) x + 1 >= 0 

sage: cube.Hrepresentation(2) 

An inequality (-1, 0, 0) x + 1 >= 0 

sage: rectangular_box_points([0]*3, [1]*3, cube, return_saturated=True) 

(((0, 0, 0), frozenset()), 

((0, 0, 1), frozenset({0})), 

((0, 1, 0), frozenset({1})), 

((0, 1, 1), frozenset({0, 1})), 

((1, 0, 0), frozenset({2})), 

((1, 0, 1), frozenset({0, 2})), 

((1, 1, 0), frozenset({1, 2})), 

((1, 1, 1), frozenset({0, 1, 2}))) 

  

TESTS: 

  

Check that this can be interrupted, see :trac:`20781`:: 

  

sage: ieqs = [(-1, -1, -1, -1, -1, -1, -1, -1, -1), 

....: (0, -1, 0, 0, 0, 0, 0, 0, 0), 

....: (0, -1, 0, 2, -1, 0, 0, 0, 0), 

....: (0, 0, -1, -1, 2, -1, 0, 0, 0), 

....: (0, 2, 0, -1, 0, 0, 0, 0, 0), 

....: (0, 0, 0, 0, 0, 0, 0, -1, 2), 

....: (1, 0, 2, 0, -1, 0, 0, 0, 0), 

....: (0, 0, 0, 0, -1, 2, -1, 0, 0), 

....: (0, 0, 0, 0, 0, 0, 0, 0, -1), 

....: (0, 0, 0, 0, 0, -1, 2, -1, 0), 

....: (0, 0, 0, 0, 0, 0, -1, 2, -1)] 

sage: P = Polyhedron(ieqs=ieqs) 

sage: alarm(0.5); P.integral_points() 

Traceback (most recent call last): 

... 

AlarmInterrupt 

""" 

assert len(box_min) == len(box_max) 

assert not (count_only and return_saturated) 

cdef int d = len(box_min) 

cdef int i, j 

cdef list diameter = sorted([ (box_max[i]-box_min[i], i) for i in range(0,d) ], 

reverse=True) 

cdef list diameter_value = [x[0] for x in diameter] 

cdef list diameter_index = [x[1] for x in diameter] 

  

# Construct the inverse permutation 

cdef list orig_perm = list(xrange(len(diameter_index))) 

for i, j in enumerate(diameter_index): 

orig_perm[j] = i 

  

box_min = perm_action(diameter_index, box_min) 

box_max = perm_action(diameter_index, box_max) 

cdef InequalityCollection inequalities = InequalityCollection(polyhedron, diameter_index, box_min, box_max) 

  

if count_only: 

return loop_over_rectangular_box_points(box_min, box_max, inequalities, d, count_only) 

  

cdef list points = [] 

cdef Vector_integer_dense v = vector(ZZ, d) 

if not return_saturated: 

for p in loop_over_rectangular_box_points(box_min, box_max, inequalities, d, count_only): 

# v = vector(ZZ, perm_action(orig_perm, p)) # too slow 

for i in range(d): 

v.set_unsafe(i, Integer(p[orig_perm[i]])) 

v_copy = copy.copy(v) 

v_copy.set_immutable() 

points.append(v_copy) 

else: 

for p, saturated in loop_over_rectangular_box_points_saturated(box_min, box_max, 

inequalities, d): 

for i in range(d): 

v.set_unsafe(i, Integer(p[orig_perm[i]])) 

v_copy = copy.copy(v) 

v_copy.set_immutable() 

points.append( (v_copy, saturated) ) 

  

return tuple(points) 

  

cdef list perm_action(list p, list lst): 

""" 

Return the action of a permutation ``p`` of `(0, ..., n-1)` 

on a list of length `n`. 

""" 

return [lst[i] for i in p] 

  

cdef loop_over_rectangular_box_points(list box_min, list box_max, 

InequalityCollection inequalities, 

int d, bint count_only): 

""" 

The inner loop of :func:`rectangular_box_points`. 

  

INPUT: 

  

- ``box_min``, ``box_max`` -- the bounding box. 

  

- ``inequalities`` -- a :class:`InequalityCollection` containing 

the inequalities defining the polyhedron. 

  

- ``d`` -- the ambient space dimension. 

  

- ``count_only`` -- whether to only return the total number of 

lattice points. 

  

OUTPUT: 

  

The integral points in the bounding box satisfying all 

inequalities. 

""" 

cdef int inc 

cdef Integer i_min, i_max 

if count_only: 

points = 0 

else: 

points = [] 

cdef list p = list(box_min) 

inequalities.prepare_next_to_inner_loop(p) 

while True: 

sig_check() 

inequalities.prepare_inner_loop(p) 

i_min = box_min[0] 

i_max = box_max[0] 

# Find the lower bound for the allowed region 

while i_min <= i_max: 

if inequalities.are_satisfied(i_min): 

break 

i_min += 1 

# Find the upper bound for the allowed region 

while i_min <= i_max: 

if inequalities.are_satisfied(i_max): 

break 

i_max -= 1 

# The points i_min .. i_max are contained in the polyhedron 

if count_only: 

if i_max >= i_min: 

points += i_max - i_min + 1 

else: 

i = i_min 

while i <= i_max: 

p[0] = i 

points.append(tuple(p)) 

i += 1 

# finally increment the other entries in p to move on to next inner loop 

inc = 1 

if d == 1: 

return points 

while True: 

if p[inc] == box_max[inc]: 

p[inc] = box_min[inc] 

inc += 1 

if inc == d: 

return points 

else: 

p[inc] += 1 

break 

if inc > 1: 

inequalities.prepare_next_to_inner_loop(p) 

  

  

  

cdef loop_over_rectangular_box_points_saturated(list box_min, list box_max, 

InequalityCollection inequalities, 

int d): 

""" 

The analog of :func:`rectangular_box_points` except that it keeps 

track of which inequalities are saturated. 

  

INPUT: 

  

See :func:`rectangular_box_points`. 

  

OUTPUT: 

  

The integral points in the bounding box satisfying all 

inequalities, each point being returned by a coordinate vector and 

a set of H-representation object indices. 

""" 

cdef int inc 

cdef list points = [] 

cdef list p = list(box_min) 

inequalities.prepare_next_to_inner_loop(p) 

while True: 

inequalities.prepare_inner_loop(p) 

i_min = box_min[0] 

i_max = box_max[0] 

# Find the lower bound for the allowed region 

while i_min <= i_max: 

if inequalities.are_satisfied(i_min): 

break 

i_min += 1 

# Find the upper bound for the allowed region 

while i_min <= i_max: 

if inequalities.are_satisfied(i_max): 

break 

i_max -= 1 

# The points i_min .. i_max are contained in the polyhedron 

i = i_min 

while i <= i_max: 

p[0] = i 

saturated = inequalities.satisfied_as_equalities(i) 

points.append( (tuple(p), saturated) ) 

i += 1 

# finally increment the other entries in p to move on to next inner loop 

inc = 1 

if d == 1: 

return points 

while True: 

if p[inc] == box_max[inc]: 

p[inc] = box_min[inc] 

inc += 1 

if inc == d: 

return points 

else: 

p[inc] += 1 

break 

if inc > 1: 

inequalities.prepare_next_to_inner_loop(p) 

  

  

cdef class Inequality_generic: 

""" 

An inequality whose coefficients are arbitrary Python/Sage objects 

  

INPUT: 

  

- ``A`` -- list of coefficients 

  

- ``b`` -- element 

  

OUTPUT: 

  

Inequality `A x + b \geq 0`. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import Inequality_generic 

sage: Inequality_generic([2*pi,sqrt(3),7/2], -5.5) 

generic: (2*pi, sqrt(3), 7/2) x + -5.50000000000000 >= 0 

""" 

  

cdef list A 

cdef object b 

cdef object coeff 

cdef object cache 

# The index of the inequality in the polyhedron H-representation 

cdef int index 

  

def __cinit__(self, list A, b, int index=-1): 

""" 

The Cython constructor 

  

INPUT: 

  

See :class:`Inequality_generic`. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import Inequality_generic 

sage: Inequality_generic([2*pi,sqrt(3),7/2], -5.5) 

generic: (2*pi, sqrt(3), 7/2) x + -5.50000000000000 >= 0 

""" 

self.A = A 

self.b = b 

self.coeff = 0 

self.cache = 0 

self.index = int(index) 

  

def __repr__(self): 

""" 

Return a string representation. 

  

OUTPUT: 

  

String. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import Inequality_generic 

sage: Inequality_generic([2,3,7], -5).__repr__() 

'generic: (2, 3, 7) x + -5 >= 0' 

""" 

s = 'generic: (' 

s += ', '.join([str(self.A[i]) for i in range(0,len(self.A))]) 

s += ') x + ' + str(self.b) + ' >= 0' 

return s 

  

cdef prepare_next_to_inner_loop(self, p): 

""" 

In :class:`Inequality_int` this method is used to peel of the 

next-to-inner loop. 

  

See :meth:`InequalityCollection.prepare_inner_loop` for more details. 

""" 

pass 

  

cdef prepare_inner_loop(self, p): 

""" 

Peel off the inner loop. 

  

See :meth:`InequalityCollection.prepare_inner_loop` for more details. 

""" 

cdef int j 

self.coeff = self.A[0] 

self.cache = self.b 

for j in range(1, len(self.A)): 

self.cache += self.A[j] * p[j] 

  

cdef bint is_not_satisfied(self, inner_loop_variable) except -1: 

r""" 

Test the inequality, using the cached value from :meth:`prepare_inner_loop` 

  

OUTPUT: 

  

Boolean. Whether the inequality is not satisfied. 

""" 

return inner_loop_variable * self.coeff + self.cache < 0 

  

cdef bint is_equality(Inequality_generic self, int inner_loop_variable) except -1: 

r""" 

Test the inequality, using the cached value from :meth:`prepare_inner_loop` 

  

OUTPUT: 

  

Boolean. Given the inequality `Ax + b \geq 0`, this method 

returns whether the equality `Ax + b = 0` is satisfied. 

""" 

return inner_loop_variable * self.coeff + self.cache == 0 

  

  

# if dim>20 then we always use the generic inequalities (Inequality_generic) 

DEF INEQ_INT_MAX_DIM = 20 

  

cdef class Inequality_int: 

""" 

Fast version of inequality in the case that all coefficients fit 

into machine ints. 

  

INPUT: 

  

- ``A`` -- list of integers 

  

- ``b`` -- integer 

  

- ``max_abs_coordinates`` -- the maximum of the coordinates that 

one wants to evaluate the coordinates on; used for overflow 

checking 

  

OUTPUT: 

  

Inequality `A x + b \geq 0`. A ``OverflowError`` is raised if a 

machine integer is not long enough to hold the results. A 

``ValueError`` is raised if some of the input is not integral. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import Inequality_int 

sage: Inequality_int([2,3,7], -5, [10]*3) 

integer: (2, 3, 7) x + -5 >= 0 

  

sage: Inequality_int([1]*21, -5, [10]*21) 

Traceback (most recent call last): 

... 

OverflowError: Dimension limit exceeded. 

  

sage: Inequality_int([2,3/2,7], -5, [10]*3) 

Traceback (most recent call last): 

... 

ValueError: Not integral. 

  

sage: Inequality_int([2,3,7], -5.2, [10]*3) 

Traceback (most recent call last): 

... 

ValueError: Not integral. 

  

sage: Inequality_int([2,3,7], -5*10^50, [10]*3) # actual error message can differ between 32 and 64 bit 

Traceback (most recent call last): 

... 

OverflowError: ... 

  

TESTS: 

  

Check that :trac:`21993` is fixed:: 

  

sage: Inequality_int([18560500, -89466500], 108027, [178933, 37121]) 

Traceback (most recent call last): 

... 

OverflowError: ... 

  

""" 

cdef int A[INEQ_INT_MAX_DIM] 

cdef int b 

cdef int dim 

# the innermost coefficient 

cdef int coeff 

cdef int cache 

# the next-to-innermost coefficient 

cdef int coeff_next 

cdef int cache_next 

# The index of the inequality in the polyhedron H-representation 

cdef int index 

  

cdef int _to_int(self, x) except? -999: 

if not x in ZZ: 

raise ValueError('Not integral.') 

return int(x) # raises OverflowError in Cython if necessary 

  

def __cinit__(self, list A, b, list max_abs_coordinates, int index=-1): 

""" 

The Cython constructor 

  

See :class:`Inequality_int` for input. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import Inequality_int 

sage: Inequality_int([2,3,7], -5, [10]*3) 

integer: (2, 3, 7) x + -5 >= 0 

""" 

cdef int i 

self.dim = self._to_int(len(A)) 

self.index = int(index) 

if self.dim < 1 or self.dim > INEQ_INT_MAX_DIM: 

raise OverflowError('Dimension limit exceeded.') 

for i in range(self.dim): 

self.A[i] = self._to_int(A[i]) 

self.b = self._to_int(b) 

self.coeff = self.A[0] 

if self.dim > 0: 

self.coeff_next = self.A[1] 

# finally, make sure that there cannot be any overflow during the enumeration 

self._to_int(abs(ZZ(b)) + sum( abs(ZZ(A[i])) * ZZ(max_abs_coordinates[i]) 

for i in range(self.dim) )) 

  

def __repr__(self): 

""" 

Return a string representation. 

  

OUTPUT: 

  

String. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import Inequality_int 

sage: Inequality_int([2,3,7], -5, [10]*3).__repr__() 

'integer: (2, 3, 7) x + -5 >= 0' 

""" 

s = 'integer: (' 

s += ', '.join([str(self.A[i]) for i in range(self.dim)]) 

s += ') x + ' + str(self.b) + ' >= 0' 

return s 

  

cdef prepare_next_to_inner_loop(Inequality_int self, p): 

""" 

Peel off the next-to-inner loop. 

  

See :meth:`InequalityCollection.prepare_inner_loop` for more details. 

""" 

cdef int j 

self.cache_next = self.b 

for j in range(2, self.dim): 

self.cache_next += self.A[j] * p[j] 

  

cdef prepare_inner_loop(Inequality_int self, p): 

""" 

Peel off the inner loop. 

  

See :meth:`InequalityCollection.prepare_inner_loop` for more details. 

""" 

cdef int j 

if self.dim > 1: 

self.cache = self.cache_next + self.coeff_next * p[1] 

else: 

self.cache = self.cache_next 

  

cdef bint is_not_satisfied(Inequality_int self, int inner_loop_variable): 

return inner_loop_variable * self.coeff + self.cache < 0 

  

cdef bint is_equality(Inequality_int self, int inner_loop_variable): 

return inner_loop_variable * self.coeff + self.cache == 0 

  

  

  

cdef class InequalityCollection: 

""" 

A collection of inequalities. 

  

INPUT: 

  

- ``polyhedron`` -- a polyhedron defining the inequalities. 

  

- ``permutation`` -- list; a 0-based permutation of the coordinates. 

Will be used to permute the coordinates of the inequality. 

  

- ``box_min``, ``box_max`` -- the (not permuted) minimal and maximal 

coordinates of the bounding box. Used for bounds checking. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: P_QQ = Polyhedron(identity_matrix(3).columns() + [(-2, -1,-1)], base_ring=QQ) 

sage: ieq = InequalityCollection(P_QQ, [0,1,2], [0]*3,[1]*3); ieq 

The collection of inequalities 

integer: (3, -2, -2) x + 2 >= 0 

integer: (-1, 4, -1) x + 1 >= 0 

integer: (-1, -1, 4) x + 1 >= 0 

integer: (-1, -1, -1) x + 1 >= 0 

  

sage: P_RR = Polyhedron(identity_matrix(2).columns() + [(-2.7, -1)], base_ring=RDF) 

sage: InequalityCollection(P_RR, [0,1], [0]*2, [1]*2) 

The collection of inequalities 

integer: (-1, -1) x + 1 >= 0 

generic: (-1.0, 3.7) x + 1.0 >= 0 

generic: (1.0, -1.35) x + 1.35 >= 0 

  

sage: line = Polyhedron(eqns=[(2,3,7)]) 

sage: InequalityCollection(line, [0,1], [0]*2, [1]*2 ) 

The collection of inequalities 

integer: (3, 7) x + 2 >= 0 

integer: (-3, -7) x + -2 >= 0 

  

TESTS:: 

  

sage: TestSuite(ieq).run(skip='_test_pickling') 

""" 

cdef list ineqs_int 

cdef list ineqs_generic 

  

def __repr__(self): 

r""" 

Return a string representation. 

  

OUTPUT: 

  

String. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: line = Polyhedron(eqns=[(2,3,7)]) 

sage: InequalityCollection(line, [0,1], [0]*2, [1]*2 ).__repr__() 

'The collection of inequalities\ninteger: (3, 7) x + 2 >= 0\ninteger: (-3, -7) x + -2 >= 0' 

""" 

s = 'The collection of inequalities\n' 

for ineq in self.ineqs_int: 

s += str(<Inequality_int>ineq) + '\n' 

for ineq in self.ineqs_generic: 

s += str(<Inequality_generic>ineq) + '\n' 

return s.strip() 

  

cpdef tuple _make_A_b(self, Hrep_obj, list permutation): 

r""" 

Return the coefficients and constant of the H-representation 

object. 

  

INPUT: 

  

- ``Hrep_obj`` -- a H-representation object of the polyhedron 

- ``permutation`` -- the permutation of the coordinates to 

apply to ``A`` 

  

OUTPUT: 

  

A pair ``(A,b)``. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: line = Polyhedron(eqns=[(2,3,7)]) 

sage: ieq = InequalityCollection(line, [0,1], [0]*2, [1]*2 ) 

sage: ieq._make_A_b(line.Hrepresentation(0), [0,1]) 

([3, 7], 2) 

sage: ieq._make_A_b(line.Hrepresentation(0), [1,0]) 

([7, 3], 2) 

""" 

cdef list v = list(Hrep_obj) 

cdef list A = perm_action(permutation, v[1:]) 

b = v[0] 

try: 

x = lcm([a.denominator() for a in A] + [b.denominator()]) 

A = [a * x for a in A] 

b = b * x 

except AttributeError: 

pass 

return (A, b) 

  

def __cinit__(self, polyhedron, list permutation, box_min, box_max): 

""" 

The Cython constructor 

  

See the class documentation for the description of the arguments. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: line = Polyhedron(eqns=[(2,3,7)]) 

sage: InequalityCollection(line, [0,1], [0]*2, [1]*2 ) 

The collection of inequalities 

integer: (3, 7) x + 2 >= 0 

integer: (-3, -7) x + -2 >= 0 

""" 

cdef list max_abs_coordinates = [max(abs(c_min), abs(c_max)) 

for c_min, c_max in zip(box_min, box_max)] 

max_abs_coordinates = perm_action(permutation, max_abs_coordinates) 

self.ineqs_int = [] 

self.ineqs_generic = [] 

if polyhedron is None: 

return 

  

try: 

# polyhedron is a PPL C_Polyhedron class? 

self._cinit_from_PPL(max_abs_coordinates, permutation, polyhedron) 

except AttributeError: 

try: 

# polyhedron is a Polyhedron class? 

self._cinit_from_Polyhedron(max_abs_coordinates, permutation, polyhedron) 

except AttributeError: 

raise TypeError('Cannot extract Hrepresentation data from polyhedron.') 

  

cdef _cinit_from_PPL(self, list max_abs_coordinates, list permutation, 

polyhedron): 

""" 

Initialize the inequalities from a PPL C_Polyhedron 

  

See __cinit__ for a description of the arguments. 

  

EXAMPLES:: 

  

sage: from sage.libs.ppl import Variable, Generator_System, C_Polyhedron, point 

sage: gs = Generator_System() 

sage: x = Variable(0); y = Variable(1); z = Variable(2) 

sage: gs.insert(point(0*x + 0*y + 1*z)) 

sage: gs.insert(point(0*x + 3*y + 1*z)) 

sage: gs.insert(point(3*x + 0*y + 1*z)) 

sage: gs.insert(point(3*x + 3*y + 1*z)) 

sage: poly = C_Polyhedron(gs) 

sage: from sage.geometry.integral_points import InequalityCollection 

sage: InequalityCollection(poly, [0,2,1], [0]*3, [3]*3 ) 

The collection of inequalities 

integer: (0, 1, 0) x + -1 >= 0 

integer: (0, -1, 0) x + 1 >= 0 

integer: (1, 0, 0) x + 0 >= 0 

integer: (0, 0, 1) x + 0 >= 0 

integer: (-1, 0, 0) x + 3 >= 0 

integer: (0, 0, -1) x + 3 >= 0 

""" 

cdef list A 

cdef int index 

for index,c in enumerate(polyhedron.minimized_constraints()): 

A = perm_action(permutation, list(c.coefficients())) 

b = c.inhomogeneous_term() 

try: 

H = Inequality_int(A, b, max_abs_coordinates, index) 

self.ineqs_int.append(H) 

except (OverflowError, ValueError): 

H = Inequality_generic(A, b, index) 

self.ineqs_generic.append(H) 

if c.is_equality(): 

A = [ -a for a in A ] 

b = -b 

try: 

H = Inequality_int(A, b, max_abs_coordinates, index) 

self.ineqs_int.append(H) 

except (OverflowError, ValueError): 

H = Inequality_generic(A, b, index) 

self.ineqs_generic.append(H) 

  

cdef _cinit_from_Polyhedron(self, list max_abs_coordinates, 

list permutation, polyhedron): 

""" 

Initialize the inequalities from a Sage Polyhedron 

  

See __cinit__ for a description of the arguments. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: line = Polyhedron(eqns=[(2,3,7)]) 

sage: InequalityCollection(line, [0,1], [0]*2, [1]*2 ) 

The collection of inequalities 

integer: (3, 7) x + 2 >= 0 

integer: (-3, -7) x + -2 >= 0 

  

TESTS: 

  

Check that :trac:`21037` is fixed:: 

  

sage: P = Polyhedron(vertices=((0, 0), (17,3))) 

sage: P += 1/1000*polytopes.regular_polygon(5) 

sage: P.integral_points() 

((0, 0), (17, 3)) 

""" 

cdef list A 

for Hrep_obj in polyhedron.inequality_generator(): 

A, b = self._make_A_b(Hrep_obj, permutation) 

try: 

H = Inequality_int(A, b, max_abs_coordinates, Hrep_obj.index()) 

self.ineqs_int.append(H) 

except (OverflowError, ValueError, TypeError): 

H = Inequality_generic(A, b, Hrep_obj.index()) 

self.ineqs_generic.append(H) 

for Hrep_obj in polyhedron.equation_generator(): 

A, b = self._make_A_b(Hrep_obj, permutation) 

# add inequality 

try: 

H = Inequality_int(A, b, max_abs_coordinates, Hrep_obj.index()) 

self.ineqs_int.append(H) 

except (OverflowError, ValueError, TypeError): 

H = Inequality_generic(A, b, Hrep_obj.index()) 

self.ineqs_generic.append(H) 

# add sign-reversed inequality 

A = [ -a for a in A ] 

b = -b 

try: 

H = Inequality_int(A, b, max_abs_coordinates, Hrep_obj.index()) 

self.ineqs_int.append(H) 

except (OverflowError, ValueError, TypeError): 

H = Inequality_generic(A, b, Hrep_obj.index()) 

self.ineqs_generic.append(H) 

  

cpdef prepare_next_to_inner_loop(self, p): 

r""" 

Peel off the next-to-inner loop. 

  

In the next-to-inner loop of :func:`rectangular_box_points`, 

we have to repeatedly evaluate `A x-A_0 x_0+b`. To speed up 

computation, we pre-evaluate 

  

.. MATH:: 

  

c = b + \sum_{i=2} A_i x_i 

  

and only compute `A x-A_0 x_0+b = A_1 x_1 +c \geq 0` in the 

next-to-inner loop. 

  

INPUT: 

  

- ``p`` -- the point coordinates. Only ``p[2:]`` coordinates 

are potentially used by this method. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection, print_cache 

sage: P = Polyhedron(ieqs=[(2,3,7,11)]) 

sage: ieq = InequalityCollection(P, [0,1,2], [0]*3,[1]*3); ieq 

The collection of inequalities 

integer: (3, 7, 11) x + 2 >= 0 

sage: ieq.prepare_next_to_inner_loop([2,1,3]) 

sage: ieq.prepare_inner_loop([2,1,3]) 

sage: print_cache(ieq) 

Cached inner loop: 3 * x_0 + 42 >= 0 

Cached next-to-inner loop: 3 * x_0 + 7 * x_1 + 35 >= 0 

""" 

for ineq in self.ineqs_int: 

(<Inequality_int>ineq).prepare_next_to_inner_loop(p) 

for ineq in self.ineqs_generic: 

(<Inequality_generic>ineq).prepare_next_to_inner_loop(p) 

  

cpdef prepare_inner_loop(self, p): 

r""" 

Peel off the inner loop. 

  

In the inner loop of :func:`rectangular_box_points`, we have 

to repeatedly evaluate `A x+b\geq 0`. To speed up computation, we pre-evaluate 

  

.. MATH:: 

  

c = A x - A_0 x_0 +b = b + \sum_{i=1} A_i x_i 

  

and only test `A_0 x_0 +c \geq 0` in the inner loop. 

  

You must call :meth:`prepare_next_to_inner_loop` before 

calling this method. 

  

INPUT: 

  

- ``p`` -- the coordinates of the point to loop over. Only the 

``p[1:]`` entries are used. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection, print_cache 

sage: P = Polyhedron(ieqs=[(2,3,7,11)]) 

sage: ieq = InequalityCollection(P, [0,1,2], [0]*3,[1]*3); ieq 

The collection of inequalities 

integer: (3, 7, 11) x + 2 >= 0 

sage: ieq.prepare_next_to_inner_loop([2,1,3]) 

sage: ieq.prepare_inner_loop([2,1,3]) 

sage: print_cache(ieq) 

Cached inner loop: 3 * x_0 + 42 >= 0 

Cached next-to-inner loop: 3 * x_0 + 7 * x_1 + 35 >= 0 

""" 

for ineq in self.ineqs_int: 

(<Inequality_int>ineq).prepare_inner_loop(p) 

for ineq in self.ineqs_generic: 

(<Inequality_generic>ineq).prepare_inner_loop(p) 

  

cpdef swap_ineq_to_front(self, int i): 

r""" 

Swap the ``i``-th entry of the list to the front of the list of inequalities. 

  

INPUT: 

  

- ``i`` -- Integer. The :class:`Inequality_int` to swap to the 

beginning of the list of integral inequalities. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: P_QQ = Polyhedron(identity_matrix(3).columns() + [(-2, -1,-1)], base_ring=QQ) 

sage: iec = InequalityCollection(P_QQ, [0,1,2], [0]*3,[1]*3) 

sage: iec 

The collection of inequalities 

integer: (3, -2, -2) x + 2 >= 0 

integer: (-1, 4, -1) x + 1 >= 0 

integer: (-1, -1, 4) x + 1 >= 0 

integer: (-1, -1, -1) x + 1 >= 0 

sage: iec.swap_ineq_to_front(3) 

sage: iec 

The collection of inequalities 

integer: (-1, -1, -1) x + 1 >= 0 

integer: (3, -2, -2) x + 2 >= 0 

integer: (-1, 4, -1) x + 1 >= 0 

integer: (-1, -1, 4) x + 1 >= 0 

""" 

i_th_entry = self.ineqs_int[i] 

cdef int j 

for j in range(i-1,-1,-1): 

self.ineqs_int[j+1] = self.ineqs_int[j] 

self.ineqs_int[0] = i_th_entry 

  

cpdef bint are_satisfied(self, inner_loop_variable) except -1: 

r""" 

Return whether all inequalities are satisfied. 

  

You must call :meth:`prepare_inner_loop` before calling this 

method. 

  

INPUT: 

  

- ``inner_loop_variable`` -- Integer. the 0-th coordinate of 

the lattice point. 

  

OUTPUT: 

  

Boolean. Whether the lattice point is in the polyhedron. 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: line = Polyhedron(eqns=[(2,3,7)]) 

sage: ieq = InequalityCollection(line, [0,1], [0]*2, [1]*2 ) 

sage: ieq.prepare_next_to_inner_loop([3,4]) 

sage: ieq.prepare_inner_loop([3,4]) 

sage: ieq.are_satisfied(3) 

False 

""" 

cdef int i 

for i in range(len(self.ineqs_int)): 

sig_check() 

ineq = self.ineqs_int[i] 

if (<Inequality_int>ineq).is_not_satisfied(inner_loop_variable): 

if i > 0: 

self.swap_ineq_to_front(i) 

return False 

for i in range(len(self.ineqs_generic)): 

sig_check() 

ineq = self.ineqs_generic[i] 

if (<Inequality_generic>ineq).is_not_satisfied(inner_loop_variable): 

return False 

return True 

  

cpdef frozenset satisfied_as_equalities(self, inner_loop_variable): 

""" 

Return the inequalities (by their index) that are satisfied as 

equalities. 

  

INPUT: 

  

- ``inner_loop_variable`` -- Integer. the 0-th coordinate of 

the lattice point. 

  

OUTPUT: 

  

A set of integers in ascending order. Each integer is the 

index of a H-representation object of the polyhedron (either a 

inequality or an equation). 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection 

sage: quadrant = Polyhedron(rays=[(1,0), (0,1)]) 

sage: ieqs = InequalityCollection(quadrant, [0,1], [-1]*2, [1]*2 ) 

sage: ieqs.prepare_next_to_inner_loop([-1,0]) 

sage: ieqs.prepare_inner_loop([-1,0]) 

sage: ieqs.satisfied_as_equalities(-1) 

frozenset({1}) 

sage: ieqs.satisfied_as_equalities(0) 

frozenset({0, 1}) 

sage: ieqs.satisfied_as_equalities(1) 

frozenset({1}) 

""" 

cdef int i 

cdef list result = [] 

for i in range(len(self.ineqs_int)): 

sig_check() 

ineq = self.ineqs_int[i] 

if (<Inequality_int>ineq).is_equality(inner_loop_variable): 

result.append( (<Inequality_int>ineq).index ) 

for i in range(len(self.ineqs_generic)): 

sig_check() 

ineq = self.ineqs_generic[i] 

if (<Inequality_generic>ineq).is_equality(inner_loop_variable): 

result.append( (<Inequality_generic>ineq).index ) 

return frozenset(result) 

  

  

  

cpdef print_cache(InequalityCollection inequality_collection): 

r""" 

Print the cached values in :class:`Inequality_int` (for 

debugging/doctesting only). 

  

EXAMPLES:: 

  

sage: from sage.geometry.integral_points import InequalityCollection, print_cache 

sage: P = Polyhedron(ieqs=[(2,3,7)]) 

sage: ieq = InequalityCollection(P, [0,1], [0]*2,[1]*2); ieq 

The collection of inequalities 

integer: (3, 7) x + 2 >= 0 

sage: ieq.prepare_next_to_inner_loop([3,5]) 

sage: ieq.prepare_inner_loop([3,5]) 

sage: print_cache(ieq) 

Cached inner loop: 3 * x_0 + 37 >= 0 

Cached next-to-inner loop: 3 * x_0 + 7 * x_1 + 2 >= 0 

""" 

cdef Inequality_int ieq = <Inequality_int>(inequality_collection.ineqs_int[0]) 

print('Cached inner loop: ' + 

str(ieq.coeff) + ' * x_0 + ' + str(ieq.cache) + ' >= 0') 

print('Cached next-to-inner loop: ' + 

str(ieq.coeff) + ' * x_0 + ' + 

str(ieq.coeff_next) + ' * x_1 + ' + str(ieq.cache_next) + ' >= 0')