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

r""" 

Dense univariate polynomials over `\RR`, implemented using MPFR 

  

TESTS: 

  

Check that operations with numpy elements work well (see :trac:`18076` and 

:trac:`8426`):: 

  

sage: import numpy 

sage: x = polygen(RR) 

sage: x * numpy.int32('1') 

x 

sage: numpy.int32('1') * x 

x 

sage: x * numpy.int64('1') 

x 

sage: numpy.int64('1') * x 

x 

sage: x * numpy.float32('1.5') 

1.50000000000000*x 

sage: numpy.float32('1.5') * x 

1.50000000000000*x 

""" 

from __future__ import absolute_import 

  

from cysignals.memory cimport check_allocarray, check_reallocarray, sig_free 

from cysignals.signals cimport sig_on, sig_off 

  

from cpython cimport PyInt_AS_LONG, PyFloat_AS_DOUBLE 

  

from sage.structure.parent cimport Parent 

from .polynomial_element cimport Polynomial, _dict_to_list 

from sage.rings.real_mpfr cimport RealField_class, RealNumber 

from sage.rings.integer cimport Integer, smallInteger 

from sage.rings.rational cimport Rational 

  

from sage.structure.element cimport Element, ModuleElement, RingElement 

from sage.structure.element cimport parent 

from sage.structure.element import coerce_binop 

from sage.libs.mpfr cimport * 

  

from sage.libs.all import pari_gen 

  

cdef class PolynomialRealDense(Polynomial): 

r""" 

  

TESTS:: 

  

sage: f = RR['x'].random_element() 

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: isinstance(f, PolynomialRealDense) 

True 

  

""" 

  

cdef Py_ssize_t _degree 

cdef mpfr_t* _coeffs 

cdef RealField_class _base_ring 

  

def __cinit__(self): 

""" 

TESTS:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: PolynomialRealDense(RR['x']) 

0 

""" 

self._coeffs = NULL 

  

def __init__(self, Parent parent, x=0, check=None, bint is_gen=False, construct=None): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: PolynomialRealDense(RR['x'], [1, int(2), RR(3), 4/1, pi]) 

3.14159265358979*x^4 + 4.00000000000000*x^3 + 3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 

sage: PolynomialRealDense(RR['x'], None) 

0 

  

TESTS: 

  

Check that errors and interrupts are handled properly (see :trac:`10100`):: 

  

sage: a = var('a') 

sage: PolynomialRealDense(RR['x'], [1,a]) 

Traceback (most recent call last): 

... 

TypeError: Cannot evaluate symbolic expression to a numeric value. 

sage: R.<x> = SR[] 

sage: (x-a).change_ring(RR) 

Traceback (most recent call last): 

... 

TypeError: Cannot evaluate symbolic expression to a numeric value. 

sage: sig_on_count() 

0 

  

Test that we don't clean up uninitialized coefficients (:trac:`9826`):: 

  

sage: k.<a> = GF(7^3) 

sage: P.<x> = PolynomialRing(k) 

sage: (a*x).complex_roots() 

Traceback (most recent call last): 

... 

TypeError: unable to convert 'a' to a real number 

  

Check that :trac:`17190` is fixed:: 

  

sage: RR['x']({}) 

0 

""" 

Polynomial.__init__(self, parent, is_gen=is_gen) 

self._base_ring = parent._base 

cdef Py_ssize_t i, degree 

cdef int prec = self._base_ring.__prec 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

if x is None: 

self._coeffs = <mpfr_t*>check_allocarray(1, sizeof(mpfr_t)) # degree zero 

mpfr_init2(self._coeffs[0], prec) 

mpfr_set_si(self._coeffs[0], 0, rnd) 

self._normalize() 

return 

if is_gen: 

x = [0, 1] 

elif isinstance(x, (int, float, Integer, Rational, RealNumber)): 

x = [x] 

elif isinstance(x, dict): 

x = _dict_to_list(x, self._base_ring.zero()) 

elif isinstance(x, pari_gen): 

x = [self._base_ring(w) for w in x.list()] 

elif not isinstance(x, list): 

try: 

x = list(x) 

except TypeError: # x is not iterable 

x = [self._base_ring(x)] 

  

sig_on() 

degree = len(x) - 1 

self._degree = -1 

cdef mpfr_t* coeffs 

coeffs = <mpfr_t*>check_allocarray(degree+1, sizeof(mpfr_t)) 

try: # We might get Python exceptions here 

for i from 0 <= i <= degree: 

mpfr_init2(coeffs[i], prec) 

self._degree += 1 

a = x[i] 

if type(a) is RealNumber: 

mpfr_set(coeffs[i], (<RealNumber>a).value, rnd) 

elif type(a) is int: 

mpfr_set_si(coeffs[i], PyInt_AS_LONG(a), rnd) 

elif type(a) is float: 

mpfr_set_d(coeffs[i], PyFloat_AS_DOUBLE(a), rnd) 

elif type(a) is Integer: 

mpfr_set_z(coeffs[i], (<Integer>a).value, rnd) 

elif type(a) is Rational: 

mpfr_set_q(coeffs[i], (<Rational>a).value, rnd) 

else: 

a = self._base_ring(a) 

mpfr_set(coeffs[i], (<RealNumber>a).value, rnd) 

finally: 

sig_off() 

self._coeffs = coeffs 

self._normalize() 

  

def __dealloc__(self): 

cdef Py_ssize_t i 

if self._coeffs != NULL: 

for i from 0 <= i <= self._degree: 

mpfr_clear(self._coeffs[i]) 

sig_free(self._coeffs) 

  

def __reduce__(self): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-2, 0, 1]) 

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

True 

""" 

return make_PolynomialRealDense, (self._parent, self.list()) 

  

cdef _normalize(self): 

""" 

Remove all leading 0's. 

""" 

cdef Py_ssize_t i 

if self._degree >= 0 and mpfr_zero_p(self._coeffs[self._degree]): 

i = self._degree 

while i >= 0 and mpfr_zero_p(self._coeffs[i]): 

mpfr_clear(self._coeffs[i]) 

i -= 1 

self._coeffs = <mpfr_t*>check_reallocarray(self._coeffs, i+1, sizeof(mpfr_t)) 

self._degree = i 

  

cdef get_unsafe(self, Py_ssize_t i): 

""" 

Return the `i`-th coefficient of ``self``. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], range(5)); f 

4.00000000000000*x^4 + 3.00000000000000*x^3 + 2.00000000000000*x^2 + x 

sage: f[0] 

0.000000000000000 

sage: f[3] 

3.00000000000000 

sage: f[5] 

0.000000000000000 

  

Test slices:: 

  

sage: R.<x> = RealField(10)[] 

sage: f = (x+1)^5; f 

x^5 + 5.0*x^4 + 10.*x^3 + 10.*x^2 + 5.0*x + 1.0 

sage: f[:3] 

10.*x^2 + 5.0*x + 1.0 

""" 

cdef RealNumber r = <RealNumber>RealNumber.__new__(RealNumber, self._base_ring) 

mpfr_set(r.value, self._coeffs[i], self._base_ring.rnd) 

return r 

  

cdef PolynomialRealDense _new(self, Py_ssize_t degree): 

cdef Py_ssize_t i 

cdef int prec = self._base_ring.__prec 

cdef PolynomialRealDense f = <PolynomialRealDense>PolynomialRealDense.__new__(PolynomialRealDense) 

f._parent = self._parent 

f._base_ring = self._base_ring 

f._degree = degree 

if degree >= 0: 

f._coeffs = <mpfr_t*>check_allocarray(degree+1, sizeof(mpfr_t)) 

for i from 0 <= i <= degree: 

mpfr_init2(f._coeffs[i], prec) 

return f 

  

def degree(self): 

""" 

Return the degree of the polynomial. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [1, 2, 3]); f 

3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 

sage: f.degree() 

2 

  

TESTS:: 

  

sage: type(f.degree()) 

<type 'sage.rings.integer.Integer'> 

""" 

return smallInteger(self._degree) 

  

cpdef Polynomial truncate(self, long n): 

r""" 

Returns the polynomial of degree `< n` which is equivalent to self 

modulo `x^n`. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RealField(10)['x'], [1, 2, 4, 8]) 

sage: f.truncate(3) 

4.0*x^2 + 2.0*x + 1.0 

sage: f.truncate(100) 

8.0*x^3 + 4.0*x^2 + 2.0*x + 1.0 

sage: f.truncate(1) 

1.0 

sage: f.truncate(0) 

0 

""" 

if n <= 0: 

return self._new(-1) 

if n > self._degree: 

return self 

cdef PolynomialRealDense f = self._new(n-1) 

cdef Py_ssize_t i 

for i from 0 <= i < n: 

mpfr_set(f._coeffs[i], self._coeffs[i], self._base_ring.rnd) 

f._normalize() 

return f 

  

def truncate_abs(self, RealNumber bound): 

""" 

Truncate all high order coefficients below bound. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RealField(10)['x'], [10^-k for k in range(10)]) 

sage: f 

1.0e-9*x^9 + 1.0e-8*x^8 + 1.0e-7*x^7 + 1.0e-6*x^6 + 0.000010*x^5 + 0.00010*x^4 + 0.0010*x^3 + 0.010*x^2 + 0.10*x + 1.0 

sage: f.truncate_abs(0.5e-6) 

1.0e-6*x^6 + 0.000010*x^5 + 0.00010*x^4 + 0.0010*x^3 + 0.010*x^2 + 0.10*x + 1.0 

sage: f.truncate_abs(10.0) 

0 

sage: f.truncate_abs(1e-100) == f 

True 

""" 

cdef Py_ssize_t i 

for i from self._degree >= i >= 0: 

if mpfr_cmpabs(self._coeffs[i], bound.value) >= 0: 

return self.truncate(i+1) 

return self._new(-1) 

  

cpdef shift(self, Py_ssize_t n): 

r""" 

Returns this polynomial multiplied by the power `x^n`. If `n` 

is negative, terms below `x^n` will be discarded. Does not 

change this polynomial. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [1, 2, 3]); f 

3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 

sage: f.shift(10) 

3.00000000000000*x^12 + 2.00000000000000*x^11 + x^10 

sage: f.shift(-1) 

3.00000000000000*x + 2.00000000000000 

sage: f.shift(-10) 

0 

  

TESTS:: 

  

sage: f = RR['x'](0) 

sage: f.shift(3).is_zero() 

True 

sage: f.shift(-3).is_zero() 

True 

""" 

cdef Py_ssize_t i 

cdef Py_ssize_t nn = 0 if n < 0 else n 

cdef PolynomialRealDense f 

if n == 0 or self._degree < 0: 

return self 

elif self._degree < -n: 

return self._new(-1) 

else: 

f = self._new(self._degree + n) 

for i from 0 <= i < n: 

mpfr_set_ui(f._coeffs[i], 0, self._base_ring.rnd) 

for i from nn <= i <= self._degree + n: 

mpfr_set(f._coeffs[i], self._coeffs[i-n], self._base_ring.rnd) 

return f 

  

cpdef list list(self, bint copy=True): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [1, 0, -2]); f 

-2.00000000000000*x^2 + 1.00000000000000 

sage: f.list() 

[1.00000000000000, 0.000000000000000, -2.00000000000000] 

""" 

cdef RealNumber r 

cdef Py_ssize_t i 

cdef list L = [] 

for i from 0 <= i <= self._degree: 

r = <RealNumber>RealNumber(self._base_ring) 

mpfr_set(r.value, self._coeffs[i], self._base_ring.rnd) 

L.append(r) 

return L 

  

def __neg__(self): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-2,0,1]) 

sage: -f 

-x^2 + 2.00000000000000 

""" 

cdef Py_ssize_t i 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense f = self._new(self._degree) 

for i from 0 <= i <= f._degree: 

mpfr_neg(f._coeffs[i], self._coeffs[i], rnd) 

return f 

  

cpdef _add_(left, _right): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-2,0,1]); f 

x^2 - 2.00000000000000 

sage: g = PolynomialRealDense(RR['x'], range(5)); g 

4.00000000000000*x^4 + 3.00000000000000*x^3 + 2.00000000000000*x^2 + x 

sage: f+g 

4.00000000000000*x^4 + 3.00000000000000*x^3 + 3.00000000000000*x^2 + x - 2.00000000000000 

sage: g + f == f + g 

True 

sage: f + (-f) 

0 

""" 

cdef Py_ssize_t i 

cdef mpfr_rnd_t rnd = left._base_ring.rnd 

cdef PolynomialRealDense right = _right 

cdef Py_ssize_t min = left._degree if left._degree < right._degree else right._degree 

cdef Py_ssize_t max = left._degree if left._degree > right._degree else right._degree 

cdef PolynomialRealDense f = left._new(max) 

for i from 0 <= i <= min: 

mpfr_add(f._coeffs[i], left._coeffs[i], right._coeffs[i], rnd) 

if left._degree < right._degree: 

for i from min < i <= max: 

mpfr_set(f._coeffs[i], right._coeffs[i], rnd) 

else: 

for i from min < i <= max: 

mpfr_set(f._coeffs[i], left._coeffs[i], rnd) 

f._normalize() 

return f 

  

cpdef _sub_(left, _right): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-3,0,1]); f 

x^2 - 3.00000000000000 

sage: g = PolynomialRealDense(RR['x'], range(4)); g 

3.00000000000000*x^3 + 2.00000000000000*x^2 + x 

sage: f-g 

-3.00000000000000*x^3 - x^2 - x - 3.00000000000000 

sage: (f-g) == -(g-f) 

True 

""" 

cdef Py_ssize_t i 

cdef mpfr_rnd_t rnd = left._base_ring.rnd 

cdef PolynomialRealDense right = _right 

cdef Py_ssize_t min = left._degree if left._degree < right._degree else right._degree 

cdef Py_ssize_t max = left._degree if left._degree > right._degree else right._degree 

cdef PolynomialRealDense f = left._new(max) 

for i from 0 <= i <= min: 

mpfr_sub(f._coeffs[i], left._coeffs[i], right._coeffs[i], rnd) 

if left._degree < right._degree: 

for i from min < i <= max: 

mpfr_neg(f._coeffs[i], right._coeffs[i], rnd) 

else: 

for i from min < i <= max: 

mpfr_set(f._coeffs[i], left._coeffs[i], rnd) 

f._normalize() 

return f 

  

cpdef _lmul_(self, Element c): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-5,0,0,1]); f 

x^3 - 5.00000000000000 

sage: 4.0 * f 

4.00000000000000*x^3 - 20.0000000000000 

sage: f * -0.2 

-0.200000000000000*x^3 + 1.00000000000000 

""" 

cdef Py_ssize_t i 

cdef RealNumber a = c 

if mpfr_zero_p(a.value): 

return self._new(-1) 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense f = self._new(self._degree) 

for i from 0 <= i <= self._degree: 

mpfr_mul(f._coeffs[i], self._coeffs[i], a.value, rnd) 

return f 

  

cpdef _mul_(left, _right): 

""" 

Here we use the naive `O(n^2)` algorithm, as asymptotically faster algorithms such 

as Karatsuba can have very inaccurate results due to intermediate rounding errors. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [1e20, 1]) 

sage: g = PolynomialRealDense(RR['x'], [1e30, 1]) 

sage: f*g 

x^2 + 1.00000000010000e30*x + 1.00000000000000e50 

sage: f._mul_karatsuba(g,0) 

x^2 + 1.00000000000000e50 

sage: f = PolynomialRealDense(RR['x'], range(5)) 

sage: g = PolynomialRealDense(RR['x'], range(3)) 

sage: f*g 

8.00000000000000*x^6 + 10.0000000000000*x^5 + 7.00000000000000*x^4 + 4.00000000000000*x^3 + x^2 

""" 

cdef Py_ssize_t i, j 

cdef mpfr_rnd_t rnd = left._base_ring.rnd 

cdef PolynomialRealDense right = _right 

cdef PolynomialRealDense f 

cdef mpfr_t tmp 

if left._degree < 0 or right._degree < 0: 

f = left._new(-1) 

else: 

f = left._new(left._degree + right._degree) 

sig_on() 

mpfr_init2(tmp, left._base_ring.__prec) 

for i from 0 <= i <= f._degree: 

# Yes, we could make this more efficient by initializing with 

# a multiple of left rather than all zeros... 

mpfr_set_ui(f._coeffs[i], 0, rnd) 

for i from 0 <= i <= left._degree: 

for j from 0 <= j <= right._degree: 

mpfr_mul(tmp, left._coeffs[i], right._coeffs[j], rnd) 

mpfr_add(f._coeffs[i+j], f._coeffs[i+j], tmp, rnd) 

mpfr_clear(tmp) 

sig_off() 

return f 

  

def _derivative(self, var=None): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [pi, 0, 2, 1]); 

sage: f.derivative() 

3.00000000000000*x^2 + 4.00000000000000*x 

""" 

if var is not None and var != self._parent.gen(): 

return self._new(-1) 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense f = self._new(self._degree-1) 

for i from 0 <= i < self._degree: 

mpfr_mul_ui(f._coeffs[i], self._coeffs[i+1], i+1, rnd) 

return f 

  

def integral(self): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [3, pi, 1]) 

sage: f.integral() 

0.333333333333333*x^3 + 1.57079632679490*x^2 + 3.00000000000000*x 

""" 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense f = self._new(self._degree+1) 

mpfr_set_ui(f._coeffs[0], 0, rnd) 

for i from 0 <= i <= self._degree: 

mpfr_div_ui(f._coeffs[i+1], self._coeffs[i], i+1, rnd) 

return f 

  

def reverse(self): 

""" 

Returns `x^d f(1/x)` where `d` is the degree of `f`. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-3, pi, 0, 1]) 

sage: f.reverse() 

-3.00000000000000*x^3 + 3.14159265358979*x^2 + 1.00000000000000 

""" 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense f = self._new(self._degree) 

for i from 0 <= i <= self._degree: 

mpfr_set(f._coeffs[self._degree-i], self._coeffs[i], rnd) 

f._normalize() 

return f 

  

@coerce_binop 

def quo_rem(self, PolynomialRealDense other): 

""" 

Return the quotient with remainder of ``self`` by ``other``. 

  

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-2, 0, 1]) 

sage: g = PolynomialRealDense(RR['x'], [5, 1]) 

sage: q, r = f.quo_rem(g) 

sage: q 

x - 5.00000000000000 

sage: r 

23.0000000000000 

sage: q*g + r == f 

True 

sage: fg = f*g 

sage: fg.quo_rem(f) 

(x + 5.00000000000000, 0) 

sage: fg.quo_rem(g) 

(x^2 - 2.00000000000000, 0) 

  

sage: f = PolynomialRealDense(RR['x'], range(5)) 

sage: g = PolynomialRealDense(RR['x'], [pi,3000,4]) 

sage: q, r = f.quo_rem(g) 

sage: g*q + r == f 

True 

  

TESTS: 

  

Check that :trac:`18467` is fixed:: 

  

sage: S.<x> = RR[] 

sage: z = S.zero() 

sage: z.degree() 

-1 

sage: q, r = z.quo_rem(x) 

sage: q.degree() 

-1 

""" 

if other._degree < 0: 

raise ZeroDivisionError("other must be nonzero") 

elif other._degree == 0: 

return self * ~other[0], self._parent.zero() 

elif other._degree > self._degree: 

return self._parent.zero(), self 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense q, r 

cdef Py_ssize_t i, j 

cdef mpfr_t tmp 

# Make divisor monic for simplicity 

leading = other[other._degree] 

other = other * ~leading 

r = self * ~leading 

q = self._new(self._degree - other._degree) 

# This is the standard division algorithm 

sig_on() 

mpfr_init2(tmp, self._base_ring.__prec) 

for i from self._degree >= i >= other._degree: 

mpfr_set(q._coeffs[i-other._degree], r._coeffs[i], rnd) 

for j from 0 <= j < other._degree: 

mpfr_mul(tmp, r._coeffs[i], other._coeffs[j], rnd) 

mpfr_sub(r._coeffs[i-other._degree+j], r._coeffs[i-other._degree+j], tmp, rnd) 

r._degree -= 1 

mpfr_clear(r._coeffs[i]) 

mpfr_clear(tmp) 

sig_off() 

r._normalize() 

return q, r * leading 

  

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

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-2, 0, 1]) 

sage: f(10) 

98.0000000000000 

sage: f(CC.0) 

-3.00000000000000 

sage: f(2.0000000000000000000000000000000000000000000) 

2.00000000000000 

sage: f(RealField(10)(2)) 

2.0 

sage: f(pi) 

1.00000000000000*pi^2 - 2.00000000000000 

  

  

sage: f = PolynomialRealDense(RR['x'], range(5)) 

sage: f(1) 

10.0000000000000 

sage: f(-1) 

2.00000000000000 

sage: f(0) 

0.000000000000000 

sage: f = PolynomialRealDense(RR['x']) 

sage: f(12) 

0.000000000000000 

  

TESTS:: 

  

sage: R.<x> = RR[] # trac #17311 

sage: (x^2+1)(x=5) 

26.0000000000000 

""" 

if len(args) == 1: 

xx = args[0] 

else: 

return Polynomial.__call__(self, *args, **kwds) 

  

if not isinstance(xx, RealNumber): 

if self._base_ring.has_coerce_map_from(parent(xx)): 

xx = self._base_ring(xx) 

else: 

return Polynomial.__call__(self, xx) 

  

cdef Py_ssize_t i 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef RealNumber x = <RealNumber>xx 

cdef RealNumber res 

  

if (<RealField_class>x._parent).__prec < self._base_ring.__prec: 

res = RealNumber(x._parent) 

else: 

res = RealNumber(self._base_ring) 

# Optimize some very useful and common cases: 

if self._degree < 0: 

mpfr_set_ui(res.value, 0, rnd) 

elif mpfr_zero_p(x.value): 

mpfr_set(res.value, self._coeffs[0], rnd) 

elif mpfr_cmp_ui(x.value, 1) == 0: 

mpfr_set(res.value, self._coeffs[0], rnd) 

for i from 0 < i <= self._degree: 

mpfr_add(res.value, res.value, self._coeffs[i], rnd) 

elif mpfr_cmp_si(x.value, -1) == 0: 

mpfr_set(res.value, self._coeffs[0], rnd) 

for i from 2 <= i <= self._degree by 2: 

mpfr_add(res.value, res.value, self._coeffs[i], rnd) 

for i from 1 <= i <= self._degree by 2: 

mpfr_sub(res.value, res.value, self._coeffs[i], rnd) 

else: 

mpfr_set(res.value, self._coeffs[self._degree], rnd) 

for i from self._degree > i >= 0: 

mpfr_mul(res.value, res.value, x.value, rnd) 

mpfr_add(res.value, res.value, self._coeffs[i], rnd) 

return res 

  

def change_ring(self, R): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import PolynomialRealDense 

sage: f = PolynomialRealDense(RR['x'], [-2, 0, 1.5]) 

sage: f.change_ring(QQ) 

3/2*x^2 - 2 

sage: f.change_ring(RealField(10)) 

1.5*x^2 - 2.0 

sage: f.change_ring(RealField(100)) 

1.5000000000000000000000000000*x^2 - 2.0000000000000000000000000000 

""" 

cdef Py_ssize_t i 

cdef mpfr_rnd_t rnd = self._base_ring.rnd 

cdef PolynomialRealDense f 

if isinstance(R, RealField_class): 

f = PolynomialRealDense(R[self.variable_name()]) 

f = f._new(self._degree) 

for i from 0 <= i <= self._degree: 

mpfr_set(f._coeffs[i], self._coeffs[i], rnd) 

return f 

else: 

return Polynomial.change_ring(self, R) 

  

  

def make_PolynomialRealDense(parent, data): 

""" 

EXAMPLES:: 

  

sage: from sage.rings.polynomial.polynomial_real_mpfr_dense import make_PolynomialRealDense 

sage: make_PolynomialRealDense(RR['x'], [1,2,3]) 

3.00000000000000*x^2 + 2.00000000000000*x + 1.00000000000000 

""" 

return PolynomialRealDense(parent, data)