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

""" 

A class to keep information about faces of a polyhedron 

 

This module gives you a tool to work with the faces of a polyhedron 

and their relative position. First, you need to find the faces. To get 

the faces in a particular dimension, use the 

:meth:`~sage.geometry.polyhedron.base.face` method:: 

 

sage: P = polytopes.cross_polytope(3) 

sage: P.faces(3) 

(<0,1,2,3,4,5>,) 

sage: P.faces(2) 

(<0,1,2>, <0,1,3>, <0,2,4>, <0,3,4>, <3,4,5>, <2,4,5>, <1,3,5>, <1,2,5>) 

sage: P.faces(1) 

(<0,1>, <0,2>, <1,2>, <0,3>, <1,3>, <0,4>, <2,4>, <3,4>, <2,5>, <3,5>, <4,5>, <1,5>) 

 

or :meth:`~sage.geometry.polyhedron.base.face_lattice` to get the 

whole face lattice as a poset:: 

 

sage: P.face_lattice() 

Finite poset containing 28 elements with distinguished linear extension 

 

The faces are printed in shorthand notation where each integer is the 

index of a vertex/ray/line in the same order as the containing 

Polyhedron's :meth:`~sage.geometry.polyhedron.base.Vrepresentation` :: 

 

sage: face = P.faces(1)[3]; face 

<0,3> 

sage: P.Vrepresentation(0) 

A vertex at (-1, 0, 0) 

sage: P.Vrepresentation(3) 

A vertex at (0, 0, 1) 

sage: face.vertices() 

(A vertex at (-1, 0, 0), A vertex at (0, 0, 1)) 

 

The face itself is not represented by Sage's 

:func:`sage.geometry.polyhedron.constructor.Polyhedron` class, but by 

an auxiliary class to keep the information. You can get the face as a 

polyhedron with the :meth:`PolyhedronFace.as_polyhedron` method:: 

 

sage: face.as_polyhedron() 

A 1-dimensional polyhedron in ZZ^3 defined as the convex hull of 2 vertices 

sage: _.equations() 

(An equation (0, 1, 0) x + 0 == 0, 

An equation (1, 0, -1) x + 1 == 0) 

""" 

 

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

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

# 

# Distributed under the terms of the GNU General Public License (GPL) 

# 

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

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

from __future__ import print_function 

 

from sage.structure.sage_object import SageObject 

from sage.structure.richcmp import richcmp_method, richcmp 

from sage.misc.all import cached_method 

from sage.modules.free_module_element import vector 

from sage.matrix.constructor import matrix 

 

 

 

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

@richcmp_method 

class PolyhedronFace(SageObject): 

r""" 

A face of a polyhedron. 

 

This class is for use in 

:meth:`~sage.geometry.polyhedron.base.Polyhedron_base.face_lattice`. 

 

INPUT: 

 

No checking is performed whether the H/V-representation indices 

actually determine a face of the polyhedron. You should not 

manually create :class:`PolyhedronFace` objects unless you know 

what you are doing. 

 

OUTPUT: 

 

A :class:`PolyhedronFace`. 

 

EXAMPLES:: 

 

sage: octahedron = polytopes.cross_polytope(3) 

sage: inequality = octahedron.Hrepresentation(2) 

sage: face_h = tuple([ inequality ]) 

sage: face_v = tuple( inequality.incident() ) 

sage: face_h_indices = [ h.index() for h in face_h ] 

sage: face_v_indices = [ v.index() for v in face_v ] 

sage: from sage.geometry.polyhedron.face import PolyhedronFace 

sage: face = PolyhedronFace(octahedron, face_v_indices, face_h_indices) 

sage: face 

<0,1,2> 

sage: face.dim() 

2 

sage: face.ambient_Hrepresentation() 

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

sage: face.ambient_Vrepresentation() 

(A vertex at (-1, 0, 0), A vertex at (0, -1, 0), A vertex at (0, 0, -1)) 

""" 

 

def __init__(self, polyhedron, V_indices, H_indices): 

r""" 

The constructor. 

 

See :class:`PolyhedronFace` for more information. 

 

INPUT: 

 

- ``polyhedron`` -- a :class:`Polyhedron`. The ambient 

polyhedron. 

 

- ``V_indices`` -- list of sorted integers. The indices of the 

face-spanning V-representation objects in the ambient 

polyhedron. 

 

- ``H_indices`` -- list of sorted integers. The indices of the 

H-representation objects of the ambient polyhedron that are 

saturated on the face. 

 

TESTS:: 

 

sage: from sage.geometry.polyhedron.face import PolyhedronFace 

sage: PolyhedronFace(Polyhedron(), [], []) # indirect doctest 

<> 

""" 

self._polyhedron = polyhedron 

self._ambient_Vrepresentation_indices = tuple(V_indices) 

self._ambient_Hrepresentation_indices = tuple(H_indices) 

self._ambient_Vrepresentation = tuple( polyhedron.Vrepresentation(i) for i in V_indices ) 

self._ambient_Hrepresentation = tuple( polyhedron.Hrepresentation(i) for i in H_indices ) 

 

def __hash__(self): 

r""" 

TESTS:: 

 

sage: P = Polyhedron([[0,0],[0,1],[23,3],[9,12]]) 

sage: list(map(hash, P.faces(1))) # random 

[2377119663630407734, 

2377136578164722109, 

5966674064902575359, 

4795242501625591634] 

""" 

return hash((self._polyhedron, self._ambient_Vrepresentation_indices)) 

 

def vertex_generator(self): 

""" 

Return a generator for the vertices of the face. 

 

EXAMPLES:: 

 

sage: triangle = Polyhedron(vertices=[[1,0],[0,1],[1,1]]) 

sage: face = triangle.faces(1)[0] 

sage: for v in face.vertex_generator(): print(v) 

A vertex at (0, 1) 

A vertex at (1, 0) 

sage: type(face.vertex_generator()) 

<... 'generator'> 

""" 

for V in self.ambient_Vrepresentation(): 

if V.is_vertex(): 

yield V 

 

@cached_method 

def vertices(self): 

""" 

Return all vertices of the face. 

 

OUTPUT: 

 

A tuple of vertices. 

 

EXAMPLES:: 

 

sage: triangle = Polyhedron(vertices=[[1,0],[0,1],[1,1]]) 

sage: face = triangle.faces(1)[0] 

sage: face.vertices() 

(A vertex at (0, 1), A vertex at (1, 0)) 

""" 

return tuple(self.vertex_generator()) 

 

def ray_generator(self): 

""" 

Return a generator for the rays of the face. 

 

EXAMPLES:: 

 

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

sage: face = pi.faces(1)[0] 

sage: next(face.ray_generator()) 

A ray in the direction (1, 0) 

""" 

for V in self.ambient_Vrepresentation(): 

if V.is_ray(): 

yield V 

 

@cached_method 

def rays(self): 

""" 

Return the rays of the face. 

 

OUTPUT: 

 

A tuple of rays. 

 

EXAMPLES:: 

 

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

sage: face = p.faces(2)[0] 

sage: face.rays() 

(A ray in the direction (1, 0, 0), A ray in the direction (0, 1, 0)) 

""" 

return tuple(self.ray_generator()) 

 

def line_generator(self): 

""" 

Return a generator for the lines of the face. 

 

EXAMPLES:: 

 

sage: pr = Polyhedron(rays = [[1,0],[-1,0],[0,1]], vertices = [[-1,-1]]) 

sage: face = pr.faces(1)[0] 

sage: next(face.line_generator()) 

A line in the direction (1, 0) 

""" 

for V in self.ambient_Vrepresentation(): 

if V.is_line(): 

yield V 

 

@cached_method 

def lines(self): 

""" 

Return all lines of the face. 

 

OUTPUT: 

 

A tuple of lines. 

 

EXAMPLES:: 

 

sage: p = Polyhedron(rays = [[1,0],[-1,0],[0,1],[1,1]], vertices = [[-2,-2],[2,3]]) 

sage: p.lines() 

(A line in the direction (1, 0),) 

""" 

return tuple(self.line_generator()) 

 

def __richcmp__(self, other, op): 

""" 

Compare ``self`` and ``other``. 

 

INPUT: 

 

- ``other`` -- anything. 

 

OUTPUT: 

 

Two faces test equal if and only if they are faces of the same 

(not just isomorphic) polyhedron and their generators have the 

same indices. 

 

EXAMPLES:: 

 

sage: square = polytopes.hypercube(2) 

sage: f = square.faces(1) 

sage: matrix(4,4, lambda i,j: ZZ(f[i] <= f[j])) 

[1 1 1 1] 

[0 1 1 1] 

[0 0 1 1] 

[0 0 0 1] 

sage: matrix(4,4, lambda i,j: ZZ(f[i] == f[j])) == 1 

True 

""" 

if not isinstance(other, PolyhedronFace): 

return NotImplemented 

if self._polyhedron is not other._polyhedron: 

return NotImplemented 

return richcmp(self._ambient_Vrepresentation_indices, 

other._ambient_Vrepresentation_indices, op) 

 

def ambient_Hrepresentation(self, index=None): 

r""" 

Return the H-representation objects of the ambient polytope 

defining the face. 

 

INPUT: 

 

- ``index`` -- optional. Either an integer or ``None`` 

(default). 

 

OUTPUT: 

 

If the optional argument is not present, a tuple of 

H-representation objects. Each entry is either an inequality 

or an equation. 

 

If the optional integer ``index`` is specified, the 

``index``-th element of the tuple is returned. 

 

EXAMPLES:: 

 

sage: square = polytopes.hypercube(2) 

sage: for face in square.face_lattice(): 

....: print(face.ambient_Hrepresentation()) 

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

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

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

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

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

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

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

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

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

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

() 

""" 

if index is None: 

return self._ambient_Hrepresentation 

else: 

return self._ambient_Hrepresentation[index] 

 

def ambient_Vrepresentation(self, index=None): 

r""" 

Return the V-representation objects of the ambient polytope 

defining the face. 

 

INPUT: 

 

- ``index`` -- optional. Either an integer or ``None`` 

(default). 

 

OUTPUT: 

 

If the optional argument is not present, a tuple of 

V-representation objects. Each entry is either a vertex, a 

ray, or a line. 

 

If the optional integer ``index`` is specified, the 

``index``-th element of the tuple is returned. 

 

EXAMPLES:: 

 

sage: square = polytopes.hypercube(2) 

sage: for fl in square.face_lattice(): 

....: print(fl.ambient_Vrepresentation()) 

() 

(A vertex at (-1, -1),) 

(A vertex at (-1, 1),) 

(A vertex at (1, -1),) 

(A vertex at (1, 1),) 

(A vertex at (-1, -1), A vertex at (-1, 1)) 

(A vertex at (-1, -1), A vertex at (1, -1)) 

(A vertex at (1, -1), A vertex at (1, 1)) 

(A vertex at (-1, 1), A vertex at (1, 1)) 

(A vertex at (-1, -1), A vertex at (-1, 1), 

A vertex at (1, -1), A vertex at (1, 1)) 

""" 

if index is None: 

return self._ambient_Vrepresentation 

else: 

return self._ambient_Vrepresentation[index] 

 

def n_ambient_Hrepresentation(self): 

""" 

Return the number of objects that make up the ambient 

H-representation of the polyhedron. 

 

See also :meth:`ambient_Hrepresentation`. 

 

OUTPUT: 

 

Integer. 

 

EXAMPLES:: 

 

sage: p = polytopes.cross_polytope(4) 

sage: face = p.face_lattice()[10] 

sage: face 

<0,2> 

sage: face.ambient_Hrepresentation() 

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

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

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

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

sage: face.n_ambient_Hrepresentation() 

4 

""" 

return len(self.ambient_Hrepresentation()) 

 

def n_ambient_Vrepresentation(self): 

""" 

Return the number of objects that make up the ambient 

V-representation of the polyhedron. 

 

See also :meth:`ambient_Vrepresentation`. 

 

OUTPUT: 

 

Integer. 

 

EXAMPLES:: 

 

sage: p = polytopes.cross_polytope(4) 

sage: face = p.face_lattice()[10] 

sage: face 

<0,2> 

sage: face.ambient_Vrepresentation() 

(A vertex at (-1, 0, 0, 0), A vertex at (0, 0, -1, 0)) 

sage: face.n_ambient_Vrepresentation() 

2 

""" 

return len(self.ambient_Vrepresentation()) 

 

def ambient_dim(self): 

r""" 

Return the dimension of the containing polyhedron. 

 

EXAMPLES:: 

 

sage: P = Polyhedron(vertices = [[1,0,0,0],[0,1,0,0]]) 

sage: face = P.faces(1)[0] 

sage: face.ambient_dim() 

4 

""" 

return self._polyhedron.ambient_dim() 

 

@cached_method 

def dim(self): 

""" 

Return the dimension of the face. 

 

OUTPUT: 

 

Integer. 

 

EXAMPLES:: 

 

sage: fl = polytopes.dodecahedron().face_lattice() 

sage: [ x.dim() for x in fl ] 

[-1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 

1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 

1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 3] 

""" 

if self.n_ambient_Vrepresentation()==0: 

return -1 

else: 

origin = vector(self.ambient_Vrepresentation(0)) 

v_list = [ vector(v)-origin for v in self.ambient_Vrepresentation() ] 

return matrix(v_list).rank() 

 

def _repr_(self): 

r""" 

Return a string representation. 

 

OUTPUT: 

 

A string listing the V-representation indices of the face. 

 

EXAMPLES:: 

 

sage: square = polytopes.hypercube(2) 

sage: a_face = list( square.face_lattice() )[8] 

sage: a_face.__repr__() 

'<1,3>' 

""" 

s = '<' 

s += ','.join([ str(v.index()) for v in self.ambient_Vrepresentation() ]) 

s += '>' 

return s 

 

def polyhedron(self): 

""" 

Return the containing polyhedron. 

 

EXAMPLES:: 

 

sage: P = polytopes.cross_polytope(3); P 

A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 6 vertices 

sage: face = P.faces(2)[3] 

sage: face 

<0,3,4> 

sage: face.polyhedron() 

A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 6 vertices 

""" 

return self._polyhedron 

 

@cached_method 

def as_polyhedron(self): 

""" 

Return the face as an independent polyhedron. 

 

OUTPUT: 

 

A polyhedron. 

 

EXAMPLES:: 

 

sage: P = polytopes.cross_polytope(3); P 

A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 6 vertices 

sage: face = P.faces(2)[3] 

sage: face 

<0,3,4> 

sage: face.as_polyhedron() 

A 2-dimensional polyhedron in ZZ^3 defined as the convex hull of 3 vertices 

 

sage: P.intersection(face.as_polyhedron()) == face.as_polyhedron() 

True 

""" 

P = self._polyhedron 

parent = P.parent() 

Vrep = (self.vertices(), self.rays(), self.lines()) 

return P.__class__(parent, Vrep, None)