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

r""" 

Linear code constructors that do not preserve the structural information 

 

This file contains a variety of constructions which builds the generator matrix 

of special (or random) linear codes and wraps them in a 

:class:`sage.coding.linear_code.LinearCode` object. These constructions are 

therefore not rich objects such as 

:class:`sage.coding.grs.GeneralizedReedSolomonCode`. 

 

For deprecation reasons, this file also contains some constructions for which 

Sage now does have rich representations. 

 

All codes available here can be accessed through the ``codes`` object:: 

 

sage: codes.GolayCode(GF(2),extended=False) 

[23, 12, 7] Golay code over GF(2) 

 

REFERENCES: 

 

- [HP2003]_ 

 

AUTHOR: 

 

- David Joyner (2007-05): initial version 

 

- " (2008-02): added cyclic codes, Hamming codes 

 

- " (2008-03): added BCH code, LinearCodeFromCheckmatrix, ReedSolomonCode, WalshCode, 

DuadicCodeEvenPair, DuadicCodeOddPair, QR codes (even and odd) 

 

- " (2008-09) fix for bug in BCHCode reported by F. Voloch 

 

- " (2008-10) small docstring changes to WalshCode and walsh_matrix 

 

""" 

 

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

# Copyright (C) 2007 David Joyner <wdjoyner@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 

from __future__ import absolute_import 

 

from sage.matrix.matrix_space import MatrixSpace 

from sage.matrix.constructor import matrix 

from sage.matrix.special import random_matrix 

from sage.rings.finite_rings.finite_field_constructor import FiniteField as GF 

from sage.groups.perm_gps.permgroup_named import SymmetricGroup 

from sage.misc.all import prod 

from .linear_code import LinearCode 

from sage.modules.free_module import span 

from sage.schemes.projective.projective_space import ProjectiveSpace 

from sage.structure.sequence import Sequence, Sequence_generic 

from sage.arith.all import GCD, LCM, divisors, quadratic_residues, gcd 

from sage.rings.finite_rings.integer_mod_ring import IntegerModRing 

from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing 

from sage.rings.integer import Integer 

from sage.sets.set import Set 

from sage.rings.finite_rings.integer_mod import Mod 

 

from sage.misc.superseded import deprecation, deprecated_function_alias 

 

############### utility functions ################ 

 

 

def _is_a_splitting(S1, S2, n, return_automorphism=False): 

""" 

Check wether ``(S1,S2)`` is a splitting of `\ZZ/n\ZZ`. 

 

A splitting of `R = \ZZ/n\ZZ` is a pair of subsets of `R` which is a 

partition of `R \\backslash \{0\}` and such that there exists an element `r` 

of `R` such that `r S_1 = S_2` and `r S_2 = S_1` (where `r S` is the 

point-wise multiplication of the elements of `S` by `r`). 

 

Splittings are useful for computing idempotents in the quotient 

ring `Q = GF(q)[x]/(x^n-1)`. 

 

INPUT: 

 

- ``S1, S2`` -- disjoint sublists partitioning ``[1, 2, ..., n-1]`` 

 

- ``n`` (integer) 

 

- ``return_automorphism`` (boolean) -- whether to return the automorphism 

exchanging `S_1` and `S_2`. 

 

OUTPUT: 

 

If ``return_automorphism is False`` (default) the function returns boolean values. 

 

Otherwise, it returns a pair ``(b, r)`` where ``b`` is a boolean indicating 

whether `S1`, `S2` is a splitting of `n`, and `r` is such that `r S_1 = S_2` 

and `r S_2 = S_1` (if `b` is ``False``, `r` is equal to ``None``). 

 

EXAMPLES:: 

 

sage: from sage.coding.code_constructions import _is_a_splitting 

sage: _is_a_splitting([1,2],[3,4],5) 

True 

sage: _is_a_splitting([1,2],[3,4],5,return_automorphism=True) 

(True, 4) 

 

sage: _is_a_splitting([1,3],[2,4,5,6],7) 

False 

sage: _is_a_splitting([1,3,4],[2,5,6],7) 

False 

 

sage: for P in SetPartitions(6,[3,3]): 

....: res,aut= _is_a_splitting(P[0],P[1],7,return_automorphism=True) 

....: if res: 

....: print((aut, P[0], P[1])) 

(6, {1, 2, 3}, {4, 5, 6}) 

(3, {1, 2, 4}, {3, 5, 6}) 

(6, {1, 3, 5}, {2, 4, 6}) 

(6, {1, 4, 5}, {2, 3, 6}) 

 

We illustrate now how to find idempotents in quotient rings:: 

 

sage: n = 11; q = 3 

sage: C = Zmod(n).cyclotomic_cosets(q); C 

[[0], [1, 3, 4, 5, 9], [2, 6, 7, 8, 10]] 

sage: S1 = C[1] 

sage: S2 = C[2] 

sage: _is_a_splitting(S1,S2,11) 

True 

sage: F = GF(q) 

sage: P.<x> = PolynomialRing(F,"x") 

sage: I = Ideal(P,[x^n-1]) 

sage: Q.<x> = QuotientRing(P,I) 

sage: i1 = -sum([x^i for i in S1]); i1 

2*x^9 + 2*x^5 + 2*x^4 + 2*x^3 + 2*x 

sage: i2 = -sum([x^i for i in S2]); i2 

2*x^10 + 2*x^8 + 2*x^7 + 2*x^6 + 2*x^2 

sage: i1^2 == i1 

True 

sage: i2^2 == i2 

True 

sage: (1-i1)^2 == 1-i1 

True 

sage: (1-i2)^2 == 1-i2 

True 

 

We return to dealing with polynomials (rather than elements of 

quotient rings), so we can construct cyclic codes:: 

 

sage: P.<x> = PolynomialRing(F,"x") 

sage: i1 = -sum([x^i for i in S1]) 

sage: i2 = -sum([x^i for i in S2]) 

sage: i1_sqrd = (i1^2).quo_rem(x^n-1)[1] 

sage: i1_sqrd == i1 

True 

sage: i2_sqrd = (i2^2).quo_rem(x^n-1)[1] 

sage: i2_sqrd == i2 

True 

sage: C1 = codes.CyclicCode(length = n, generator_pol = gcd(i1, x^n - 1)) 

sage: C2 = codes.CyclicCode(length = n, generator_pol = gcd(1-i2, x^n - 1)) 

sage: C1.dual_code().systematic_generator_matrix() == C2.systematic_generator_matrix() 

True 

 

This is a special case of Theorem 6.4.3 in [HP2003]_. 

""" 

R = IntegerModRing(n) 

S1 = set(R(x) for x in S1) 

S2 = set(R(x) for x in S2) 

 

# we first check whether (S1,S2) is a partition of R - {0} 

if (len(S1) + len(S2) != n-1 or len(S1) != len(S2) or 

R.zero() in S1 or R.zero() in S2 or not S1.isdisjoint(S2)): 

if return_automorphism: 

return False, None 

else: 

return False 

 

# now that we know that (S1,S2) is a partition, we look for an invertible 

# element b that maps S1 to S2 by multiplication 

for b in range(2,n): 

if GCD(b,n) == 1 and all(b*x in S2 for x in S1): 

if return_automorphism: 

return True, b 

else: 

return True 

if return_automorphism: 

return False, None 

else: 

return False 

 

is_a_splitting = deprecated_function_alias(21165, _is_a_splitting) 

 

def _lift2smallest_field(a): 

""" 

INPUT: a is an element of a finite field GF(q) 

 

OUTPUT: the element b of the smallest subfield F of GF(q) for 

which F(b)=a. 

 

EXAMPLES:: 

 

sage: from sage.coding.code_constructions import lift2smallest_field 

sage: FF.<z> = GF(3^4,"z") 

sage: a = z^10 

sage: lift2smallest_field(a) 

doctest:...: DeprecationWarning: lift2smallest_field is deprecated. Please use sage.coding.code_constructions._lift2smallest_field instead. 

See http://trac.sagemath.org/21165 for details. 

(2*z + 1, Finite Field in z of size 3^2) 

sage: a = z^40 

sage: lift2smallest_field(a) 

(2, Finite Field of size 3) 

 

AUTHORS: 

 

- John Cremona 

""" 

FF = a.parent() 

k = FF.degree() 

if k == 1: 

return a, FF 

pol = a.minimal_polynomial() 

d = pol.degree() 

if d == k: 

return a, FF 

p = FF.characteristic() 

F = GF(p**d,"z") 

b = pol.roots(F,multiplicities=False)[0] 

return b, F 

 

 

lift2smallest_field = deprecated_function_alias(21165, _lift2smallest_field) 

 

def lift2smallest_field2(a): 

""" 

INPUT: a is an element of a finite field GF(q) 

 

OUTPUT: the element b of the smallest subfield F of GF(q) for which F(b)=a. 

 

EXAMPLES:: 

 

sage: from sage.coding.code_constructions import lift2smallest_field2 

sage: FF.<z> = GF(3^4,"z") 

sage: a = z^40 

sage: lift2smallest_field2(a) 

doctest:...: DeprecationWarning: lift2smallest_field2 will be removed in a future release of Sage. Consider using sage.coding.code_constructions._lift2smallest_field instead, though this is private and may be removed in the future without deprecation warning. If you care about this functionality being in Sage, consider opening a Trac ticket for promoting the function to public. 

See http://trac.sagemath.org/21165 for details. 

(2, Finite Field of size 3) 

sage: FF.<z> = GF(2^4,"z") 

sage: a = z^15 

sage: lift2smallest_field2(a) 

(1, Finite Field of size 2) 

 

.. warning:: 

 

Since coercion (the FF(b) step) has a bug in it, this 

*only works* in the case when you *know* F is a prime field. 

 

AUTHORS: 

 

- David Joyner 

""" 

deprecation(21165, "lift2smallest_field2 will be removed in a future release of Sage. Consider using sage.coding.code_constructions._lift2smallest_field instead, though this is private and may be removed in the future without deprecation warning. If you care about this functionality being in Sage, consider opening a Trac ticket for promoting the function to public.") 

FF = a.parent() 

q = FF.order() 

if q.is_prime(): 

return a,FF 

p = q.factor()[0][0] 

k = q.factor()[0][1] 

for d in divisors(k): 

F = GF(p**d,"zz") 

for b in F: 

if FF(b) == a: 

return b, F 

 

 

def permutation_action(g,v): 

""" 

Returns permutation of rows g\*v. Works on lists, matrices, 

sequences and vectors (by permuting coordinates). The code requires 

switching from i to i+1 (and back again) since the SymmetricGroup 

is, by convention, the symmetric group on the "letters" 1, 2, ..., 

n (not 0, 1, ..., n-1). 

 

EXAMPLES:: 

 

sage: V = VectorSpace(GF(3),5) 

sage: v = V([0,1,2,0,1]) 

sage: G = SymmetricGroup(5) 

sage: g = G([(1,2,3)]) 

sage: permutation_action(g,v) 

(1, 2, 0, 0, 1) 

sage: g = G([()]) 

sage: permutation_action(g,v) 

(0, 1, 2, 0, 1) 

sage: g = G([(1,2,3,4,5)]) 

sage: permutation_action(g,v) 

(1, 2, 0, 1, 0) 

sage: L = Sequence([1,2,3,4,5]) 

sage: permutation_action(g,L) 

[2, 3, 4, 5, 1] 

sage: MS = MatrixSpace(GF(3),3,7) 

sage: A = MS([[1,0,0,0,1,1,0],[0,1,0,1,0,1,0],[0,0,0,0,0,0,1]]) 

sage: S5 = SymmetricGroup(5) 

sage: g = S5([(1,2,3)]) 

sage: A 

[1 0 0 0 1 1 0] 

[0 1 0 1 0 1 0] 

[0 0 0 0 0 0 1] 

sage: permutation_action(g,A) 

[0 1 0 1 0 1 0] 

[0 0 0 0 0 0 1] 

[1 0 0 0 1 1 0] 

 

It also works on lists and is a "left action":: 

 

sage: v = [0,1,2,0,1] 

sage: G = SymmetricGroup(5) 

sage: g = G([(1,2,3)]) 

sage: gv = permutation_action(g,v); gv 

[1, 2, 0, 0, 1] 

sage: permutation_action(g,v) == g(v) 

True 

sage: h = G([(3,4)]) 

sage: gv = permutation_action(g,v) 

sage: hgv = permutation_action(h,gv) 

sage: hgv == permutation_action(h*g,v) 

True 

 

AUTHORS: 

 

- David Joyner, licensed under the GPL v2 or greater. 

""" 

v_type_list = False 

if isinstance(v, list): 

v_type_list = True 

v = Sequence(v) 

if isinstance(v, Sequence_generic): 

V = v.universe() 

else: 

V = v.parent() 

n = len(list(v)) 

gv = [] 

for i in range(n): 

gv.append(v[g(i+1)-1]) 

if v_type_list: 

return gv 

return V(gv) 

 

def walsh_matrix(m0): 

""" 

This is the generator matrix of a Walsh code. The matrix of 

codewords correspond to a Hadamard matrix. 

 

EXAMPLES:: 

 

sage: walsh_matrix(2) 

[0 0 1 1] 

[0 1 0 1] 

sage: walsh_matrix(3) 

[0 0 0 0 1 1 1 1] 

[0 0 1 1 0 0 1 1] 

[0 1 0 1 0 1 0 1] 

sage: C = LinearCode(walsh_matrix(4)); C 

[16, 4] linear code over GF(2) 

sage: C.spectrum() 

[1, 0, 0, 0, 0, 0, 0, 0, 15, 0, 0, 0, 0, 0, 0, 0, 0] 

 

This last code has minimum distance 8. 

 

REFERENCES: 

 

- :wikipedia:`Hadamard_matrix` 

""" 

m = int(m0) 

if m == 1: 

return matrix(GF(2), 1, 2, [ 0, 1]) 

if m > 1: 

row2 = [x.list() for x in walsh_matrix(m-1).augment(walsh_matrix(m-1)).rows()] 

return matrix(GF(2), m, 2**m, [[0]*2**(m-1) + [1]*2**(m-1)] + row2) 

raise ValueError("%s must be an integer > 0."%m0) 

 

 

##################### main constructions ##################### 

 

 

def BinaryGolayCode(): 

""" 

This method is now deprecated. 

Please use :class:`sage.coding.golay_code.GolayCode` instead. 

""" 

from sage.misc.superseded import deprecation 

from .golay_code import GolayCode 

deprecation(20787, "codes.BinaryGolayCode is now deprecated. Please use codes.GolayCode instead.") 

return GolayCode(GF(2), False) 

 

def CyclicCodeFromGeneratingPolynomial(n,g,ignore=True): 

r""" 

If g is a polynomial over GF(q) which divides `x^n-1` then 

this constructs the code "generated by g" (ie, the code associated 

with the principle ideal `gR` in the ring 

`R = GF(q)[x]/(x^n-1)` in the usual way). 

 

The option "ignore" says to ignore the condition that (a) the 

characteristic of the base field does not divide the length (the 

usual assumption in the theory of cyclic codes), and (b) `g` 

must divide `x^n-1`. If ignore=True, instead of returning 

an error, a code generated by `gcd(x^n-1,g)` is created. 

 

EXAMPLES:: 

 

sage: P.<x> = PolynomialRing(GF(3),"x") 

sage: g = x-1 

sage: C = codes.CyclicCodeFromGeneratingPolynomial(4,g); C 

doctest:... 

DeprecationWarning: codes.CyclicCodeFromGeneratingPolynomial is now deprecated. Please use codes.CyclicCode instead. 

See http://trac.sagemath.org/20100 for details. 

[4, 3] Cyclic Code over GF(3) 

sage: P.<x> = PolynomialRing(GF(4,"a"),"x") 

sage: g = x^3+1 

sage: C = codes.CyclicCodeFromGeneratingPolynomial(9,g); C 

[9, 6] Cyclic Code over GF(4) 

sage: P.<x> = PolynomialRing(GF(2),"x") 

sage: g = x^3+x+1 

sage: C = codes.CyclicCodeFromGeneratingPolynomial(7,g); C 

[7, 4] Cyclic Code over GF(2) 

sage: C.generator_matrix() 

[1 1 0 1 0 0 0] 

[0 1 1 0 1 0 0] 

[0 0 1 1 0 1 0] 

[0 0 0 1 1 0 1] 

sage: g = x+1 

sage: C = codes.CyclicCodeFromGeneratingPolynomial(4,g); C 

Traceback (most recent call last): 

... 

ValueError: Only cyclic codes whose length and field order are coprimes are implemented. 

""" 

from sage.misc.superseded import deprecation 

from sage.coding.cyclic_code import CyclicCode 

deprecation(20100, "codes.CyclicCodeFromGeneratingPolynomial is now deprecated. Please use codes.CyclicCode instead.") 

return CyclicCode(length = n, generator_pol = g) 

 

 

def CyclicCodeFromCheckPolynomial(n,h,ignore=True): 

r""" 

If h is a polynomial over GF(q) which divides `x^n-1` then 

this constructs the code "generated by `g = (x^n-1)/h`" 

(ie, the code associated with the principle ideal `gR` in 

the ring `R = GF(q)[x]/(x^n-1)` in the usual way). The 

option "ignore" says to ignore the condition that the 

characteristic of the base field does not divide the length (the 

usual assumption in the theory of cyclic codes). 

 

EXAMPLES:: 

 

sage: P.<x> = PolynomialRing(GF(3),"x") 

sage: C = codes.CyclicCodeFromCheckPolynomial(4,x + 1); C 

doctest:... 

DeprecationWarning: codes.CyclicCodeFromCheckPolynomial is now deprecated. Please use codes.CyclicCode instead. 

See http://trac.sagemath.org/20100 for details. 

[4, 1] Cyclic Code over GF(3) 

sage: C = codes.CyclicCodeFromCheckPolynomial(4,x^3 + x^2 + x + 1); C 

[4, 3] Cyclic Code over GF(3) 

sage: C.generator_matrix() 

[2 1 0 0] 

[0 2 1 0] 

[0 0 2 1] 

""" 

from sage.misc.superseded import deprecation 

from sage.coding.cyclic_code import CyclicCode 

deprecation(20100, "codes.CyclicCodeFromCheckPolynomial is now deprecated. Please use codes.CyclicCode instead.") 

P = h.parent() 

x = P.gen() 

g = P((x**n-1)/h) 

return CyclicCode(length = n, generator_pol = g) 

 

def DuadicCodeEvenPair(F,S1,S2): 

r""" 

Constructs the "even pair" of duadic codes associated to the 

"splitting" (see the docstring for ``_is_a_splitting`` 

for the definition) S1, S2 of n. 

 

.. warning:: 

 

Maybe the splitting should be associated to a sum of 

q-cyclotomic cosets mod n, where q is a *prime*. 

 

EXAMPLES:: 

 

sage: from sage.coding.code_constructions import _is_a_splitting 

sage: n = 11; q = 3 

sage: C = Zmod(n).cyclotomic_cosets(q); C 

[[0], [1, 3, 4, 5, 9], [2, 6, 7, 8, 10]] 

sage: S1 = C[1] 

sage: S2 = C[2] 

sage: _is_a_splitting(S1,S2,11) 

True 

sage: codes.DuadicCodeEvenPair(GF(q),S1,S2) 

([11, 5] Cyclic Code over GF(3), 

[11, 5] Cyclic Code over GF(3)) 

""" 

from .cyclic_code import CyclicCode 

n = len(S1) + len(S2) + 1 

if not _is_a_splitting(S1,S2,n): 

raise TypeError("%s, %s must be a splitting of %s."%(S1,S2,n)) 

q = F.order() 

k = Mod(q,n).multiplicative_order() 

FF = GF(q**k,"z") 

z = FF.gen() 

zeta = z**((q**k-1)/n) 

P1 = PolynomialRing(FF,"x") 

x = P1.gen() 

g1 = prod([x-zeta**i for i in S1+[0]]) 

g2 = prod([x-zeta**i for i in S2+[0]]) 

P2 = PolynomialRing(F,"x") 

x = P2.gen() 

gg1 = P2([_lift2smallest_field(c)[0] for c in g1.coefficients(sparse=False)]) 

gg2 = P2([_lift2smallest_field(c)[0] for c in g2.coefficients(sparse=False)]) 

C1 = CyclicCode(length = n, generator_pol = gg1) 

C2 = CyclicCode(length = n, generator_pol = gg2) 

return C1,C2 

 

def DuadicCodeOddPair(F,S1,S2): 

""" 

Constructs the "odd pair" of duadic codes associated to the 

"splitting" S1, S2 of n. 

 

.. warning:: 

 

Maybe the splitting should be associated to a sum of 

q-cyclotomic cosets mod n, where q is a *prime*. 

 

EXAMPLES:: 

 

sage: from sage.coding.code_constructions import _is_a_splitting 

sage: n = 11; q = 3 

sage: C = Zmod(n).cyclotomic_cosets(q); C 

[[0], [1, 3, 4, 5, 9], [2, 6, 7, 8, 10]] 

sage: S1 = C[1] 

sage: S2 = C[2] 

sage: _is_a_splitting(S1,S2,11) 

True 

sage: codes.DuadicCodeOddPair(GF(q),S1,S2) 

([11, 6] Cyclic Code over GF(3), 

[11, 6] Cyclic Code over GF(3)) 

 

This is consistent with Theorem 6.1.3 in [HP2003]_. 

""" 

from .cyclic_code import CyclicCode 

n = len(S1) + len(S2) + 1 

if not _is_a_splitting(S1,S2,n): 

raise TypeError("%s, %s must be a splitting of %s."%(S1,S2,n)) 

q = F.order() 

k = Mod(q,n).multiplicative_order() 

FF = GF(q**k,"z") 

z = FF.gen() 

zeta = z**((q**k-1)/n) 

P1 = PolynomialRing(FF,"x") 

x = P1.gen() 

g1 = prod([x-zeta**i for i in S1+[0]]) 

g2 = prod([x-zeta**i for i in S2+[0]]) 

j = sum([x**i/n for i in range(n)]) 

P2 = PolynomialRing(F,"x") 

x = P2.gen() 

coeffs1 = [_lift2smallest_field(c)[0] for c in (g1+j).coefficients(sparse=False)] 

coeffs2 = [_lift2smallest_field(c)[0] for c in (g2+j).coefficients(sparse=False)] 

gg1 = P2(coeffs1) 

gg2 = P2(coeffs2) 

gg1 = gcd(gg1, x**n - 1) 

gg2 = gcd(gg2, x**n - 1) 

C1 = CyclicCode(length = n, generator_pol = gg1) 

C2 = CyclicCode(length = n, generator_pol = gg2) 

return C1,C2 

 

 

def ExtendedBinaryGolayCode(): 

""" 

This method is now deprecated. 

Please use :class:`sage.coding.golay_code.GolayCode` instead. 

""" 

from sage.misc.superseded import deprecation 

from .golay_code import GolayCode 

deprecation(20787, "codes.ExtendedBinaryGolayCode is now deprecated. Please use codes.GolayCode instead.") 

return GolayCode(GF(2)) 

 

 

def ExtendedQuadraticResidueCode(n,F): 

r""" 

The extended quadratic residue code (or XQR code) is obtained from 

a QR code by adding a check bit to the last coordinate. (These 

codes have very remarkable properties such as large automorphism 

groups and duality properties - see [HP2003]_, Section 6.6.3-6.6.4.) 

 

INPUT: 

 

 

- ``n`` - an odd prime 

 

- ``F`` - a finite prime field F whose order must be a 

quadratic residue modulo n. 

 

 

OUTPUT: Returns an extended quadratic residue code. 

 

EXAMPLES:: 

 

sage: C1 = codes.QuadraticResidueCode(7,GF(2)) 

sage: C2 = C1.extended_code() 

sage: C3 = codes.ExtendedQuadraticResidueCode(7,GF(2)); C3 

Extension of [7, 4] Cyclic Code over GF(2) 

sage: C2 == C3 

True 

sage: C = codes.ExtendedQuadraticResidueCode(17,GF(2)) 

sage: C 

Extension of [17, 9] Cyclic Code over GF(2) 

sage: C3 = codes.QuadraticResidueCodeOddPair(7,GF(2))[0] 

sage: C3x = C3.extended_code() 

sage: C4 = codes.ExtendedQuadraticResidueCode(7,GF(2)) 

sage: C3x == C4 

True 

 

AUTHORS: 

 

- David Joyner (07-2006) 

""" 

C = QuadraticResidueCodeOddPair(n,F)[0] 

return C.extended_code() 

 

def ExtendedTernaryGolayCode(): 

""" 

This method is now deprecated. 

Please use :class:`sage.coding.golay_code.GolayCode` instead. 

""" 

from sage.misc.superseded import deprecation 

from .golay_code import GolayCode 

deprecation(20787, "codes.ExtendedTernaryGolayCode is now deprecated. Please use codes.GolayCode instead.") 

return GolayCode(GF(3)) 

 

def from_parity_check_matrix(H): 

r""" 

Return the linear code that has ``H`` as a parity check matrix. 

 

If ``H`` has dimensions `h \times n` then the linear code will have 

dimension `n-h` and length `n`. 

 

EXAMPLES:: 

 

sage: C = codes.HammingCode(GF(2), 3); C 

[7, 4] Hamming Code over GF(2) 

sage: H = C.parity_check_matrix(); H 

[1 0 1 0 1 0 1] 

[0 1 1 0 0 1 1] 

[0 0 0 1 1 1 1] 

sage: C2 = codes.from_parity_check_matrix(H); C2 

[7, 4] linear code over GF(2) 

sage: C2.systematic_generator_matrix() == C.systematic_generator_matrix() 

True 

""" 

Cd = LinearCode(H) 

return Cd.dual_code() 

 

LinearCodeFromCheckMatrix = deprecated_function_alias(21165, from_parity_check_matrix) 

 

 

def QuadraticResidueCode(n,F): 

r""" 

A quadratic residue code (or QR code) is a cyclic code whose 

generator polynomial is the product of the polynomials 

`x-\alpha^i` (`\alpha` is a primitive 

`n^{th}` root of unity; `i` ranges over the set of 

quadratic residues modulo `n`). 

 

See QuadraticResidueCodeEvenPair and QuadraticResidueCodeOddPair 

for a more general construction. 

 

INPUT: 

 

 

- ``n`` - an odd prime 

 

- ``F`` - a finite prime field F whose order must be a 

quadratic residue modulo n. 

 

 

OUTPUT: Returns a quadratic residue code. 

 

EXAMPLES:: 

 

sage: C = codes.QuadraticResidueCode(7,GF(2)) 

sage: C 

[7, 4] Cyclic Code over GF(2) 

sage: C = codes.QuadraticResidueCode(17,GF(2)) 

sage: C 

[17, 9] Cyclic Code over GF(2) 

sage: C1 = codes.QuadraticResidueCodeOddPair(7,GF(2))[0] 

sage: C2 = codes.QuadraticResidueCode(7,GF(2)) 

sage: C1 == C2 

True 

sage: C1 = codes.QuadraticResidueCodeOddPair(17,GF(2))[0] 

sage: C2 = codes.QuadraticResidueCode(17,GF(2)) 

sage: C1 == C2 

True 

 

AUTHORS: 

 

- David Joyner (11-2005) 

""" 

return QuadraticResidueCodeOddPair(n,F)[0] 

 

def QuadraticResidueCodeEvenPair(n,F): 

""" 

Quadratic residue codes of a given odd prime length and base ring 

either don't exist at all or occur as 4-tuples - a pair of 

"odd-like" codes and a pair of "even-like" codes. If `n > 2` is prime 

then (Theorem 6.6.2 in [HP2003]_) a QR code exists over `GF(q)` iff q is a 

quadratic residue mod `n`. 

 

They are constructed as "even-like" duadic codes associated the 

splitting (Q,N) mod n, where Q is the set of non-zero quadratic 

residues and N is the non-residues. 

 

EXAMPLES:: 

 

sage: codes.QuadraticResidueCodeEvenPair(17, GF(13)) 

([17, 8] Cyclic Code over GF(13), 

[17, 8] Cyclic Code over GF(13)) 

sage: codes.QuadraticResidueCodeEvenPair(17, GF(2)) 

([17, 8] Cyclic Code over GF(2), 

[17, 8] Cyclic Code over GF(2)) 

sage: codes.QuadraticResidueCodeEvenPair(13,GF(9,"z")) 

([13, 6] Cyclic Code over GF(9), 

[13, 6] Cyclic Code over GF(9)) 

sage: C1,C2 = codes.QuadraticResidueCodeEvenPair(7,GF(2)) 

sage: C1.is_self_orthogonal() 

True 

sage: C2.is_self_orthogonal() 

True 

sage: C3 = codes.QuadraticResidueCodeOddPair(17,GF(2))[0] 

sage: C4 = codes.QuadraticResidueCodeEvenPair(17,GF(2))[1] 

sage: C3.systematic_generator_matrix() == C4.dual_code().systematic_generator_matrix() 

True 

 

This is consistent with Theorem 6.6.9 and Exercise 365 in [HP2003]_. 

 

TESTS:: 

 

sage: codes.QuadraticResidueCodeEvenPair(14,Zmod(4)) 

Traceback (most recent call last): 

... 

ValueError: the argument F must be a finite field 

sage: codes.QuadraticResidueCodeEvenPair(14,GF(2)) 

Traceback (most recent call last): 

... 

ValueError: the argument n must be an odd prime 

sage: codes.QuadraticResidueCodeEvenPair(5,GF(2)) 

Traceback (most recent call last): 

... 

ValueError: the order of the finite field must be a quadratic residue modulo n 

""" 

from sage.arith.srange import srange 

from sage.categories.finite_fields import FiniteFields 

if F not in FiniteFields(): 

raise ValueError("the argument F must be a finite field") 

q = F.order() 

n = Integer(n) 

if n <= 2 or not n.is_prime(): 

raise ValueError("the argument n must be an odd prime") 

Q = quadratic_residues(n); Q.remove(0) # non-zero quad residues 

N = [x for x in srange(1,n) if x not in Q] # non-zero quad non-residues 

if q not in Q: 

raise ValueError("the order of the finite field must be a quadratic residue modulo n") 

return DuadicCodeEvenPair(F,Q,N) 

 

def QuadraticResidueCodeOddPair(n,F): 

""" 

Quadratic residue codes of a given odd prime length and base ring 

either don't exist at all or occur as 4-tuples - a pair of 

"odd-like" codes and a pair of "even-like" codes. If n 2 is prime 

then (Theorem 6.6.2 in [HP2003]_) a QR code exists over GF(q) iff q is a 

quadratic residue mod n. 

 

They are constructed as "odd-like" duadic codes associated the 

splitting (Q,N) mod n, where Q is the set of non-zero quadratic 

residues and N is the non-residues. 

 

EXAMPLES:: 

 

sage: codes.QuadraticResidueCodeOddPair(17, GF(13)) 

([17, 9] Cyclic Code over GF(13), 

[17, 9] Cyclic Code over GF(13)) 

sage: codes.QuadraticResidueCodeOddPair(17, GF(2)) 

([17, 9] Cyclic Code over GF(2), 

[17, 9] Cyclic Code over GF(2)) 

sage: codes.QuadraticResidueCodeOddPair(13, GF(9,"z")) 

([13, 7] Cyclic Code over GF(9), 

[13, 7] Cyclic Code over GF(9)) 

sage: C1 = codes.QuadraticResidueCodeOddPair(17, GF(2))[1] 

sage: C1x = C1.extended_code() 

sage: C2 = codes.QuadraticResidueCodeOddPair(17, GF(2))[0] 

sage: C2x = C2.extended_code() 

sage: C2x.spectrum(); C1x.spectrum() 

[1, 0, 0, 0, 0, 0, 102, 0, 153, 0, 153, 0, 102, 0, 0, 0, 0, 0, 1] 

[1, 0, 0, 0, 0, 0, 102, 0, 153, 0, 153, 0, 102, 0, 0, 0, 0, 0, 1] 

sage: C3 = codes.QuadraticResidueCodeOddPair(7, GF(2))[0] 

sage: C3x = C3.extended_code() 

sage: C3x.spectrum() 

[1, 0, 0, 0, 14, 0, 0, 0, 1] 

 

This is consistent with Theorem 6.6.14 in [HP2003]_. 

 

TESTS:: 

 

sage: codes.QuadraticResidueCodeOddPair(9,GF(2)) 

Traceback (most recent call last): 

... 

ValueError: the argument n must be an odd prime 

""" 

from sage.arith.srange import srange 

from sage.categories.finite_fields import FiniteFields 

if F not in FiniteFields(): 

raise ValueError("the argument F must be a finite field") 

q = F.order() 

n = Integer(n) 

if n <= 2 or not n.is_prime(): 

raise ValueError("the argument n must be an odd prime") 

Q = quadratic_residues(n); Q.remove(0) # non-zero quad residues 

N = [x for x in srange(1,n) if x not in Q] # non-zero quad non-residues 

if q not in Q: 

raise ValueError("the order of the finite field must be a quadratic residue modulo n") 

return DuadicCodeOddPair(F,Q,N) 

 

def random_linear_code(F, length, dimension): 

r""" 

Generate a random linear code of length ``length``, dimension ``dimension`` 

and over the field ``F``. 

 

This function is Las Vegas probabilistic: always correct, usually fast. 

Random matrices over the ``F`` are drawn until one with full rank is hit. 

 

If ``F`` is infinite, the distribution of the elements in the random 

generator matrix will be random according to the distribution of 

``F.random_element()``. 

 

EXAMPLES:: 

 

sage: C = codes.random_linear_code(GF(2), 10, 3) 

sage: C 

[10, 3] linear code over GF(2) 

sage: C.generator_matrix().rank() 

3 

""" 

while True: 

G = random_matrix(F, dimension, length) 

if G.rank() == dimension: 

return LinearCode(G) 

 

def RandomLinearCode(n, k, F): 

r""" 

Deprecated alias of :func:`random_linear_code`. 

 

EXAMPLES:: 

 

sage: C = codes.RandomLinearCode(10, 3, GF(2)) 

doctest:...: DeprecationWarning: codes.RandomLinearCode(n, k, F) is deprecated. Please use codes.random_linear_code(F, n, k) instead 

See http://trac.sagemath.org/21165 for details. 

sage: C 

[10, 3] linear code over GF(2) 

sage: C.generator_matrix().rank() 

3 

""" 

deprecation(21165, "codes.RandomLinearCode(n, k, F) is deprecated. Please use codes.random_linear_code(F, n, k) instead") 

return random_linear_code(F, n, k) 

 

def ReedSolomonCode(n,k,F,pts = None): 

from sage.coding.grs import GeneralizedReedSolomonCode 

deprecation(18928, "codes.ReedSolomonCode is now deprecated. Please use codes.GeneralizedReedSolomonCode instead.") 

q = F.order() 

if n>q or k>n or k>q: 

raise ValueError("RS codes does not exist with the given input.") 

if pts is not None and len(pts) != n: 

raise ValueError("You must provide exactly %s distinct points of %s"%(n,F)) 

if (pts is None): 

pts = [] 

i = 0 

for x in F: 

if i<n: 

pts.append(x) 

i = i+1 

return GeneralizedReedSolomonCode(pts, k) 

 

 

def TernaryGolayCode(): 

""" 

This method is now deprecated. 

Please use :class:`sage.coding.golay_code.GolayCode` instead. 

""" 

from sage.misc.superseded import deprecation 

from .golay_code import GolayCode 

deprecation(20787, "codes.TernaryGolayCode is now deprecated. Please use codes.GolayCode instead.") 

return GolayCode(GF(3), False) 

 

def ToricCode(P,F): 

r""" 

Let `P` denote a list of lattice points in 

`\ZZ^d` and let `T` denote the set of all 

points in `(F^x)^d` (ordered in some fixed way). Put 

`n=|T|` and let `k` denote the dimension of the 

vector space of functions `V = \mathrm{Span}\{x^e \ |\ e \in P\}`. 

The associated toric code `C` is the evaluation code which 

is the image of the evaluation map 

 

.. MATH:: 

 

\mathrm{eval_T} : V \rightarrow F^n, 

 

 

where `x^e` is the multi-index notation 

(`x=(x_1,...,x_d)`, `e=(e_1,...,e_d)`, and 

`x^e = x_1^{e_1}...x_d^{e_d}`), where 

`eval_T (f(x)) = (f(t_1),...,f(t_n))`, and where 

`T=\{t_1,...,t_n\}`. This function returns the toric 

codes discussed in [Joy2004]_. 

 

INPUT: 

 

 

- ``P`` - all the integer lattice points in a polytope 

defining the toric variety. 

 

- ``F`` - a finite field. 

 

 

OUTPUT: Returns toric code with length n = , dimension k over field 

F. 

 

EXAMPLES:: 

 

sage: C = codes.ToricCode([[0,0],[1,0],[2,0],[0,1],[1,1]],GF(7)) 

sage: C 

[36, 5] linear code over GF(7) 

sage: C.minimum_distance() 

24 

sage: C = codes.ToricCode([[-2,-2],[-1,-2],[-1,-1],[-1,0],[0,-1],[0,0],[0,1],[1,-1],[1,0]],GF(5)) 

sage: C 

[16, 9] linear code over GF(5) 

sage: C.minimum_distance() 

6 

sage: C = codes.ToricCode([ [0,0],[1,1],[1,2],[1,3],[1,4],[2,1],[2,2],[2,3],[3,1],[3,2],[4,1]],GF(8,"a")) 

sage: C 

[49, 11] linear code over GF(8) 

 

This is in fact a [49,11,28] code over GF(8). If you type next 

``C.minimum_distance()`` and wait overnight (!), you 

should get 28. 

 

AUTHOR: 

 

- David Joyner (07-2006) 

""" 

from sage.combinat.all import Tuples 

mset = [x for x in F if x!=0] 

d = len(P[0]) 

pts = Tuples(mset,d).list() 

n = len(pts) # (q-1)^d 

k = len(P) 

e = P[0] 

B = [] 

for e in P: 

tmpvar = [prod([t[i]**e[i] for i in range(d)]) for t in pts] 

B.append(tmpvar) 

# now B0 *should* be a full rank matrix 

MS = MatrixSpace(F,k,n) 

return LinearCode(MS(B)) 

 

 

def WalshCode(m): 

r""" 

Return the binary Walsh code of length `2^m`. 

 

The matrix 

of codewords correspond to a Hadamard matrix. This is a (constant 

rate) binary linear `[2^m,m,2^{m-1}]` code. 

 

EXAMPLES:: 

 

sage: C = codes.WalshCode(4); C 

[16, 4] linear code over GF(2) 

sage: C = codes.WalshCode(3); C 

[8, 3] linear code over GF(2) 

sage: C.spectrum() 

[1, 0, 0, 0, 7, 0, 0, 0, 0] 

sage: C.minimum_distance() 

4 

sage: C.minimum_distance(algorithm='gap') # check d=2^(m-1) 

4 

 

REFERENCES: 

 

- :wikipedia:`Hadamard_matrix` 

 

- :wikipedia:`Walsh_code` 

""" 

return LinearCode(walsh_matrix(m), d=2**(m - 1))