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

# distutils: sources = sage/modular/arithgroup/sl2z.cpp sage/modular/arithgroup/farey.cpp 

r""" 

Farey Symbol for arithmetic subgroups of `{\rm PSL}_2(\ZZ)` 

  

AUTHORS: 

  

- Hartmut Monien (08 - 2011) 

  

based on the *KFarey* package by Chris Kurth. Implemented as C++ module 

for speed. 

""" 

  

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

# Copyright (C) 2011 Hartmut Monien <monien@th.physik.uni-bonn.de> 

# 

# This program is free software: you can redistribute it and/or modify 

# it under the terms of the GNU General Public License as published by 

# the Free Software Foundation, either version 2 of the License, or 

# (at your option) any later version. 

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

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

  

from __future__ import absolute_import, division 

  

from cpython.object cimport PyObject_RichCompare 

from itertools import groupby 

from cysignals.signals cimport sig_on, sig_off 

  

from sage.libs.gmpxx cimport * 

  

from sage.rings.all import CC, RR 

from sage.rings.integer cimport Integer 

from sage.rings.infinity import infinity 

from .congroup_gammaH import is_GammaH 

from .congroup_gamma1 import is_Gamma1 

from .congroup_gamma0 import is_Gamma0 

from .congroup_gamma import is_Gamma 

from .congroup_sl2z import SL2Z 

from sage.modular.cusps import Cusp 

  

from sage.plot.all import Graphics 

from sage.plot.colors import to_mpl_color 

from sage.plot.misc import options, rename_keyword 

from sage.plot.all import hyperbolic_arc, hyperbolic_triangle, text 

  

from sage.misc.latex import latex 

from sage.misc.lazy_attribute import lazy_attribute 

from sage.misc.cachefunc import cached_method 

from sage.structure.richcmp cimport richcmp_not_equal 

  

  

cdef extern from "sage/modular/arithgroup/sl2z.hpp": 

cppclass cpp_SL2Z "SL2Z": 

mpz_class a, b, c, d 

cpp_SL2Z(int, int, int, int) 

cpp_SL2Z(mpz_class, mpz_class, mpz_class, mpz_class) 

mpz_class a() 

mpz_class b() 

mpz_class c() 

mpz_class d() 

  

cdef extern from "sage/modular/arithgroup/farey.hpp": 

cppclass is_element_Gamma0: 

is_element_Gamma0(int) 

cppclass is_element_Gamma1: 

is_element_Gamma1(int) 

cppclass is_element_Gamma: 

is_element_Gamma(int) 

cppclass is_element_GammaH: 

is_element_GammaH(int, object) 

cppclass cpp_farey "FareySymbol": 

cpp_farey() 

cpp_farey(object) 

cpp_farey(object, is_element_Gamma*) 

cpp_farey(object, is_element_Gamma0*) 

cpp_farey(object, is_element_Gamma1*) 

cpp_farey(object, is_element_GammaH*) 

size_t genus() 

size_t index() 

size_t level() 

size_t nu2() 

size_t nu3() 

object is_element(mpz_t, mpz_t, mpz_t, mpz_t) 

object word_problem(mpz_t, mpz_t, mpz_t, mpz_t, cpp_SL2Z *) 

size_t get_cusp_class(mpz_t, mpz_t) 

object get_cusps() 

object get_cusp_widths() 

object get_transformation_to_cusp(mpz_t, mpz_t) 

object get_fractions() 

object get_coset() 

object get_generators() 

object get_pairings() 

object get_paired_sides() 

object get_pairing_matrices() 

object dumps() 

  

  

cdef class Farey: 

r""" 

A class for calculating Farey symbols of arithmetics subgroups of 

`{\rm PSL}_2(\ZZ)`. 

  

The arithmetic subgroup can be either any of 

the congruence subgroups implemented in Sage, i.e. Gamma, Gamma0, 

Gamma1 and GammaH or a subgroup of `{\rm PSL}_2(\ZZ)` which is 

given by a user written helper class defining membership in that 

group. 

  

REFERENCES: 

  

- Ravi S. Kulkarni, ''An arithmetic-geometric method in the study of the 

subgroups of the modular group'', `Amer. J. Math., 113(6):1053--1133, 

1991. <http://www.jstor.org/stable/2374900>`_ 

  

INPUT: 

  

- `G` - an arithmetic subgroup of `{\rm PSL}_2(\ZZ)` 

  

EXAMPLES: 

  

Create a Farey symbol for the group `\Gamma_0(11)`:: 

  

sage: f = FareySymbol(Gamma0(11)); f 

FareySymbol(Congruence Subgroup Gamma0(11)) 

  

Calculate the generators:: 

  

sage: f.generators() 

[ 

[1 1] [ 7 -2] [ 8 -3] [-1 0] 

[0 1], [11 -3], [11 -4], [ 0 -1] 

] 

  

Pickling the FareySymbol and recovering it:: 

  

sage: f == loads(dumps(f)) 

True 

  

Calculate the index of `\Gamma_H(33, [2, 5])` in 

`{\rm PSL}_2(\ZZ)` via FareySymbol:: 

  

sage: FareySymbol(GammaH(33, [2, 5])).index() 

48 

  

Calculate the generators of `\Gamma_1(4)`:: 

  

sage: FareySymbol(Gamma1(4)).generators() 

[ 

[1 1] [-3 1] 

[0 1], [-4 1] 

] 

  

Calculate the generators of the :meth:`example 

<sage.modular.arithgroup.arithgroup_perm.HsuExample10>` of an 

index 10 arithmetic subgroup given by Tim Hsu:: 

  

sage: from sage.modular.arithgroup.arithgroup_perm import HsuExample10 

sage: FareySymbol(HsuExample10()).generators() 

[ 

[1 2] [-2 1] [ 4 -3] 

[0 1], [-7 3], [ 3 -2] 

] 

  

Calculate the generators of the group `\Gamma' = 

\Gamma_0(8)\cap\Gamma_1(4)` using a helper class to define group membership:: 

  

sage: class GPrime: 

....: def __contains__(self, M): 

....: return M in Gamma0(8) and M in Gamma1(4) 

....: 

  

sage: FareySymbol(GPrime()).generators() 

[ 

[1 1] [ 5 -1] [ 5 -2] 

[0 1], [16 -3], [ 8 -3] 

] 

  

Calculate cusps of arithmetic subgroup defined via permutation group:: 

  

sage: L = SymmetricGroup(4)('(1, 2, 3)') 

  

sage: R = SymmetricGroup(4)('(1, 2, 4)') 

  

sage: FareySymbol(ArithmeticSubgroup_Permutation(L, R)).cusps() 

[-1, Infinity] 

  

Calculate the left coset representation of `\Gamma_H(8, [3])`:: 

  

sage: FareySymbol(GammaH(8, [3])).coset_reps() 

[ 

[1 0] [ 4 -1] [ 3 -1] [ 2 -1] [ 1 -1] [ 3 -1] [ 2 -1] [-1 0] 

[0 1], [ 1 0], [ 1 0], [ 1 0], [ 1 0], [ 4 -1], [ 3 -1], [ 3 -1], 

[ 1 -1] [-1 0] [ 0 -1] [-1 0] 

[ 2 -1], [ 2 -1], [ 1 -1], [ 1 -1] 

] 

""" 

cdef cpp_farey *this_ptr 

cdef object group 

  

def __cinit__(self, group, data=None): 

r""" 

Initialize FareySymbol:: 

  

sage: FareySymbol(Gamma0(23)) 

FareySymbol(Congruence Subgroup Gamma0(23)) 

""" 

self.group = group 

# if data is present we want to restore 

if data is not None: 

sig_on() 

self.this_ptr = new cpp_farey(data) 

sig_off() 

return 

## to accelerate the calculation of the FareySymbol 

## we implement the tests for the standard congruence groups 

## in the c++ module. For a general group the test if an element 

## of SL2Z is in the group the python __contains__ attribute 

## of the group is called 

cdef int p 

if hasattr(group, "level"): 

p = group.level() 

if group == SL2Z: 

sig_on() 

self.this_ptr = new cpp_farey() 

sig_off() 

elif is_Gamma0(group): 

sig_on() 

self.this_ptr = new cpp_farey(group, new is_element_Gamma0(p)) 

sig_off() 

elif is_Gamma1(group): 

sig_on() 

self.this_ptr = new cpp_farey(group, new is_element_Gamma1(p)) 

sig_off() 

elif is_Gamma(group): 

sig_on() 

self.this_ptr = new cpp_farey(group, new is_element_Gamma(p)) 

sig_off() 

elif is_GammaH(group): 

sig_on() 

l = group._GammaH_class__H 

self.this_ptr = new cpp_farey(group, new is_element_GammaH(p, l)) 

sig_off() 

else: 

sig_on() 

self.this_ptr = new cpp_farey(group) 

sig_off() 

  

def __dealloc__(self): 

r""" 

Remove reference to FareySymbol. 

  

TESTS:: 

  

sage: F = FareySymbol(Gamma0(23)) 

sage: del F 

""" 

del self.this_ptr 

  

@cached_method 

def pairing_matrices_to_tietze_index(self): 

r""" 

Obtain the translation table from pairing matrices 

to generators. 

  

The result is cached. 

  

OUTPUT: 

  

a list where the `i`-th entry is a nonzero integer `k`, 

such that if `k > 0` then the `i`-th pairing matrix is (up to sign) 

the `(k-1)`-th generator and, if `k < 0`, then the `i`-th pairing 

matrix is (up to sign) the inverse of the `(-k-1)`-th generator. 

  

EXAMPLES:: 

  

sage: F = Gamma0(40).farey_symbol() 

sage: table = F.pairing_matrices_to_tietze_index() 

sage: table[12] 

(-2, -1) 

sage: F.pairing_matrices()[12] 

[ 3 -1] 

[ 40 -13] 

sage: F.generators()[1]**-1 

[ -3 1] 

[-40 13] 

""" 

gens_dict = {g:i+1 for i,g in enumerate(self.generators())} 

ans = [] 

for pm in self.pairing_matrices(): 

a, b, c, d = pm.matrix().list() 

newval = gens_dict.get(SL2Z([a, b, c, d])) 

if newval is not None: 

ans.append((newval,1)) 

continue 

newval = gens_dict.get(SL2Z([-a, -b, -c, -d])) 

if newval is not None: 

ans.append((newval,-1)) 

continue 

newval = gens_dict.get(SL2Z([d, -b, -c, a])) 

if newval is not None: 

ans.append((-newval,1)) 

continue 

newval = gens_dict.get(SL2Z([-d, b, c, -a])) 

if newval is not None: 

ans.append((-newval,-1)) 

continue 

raise RuntimeError("This should have not happened") 

return ans 

  

@cached_method 

def _get_minus_one(self): 

r""" 

If -I belongs to self, return a Tietze word representing it. 

  

OUTPUT: 

  

A Tietze word representing the element -I if it belongs to self. 

Otherwise return [] 

  

EXAMPLES:: 

  

sage: Gamma1(30).farey_symbol()._get_minus_one() 

[] 

sage: G = Gamma0(30).farey_symbol() 

sage: G._get_minus_one() 

[14] 

sage: g = G.generators()[13] 

sage: (-g.matrix()).is_one() 

True 

sage: G = Gamma(1).farey_symbol() 

sage: G._get_minus_one() 

[1, 1] 

sage: g = G.generators()[0]**2 

sage: (-g.matrix()).is_one() 

True 

sage: G = Gamma0(3).farey_symbol() 

sage: G._get_minus_one() 

[2, 2, 2] 

sage: g = G.generators()[1]**3 

sage: (-g.matrix()).is_one() 

True 

""" 

for i,g in enumerate(self.generators()): 

m = g.matrix() 

if (-m).is_one(): 

return [i+1] 

t = m.trace() 

if t == 0: # order = 4 

return 2 * [i+1] 

elif t == 1: # order = 6 

return 3 * [i+1] 

return [] 

  

def word_problem(self, M, output = 'standard', check = True): 

r""" 

Solve the word problem (up to sign) using this Farey symbol. 

  

INPUT: 

  

- ``M`` -- An element `M` of `{\rm SL}_2(\ZZ)`. 

- ``output`` -- (default: ``'standard'``) Should be one of ``'standard'``, 

``'syllables'``, ``'gens'``. 

- ``check`` -- (default: ``True``) Whether to check for correct input and output. 

  

OUTPUT: 

  

A solution to the word problem for the matrix `M`. 

The format depends on the ``output`` parameter, as follows. 

  

- ``standard`` returns the so called the Tietze representation, 

consists of a tuple of nonzero integers `i`, where if `i` > 0 

then it indicates the `i`th generator (that is, ``self.generators()[0]`` 

would correspond to `i` = 1), and if `i` < 0 then it indicates 

the inverse of the `i`-th generator. 

- ``syllables`` returns a tuple of tuples of the form `(i,n)`, where 

`(i,n)` represents ``self.generators()[i] ^ n``, 

whose product equals `M` up to sign. 

- ``gens`` returns tuple of tuples of the form `(g,n)`, 

`(g,n)` such that the product of the matrices `g^n` 

equals `M` up to sign. 

  

EXAMPLES:: 

  

sage: F = Gamma0(30).farey_symbol() 

sage: gens = F.generators() 

sage: g = gens[3] * gens[10] * gens[8]^-1 * gens[5] 

sage: g 

[-628597 73008] 

[-692130 80387] 

sage: F.word_problem(g) 

(4, 11, -9, 6) 

sage: g = gens[3] * gens[10]^2 * gens[8]^-1 * gens[5] 

sage: g 

[-5048053 586303] 

[-5558280 645563] 

sage: F.word_problem(g, output = 'gens') 

(( 

[109 -10] 

[120 -11], 1 

), 

( 

[ 19 -7] 

[ 30 -11], 2 

), 

( 

[ 49 -9] 

[ 60 -11], -1 

), 

( 

[17 -2] 

[60 -7], 1 

)) 

sage: F.word_problem(g, output = 'syllables') 

((3, 1), (10, 2), (8, -1), (5, 1)) 

  

TESTS: 

  

Check that problem with forgotten generator is fixed:: 

  

sage: from sage.misc.misc_c import prod 

sage: G = Gamma0(10) 

sage: F = G.farey_symbol() 

sage: g = G([-701,-137,4600,899]) 

sage: g1 = prod(F.generators()[i]**a for i,a in F.word_problem(g, output = 'syllables')) 

sage: g == g1 

True 

  

Check that it works for GammaH as well (:trac:`19660`):: 

  

sage: G = GammaH(147, [8]) 

sage: G.farey_symbol().word_problem(G([1,1,0,1])) 

(1,) 

  

Check that :trac:`20347` is solved:: 

  

sage: from sage.misc.misc_c import prod 

sage: G = ArithmeticSubgroup_Permutation(S2="(1,2)(3,4)",S3="(1,2,3)") 

sage: S = G.farey_symbol() 

sage: g1,g2 = S.generators() 

sage: g = g1^3 * g2^-2 * g1 * g2 

sage: S.word_problem(g) 

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

sage: h = prod(S.generators()[i]**a for i,a in S.word_problem(g, output = 'syllables')) 

sage: g == h 

True 

""" 

if output not in ['standard', 'syllables', 'gens']: 

raise ValueError('Unrecognized output format') 

if check: 

if M not in self.group: 

raise ValueError("Matrix ( %s ) is not in group ( %s )"%(M,self.group)) 

cdef Integer a = M.d() 

cdef Integer b = -M.b() 

cdef Integer c = -M.c() 

cdef Integer d = M.a() 

cdef cpp_SL2Z *cpp_beta = new cpp_SL2Z(1,0,0,1) 

E = SL2Z([-1,0,0,-1]) 

sig_on() 

result = self.this_ptr.word_problem(a.value, b.value, c.value, d.value, cpp_beta) 

sig_off() 

beta = convert_to_SL2Z(cpp_beta[0])**-1 

mbeta = SL2Z([-beta.a(),-beta.b(),-beta.c(),-beta.d()]) 

V = self.pairing_matrices_to_tietze_index() 

sgn = 1 

tietze = [] 

for o in result: 

if o > 0: 

tietze.append(V[o-1][0]) 

sgn *= V[o-1][1] 

else: 

tietze.append(-V[-o-1][0]) 

sgn *= V[-o-1][1] 

if sgn == -1: 

beta, mbeta = mbeta, beta 

  

gens_dict = {g:i+1 for i,g in enumerate(self.generators())} 

extra_tietze = [] 

if beta.is_one(): 

found = True 

elif mbeta.is_one(): 

found = True 

extra_tietze = self._get_minus_one() 

else: 

found = False 

if not found: 

newval = gens_dict.get(beta) 

if newval is not None: 

found = True 

extra_tietze = [newval] 

if not found: 

newval = gens_dict.get(beta**-1) 

if newval is not None: 

found = True 

extra_tietze = [-newval] 

if not found: 

newval = gens_dict.get(mbeta) 

if newval is not None: 

found = True 

extra_tietze = [newval] + self._get_minus_one() 

if not found: 

newval = gens_dict.get(mbeta**-1) 

if newval is not None: 

found = True 

extra_tietze = [-newval] + self._get_minus_one() 

tietze.extend(extra_tietze) 

tietze.reverse() 

gens = self.generators() 

if check: 

tmp = SL2Z([1,0,0,1]) 

for i in range(len(tietze)): 

t = tietze[i] 

tmp = tmp * gens[t-1] if t > 0 else tmp * gens[-t-1]**-1 

assert tmp.matrix() == M.matrix(),'%s %s %s'%(tietze, tmp.matrix(),M.matrix()) 

if output == 'standard': 

return tuple(tietze) 

if output == 'syllables': 

return tuple((a-1,len(list(g))) if a > 0 else (-a-1,-len(list(g))) for a,g in groupby(tietze)) 

else: # output == 'gens' 

return tuple((gens[a-1],len(list(g))) if a > 0 else (gens[-a-1],-len(list(g))) for a,g in groupby(tietze)) 

  

def __contains__(self, M): 

r""" 

Tests if element is in the arithmetic group of the Farey symbol 

via LLT algorithm. 

  

EXAMPLES:: 

  

sage: SL2Z([0, -1, 1, 0]) in FareySymbol(Gamma0(6)) 

False 

  

sage: SL2Z([1, 1, 0, 1]) in FareySymbol(Gamma0(6)) 

True 

""" 

cdef Integer a = M.a() 

cdef Integer b = M.b() 

cdef Integer c = M.c() 

cdef Integer d = M.d() 

sig_on() 

result = self.this_ptr.is_element(a.value, b.value, c.value, d.value) 

sig_off() 

return result 

  

def __richcmp__(self, other, op): 

r""" 

Compare self to others. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma(2)) == FareySymbol(Gamma0(7)) 

False 

  

sage: FareySymbol(Gamma0(23)) == loads(dumps(FareySymbol(Gamma0(23)))) 

True 

""" 

if not isinstance(other, Farey): 

return NotImplemented 

  

cosetA = self.coset_reps() 

cosetB = other.coset_reps() 

if cosetA != cosetB: 

return richcmp_not_equal(cosetA, cosetB, op) 

  

cuspA = self.cusps() 

cuspB = other.cusps() 

if cuspA != cuspB: 

return richcmp_not_equal(cuspA, cuspB, op) 

  

return PyObject_RichCompare(self.fractions(), other.fractions(), op) 

  

def __reduce__(self): 

r""" 

Serialization for pickling:: 

  

sage: FareySymbol(Gamma0(4)).__reduce__() 

(<type 'sage.modular.arithgroup.farey_symbol.Farey'>, ...)) 

  

""" 

return Farey, (self.group, self.this_ptr.dumps()) 

  

def __repr__(self): 

r""" 

Return the string representation of ``self``. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(23)).__repr__() 

'FareySymbol(Congruence Subgroup Gamma0(23))' 

""" 

if hasattr(self.group, "_repr_"): 

return "FareySymbol(%s)" % self.group._repr_() 

elif hasattr(self.group, "__repr__"): 

return "FareySymbol(%r)" % self.group 

else: 

return "FareySymbol(?)" 

  

def _latex_(self, forced_format=None): 

r""" 

Return the LaTeX representation of ``self``. 

  

INPUT: 

  

- ``forced_format`` -- A format string ('plain' or 'xymatrix') 

or ``None``. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(11))._latex_(forced_format = 'plain') 

'\\left( -\\infty\\underbrace{\\quad}_{1} 0\\underbrace{\\quad}_{2} \\frac{1}{3}\\underbrace{\\quad}_{3} \\frac{1}{2}\\underbrace{\\quad}_{2} \\frac{2}{3}\\underbrace{\\quad}_{3} 1\\underbrace{\\quad}_{1} \\infty\\right)' 

sage: FareySymbol(Gamma0(11))._latex_(forced_format = 'xymatrix') 

'\\begin{xy}\\xymatrix{& -\\infty \\ar@{-}@/_1pc/[r]_{1}& 0 \\ar@{-}@/_1pc/[r]_{2}& \\frac{1}{3} \\ar@{-}@/_1pc/[r]_{3}& \\frac{1}{2} \\ar@{-}@/_1pc/[r]_{2}& \\frac{2}{3} \\ar@{-}@/_1pc/[r]_{3}& 1 \\ar@{-}@/_1pc/[r]_{1}& \\infty }\\end{xy}' 

  

sage: if '\\xymatrix' in sage.misc.latex.latex.mathjax_avoid_list(): 

....: 'xymatrix' not in FareySymbol(Gamma0(11))._latex_() 

....: else: 

....: 'xymatrix' in FareySymbol(Gamma0(11))._latex_() 

True 

""" 

if forced_format == 'plain' or \ 

(forced_format is None and '\\xymatrix' in latex.mathjax_avoid_list()): 

# output not using xymatrix 

s = r'\left( -\infty' 

a = [x._latex_() for x in self.fractions()] + ['\infty'] 

b = self.pairings() 

for i in xrange(len(a)): 

u = b[i] 

if u == -3: 

u = r'\bullet' 

elif u == -2: 

u = r'\circ' 

s += r'\underbrace{\quad}_{%s} %s' % (u, a[i]) 

return s + r'\right)' 

else: 

# output using xymatrix 

s = r'\begin{xy}\xymatrix{& -\infty ' 

f = [x._latex_() for x in self.fractions()]+[r'\infty'] 

f.reverse() 

for p in self.pairings(): 

if p >= 0: 

s += r'\ar@{-}@/_1pc/[r]_{%s}' % p 

elif p == -2: 

s += r'\ar@{-}@/_1pc/[r]_{\circ}' 

elif p == -3: 

s += r'\ar@{-}@/_1pc/[r]_{\bullet}' 

s += r'& %s ' % f.pop() 

s += r'}\end{xy}' 

return s 

  

def index(self): 

r""" 

Return the index of the arithmetic group of the FareySymbol 

in `{\rm PSL}_2(\ZZ)`. 

  

EXAMPLES:: 

  

sage: [FareySymbol(Gamma0(n)).index() for n in range(1, 16)] 

[1, 3, 4, 6, 6, 12, 8, 12, 12, 18, 12, 24, 14, 24, 24] 

""" 

return self.this_ptr.index() 

  

def genus(self): 

r""" 

Return the genus of the arithmetic group of the FareySymbol. 

  

EXAMPLES:: 

  

sage: [FareySymbol(Gamma0(n)).genus() for n in range(16, 32)] 

[0, 1, 0, 1, 1, 1, 2, 2, 1, 0, 2, 1, 2, 2, 3, 2] 

""" 

return self.this_ptr.genus() 

  

def level(self): 

r""" 

Return the level of the arithmetic group of the FareySymbol. 

  

EXAMPLES:: 

  

sage: [FareySymbol(Gamma0(n)).level() for n in range(1, 16)] 

[1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15] 

""" 

return self.this_ptr.level() 

  

def nu2(self): 

r""" 

Return the number of elliptic points of order two. 

  

EXAMPLES:: 

  

sage: [FareySymbol(Gamma0(n)).nu2() for n in range(1, 16)] 

[1, 1, 0, 0, 2, 0, 0, 0, 0, 2, 0, 0, 2, 0, 0] 

""" 

return self.this_ptr.nu2() 

  

def nu3(self): 

r""" 

Return the number of elliptic points of order three. 

  

EXAMPLES:: 

  

sage: [FareySymbol(Gamma0(n)).nu3() for n in range(1, 16)] 

[1, 0, 1, 0, 0, 0, 2, 0, 0, 0, 0, 0, 2, 0, 0] 

""" 

return self.this_ptr.nu3() 

  

def coset_reps(self): 

r""" 

Left coset of the arithmetic group of the FareySymbol. 

  

EXAMPLES: 

  

Calculate the left coset of `\Gamma_0(6)`:: 

  

sage: FareySymbol(Gamma0(6)).coset_reps() 

[ 

[1 0] [ 3 -1] [ 2 -1] [ 1 -1] [ 2 -1] [ 3 -2] [ 1 -1] [-1 0] 

[0 1], [ 1 0], [ 1 0], [ 1 0], [ 3 -1], [ 2 -1], [ 2 -1], [ 2 -1], 

[ 1 -1] [ 0 -1] [-1 0] [-2 1] 

[ 3 -2], [ 1 -1], [ 1 -1], [ 1 -1] 

] 

""" 

return self.this_ptr.get_coset() 

  

def generators(self): 

r""" 

Minimal set of generators of the group of the FareySymbol. 

  

EXAMPLES: 

  

Calculate the generators of `\Gamma_0(6)`:: 

  

sage: FareySymbol(Gamma0(6)).generators() 

[ 

[1 1] [ 5 -1] [ 7 -3] [-1 0] 

[0 1], [ 6 -1], [12 -5], [ 0 -1] 

] 

  

Calculate the generators of `{\rm SL}_2(\ZZ)`:: 

  

sage: FareySymbol(SL2Z).generators() 

[ 

[ 0 -1] [ 0 -1] 

[ 1 0], [ 1 -1] 

] 

  

The unique index 2 even subgroup and index 4 odd subgroup each get handled correctly:: 

  

sage: FareySymbol(ArithmeticSubgroup_Permutation(S2="(1,2)", S3="()")).generators() 

[ 

[ 0 1] [-1 1] 

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

] 

sage: FareySymbol(ArithmeticSubgroup_Permutation(S2="(1,2, 3, 4)", S3="(1,3)(2,4)")).generators() 

[ 

[ 0 1] [-1 1] 

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

] 

""" 

return self.this_ptr.get_generators() 

  

def fractions(self): 

r""" 

Fractions of the FareySymbol. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma(4)).fractions() 

[0, 1/2, 1, 3/2, 2, 5/2, 3, 7/2, 4] 

""" 

return self.this_ptr.get_fractions() 

  

def pairings(self): 

r""" 

Pairings of the sides of the fundamental domain of the Farey symbol 

of the arithmetic group. 

  

The sides of the hyperbolic polygon are 

numbered 0, 1, ... from left to right. Conventions: even pairings are 

denoted by -2, odd pairings by -3 while free pairings are denoted by 

an integer number greater than zero. 

  

EXAMPLES: 

  

Odd pairings:: 

  

sage: FareySymbol(Gamma0(7)).pairings() 

[1, -3, -3, 1] 

  

Even and odd pairings:: 

  

FareySymbol(Gamma0(13)).pairings() 

[1, -3, -2, -2, -3, 1] 

  

Only free pairings:: 

  

sage: FareySymbol(Gamma0(23)).pairings() 

[1, 2, 3, 5, 3, 4, 2, 4, 5, 1] 

""" 

return self.this_ptr.get_pairings() 

  

def paired_sides(self): 

r""" 

Pairs of index of the sides of the fundamental domain of the 

Farey symbol of the arithmetic group. The sides of the 

hyperbolic polygon are numbered 0, 1, ... from left to right. 

  

.. image:: ../../../media/pairing.png 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(11)).paired_sides() 

[(0, 5), (1, 3), (2, 4)] 

  

indicating that the side 0 is paired with 5, 1 with 3 and 2 with 4. 

""" 

return self.this_ptr.get_paired_sides() 

  

def pairing_matrices(self): 

r""" 

Pairing matrices of the sides of the fundamental domain. The sides 

of the hyperbolic polygon are numbered 0, 1, ... from left to right. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(6)).pairing_matrices() 

[ 

[1 1] [ 5 -1] [ 7 -3] [ 5 -3] [ 1 -1] [-1 1] 

[0 1], [ 6 -1], [12 -5], [12 -7], [ 6 -5], [ 0 -1] 

] 

""" 

  

return self.this_ptr.get_pairing_matrices() 

  

def cusps(self): 

r""" 

Cusps of the FareySymbol. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(6)).cusps() 

[0, 1/3, 1/2, Infinity] 

""" 

return self.this_ptr.get_cusps()+[Cusp(infinity)] 

  

def cusp_widths(self): 

r""" 

Cusps widths of the FareySymbol. 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(6)).cusp_widths() 

[6, 2, 3, 1] 

""" 

return self.this_ptr.get_cusp_widths() 

  

def cusp_class(self, c): 

r""" 

Cusp class of a cusp in the FareySymbol. 

  

INPUT: 

  

``c`` -- a cusp 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(12)).cusp_class(Cusp(1, 12)) 

5 

""" 

cusp = Cusp(c) 

cdef Integer p = cusp.numerator() 

cdef Integer q = cusp.denominator() 

sig_on() 

result = self.this_ptr.get_cusp_class(p.value, q.value) 

sig_off() 

return result 

  

def reduce_to_cusp(self, r): 

r""" 

Transformation of a rational number to cusp representative. 

  

INPUT: 

  

``r`` -- a rational number 

  

EXAMPLES:: 

  

sage: FareySymbol(Gamma0(12)).reduce_to_cusp(5/8) 

[ 5 -3] 

[12 -7] 

  

Reduce 11/17 to a cusp of for HsuExample10():: 

  

sage: from sage.modular.arithgroup.arithgroup_perm import HsuExample10 

sage: f = FareySymbol(HsuExample10()) 

sage: f.reduce_to_cusp(11/17) 

[14 -9] 

[-3 2] 

sage: _.acton(11/17) 

1 

sage: f.cusps()[f.cusp_class(11/17)] 

1 

""" 

cdef Integer p = r.numerator() 

cdef Integer q = r.denominator() 

sig_on() 

result = self.this_ptr.get_transformation_to_cusp(p.value, q.value) 

sig_off() 

return result 

  

@rename_keyword(rgbcolor='color') 

@options(alpha=1, fill=True, thickness=1, color='lightgray', 

color_even='white', 

zorder=2, linestyle='solid', show_pairing=True, 

tesselation='Dedekind', ymax=1) 

def fundamental_domain(self, **options): 

r""" 

Plot a fundamental domain of an arithmetic subgroup of 

`{\rm PSL}_2(\ZZ)` corresponding to the Farey symbol. 

  

OPTIONS: 

  

- ``fill`` -- boolean (default ``True``) fill the fundamental domain 

  

- ``linestyle`` -- string (default: 'solid') The style of the line, 

which is one of 'dashed', 'dotted', 'solid', 'dashdot', or '--', 

':', '-', '-.', respectively 

  

- ``color`` -- (default: 'lightgray') fill color; fill 

color for odd part of Dedekind tesselation. 

  

- ``show_pairing`` -- boolean (default: ``True``) flag for pairing 

  

- ``tesselation`` -- (default: 'Dedekind') The type of 

hyperbolic tesselation which is one of 

'coset', 'Dedekind' or ``None`` respectively 

  

- ``color_even`` -- fill color for even parts of Dedekind 

tesselation (default 'white'); ignored for other tesselations 

  

- ``thickness`` -- float (default: `1`) the thickness of the line 

  

- ``ymax`` -- float (default: `1`) maximal height 

  

EXAMPLES: 

  

For example, to plot the fundamental domain of `\Gamma_0(11)` 

with pairings use the following command:: 

  

sage: FareySymbol(Gamma0(11)).fundamental_domain() 

Graphics object consisting of 54 graphics primitives 

  

indicating that side 1 is paired with side 3 and side 2 is 

paired with side 4, see also :meth:`.paired_sides`. 

  

To plot the fundamental domain of `\Gamma(3)` without pairings 

use the following command:: 

  

sage: FareySymbol(Gamma(3)).fundamental_domain(show_pairing=False) 

Graphics object consisting of 48 graphics primitives 

  

Plot the fundamental domain of `\Gamma_0(23)` showing the left 

coset representatives:: 

  

sage: FareySymbol(Gamma0(23)).fundamental_domain(tesselation='coset') 

Graphics object consisting of 58 graphics primitives 

  

The same as above but with a custom linestyle:: 

  

sage: FareySymbol(Gamma0(23)).fundamental_domain(tesselation='coset', linestyle=':', thickness='2') 

Graphics object consisting of 58 graphics primitives 

  

""" 

from sage.plot.colors import rainbow 

I = CC(0, 1) 

w = RR(3).sqrt() 

L = 1000 

g = Graphics() 

## show coset 

for x in self.coset_reps(): 

a, b, c, d = x[1, 1], -x[0, 1], -x[1, 0], x[0, 0] 

A, B = CC(0, L), CC(0, L) 

if d != 0: 

A = b / d 

if c != 0: 

B = a / c 

C = (a*c+b*d+(a*d+b*c)/2+I*w/2)/(c*c+c*d+d*d) 

D = (a*c+b*d + CC(0, 1))/(c*c+d*d) 

if options['tesselation'] == 'Dedekind': 

g += hyperbolic_triangle(A, D, C, 

alpha=options['alpha'], 

color=options['color'], 

fill=options['fill'], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

g += hyperbolic_triangle(D, C, B, 

alpha=options['alpha'], 

color=options['color_even'], 

fill=options['fill'], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

g += hyperbolic_triangle(A, D, C, color='gray') 

g += hyperbolic_triangle(D, C, B, color='gray') 

elif options['tesselation'] == 'coset': 

g += hyperbolic_triangle(A, B, C, 

alpha=options['alpha'], 

color=options['color'], 

fill=options['fill'], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

g += hyperbolic_triangle(A, B, C, color='gray', 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

else: 

g += hyperbolic_triangle(A, B, C, 

alpha=options['alpha'], 

color=options['color'], 

fill=options['fill'], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

## show pairings 

p = self.pairings() 

x = self.fractions() 

if options['show_pairing']: 

rc = rainbow(max(p)-min(p)+1) 

if p[0] > 0: 

g += hyperbolic_arc(CC(0, L), x[0], color=rc[p[0]-min(p)], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

if p[-1] > 0: 

g += hyperbolic_arc(CC(0, L), x[-1], color=rc[p[-1]-min(p)], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

for i in range(len(x)-1): 

if p[i+1] > 0: 

g += hyperbolic_arc(x[i], x[i+1], color=rc[p[i+1]-min(p)], 

linestyle=options['linestyle'], 

thickness=options['thickness']) 

d = g.get_minmax_data() 

g.set_axes_range(d['xmin'], d['xmax'], 0, options['ymax']) 

return g 

  

  

#--- conversions ------------------------------------------------------------ 

  

cdef public long convert_to_long(n): 

cdef long m = n 

return m 

  

cdef public object convert_to_Integer(mpz_class a): 

A = Integer() 

A.set_from_mpz(a.get_mpz_t()) 

return A 

  

cdef public object convert_to_rational(mpq_class r): 

a = Integer() 

a.set_from_mpz(r.get_num_mpz_t()) 

b = Integer() 

b.set_from_mpz(r.get_den_mpz_t()) 

return a/b 

  

cdef public object convert_to_cusp(mpq_class r): 

a = Integer() 

a.set_from_mpz(r.get_num_mpz_t()) 

b = Integer() 

b.set_from_mpz(r.get_den_mpz_t()) 

return Cusp(a/b) 

  

cdef public object convert_to_SL2Z(cpp_SL2Z M): 

a = convert_to_Integer(M.a()) 

b = convert_to_Integer(M.b()) 

c = convert_to_Integer(M.c()) 

d = convert_to_Integer(M.d()) 

return SL2Z([a, b, c, d])