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

1453

1454

1455

1456

1457

1458

1459

1460

1461

1462

1463

1464

1465

1466

1467

1468

1469

1470

1471

1472

1473

1474

1475

1476

1477

1478

1479

1480

1481

1482

1483

1484

1485

1486

1487

1488

1489

1490

1491

1492

1493

1494

1495

1496

1497

1498

1499

1500

1501

1502

1503

1504

1505

1506

1507

1508

1509

1510

1511

1512

1513

1514

1515

1516

1517

1518

1519

1520

1521

1522

1523

1524

1525

1526

1527

1528

1529

1530

1531

1532

1533

1534

1535

1536

1537

1538

1539

1540

1541

1542

1543

1544

1545

1546

1547

1548

1549

1550

1551

1552

1553

1554

1555

1556

1557

1558

1559

1560

1561

1562

1563

1564

1565

1566

1567

1568

1569

1570

1571

1572

1573

1574

1575

1576

1577

1578

1579

1580

1581

1582

1583

1584

1585

1586

1587

1588

1589

1590

1591

1592

1593

1594

1595

1596

1597

1598

1599

1600

1601

1602

1603

1604

1605

1606

1607

1608

1609

1610

1611

1612

1613

1614

1615

1616

1617

1618

1619

1620

1621

1622

1623

1624

1625

1626

1627

1628

1629

1630

1631

1632

1633

1634

1635

1636

1637

1638

1639

1640

1641

1642

1643

1644

1645

1646

1647

1648

1649

1650

1651

1652

1653

1654

1655

1656

1657

1658

1659

1660

1661

1662

1663

1664

1665

1666

1667

1668

1669

1670

1671

1672

1673

1674

1675

1676

1677

1678

1679

1680

1681

1682

1683

1684

1685

1686

1687

1688

1689

1690

1691

1692

1693

1694

1695

1696

1697

1698

1699

1700

1701

1702

1703

1704

1705

1706

1707

1708

1709

1710

1711

1712

1713

1714

1715

1716

1717

1718

1719

1720

1721

1722

1723

1724

1725

1726

1727

""" 

Affine curves. 

 

EXAMPLES: 

 

We can construct curves in either an affine plane:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y - x^2], A); C 

Affine Plane Curve over Rational Field defined by -x^2 + y 

 

or in higher dimensional affine space:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = Curve([y - x^2, z - w^3, w - y^4], A); C 

Affine Curve over Rational Field defined by -x^2 + y, -w^3 + z, -y^4 + w 

 

AUTHORS: 

 

- William Stein (2005-11-13) 

 

- David Joyner (2005-11-13) 

 

- David Kohel (2006-01) 

 

- Grayson Jorgenson (2016-8) 

""" 

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

# Copyright (C) 2005 William Stein <wstein@gmail.com> 

# 

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

# 

# The full text of the GPL is available at: 

# 

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

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

from __future__ import absolute_import 

 

from sage.arith.misc import binomial 

from sage.categories.fields import Fields 

from sage.categories.finite_fields import FiniteFields 

from copy import copy 

from sage.categories.homset import Hom, End 

from sage.categories.number_fields import NumberFields 

from sage.interfaces.all import singular 

import sage.libs.singular 

 

from sage.misc.all import add 

 

from sage.rings.all import degree_lowest_rational_function 

 

from sage.rings.number_field.number_field import NumberField 

from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing 

from sage.rings.qqbar import number_field_elements_from_algebraics, QQbar 

from sage.rings.rational_field import is_RationalField 

from sage.schemes.affine.affine_space import (AffineSpace, 

is_AffineSpace) 

from . import point 

 

from sage.schemes.affine.affine_subscheme import AlgebraicScheme_subscheme_affine 

 

from sage.schemes.affine.affine_space import AffineSpace, is_AffineSpace 

from sage.schemes.projective.projective_space import ProjectiveSpace 

 

from .curve import Curve_generic 

 

class AffineCurve(Curve_generic, AlgebraicScheme_subscheme_affine): 

 

_point = point.AffineCurvePoint_field 

 

def _repr_type(self): 

r""" 

Return a string representation of the type of this curve. 

 

EXAMPLES:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = Curve([x - y, z - w, w - x], A) 

sage: C._repr_type() 

'Affine' 

""" 

return "Affine" 

 

def __init__(self, A, X): 

r""" 

Initialization function. 

 

EXAMPLES:: 

 

sage: R.<v> = QQ[] 

sage: K.<u> = NumberField(v^2 + 3) 

sage: A.<x,y,z> = AffineSpace(K, 3) 

sage: C = Curve([z - u*x^2, y^2], A); C 

Affine Curve over Number Field in u with defining polynomial v^2 + 3 

defined by (-u)*x^2 + z, y^2 

 

:: 

 

sage: A.<x,y,z> = AffineSpace(GF(7), 3) 

sage: C = Curve([x^2 - z, z - 8*x], A); C 

Affine Curve over Finite Field of size 7 defined by x^2 - z, -x + z 

""" 

if not is_AffineSpace(A): 

raise TypeError("A (=%s) must be an affine space"%A) 

Curve_generic.__init__(self, A, X) 

d = self.dimension() 

if d != 1: 

raise ValueError("defining equations (=%s) define a scheme of dimension %s != 1"%(X,d)) 

 

def projective_closure(self, i=0, PP=None): 

r""" 

Return the projective closure of this affine curve. 

 

INPUT: 

 

- ``i`` -- (default: 0) the index of the affine coordinate chart of the projective space that the affine 

ambient space of this curve embeds into. 

 

- ``PP`` -- (default: None) ambient projective space to compute the projective closure in. This is 

constructed if it is not given. 

 

OUTPUT: 

 

- a curve in projective space. 

 

EXAMPLES:: 

 

sage: A.<x,y,z> = AffineSpace(QQ, 3) 

sage: C = Curve([y-x^2,z-x^3], A) 

sage: C.projective_closure() 

Projective Curve over Rational Field defined by x1^2 - x0*x2, 

x1*x2 - x0*x3, x2^2 - x1*x3 

 

:: 

 

sage: A.<x,y,z> = AffineSpace(QQ, 3) 

sage: C = Curve([y - x^2, z - x^3], A) 

sage: C.projective_closure() 

Projective Curve over Rational Field defined by 

x1^2 - x0*x2, x1*x2 - x0*x3, x2^2 - x1*x3 

 

:: 

 

sage: A.<x,y> = AffineSpace(CC, 2) 

sage: C = Curve(y - x^3 + x - 1, A) 

sage: C.projective_closure(1) 

Projective Plane Curve over Complex Field with 53 bits of precision defined by 

x0^3 - x0*x1^2 + x1^3 - x1^2*x2 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: P.<u,v,w> = ProjectiveSpace(QQ, 2) 

sage: C = Curve([y - x^2], A) 

sage: C.projective_closure(1, P).ambient_space() == P 

True 

""" 

from .constructor import Curve 

return Curve(AlgebraicScheme_subscheme_affine.projective_closure(self, i, PP)) 

 

def projection(self, indices, AS=None): 

r""" 

Return the projection of this curve onto the coordinates specified by ``indices``. 

 

INPUT: 

 

- ``indices`` -- a list or tuple of distinct integers specifying the indices of the coordinates to use 

in the projection. Can also be a list or tuple consisting of variables of the coordinate ring of the 

ambient space of this curve. If integers are used to specify the coordinates, 0 denotes the first 

coordinate. The length of ``indices`` must be between two and one less than the dimension of the ambient 

space of this curve, inclusive. 

 

- ``AS`` -- (default: None) the affine space the projected curve will be defined in. This space must be 

defined over the same base field as this curve, and must have dimension equal to the length of ``indices``. 

This space is constructed if not specified. 

 

OUTPUT: 

 

- a tuple consisting of two elements: a scheme morphism from this curve to affine space of dimension 

equal to the number of coordinates specified in ``indices``, and the affine subscheme that is the image 

of that morphism. If the image is a curve, the second element of the tuple will be a curve. 

 

EXAMPLES:: 

 

sage: A.<x,y,z> = AffineSpace(QQ, 3) 

sage: C = Curve([y^7 - x^2 + x^3 - 2*z, z^2 - x^7 - y^2], A) 

sage: C.projection([0,1]) 

(Scheme morphism: 

From: Affine Curve over Rational Field defined by y^7 + x^3 - x^2 - 

2*z, -x^7 - y^2 + z^2 

To: Affine Space of dimension 2 over Rational Field 

Defn: Defined on coordinates by sending (x, y, z) to 

(x, y), 

Affine Plane Curve over Rational Field defined by x1^14 + 2*x0^3*x1^7 - 

2*x0^2*x1^7 - 4*x0^7 + x0^6 - 2*x0^5 + x0^4 - 4*x1^2) 

sage: C.projection([0,1,3,4]) 

Traceback (most recent call last): 

... 

ValueError: (=[0, 1, 3, 4]) must be a list or tuple of length between 2 

and (=2), inclusive 

 

:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = Curve([x - 2, y - 3, z - 1], A) 

sage: B.<a,b,c> = AffineSpace(QQ, 3) 

sage: C.projection([0,1,2], AS=B) 

(Scheme morphism: 

From: Affine Curve over Rational Field defined by x - 2, y - 3, z - 1 

To: Affine Space of dimension 3 over Rational Field 

Defn: Defined on coordinates by sending (x, y, z, w) to 

(x, y, z), 

Closed subscheme of Affine Space of dimension 3 over Rational Field 

defined by: 

c - 1, 

b - 3, 

a - 2) 

 

:: 

 

sage: A.<x,y,z,w,u> = AffineSpace(GF(11), 5) 

sage: C = Curve([x^3 - 5*y*z + u^2, x - y^2 + 3*z^2, w^2 + 2*u^3*y, y - u^2 + z*x], A) 

sage: B.<a,b,c> = AffineSpace(GF(11), 3) 

sage: proj1 = C.projection([1,2,4], AS=B) 

sage: proj1 

(Scheme morphism: 

From: Affine Curve over Finite Field of size 11 defined by x^3 - 

5*y*z + u^2, -y^2 + 3*z^2 + x, 2*y*u^3 + w^2, x*z - u^2 + y 

To: Affine Space of dimension 3 over Finite Field of size 11 

Defn: Defined on coordinates by sending (x, y, z, w, u) to 

(y, z, u), 

Affine Curve over Finite Field of size 11 defined by a^2*b - 3*b^3 - 

c^2 + a, c^6 - 5*a*b^4 + b^3*c^2 - 3*a*c^4 + 3*a^2*c^2 - a^3, a^2*c^4 - 

3*b^2*c^4 - 2*a^3*c^2 - 5*a*b^2*c^2 + a^4 - 5*a*b^3 + 2*b^4 + b^2*c^2 - 

3*b*c^2 + 3*a*b, a^4*c^2 + 2*b^4*c^2 - a^5 - 2*a*b^4 + 5*b*c^4 + a*b*c^2 

- 5*a*b^2 + 4*b^3 + b*c^2 + 5*c^2 - 5*a, a^6 - 5*b^6 - 5*b^3*c^2 + 

5*a*b^3 + 2*c^4 - 4*a*c^2 + 2*a^2 - 5*a*b + c^2) 

sage: proj1[1].ambient_space() is B 

True 

sage: proj2 = C.projection([1,2,4]) 

sage: proj2[1].ambient_space() is B 

False 

sage: C.projection([1,2,3,5], AS=B) 

Traceback (most recent call last): 

... 

TypeError: (=Affine Space of dimension 3 over Finite Field of size 11) 

must have dimension (=4) 

 

:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = A.curve([x*y - z^3, x*z - w^3, w^2 - x^3]) 

sage: C.projection([y,z]) 

(Scheme morphism: 

From: Affine Curve over Rational Field defined by -z^3 + x*y, -w^3 + 

x*z, -x^3 + w^2 

To: Affine Space of dimension 2 over Rational Field 

Defn: Defined on coordinates by sending (x, y, z, w) to 

(y, z), 

Affine Plane Curve over Rational Field defined by x1^23 - x0^7*x1^4) 

sage: B.<x,y,z> = AffineSpace(QQ, 3) 

sage: C.projection([x,y,z], AS=B) 

(Scheme morphism: 

From: Affine Curve over Rational Field defined by -z^3 + x*y, -w^3 + 

x*z, -x^3 + w^2 

To: Affine Space of dimension 3 over Rational Field 

Defn: Defined on coordinates by sending (x, y, z, w) to 

(x, y, z), 

Affine Curve over Rational Field defined by z^3 - x*y, x^8 - x*z^2, 

x^7*z^2 - x*y*z) 

sage: C.projection([y,z,z]) 

Traceback (most recent call last): 

... 

ValueError: (=[y, z, z]) must be a list or tuple of distinct indices or 

variables 

""" 

AA = self.ambient_space() 

n = AA.dimension_relative() 

if n == 2: 

raise TypeError("this curve is already a plane curve") 

if self.base_ring() not in Fields(): 

raise TypeError("this curve must be defined over a field") 

if len(indices) < 2 or len(indices) >= n: 

raise ValueError("(=%s) must be a list or tuple of length between 2 and (=%s), inclusive"%(indices,n - 1)) 

if len(set(indices)) < len(indices): 

raise ValueError("(=%s) must be a list or tuple of distinct indices or variables"%indices) 

if not AS is None: 

if not is_AffineSpace(AS): 

raise TypeError("(=%s) must be an affine space"%AS) 

if AS.dimension_relative() != len(indices): 

raise TypeError("(=%s) must have dimension (=%s)"%(AS, len(indices))) 

if AS.base_ring() != AA.base_ring(): 

raise TypeError("(=%s) must be defined over the same base field as this curve"%AS) 

indices = list(indices) 

if all([f in AA.gens() for f in indices]): 

indices = [AA.gens().index(f) for f in indices] 

indices.sort() 

else: 

indices = [int(i) for i in indices] # type checking 

indices.sort() 

if indices[0] < 0 or indices[-1] > n - 1: 

raise ValueError("index values must be between 0 and one minus the dimension of the ambient space " \ 

"of this curve") 

# construct the projection map 

if AS is None: 

AA2 = AffineSpace(self.base_ring(), len(indices)) 

else: 

AA2 = AS 

H = Hom(self, AA2) 

psi = H([AA.gens()[i] for i in indices]) 

# compute the image via elimination 

removecoords = list(AA.gens()) 

for i in range(len(indices) - 1, -1, -1): 

removecoords.pop(indices[i]) 

J = self.defining_ideal().elimination_ideal(removecoords) 

K = Hom(AA.coordinate_ring(), AA2.coordinate_ring()) 

l = [0]*(n) 

for i in range(len(indices)): 

l[indices[i]] = AA2.gens()[i] 

phi = K(l) 

G = [phi(f) for f in J.gens()] 

try: 

C = AA2.curve(G) 

except (TypeError, ValueError): 

C = AA2.subscheme(G) 

return tuple([psi, C]) 

 

def plane_projection(self, AP=None): 

r""" 

Return a projection of this curve into an affine plane so that the image of the projection is 

a plane curve. 

 

INPUT: 

 

- ``AP`` -- (default: None) the affine plane to project this curve into. This space must be defined over 

the same base field as this curve, and must have dimension two. This space will be constructed if not 

specified. 

 

OUTPUT: 

 

- a tuple consisting of two elements: a scheme morphism from this curve into an affine plane, and the plane 

curve that defines the image of that morphism. 

 

EXAMPLES:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = Curve([x^2 - y*z*w, z^3 - w, w + x*y - 3*z^3], A) 

sage: C.plane_projection() 

(Scheme morphism: 

From: Affine Curve over Rational Field defined by -y*z*w + x^2, z^3 - 

w, -3*z^3 + x*y + w 

To: Affine Space of dimension 2 over Rational Field 

Defn: Defined on coordinates by sending (x, y, z, w) to 

(x, y), Affine Plane Curve over Rational Field defined by 

x0^2*x1^7 - 16*x0^4) 

 

:: 

 

sage: R.<a> = QQ[] 

sage: K.<b> = NumberField(a^2 + 2) 

sage: A.<x,y,z> = AffineSpace(K, 3) 

sage: C = A.curve([x - b, y - 2]) 

sage: B.<a,b> = AffineSpace(K, 2) 

sage: proj1 = C.plane_projection(AP=B) 

sage: proj1 

(Scheme morphism: 

From: Affine Curve over Number Field in b with defining polynomial 

a^2 + 2 defined by x + (-b), y - 2 

To: Affine Space of dimension 2 over Number Field in b with 

defining polynomial a^2 + 2 

Defn: Defined on coordinates by sending (x, y, z) to 

(x, z), 

Affine Plane Curve over Number Field in b with defining polynomial a^2 

+ 2 defined by a + (-b)) 

sage: proj1[1].ambient_space() is B 

True 

sage: proj2 = C.plane_projection() 

sage: proj2[1].ambient_space() is B 

False 

""" 

n = self.ambient_space().dimension_relative() 

# finds a projection that will have a plane curve as its image 

# the following iterates over all pairs (i,j) with j > i to test all 

# possible projections 

for i in range(0, n - 1): 

for j in range(i + 1, n): 

L = self.projection([i,j], AP) 

if isinstance(L[1], Curve_generic): 

return L 

 

def blowup(self, P=None): 

r""" 

Return the blow up of this affine curve at the point ``P``. 

 

The blow up is described by affine charts. This curve must be irreducible. 

 

INPUT: 

 

- ``P`` -- (default: None) a point on this curve at which to blow up. If ``None``, then ``P`` is 

taken to be the origin. 

 

OUTPUT: 

 

- a tuple consisting of three elements. The first is a tuple of curves in affine space of 

the same dimension as the ambient space of this curve, which define the blow up in 

each affine chart. The second is a tuple of tuples such that the jth element of the ith 

tuple is the transition map from the ith affine patch to the jth affine patch. Lastly, 

the third element is a tuple consisting of the restrictions of the projection map from 

the blow up back to the original curve, restricted to each affine patch. There the 

ith element will be the projection from the ith affine patch. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y^2 - x^3], A) 

sage: C.blowup() 

((Affine Plane Curve over Rational Field defined by s1^2 - x, 

Affine Plane Curve over Rational Field defined by y*s0^3 - 1), 

([Scheme endomorphism of Affine Plane Curve over Rational Field defined by s1^2 - x 

Defn: Defined on coordinates by sending (x, s1) to 

(x, s1), Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by s1^2 - x 

To: Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

Defn: Defined on coordinates by sending (x, s1) to 

(x*s1, 1/s1)], [Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

To: Affine Plane Curve over Rational Field defined by s1^2 - x 

Defn: Defined on coordinates by sending (y, s0) to 

(y*s0, 1/s0), 

Scheme endomorphism of Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

Defn: Defined on coordinates by sending (y, s0) to 

(y, s0)]), 

(Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by s1^2 - x 

To: Affine Plane Curve over Rational Field defined by -x^3 + y^2 

Defn: Defined on coordinates by sending (x, s1) to 

(x, x*s1), Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

To: Affine Plane Curve over Rational Field defined by -x^3 + y^2 

Defn: Defined on coordinates by sending (y, s0) to 

(y*s0, y))) 

 

:: 

 

sage: K.<a> = QuadraticField(2) 

sage: A.<x,y,z> = AffineSpace(K, 3) 

sage: C = Curve([y^2 - a*x^5, x - z], A) 

sage: B = C.blowup() 

sage: B[0] 

(Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s2 - 1, 

2*x^3 + (-a)*s1^2, 

Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s0 - s2, 

2*y^3*s2^5 + (-a), 

Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s0 - 1, 

2*z^3 + (-a)*s1^2) 

sage: B[1][0][2] 

Scheme morphism: 

From: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s2 - 1, 

2*x^3 + (-a)*s1^2 

To: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s0 - 1,  

2*z^3 + (-a)*s1^2 

Defn: Defined on coordinates by sending (x, s1, s2) to 

(x*s2, 1/s2, s1/s2) 

sage: B[1][2][0] 

Scheme morphism: 

From: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s0 - 1, 

2*z^3 + (-a)*s1^2 

To: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by s2 - 1, 

2*x^3 + (-a)*s1^2 

Defn: Defined on coordinates by sending (z, s0, s1) to 

(z*s0, s1/s0, 1/s0) 

sage: B[2] 

(Scheme morphism: 

From: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by 

s2 - 1, 2*x^3 + (-a)*s1^2 

To: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by 

(-a)*x^5 + y^2, x - z 

Defn: Defined on coordinates by sending (x, s1, s2) to 

(x, x*s1, x*s2), Scheme morphism: 

From: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by 

s0 - s2, 2*y^3*s2^5 + (-a) 

To: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by 

(-a)*x^5 + y^2, x - z 

Defn: Defined on coordinates by sending (y, s0, s2) to 

(y*s0, y, y*s2), Scheme morphism: 

From: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by 

s0 - 1, 2*z^3 + (-a)*s1^2 

To: Affine Curve over Number Field in a with defining polynomial x^2 - 2 defined by 

(-a)*x^5 + y^2, x - z 

Defn: Defined on coordinates by sending (z, s0, s1) to 

(z*s0, z*s1, z)) 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve((y - 3/2)^3 - (x + 2)^5 - (x + 2)^6) 

sage: Q = A([-2,3/2]) 

sage: C.blowup(Q) 

((Affine Plane Curve over Rational Field defined by x^3 - s1^3 + 7*x^2 + 16*x + 12, 

Affine Plane Curve over Rational Field defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5 + 

54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8), 

([Scheme endomorphism of Affine Plane Curve over Rational Field defined by x^3 - s1^3 + 7*x^2 + 

16*x + 12 

Defn: Defined on coordinates by sending (x, s1) to 

(x, s1), Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by x^3 - s1^3 + 7*x^2 + 16*x + 12 

To: Affine Plane Curve over Rational Field defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5 + 

54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8 

Defn: Defined on coordinates by sending (x, s1) to 

(x*s1 + 2*s1 + 3/2, 1/s1)], [Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5 + 

54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8 

To: Affine Plane Curve over Rational Field defined by x^3 - s1^3 + 7*x^2 + 16*x + 12 

Defn: Defined on coordinates by sending (y, s0) to 

(y*s0 - 3/2*s0 - 2, 1/s0), 

Scheme endomorphism of Affine Plane Curve over Rational Field defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 

8*y^2*s0^5 + 54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8 

Defn: Defined on coordinates by sending (y, s0) to 

(y, s0)]), 

(Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by x^3 - s1^3 + 7*x^2 + 16*x + 12 

To: Affine Plane Curve over Rational Field defined by -x^6 - 13*x^5 - 70*x^4 - 200*x^3 + y^3 - 

320*x^2 - 9/2*y^2 - 272*x + 27/4*y - 795/8 

Defn: Defined on coordinates by sending (x, s1) to 

(x, x*s1 + 2*s1 + 3/2), Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by 8*y^3*s0^6 - 36*y^2*s0^6 + 8*y^2*s0^5 + 

54*y*s0^6 - 24*y*s0^5 - 27*s0^6 + 18*s0^5 - 8 

To: Affine Plane Curve over Rational Field defined by -x^6 - 13*x^5 - 70*x^4 - 200*x^3 + y^3 - 

320*x^2 - 9/2*y^2 - 272*x + 27/4*y - 795/8 

Defn: Defined on coordinates by sending (y, s0) to 

(y*s0 - 3/2*s0 - 2, y))) 

 

:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = A.curve([((x + 1)^2 + y^2)^3 - 4*(x + 1)^2*y^2, y - z, w - 4]) 

sage: Q = C([-1,0,0,4]) 

sage: B = C.blowup(Q) 

sage: B[0] 

(Affine Curve over Rational Field defined by s3, s1 - s2, x^2*s2^6 + 

2*x*s2^6 + 3*x^2*s2^4 + s2^6 + 6*x*s2^4 + 3*x^2*s2^2 + 3*s2^4 + 6*x*s2^2 

+ x^2 - s2^2 + 2*x + 1, 

Affine Curve over Rational Field defined by s3, s2 - 1, y^2*s0^6 + 

3*y^2*s0^4 + 3*y^2*s0^2 + y^2 - 4*s0^2, 

Affine Curve over Rational Field defined by s3, s1 - 1, z^2*s0^6 + 

3*z^2*s0^4 + 3*z^2*s0^2 + z^2 - 4*s0^2, 

Closed subscheme of Affine Space of dimension 4 over Rational Field 

defined by: 

1) 

sage: Q = A([6,2,3,1]) 

sage: B = C.blowup(Q) 

Traceback (most recent call last): 

... 

TypeError: (=(6, 2, 3, 1)) must be a point on this curve 

 

:: 

 

sage: A.<x,y> = AffineSpace(QuadraticField(-1), 2) 

sage: C = A.curve([y^2 + x^2]) 

sage: C.blowup() 

Traceback (most recent call last): 

... 

TypeError: this curve must be irreducible 

""" 

A = self.ambient_space() 

n = A.dimension_relative() 

if P is None: 

P = A([0]*n) 

try: 

self(P) 

except TypeError: 

raise TypeError("(=%s) must be a point on this curve"%P) 

if not self.base_ring() in Fields(): 

raise TypeError("the base ring of this curve must be a field") 

if not self.is_irreducible(): 

raise TypeError("this curve must be irreducible") 

# attempt to make the variable names more organized 

# the convention used here is to have the homogeneous coordinates for the projective component of the 

# product space the blow up resides in be generated from the letter 's'. The following loop is in place 

# to prevent conflicts in the names from occurring 

rf = 1 

for i in range(n): 

if str(A.gens()[i])[0] == 's' and len(str(A.gens()[i])) > rf: 

rf = len(str(A.gens()[i])) 

var_names = [str(A.gens()[i]) for i in range(n)] + ['s'*rf + str(i) for i in range(n)] 

R = PolynomialRing(A.base_ring(), 2*n, var_names) 

# move the defining polynomials of this curve into R 

H = Hom(A.coordinate_ring(), R) 

psi = H([R.gens()[i] for i in range(n)]) 

n_polys = [psi(f) for f in self.defining_polynomials()] 

# the blow up ideal of A at P is the ideal generated by 

# (z_i - p_i)*s_j - (z_j - p_j)*s_i for i != j from 0,...,n-1 

# in the mixed product space of A^n and P^{n-1} where the z_i are the gens 

# of A^n, the s_i are the gens for P^{n-1}, and P = (p_1,...,p_n). We describe the 

# blow up of this curve at P in each affine chart 

patches = [] 

for i in range(n): 

# in this chart, s_i is assumed to be 1 

# substitute in z_j = (z_i - p_i)*s_j + p_j for each j != i 

coords = list(R.gens()) 

for j in range(n): 

if j != i: 

coords[j] = (R.gens()[i] - P[i])*R.gens()[j + n] + P[j] 

c_polys = [f(coords) for f in n_polys] 

var_names = list(R.gens())[n:2*n] 

var_names.pop(i) 

var_names.insert(0, R.gens()[i]) 

c_A = AffineSpace(R.base_ring(), n, var_names) 

H = Hom(R, c_A.coordinate_ring()) 

coords = [0]*(2*n) 

coords[i] = c_A.gens()[0] 

t = 1 

for j in range(n): 

if j != i: 

coords[j + n] = c_A.gens()[t] 

t = t + 1 

else: 

coords[j + n] = 1 

psi = H(coords) 

c_polys = [psi(f) for f in c_polys] 

# choose the component of the subscheme defined by these polynomials 

# that corresponds to the proper transform 

irr_comps = c_A.subscheme(c_polys).irreducible_components() 

for j in range(len(irr_comps)): 

proper_transform = True 

for f in irr_comps[j].defining_polynomials(): 

if (c_A.gens()[0] - P[i]).divides(f): 

proper_transform = False 

break 

if proper_transform: 

patches.append(c_A.curve(irr_comps[j].defining_polynomials())) 

break 

elif j + 1 == len(irr_comps): 

# patch of blowup in this chart is empty 

patches.append(c_A.subscheme(1)) 

# create the transition maps between the charts 

t_maps = [] 

for i in range(n): 

maps = [] 

for j in range(n): 

AA = patches[i].ambient_space() 

H = Hom(patches[i], patches[j]) 

vars = AA.gens() 

homvars = list(AA.gens()) 

homvars.pop(0) 

homvars.insert(i, 1) 

coords = [(vars[0] - P[i])*homvars[j] + P[j]] 

for t in range(n): 

if t != j: 

coords.append(homvars[t]/homvars[j]) 

maps.append(H(coords)) 

t_maps.append(maps) 

# create the restrictions of the projection map 

proj_maps = [] 

for i in range(n): 

p_A = patches[i].ambient_space() 

H = Hom(patches[i], self) 

homvars = list(p_A.gens())[1:n] 

homvars.insert(i, 1) 

coords = [(p_A.gens()[0] - P[i])*homvars[j] + P[j] for j in range(n)] 

proj_maps.append(H(coords)) 

return tuple([tuple(patches), tuple(t_maps), tuple(proj_maps)]) 

 

def resolution_of_singularities(self, extend=False): 

r""" 

Return a nonsingular model for this affine curve created by blowing up its singular points. 

 

The nonsingular model is given as a collection of affine patches that cover it. If ``extend`` is ``False`` 

and if the base field is a number field, or if the base field is a finite field, the model returned may have 

singularities with coordinates not contained in the base field. An error is returned if this curve is already 

nonsingular, or if it has no singular points over its base field. This curve must be irreducible, and must be 

defined over a number field or finite field. 

 

INPUT: 

 

- ``extend`` -- (default: False) specifies whether to extend the base field when necessary to find all 

singular points when this curve is defined over a number field. If ``extend`` is ``False``, then only 

singularities with coordinates in the base field of this curve will be resolved. However, setting 

``extend`` to ``True`` will slow down computations. 

 

OUTPUT: 

 

- a tuple consisting of three elements. The first is a tuple of curves in affine space of 

the same dimension as the ambient space of this curve, which represent affine patches 

of the resolution of singularities. The second is a tuple of tuples such that the jth 

element of the ith tuple is the transition map from the ith patch to the jth patch. Lastly, 

the third element is a tuple consisting of birational maps from the patches back to the 

original curve that were created by composing the projection maps generated from the blow up 

computations. There the ith element will be a map from the ith patch. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y^2 - x^3], A) 

sage: C.resolution_of_singularities() 

((Affine Plane Curve over Rational Field defined by s1^2 - x, 

Affine Plane Curve over Rational Field defined by y*s0^3 - 1), 

((Scheme endomorphism of Affine Plane Curve over Rational Field defined by s1^2 - x 

Defn: Defined on coordinates by sending (x, s1) to 

(x, s1), Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by s1^2 - x 

To: Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

Defn: Defined on coordinates by sending (x, s1) to 

(x*s1, 1/s1)), (Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

To: Affine Plane Curve over Rational Field defined by s1^2 - x 

Defn: Defined on coordinates by sending (y, s0) to 

(y*s0, 1/s0), 

Scheme endomorphism of Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

Defn: Defined on coordinates by sending (y, s0) to 

(y, s0))), 

(Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by s1^2 - x 

To: Affine Plane Curve over Rational Field defined by -x^3 + y^2 

Defn: Defined on coordinates by sending (x, s1) to 

(x, x*s1), Scheme morphism: 

From: Affine Plane Curve over Rational Field defined by y*s0^3 - 1 

To: Affine Plane Curve over Rational Field defined by -x^3 + y^2 

Defn: Defined on coordinates by sending (y, s0) to 

(y*s0, y))) 

 

:: 

 

sage: set_verbose(-1) 

sage: K.<a> = QuadraticField(3) 

sage: A.<x,y> = AffineSpace(K, 2) 

sage: C = A.curve(x^4 + 2*x^2 + a*y^3 + 1) 

sage: C.resolution_of_singularities(extend=True)[0] # long time (2 seconds) 

(Affine Plane Curve over Number Field in a0 with defining polynomial y^4 - 4*y^2 + 16 defined by 

24*x^2*ss1^3 + 24*ss1^3 + (a0^3 - 8*a0), 

Affine Plane Curve over Number Field in a0 with defining polynomial y^4 - 4*y^2 + 16 defined by 

24*s1^2*ss0 + (a0^3 - 8*a0)*ss0^2 + (-6*a0^3)*s1, 

Affine Plane Curve over Number Field in a0 with defining polynomial y^4 - 4*y^2 + 16 defined by 

8*y^2*s0^4 + (4*a0^3)*y*s0^3 - 32*s0^2 + (a0^3 - 8*a0)*y) 

 

:: 

 

sage: A.<x,y,z> = AffineSpace(GF(5), 3) 

sage: C = Curve([y - x^3, (z - 2)^2 - y^3 - x^3], A) 

sage: R = C.resolution_of_singularities() 

sage: R[0] 

(Affine Curve over Finite Field of size 5 defined by x^2 - s1, s1^4 - x*s2^2 + s1, x*s1^3 - s2^2 + x, 

Affine Curve over Finite Field of size 5 defined by y*s2^2 - y^2 - 1, s2^4 - s0^3 - y^2 - 2, y*s0^3 

- s2^2 + y, Affine Curve over Finite Field of size 5 defined by s0^3*s1 + z*s1^3 + s1^4 - 2*s1^3 - 1, 

z*s0^3 + z*s1^3 - 2*s0^3 - 2*s1^3 - 1, z^2*s1^3 + z*s1^3 - s1^3 - z + s1 + 2) 

 

:: 

 

sage: A.<x,y,z,w> = AffineSpace(QQ, 4) 

sage: C = A.curve([((x - 2)^2 + y^2)^2 - (x - 2)^2 - y^2 + (x - 2)^3, z - y - 7, w - 4]) 

sage: B = C.resolution_of_singularities() 

sage: B[0] 

(Affine Curve over Rational Field defined by s3, s1 - s2, x^2*s2^4 - 

4*x*s2^4 + 2*x^2*s2^2 + 4*s2^4 - 8*x*s2^2 + x^2 + 7*s2^2 - 3*x + 1, 

Affine Curve over Rational Field defined by s3, s2 - 1, y^2*s0^4 + 

2*y^2*s0^2 + y*s0^3 + y^2 - s0^2 - 1, 

Affine Curve over Rational Field defined by s3, s1 - 1, z^2*s0^4 - 

14*z*s0^4 + 2*z^2*s0^2 + z*s0^3 + 49*s0^4 - 28*z*s0^2 - 7*s0^3 + z^2 + 

97*s0^2 - 14*z + 48, 

Closed subscheme of Affine Space of dimension 4 over Rational Field 

defined by: 

1) 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y - x^2 + 1], A) 

sage: C.resolution_of_singularities() 

Traceback (most recent call last): 

... 

TypeError: this curve is already nonsingular 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve([(x^2 + y^2 - y - 2)*(y - x^2 + 2) + y^3]) 

sage: C.resolution_of_singularities() 

Traceback (most recent call last): 

... 

TypeError: this curve has no singular points over its base field. If 

working over a number field use extend=True 

""" 

# helper function for extending the base field (in the case of working over a number field) 

def extension(self): 

F = self.base_ring() 

pts = self.change_ring(F.embeddings(QQbar)[0]).rational_points() 

L = [t for pt in pts for t in pt] 

K = number_field_elements_from_algebraics(L)[0] 

if is_RationalField(K): 

return F.embeddings(F)[0] 

else: 

if is_RationalField(F): 

return F.embeddings(K)[0] 

else: 

# make sure the defining polynomial variable names are the same for K, N 

N = NumberField(K.defining_polynomial().parent()(F.defining_polynomial()), str(K.gen())) 

return N.composite_fields(K, both_maps=True)[0][1]*F.embeddings(N)[0] 

# find the set of singular points of this curve 

# in the case that the base field is a number field, extend it as needed (if extend == True) 

C = self 

n = C.ambient_space().dimension_relative() 

if not self.is_irreducible(): 

raise TypeError("this curve must be irreducible") 

if not (self.base_ring() in NumberFields() or self.base_ring() in FiniteFields()): 

raise NotImplementedError("this curve must be defined over either a number field or a finite field") 

if C.base_ring() in NumberFields() and extend: 

C = C.change_ring(extension(C.singular_subscheme())) 

H = End(C) 

placeholder = H(C.ambient_space().gens()) 

# the list res holds the data for the patches of the resolution of singularities 

# each element is a list consisting of the curve defining the patch, a list 

# of the transition maps from that patch to the other patches, a projection 

# map from the patch to the original curve, and the set of singular points 

# of the patch 

res = [[C, [placeholder], placeholder, C.singular_points()]] 

if len(res[0][3]) == 0: 

if C.is_smooth(): 

raise TypeError("this curve is already nonsingular") 

else: 

raise TypeError("this curve has no singular points over its base field. If working over"\ 

" a number field use extend=True") 

not_resolved = True 

t = 0 

# loop through the patches and blow up each until no patch has singular points 

while not_resolved: 

[BC, t_maps, pi, pts] = [res[t][0], res[t][1], res[t][2], res[t][3]] 

# check if there are any singular points in this patch 

if len(pts) == 0: 

t = t + 1 

if t == len(res): 

not_resolved = False 

continue 

# the identity map should be replaced for each of the charts of the blow up 

t_maps.pop(t) 

# blow up pts[0] 

B = list(BC.blowup(pts[0])) 

B = [list(B[0]), [list(B[1][i]) for i in range(len(B[1]))], list(B[2])] 

# the t-th element of res will be replaced with the new data corresponding to the charts 

# of the blow up 

res.pop(t) 

# take out the transition maps from the other resolution patches to the t-th patch 

old_maps = [res[i][1].pop(t) for i in range(len(res))] 

patches_to_add = [] 

# generate the needed data for each patch of the blow up 

for i in range(len(B[0])): 

# check if there are any singular points where this patch meets the exceptional divisor 

AA = AffineSpace(B[0][i].base_ring(), n - 1, 'x') 

coords = [pts[0][i]] 

coords.extend(list(AA.gens())) 

H = Hom(B[0][i].ambient_space().coordinate_ring(), AA.coordinate_ring()) 

poly_hom = H(coords) 

X = AA.subscheme([poly_hom(f) for f in B[0][i].defining_polynomials()]) 

# in the case of working over a number field, it might be necessary to extend the base 

# field in order to find all intersection points 

n_pts = [] 

if B[0][i].base_ring() in NumberFields() and extend: 

emb = extension(X) 

X = X.change_ring(emb) 

tmp_curve = B[0][i].change_ring(emb) 

for pt in X.rational_points(): 

tmp_pt = tmp_curve([pts[0][i]] + list(pt)) 

if tmp_curve.is_singular(tmp_pt): 

n_pts.append(tmp_pt) 

# avoid needlessly extending the base field 

if len(n_pts) > 0: 

# coerce everything to the new base field 

BC = BC.change_ring(emb) 

t_maps = [t_maps[j].change_ring(emb) for j in range(len(t_maps))] 

old_maps = [old_maps[j].change_ring(emb) for j in range(len(old_maps))] 

pi = pi.change_ring(emb) 

pts = [pt.change_ring(emb) for pt in pts] 

# coerce the current blow up data 

for j in range(len(B[0])): 

B[0][j] = B[0][j].change_ring(emb) 

for j in range(len(B[1])): 

for k in range(len(B[1])): 

B[1][j][k] = B[1][j][k].change_ring(emb) 

for j in range(len(B[2])): 

B[2][j] = B[2][j].change_ring(emb) 

# coerce the other data in res 

for j in range(len(res)): 

res[j][0] = res[j][0].change_ring(emb) 

for k in range(len(res[j][1])): 

res[j][1][k] = res[j][1][k].change_ring(emb) 

res[j][2].change_ring(emb) 

for k in range(len(res[j][3])): 

res[j][3][k] = res[j][3][k].change_ring(emb) 

else: 

for pt in X.rational_points(): 

tmp_pt = B[0][i]([pts[0][i]] + list(pt)) 

if B[0][i].is_singular(tmp_pt): 

n_pts.append(tmp_pt) 

b_data = [B[0][i]] 

# projection map and its inverse 

t_pi = B[2][i] 

coords = [(BC.ambient_space().gens()[j] - pts[0][j])/(BC.ambient_space().gens()[i] - pts[0][i]) for\ 

j in range(n)] 

coords.pop(i) 

coords.insert(0, BC.ambient_space().gens()[i]) 

H = Hom(BC, B[0][i]) 

t_pi_inv = H(coords) 

# compose the current transition maps from the original curve to the other patches 

# with the projection map 

L = list(t_maps) 

for j in range(len(t_maps)): 

L[j] = L[j]*t_pi 

for j in range(len(B[1][i])): 

L.insert(t + j, B[1][i][j]) 

b_data.append(L) 

# update transition maps of each other element of res 

for j in range(len(res)): 

new_t_map = t_pi_inv*old_maps[j] 

res[j][1].insert(t + i, new_t_map) 

# create the projection map 

b_data.append(pi*t_pi) 

# singular points 

# translate the singular points of the parent patch (other than that which was the center of the 

# blow up) by the inverse of the first projection map 

for j in range(1, len(pts)): 

# make sure this point is in this chart before attempting to map it 

try: 

n_pts.append(t_pi_inv(BC(pts[j]))) 

except (TypeError, ZeroDivisionError): 

pass 

b_data.append(n_pts) 

patches_to_add.append(b_data) 

for i in range(len(patches_to_add)): 

res.insert(t + i, patches_to_add[i]) 

t = 0 

patches = [res[i][0] for i in range(len(res))] 

t_maps = [tuple(res[i][1]) for i in range(len(res))] 

p_maps = [res[i][2] for i in range(len(res))] 

return tuple([tuple(patches), tuple(t_maps), tuple(p_maps)]) 

 

class AffinePlaneCurve(AffineCurve): 

 

_point = point.AffinePlaneCurvePoint_field 

 

def __init__(self, A, f): 

r""" 

Initialization function. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([x^3 - y^2], A); C 

Affine Plane Curve over Rational Field defined by x^3 - y^2 

 

:: 

 

sage: A.<x,y> = AffineSpace(CC, 2) 

sage: C = Curve([y^2 + x^2], A); C 

Affine Plane Curve over Complex Field with 53 bits of precision defined 

by x^2 + y^2 

""" 

if not (is_AffineSpace(A) and A.dimension != 2): 

raise TypeError("Argument A (= %s) must be an affine plane."%A) 

Curve_generic.__init__(self, A, [f]) 

 

def _repr_type(self): 

r""" 

Return a string representation of the type of this curve. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y - 7/2*x^5 + x - 3], A) 

sage: C._repr_type() 

'Affine Plane' 

""" 

return "Affine Plane" 

 

def divisor_of_function(self, r): 

""" 

Return the divisor of a function on a curve. 

 

INPUT: r is a rational function on X 

 

OUTPUT: 

 

 

- ``list`` - The divisor of r represented as a list of 

coefficients and points. (TODO: This will change to a more 

structural output in the future.) 

 

 

EXAMPLES:: 

 

sage: F = GF(5) 

sage: P2 = AffineSpace(2, F, names = 'xy') 

sage: R = P2.coordinate_ring() 

sage: x, y = R.gens() 

sage: f = y^2 - x^9 - x 

sage: C = Curve(f) 

sage: K = FractionField(R) 

sage: r = 1/x 

sage: C.divisor_of_function(r) # todo: not implemented (broken) 

[[-1, (0, 0, 1)]] 

sage: r = 1/x^3 

sage: C.divisor_of_function(r) # todo: not implemented (broken) 

[[-3, (0, 0, 1)]] 

""" 

F = self.base_ring() 

f = self.defining_polynomial() 

pts = self.places_on_curve() 

numpts = len(pts) 

R = f.parent() 

x,y = R.gens() 

R0 = PolynomialRing(F,3,names = [str(x),str(y),"t"]) 

vars0 = R0.gens() 

t = vars0[2] 

divf = [] 

for pt0 in pts: 

if pt0[2] != F(0): 

lcs = self.local_coordinates(pt0,5) 

yt = lcs[1] 

xt = lcs[0] 

ldg = degree_lowest_rational_function(r(xt,yt),t) 

if ldg[0] != 0: 

divf.append([ldg[0],pt0]) 

return divf 

 

def local_coordinates(self, pt, n): 

r""" 

Return local coordinates to precision n at the given point. 

 

Behaviour is flaky - some choices of `n` are worst that 

others. 

 

 

INPUT: 

 

 

- ``pt`` - an F-rational point on X which is not a 

point of ramification for the projection (x,y) - x. 

 

- ``n`` - the number of terms desired 

 

 

OUTPUT: x = x0 + t y = y0 + power series in t 

 

EXAMPLES:: 

 

sage: F = GF(5) 

sage: pt = (2,3) 

sage: R = PolynomialRing(F,2, names = ['x','y']) 

sage: x,y = R.gens() 

sage: f = y^2-x^9-x 

sage: C = Curve(f) 

sage: C.local_coordinates(pt, 9) 

[t + 2, -2*t^12 - 2*t^11 + 2*t^9 + t^8 - 2*t^7 - 2*t^6 - 2*t^4 + t^3 - 2*t^2 - 2] 

""" 

f = self.defining_polynomial() 

R = f.parent() 

F = self.base_ring() 

p = F.characteristic() 

x0 = F(pt[0]) 

y0 = F(pt[1]) 

astr = ["a"+str(i) for i in range(1,2*n)] 

x,y = R.gens() 

R0 = PolynomialRing(F,2*n+2,names = [str(x),str(y),"t"]+astr) 

vars0 = R0.gens() 

t = vars0[2] 

yt = y0*t**0+add([vars0[i]*t**(i-2) for i in range(3,2*n+2)]) 

xt = x0+t 

ft = f(xt,yt) 

S = singular 

S.eval('ring s = '+str(p)+','+str(R0.gens())+',lp;') 

S.eval('poly f = '+str(ft) + ';') 

c = S('coeffs(%s, t)'%ft) 

N = int(c.size()) 

b = ["%s[%s,1],"%(c.name(), i) for i in range(2,N//2-4)] 

b = ''.join(b) 

b = b[:len(b)-1] # to cut off the trailing comma 

cmd = 'ideal I = '+b 

S.eval(cmd) 

S.eval('short=0') # print using *'s and ^'s. 

c = S.eval('slimgb(I)') 

d = c.split("=") 

d = d[1:] 

d[len(d)-1] += "\n" 

e = [x[:x.index("\n")] for x in d] 

vals = [] 

for x in e: 

for y in vars0: 

if str(y) in x: 

if len(x.replace(str(y),"")) != 0: 

i = x.find("-") 

if i>0: 

vals.append([eval(x[1:i]),x[:i],F(eval(x[i+1:]))]) 

i = x.find("+") 

if i>0: 

vals.append([eval(x[1:i]),x[:i],-F(eval(x[i+1:]))]) 

else: 

vals.append([eval(str(y)[1:]),str(y),F(0)]) 

vals.sort() 

k = len(vals) 

v = [x0+t,y0+add([vals[i][2]*t**(i+1) for i in range(k)])] 

return v 

 

def plot(self, *args, **kwds): 

""" 

Plot the real points on this affine plane curve. 

 

INPUT: 

 

 

- ``self`` - an affine plane curve 

 

- ``*args`` - optional tuples (variable, minimum, maximum) for 

plotting dimensions 

 

- ``**kwds`` - optional keyword arguments passed on to 

``implicit_plot`` 

 

 

EXAMPLES: 

 

A cuspidal curve:: 

 

sage: R.<x, y> = QQ[] 

sage: C = Curve(x^3 - y^2) 

sage: C.plot() 

Graphics object consisting of 1 graphics primitive 

 

A 5-nodal curve of degree 11. This example also illustrates 

some of the optional arguments:: 

 

sage: R.<x, y> = ZZ[] 

sage: C = Curve(32*x^2 - 2097152*y^11 + 1441792*y^9 - 360448*y^7 + 39424*y^5 - 1760*y^3 + 22*y - 1) 

sage: C.plot((x, -1, 1), (y, -1, 1), plot_points=400) 

Graphics object consisting of 1 graphics primitive 

 

A line over `\mathbf{RR}`:: 

 

sage: R.<x, y> = RR[] 

sage: C = Curve(R(y - sqrt(2)*x)) 

sage: C.plot() 

Graphics object consisting of 1 graphics primitive 

""" 

I = self.defining_ideal() 

return I.plot(*args, **kwds) 

 

def is_transverse(self, C, P): 

r""" 

Return whether the intersection of this curve with the curve ``C`` at the point ``P`` is transverse. 

 

The intersection at ``P`` is transverse if ``P`` is a nonsingular point of both curves, and if the 

tangents of the curves at ``P`` are distinct. 

 

INPUT: 

 

- ``C`` -- a curve in the ambient space of this curve. 

 

- ``P`` -- a point in the intersection of both curves. 

 

OUTPUT: Boolean. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([x^2 + y^2 - 1], A) 

sage: D = Curve([x - 1], A) 

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

sage: C.is_transverse(D, Q) 

False 

 

:: 

 

sage: R.<a> = QQ[] 

sage: K.<b> = NumberField(a^3 + 2) 

sage: A.<x,y> = AffineSpace(K, 2) 

sage: C = A.curve([x*y]) 

sage: D = A.curve([y - b*x]) 

sage: Q = A([0,0]) 

sage: C.is_transverse(D, Q) 

False 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y - x^3], A) 

sage: D = Curve([y + x], A) 

sage: Q = A([0,0]) 

sage: C.is_transverse(D, Q) 

True 

""" 

if not self.intersects_at(C, P): 

raise TypeError("(=%s) must be a point in the intersection of (=%s) and this curve"%(P,C)) 

if self.is_singular(P) or C.is_singular(P): 

return False 

 

# there is only one tangent at a nonsingular point of a plane curve 

return not self.tangents(P)[0] == C.tangents(P)[0] 

 

def multiplicity(self, P): 

r""" 

Return the multiplicity of this affine plane curve at the point ``P``. 

 

In the special case of affine plane curves, the multiplicity of an affine 

plane curve at the point (0,0) can be computed as the minimum of the degrees 

of the homogeneous components of its defining polynomial. To compute the 

multiplicity of a different point, a linear change of coordinates is used. 

 

This curve must be defined over a field. An error if raised if ``P`` is 

not a point on this curve. 

 

INPUT: 

 

- ``P`` -- a point in the ambient space of this curve. 

 

OUTPUT: 

 

An integer. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y^2 - x^3], A) 

sage: Q1 = A([1,1]) 

sage: C.multiplicity(Q1) 

1 

sage: Q2 = A([0,0]) 

sage: C.multiplicity(Q2) 

2 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQbar,2) 

sage: C = Curve([-x^7 + (-7)*x^6 + y^6 + (-21)*x^5 + 12*y^5 + (-35)*x^4 + 60*y^4 +\ 

(-35)*x^3 + 160*y^3 + (-21)*x^2 + 240*y^2 + (-7)*x + 192*y + 63], A) 

sage: Q = A([-1,-2]) 

sage: C.multiplicity(Q) 

6 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve([y^3 - x^3 + x^6]) 

sage: Q = A([1,1]) 

sage: C.multiplicity(Q) 

Traceback (most recent call last): 

... 

TypeError: (=(1, 1)) is not a point on (=Affine Plane Curve over 

Rational Field defined by x^6 - x^3 + y^3) 

""" 

if not self.base_ring() in Fields(): 

raise TypeError("curve must be defined over a field") 

 

# Check whether P is a point on this curve 

try: 

P = self(P) 

except TypeError: 

raise TypeError("(=%s) is not a point on (=%s)"%(P,self)) 

 

# Apply a linear change of coordinates to self so that P becomes (0,0) 

AA = self.ambient_space() 

f = self.defining_polynomials()[0](AA.gens()[0] + P[0], AA.gens()[1] + P[1]) 

 

# Compute the multiplicity of the new curve at (0,0), which is the minimum of the degrees of its 

# nonzero terms 

return min([g.degree() for g in f.monomials()]) 

 

def tangents(self, P, factor=True): 

r""" 

Return the tangents of this affine plane curve at the point ``P``. 

 

The point ``P`` must be a point on this curve. 

 

INPUT: 

 

- ``P`` -- a point on this curve. 

 

- ``factor`` -- (default: True) whether to attempt computing the polynomials of the individual tangent 

lines over the base field of this curve, or to just return the polynomial corresponding to the union 

of the tangent lines (which requires fewer computations). 

 

OUTPUT: 

 

- a list of polynomials in the coordinate ring of the ambient space of this curve. 

 

EXAMPLES:: 

 

sage: set_verbose(-1) 

sage: A.<x,y> = AffineSpace(QQbar, 2) 

sage: C = Curve([x^5*y^3 + 2*x^4*y^4 + x^3*y^5 + 3*x^4*y^3 + 6*x^3*y^4 + 3*x^2*y^5\ 

+ 3*x^3*y^3 + 6*x^2*y^4 + 3*x*y^5 + x^5 + 10*x^4*y + 40*x^3*y^2 + 81*x^2*y^3 + 82*x*y^4\ 

+ 33*y^5], A) 

sage: Q = A([0,0]) 

sage: C.tangents(Q) 

[x + 3.425299577684700?*y, x + (1.949159013086856? + 1.179307909383728?*I)*y, 

x + (1.949159013086856? - 1.179307909383728?*I)*y, x + (1.338191198070795? + 0.2560234251008043?*I)*y, 

x + (1.338191198070795? - 0.2560234251008043?*I)*y] 

sage: C.tangents(Q, factor=False) 

[120*x^5 + 1200*x^4*y + 4800*x^3*y^2 + 9720*x^2*y^3 + 9840*x*y^4 + 3960*y^5] 

 

:: 

 

sage: R.<a> = QQ[] 

sage: K.<b> = NumberField(a^2 - 3) 

sage: A.<x,y> = AffineSpace(K, 2) 

sage: C = Curve([(x^2 + y^2 - 2*x)^2 - x^2 - y^2], A) 

sage: Q = A([0,0]) 

sage: C.tangents(Q) 

[x + (-1/3*b)*y, x + (1/3*b)*y] 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve([y^2 - x^3 - x^2]) 

sage: Q = A([0,0]) 

sage: C.tangents(Q) 

[x - y, x + y] 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve([y*x - x^4 + 2*x^2]) 

sage: Q = A([1,1]) 

sage: C.tangents(Q) 

Traceback (most recent call last): 

... 

TypeError: (=(1, 1)) is not a point on (=Affine Plane Curve over 

Rational Field defined by -x^4 + 2*x^2 + x*y) 

""" 

r = self.multiplicity(P) 

f = self.defining_polynomial() 

# move P to (0,0) 

vars = self.ambient_space().gens() 

coords = [vars[0] + P[0], vars[1] + P[1]] 

f = f(coords) 

coords = [vars[0] - P[0], vars[1] - P[1]] # coords to change back with 

deriv = [f.derivative(vars[0],i).derivative(vars[1], r-i)([0,0]) for i in range(r+1)] 

T = sum([binomial(r,i)*deriv[i]*(vars[0])**i*(vars[1])**(r-i) for i in range(r+1)]) 

if not factor: 

return [T(coords)] 

if self.base_ring() == QQbar: 

fact = [] 

# first add tangents corresponding to vars[0], vars[1] if they divide T 

t = min([e[0] for e in T.exponents()]) 

# vars[0] divides T 

if t > 0: 

fact.append(vars[0]) 

# divide T by that power of vars[0] 

T = self.ambient_space().coordinate_ring()(dict([((v[0] - t,v[1]), h) for (v,h) in T.dict().items()])) 

t = min([e[1] for e in T.exponents()]) 

# vars[1] divides T 

if t > 0: 

fact.append(vars[1]) 

# divide T by that power of vars[1] 

T = self.ambient_space().coordinate_ring()(dict([((v[0],v[1] - t), h) for (v,h) in T.dict().items()])) 

# T is homogeneous in var[0], var[1] if nonconstant, so dehomogenize 

if not T in self.base_ring(): 

if T.degree(vars[0]) > 0: 

T = T(vars[0], 1) 

roots = T.univariate_polynomial().roots() 

fact.extend([vars[0] - roots[i][0]*vars[1] for i in range(len(roots))]) 

else: 

T = T(1, vars[1]) 

roots = T.univariate_polynomial().roots() 

fact.extend([vars[1] - roots[i][0]*vars[0] for i in range(len(roots))]) 

return [f(coords) for f in fact] 

else: 

fact = T.factor() 

return [l[0](coords) for l in fact] 

 

def is_ordinary_singularity(self, P): 

r""" 

Return whether the singular point ``P`` of this affine plane curve is an ordinary singularity. 

 

The point ``P`` is an ordinary singularity of this curve if it is a singular point, and 

if the tangents of this curve at ``P`` are distinct. 

 

INPUT: 

 

- ``P`` -- a point on this curve. 

 

OUTPUT: 

 

- Boolean. True or False depending on whether ``P`` is or is not an ordinary singularity of this 

curve, respectively. An error is raised if ``P`` is not a singular point of this curve. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y^2 - x^3], A) 

sage: Q = A([0,0]) 

sage: C.is_ordinary_singularity(Q) 

False 

 

:: 

 

sage: R.<a> = QQ[] 

sage: K.<b> = NumberField(a^2 - 3) 

sage: A.<x,y> = AffineSpace(K, 2) 

sage: C = Curve([(x^2 + y^2 - 2*x)^2 - x^2 - y^2], A) 

sage: Q = A([0,0]) 

sage: C.is_ordinary_singularity(Q) 

True 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve([x^2*y - y^2*x + y^2 + x^3]) 

sage: Q = A([-1,-1]) 

sage: C.is_ordinary_singularity(Q) 

Traceback (most recent call last): 

... 

TypeError: (=(-1, -1)) is not a singular point of (=Affine Plane Curve 

over Rational Field defined by x^3 + x^2*y - x*y^2 + y^2) 

""" 

r = self.multiplicity(P) 

if r < 2: 

raise TypeError("(=%s) is not a singular point of (=%s)"%(P,self)) 

 

T = self.tangents(P, factor=False)[0] 

vars = self.ambient_space().gens() 

 

# use resultants to determine if there is a higher multiplicity tangent 

if T.degree(vars[0]) > 0: 

return T.resultant(T.derivative(vars[0]), vars[0]) != 0 

else: 

return T.resultant(T.derivative(vars[1]), vars[1]) != 0 

 

def rational_parameterization(self): 

r""" 

Return a rational parameterization of this curve. 

 

This curve must have rational coefficients and be absolutely irreducible (i.e. irreducible 

over the algebraic closure of the rational field). The curve must also be rational (have 

geometric genus zero). 

 

The rational parameterization may have coefficients in a quadratic extension of the rational 

field. 

 

OUTPUT: 

 

- a birational map between `\mathbb{A}^{1}` and this curve, given as a scheme morphism. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([y^2 - x], A) 

sage: C.rational_parameterization() 

Scheme morphism: 

From: Affine Space of dimension 1 over Rational Field 

To: Affine Plane Curve over Rational Field defined by y^2 - x 

Defn: Defined on coordinates by sending (t) to 

(t^2, t) 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([(x^2 + y^2 - 2*x)^2 - x^2 - y^2], A) 

sage: C.rational_parameterization() 

Scheme morphism: 

From: Affine Space of dimension 1 over Rational Field 

To: Affine Plane Curve over Rational Field defined by x^4 + 

2*x^2*y^2 + y^4 - 4*x^3 - 4*x*y^2 + 3*x^2 - y^2 

Defn: Defined on coordinates by sending (t) to 

((-12*t^4 + 6*t^3 + 4*t^2 - 2*t)/(-25*t^4 + 40*t^3 - 26*t^2 + 

8*t - 1), (-9*t^4 + 12*t^3 - 4*t + 1)/(-25*t^4 + 40*t^3 - 26*t^2 + 8*t - 1)) 

 

:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = Curve([x^2 + y^2 + 7], A) 

sage: C.rational_parameterization() 

Scheme morphism: 

From: Affine Space of dimension 1 over Number Field in a with defining polynomial a^2 + 7 

To: Affine Plane Curve over Number Field in a with defining 

polynomial a^2 + 7 defined by x^2 + y^2 + 7 

Defn: Defined on coordinates by sending (t) to 

((-7*t^2 + 7)/((-a)*t^2 + (-a)), 14*t/((-a)*t^2 + (-a))) 

""" 

para = self.projective_closure(i=0).rational_parameterization().defining_polynomials() 

# these polynomials are homogeneous in two indeterminants, so dehomogenize wrt one of the variables 

R = para[0].parent() 

A_line = AffineSpace(R.base_ring(), 1, 't') 

para = [A_line.coordinate_ring()(para[i].substitute({R.gens()[0]: 1})) for i in range(3)] 

C = self.change_ring(R.base_ring()) 

# because of the parameter i=0, the projective closure is constructed with respect to the 

# affine patch corresponding to the first coordinate being nonzero. Thus para[0] will not be 

# the zero polynomial, and dehomogenization won't change this 

H = Hom(A_line, C) 

return H([para[1]/para[0], para[2]/para[0]]) 

 

def fundamental_group(self): 

r""" 

Return a presentation of the fundamental group of the complement 

of ``self``. 

 

.. NOTE:: 

 

The curve must be defined over the rationals or a number field 

with an embedding over `\QQbar`. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(QQ, 2) 

sage: C = A.curve(y^2 - x^3 - x^2) 

sage: C.fundamental_group() # optional - sirocco 

Finitely presented group < x0 | > 

 

In the case of number fields, they need to have an embedding 

to the algebraic field:: 

 

sage: a = QQ[x](x^2+5).roots(QQbar)[0][0] 

sage: F = NumberField(a.minpoly(), 'a', embedding=a) 

sage: F.inject_variables() 

Defining a 

sage: A.<x,y> = AffineSpace(F, 2) 

sage: C = A.curve(y^2 - a*x^3 - x^2) 

sage: C.fundamental_group() # optional - sirocco 

Finitely presented group < x0 | > 

 

.. WARNING:: 

 

This functionality requires the sirocco package to be installed. 

""" 

from sage.schemes.curves.zariski_vankampen import fundamental_group 

F = self.base_ring() 

from sage.rings.qqbar import QQbar 

if QQbar.coerce_map_from(F) is None: 

raise NotImplementedError("the base field must have an embedding" 

" to the algebraic field") 

f = self.defining_polynomial() 

return fundamental_group(f, projective=False) 

 

def riemann_surface(self,**kwargs): 

r"""Return the complex riemann surface determined by this curve 

 

OUTPUT: 

 

- RiemannSurface object 

 

EXAMPLES:: 

 

sage: R.<x,y>=QQ[] 

sage: C=Curve(x^3+3*y^3+5) 

sage: C.riemann_surface() 

Riemann surface defined by polynomial f = x^3 + 3*y^3 + 5 = 0, with 53 bits of precision 

 

""" 

from sage.schemes.riemann_surfaces.riemann_surface import RiemannSurface 

return RiemannSurface(self.defining_polynomial(),**kwargs) 

 

 

class AffinePlaneCurve_finite_field(AffinePlaneCurve): 

 

_point = point.AffinePlaneCurvePoint_finite_field 

 

def rational_points(self, algorithm="enum"): 

r""" 

Return sorted list of all rational points on this curve. 

 

Use *very* naive point enumeration to find all rational points on 

this curve over a finite field. 

 

EXAMPLES:: 

 

sage: A.<x,y> = AffineSpace(2,GF(9,'a')) 

sage: C = Curve(x^2 + y^2 - 1) 

sage: C 

Affine Plane Curve over Finite Field in a of size 3^2 defined by x^2 + y^2 - 1 

sage: C.rational_points() 

[(0, 1), (0, 2), (1, 0), (2, 0), (a + 1, a + 1), (a + 1, 2*a + 2), (2*a + 2, a + 1), (2*a + 2, 2*a + 2)] 

""" 

f = self.defining_polynomial() 

R = f.parent() 

K = R.base_ring() 

points = [] 

for x in K: 

for y in K: 

if f(x,y) == 0: 

points.append(self((x,y))) 

points.sort() 

return points 

 

 

class AffinePlaneCurve_prime_finite_field(AffinePlaneCurve_finite_field): 

# CHECK WHAT ASSUMPTIONS ARE MADE REGARDING AFFINE VS. PROJECTIVE MODELS!!! 

# THIS IS VERY DIRTY STILL -- NO DATASTRUCTURES FOR DIVISORS. 

 

def riemann_roch_basis(self,D): 

""" 

Interfaces with Singular's BrillNoether command. 

 

INPUT: 

 

 

- ``self`` - a plane curve defined by a polynomial eqn f(x,y) 

= 0 over a prime finite field F = GF(p) in 2 variables x,y 

representing a curve X: f(x,y) = 0 having n F-rational 

points (see the Sage function places_on_curve) 

 

- ``D`` - an n-tuple of integers 

`(d1, ..., dn)` representing the divisor 

`Div = d1*P1+...+dn*Pn`, where 

`X(F) = \{P1,...,Pn\}`. 

*The ordering is that dictated by places_on_curve.* 

 

 

OUTPUT: basis of L(Div) 

 

EXAMPLES:: 

 

sage: R = PolynomialRing(GF(5),2,names = ["x","y"]) 

sage: x, y = R.gens() 

sage: f = y^2 - x^9 - x 

sage: C = Curve(f) 

sage: D = [6,0,0,0,0,0] 

sage: C.riemann_roch_basis(D) 

[1, (y^2*z^4 - x*z^5)/x^6, (y^2*z^5 - x*z^6)/x^7, (y^2*z^6 - x*z^7)/x^8] 

""" 

f = self.defining_polynomial() 

R = f.parent() 

F = self.base_ring() 

p = F.characteristic() 

Dstr = str(tuple(D)) 

G = singular(','.join([str(x) for x in D]), type='intvec') 

 

singular.LIB('brnoeth.lib') 

 

S = singular.ring(p, R.gens(), 'lp') 

fsing = singular(str(f)) 

 

X = fsing.Adj_div() 

P = singular.NSplaces(1, X) 

T = P[1][2] 

T.set_ring() 

 

LG = G.BrillNoether(P) 

 

dim = len(LG) 

basis = [(LG[i][1], LG[i][2]) for i in range(1,dim+1)] 

x, y, z = PolynomialRing(F, 3, names = ["x","y","z"]).gens() 

V = [] 

for g in basis: 

T.set_ring() # necessary... 

V.append(eval(g[0].sage_polystring())/eval(g[1].sage_polystring())) 

return V 

 

 

def rational_points(self, algorithm="enum"): 

r""" 

Return sorted list of all rational points on this curve. 

 

INPUT: 

 

 

- ``algorithm`` - string: 

 

+ ``'enum'`` - straightforward enumeration 

 

+ ``'bn'`` - via Singular's Brill-Noether package. 

 

+ ``'all'`` - use all implemented algorithms and 

verify that they give the same answer, then return it 

 

 

.. note:: 

 

The Brill-Noether package does not always work. When it 

fails a RuntimeError exception is raised. 

 

EXAMPLES:: 

 

sage: x, y = (GF(5)['x,y']).gens() 

sage: f = y^2 - x^9 - x 

sage: C = Curve(f); C 

Affine Plane Curve over Finite Field of size 5 defined by -x^9 + y^2 - x 

sage: C.rational_points(algorithm='bn') 

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

sage: C = Curve(x - y + 1) 

sage: C.rational_points() 

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

 

We compare Brill-Noether and enumeration:: 

 

sage: x, y = (GF(17)['x,y']).gens() 

sage: C = Curve(x^2 + y^5 + x*y - 19) 

sage: v = C.rational_points(algorithm='bn') 

sage: w = C.rational_points(algorithm='enum') 

sage: len(v) 

20 

sage: v == w 

True 

""" 

if algorithm == "enum": 

 

return AffinePlaneCurve_finite_field.rational_points(self, algorithm="enum") 

 

elif algorithm == "bn": 

f = self.defining_polynomial()._singular_() 

singular = f.parent() 

singular.lib('brnoeth') 

try: 

X1 = f.Adj_div() 

except (TypeError, RuntimeError) as s: 

raise RuntimeError(str(s) + "\n\n ** Unable to use the Brill-Noether Singular package to compute all points (see above).") 

 

X2 = singular.NSplaces(1, X1) 

R = X2[5][1][1] 

singular.set_ring(R) 

 

# We use sage_flattened_str_list since iterating through 

# the entire list through the sage/singular interface directly 

# would involve hundreds of calls to singular, and timing issues 

# with the expect interface could crop up. Also, this is vastly 

# faster (and more robust). 

v = singular('POINTS').sage_flattened_str_list() 

pnts = [self(int(v[3*i]), int(v[3*i+1])) for i in range(len(v)//3) if int(v[3*i+2])!=0] 

# remove multiple points 

pnts = sorted(set(pnts)) 

return pnts 

 

elif algorithm == "all": 

 

S_enum = self.rational_points(algorithm = "enum") 

S_bn = self.rational_points(algorithm = "bn") 

if S_enum != S_bn: 

raise RuntimeError("Bug in rational_points -- different algorithms give different answers for curve %s!"%self) 

return S_enum 

 

else: 

raise ValueError("No algorithm '%s' known"%algorithm)