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

r""" 

Interface to PHC. 

 

PHC computes numerical information about systems of polynomials over 

the complex numbers. 

 

PHC implements polynomial homotopy methods to exploit structure in 

order to better approximate all isolated solutions. The package also 

includes extra tools to handle positive dimensional solution 

components. 

 

AUTHORS: 

 

- PHC was written by J. Verschelde, R. Cools, and many others (?) 

 

- William Stein and Kelly ?? -- first version of interface to PHC 

 

- Marshall Hampton -- second version of interface to PHC 

 

- Marshall Hampton and Alex Jokela -- third version, path tracking 

 

""" 

 

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

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

# Copyright (C) 2008 Marshall Hampton <hamptonio@gmail.com> 

# 

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

# as published by the Free Software Foundation; either version 2 of 

# the License, or (at your option) any later version. 

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

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

from __future__ import print_function 

 

import os 

import re 

import pexpect 

import random 

 

from sage.misc.all import tmp_filename 

from sage.rings.real_mpfr import RR 

from sage.rings.all import CC 

from sage.rings.integer import Integer 

from sage.plot.plot import line, point 

 

 

def get_solution_dicts(output_file_contents, input_ring, get_failures = True): 

""" 

Returns a list of dictionaries of variable:value (key:value) 

pairs. Only used internally; see the solution_dict function in 

the PHC_Object class definition for details. 

 

INPUT: 

 

- output_file_contents -- phc solution output as a string 

- input_ring -- a PolynomialRing that variable names can be coerced into 

 

OUTPUT: 

 

a list of dictionaries of solutions 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x1,x2> = PolynomialRing(QQ,2) 

sage: test_sys = [(x1-1)^5-x2, (x2-1)^5-1] 

sage: sol = phc.blackbox(test_sys, R2) # optional -- phc 

sage: test = get_solution_dicts(sol.output_file_contents,R2) # optional -- phc 

sage: str(sum([q[x1].real() for q in test]))[0:4] # optional -- phc 

'25.0' 

""" 

output_list = output_file_contents.splitlines() 

solution_dicts = [] 

for solution_line in range(len(output_list)-1,-1,-1): 

if output_list[solution_line].find('THE SOLUTIONS') == 0: 

break 

try: 

var_number = int(output_list[solution_line+2].split(' ')[1]) 

# sol_number = int(output_list[solution_line+2].split(' ')[0]) 

except IndexError: 

var_number = int(output_list[solution_line+1].split(' ')[1]) 

# sol_number = int(output_list[solution_line+1].split(' ')[0]) 

for i in range(solution_line + 1,len(output_list)): 

if output_list[i].count('the solution for t') == 1: 

if output_list[i-3].count('success') > 0 or get_failures: 

temp_dict = {} 

for j in range(1,var_number+1): 

rawsplit = output_list[i+j].split(': ')[1].split(' ') 

for extras in range(rawsplit.count('')): 

rawsplit.remove('') 

temp_var = output_list[i+j].split(': ')[0].replace(' ','') 

temp_dict[input_ring(temp_var)] = CC(rawsplit[0],rawsplit[1]) 

solution_dicts.append(temp_dict) 

return solution_dicts 

 

def get_classified_solution_dicts(output_file_contents, input_ring, get_failures = True): 

""" 

Returns a dictionary of lists of dictionaries of variable:value (key:value) 

pairs. Only used internally; see the classified_solution_dict function in 

the PHC_Object class definition for details. 

 

INPUT: 

 

- output_file_contents -- phc solution output as a string 

- input_ring -- a PolynomialRing that variable names can be coerced into 

 

OUTPUT: 

 

- a dictionary of lists if dictionaries of solutions, classifies by type 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x1,x2> = PolynomialRing(QQ,2) 

sage: test_sys = [(x1-2)^5-x2, (x2-1)^5-1] 

sage: sol = phc.blackbox(test_sys, R2) # optional -- phc 

sage: sol_classes = get_classified_solution_dicts(sol.output_file_contents,R2) # optional -- phc 

sage: len(sol_classes['real']) # optional -- phc 

1 

""" 

output_list = output_file_contents.splitlines() 

solution_dicts = {} 

solution_types = ['complex', 'real','failure'] 

for sol_type in solution_types: 

solution_dicts[sol_type] = [] 

for solution_line in range(len(output_list)-1,-1,-1): 

if output_list[solution_line].find('THE SOLUTIONS') == 0: 

break 

var_number = int(output_list[solution_line+2].split(' ')[1]) 

# sol_number = int(output_list[solution_line+2].split(' ')[0]) 

for i in range(solution_line + 1,len(output_list)): 

if output_list[i].count('the solution for t') == 1: 

phc_type = output_list[i+var_number+1].split(' = ')[-1] 

if phc_type.find('complex') != -1: 

phc_type = 'complex' 

elif phc_type.find('real') != -1: 

phc_type = 'real' 

else: 

phc_type = 'failure' 

temp_dict = {} 

for j in range(1,var_number+1): 

rawsplit = output_list[i+j].split(': ')[1].split(' ') 

for extras in range(rawsplit.count('')): 

rawsplit.remove('') 

temp_var = output_list[i+j].split(': ')[0].replace(' ','') 

if phc_type == 'real': 

temp_dict[input_ring(temp_var)] = RR(rawsplit[0]) 

else: 

temp_dict[input_ring(temp_var)] = CC(rawsplit[0],rawsplit[1]) 

solution_dicts[phc_type].append(temp_dict) 

return solution_dicts 

 

def get_variable_list(output_file_contents): 

""" 

Returns the variables, as strings, in the order in which PHCpack has processed them. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x1,x2> = PolynomialRing(QQ,2) 

sage: test_sys = [(x1-2)^5-x2, (x2-1)^5-1] 

sage: sol = phc.blackbox(test_sys, R2) # optional -- phc 

sage: get_variable_list(sol.output_file_contents) # optional -- phc 

['x1', 'x2'] 

""" 

output_list = output_file_contents.splitlines() 

for solution_line in range(len(output_list)-1,-1,-1): 

if output_list[solution_line].find('THE SOLUTIONS') == 0: 

break 

var_number = int(output_list[solution_line+2].split(' ')[1]) 

varlist = [] 

for var_ind in range(var_number): 

var = output_list[solution_line + 8 + var_ind].split(' ')[1] 

varlist.append(var) 

return varlist 

 

class PHC_Object: 

 

def __init__(self, output_file_contents, input_ring): 

""" 

A container for data from the PHCpack program - lists of float 

solutions, etc. Currently the file contents are kept as a string; 

for really large outputs this would be bad. 

 

INPUT: 

 

- output_file_contents: the string output of PHCpack 

- input_ring: for coercion of the variables into the desired ring. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import phc 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [(x-1)^2+(y-1)-1, x^2+y^2-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: str(sum([x[0] for x in sol.solutions()]).real())[0:3] # optional -- phc 

'2.0' 

""" 

self.output_file_contents = output_file_contents 

self.input_ring = input_ring 

 

 

def save_as_start(self, start_filename = None, sol_filter = ''): 

""" 

Saves a solution as a phcpack start file. The usual output is 

just as a string, but it can be saved to a file as well. Even 

if saved to a file, it still returns the output string. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import phc 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^3-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: start_save = sol.save_as_start() # optional -- phc 

sage: end_sys = [x^7-2,y^5-x^2] # optional -- phc 

sage: sol = phc.start_from(start_save, end_sys, R2) # optional -- phc 

sage: len(sol.solutions()) # optional -- phc 

15 

""" 

start_data = '' 

output_list = self.output_file_contents.splitlines() 

for a_line in output_list: 

if a_line.find('ROOT COUNTS') != -1 or a_line.find('START SOLUTIONS') != -1: 

break 

else: 

start_data += a_line + '\n' 

for index in range(len(output_list)-1,0,-1): 

a_line = output_list[index] 

if a_line.find('THE SOLUTIONS') != -1: 

found_solutions = index 

break 

start_data += output_list[found_solutions] + '\n\n' 

try: 

var_number = int(output_list[found_solutions+1].split(' ')[1]) 

except Exception: 

# bad error handling 

var_number = int(output_list[found_solutions+2].split(' ')[1]) 

sol_count = 0 

sol_data = '' 

for i in range(found_solutions + 2, len(output_list)): 

if output_list[i].count('the solution for t') == 1 and output_list[i+1+var_number].find(sol_filter) != -1: 

phc_type = output_list[i+var_number+1].split(' = ')[-1] 

if phc_type.find('no solution') == -1: 

sol_count += 1 

for ind2 in range(i-3,i+var_number+2): 

sol_data += output_list[ind2] + '\n' 

jan_bar = '===========================================================================\n' 

sol_data += jan_bar 

start_data += str(sol_count) + ' ' + str(var_number) + '\n' 

start_data += jan_bar + sol_data 

if start_filename is not None: 

start_file = open(start_filename, 'w') 

start_file.write(start_data) 

start_file.close() 

return start_data 

 

def classified_solution_dicts(self): 

""" 

Returns a dictionary of lists of dictionaries of solutions. 

Its not as crazy as it sounds; the keys are the types of solutions as 

classified by phcpack: regular vs. singular, complex vs. real 

 

INPUT: 

 

- None 

 

OUTPUT: 

 

- A dictionary of lists of dictionaries of solutions 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import phc 

sage: R.<x,y> = PolynomialRing(CC,2) 

sage: p_sys = [x^10-y,y^2-1] 

sage: sol = phc.blackbox(p_sys,R) # optional -- phc 

sage: classifieds = sol.classified_solution_dicts() # optional -- phc 

sage: str(sum([q[y] for q in classifieds['real']]))[0:3] # optional -- phc 

'2.0' 

""" 

try: 

return self.__classified_sols 

except AttributeError: 

pass 

classified_sols = get_classified_solution_dicts(self.output_file_contents, self.input_ring) 

self.__classified_sols = classified_sols 

return classified_sols 

 

def solution_dicts(self, get_failures = False): 

""" 

Returns a list of solutions in dictionary form: variable:value. 

 

INPUT: 

 

- self -- for access to self_out_file_contents, the string 

of raw PHCpack output. 

 

- get_failures (optional) -- a boolean. The default (False) 

is to not process failed homotopies. These either lie on 

positive-dimensional components or at infinity. 

 

OUTPUT: 

 

- solution_dicts: a list of dictionaries. Each dictionary 

element is of the form variable:value, where the variable 

is an element of the input_ring, and the value is in 

ComplexField. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R.<x,y,z> = PolynomialRing(QQ,3) 

sage: fs = [x^2-1,y^2-x,z^2-y] 

sage: sol = phc.blackbox(fs,R) # optional -- phc 

sage: s_list = sol.solution_dicts() # optional -- phc 

sage: s_list.sort() # optional -- phc 

sage: s_list[0] # optional -- phc 

{y: 1.00000000000000, z: -1.00000000000000, x: 1.00000000000000} 

""" 

try: 

return self.__solution_dicts 

except AttributeError: 

pass 

solution_dicts = get_solution_dicts(self.output_file_contents, self.input_ring, get_failures = get_failures) 

self.__solution_dicts = solution_dicts 

return solution_dicts 

 

def solutions(self, get_failures = False): 

""" 

Returns a list of solutions in the ComplexField. 

 

Use the variable_list function to get the order of variables used by 

PHCpack, which is usually different than the term order of the 

input_ring. 

 

INPUT: 

 

- self -- for access to self_out_file_contents, the string 

of raw PHCpack output. 

- get_failures (optional) -- a boolean. The default (False) 

is to not process failed homotopies. These either lie on 

positive-dimensional components or at infinity. 

 

OUTPUT: 

 

- solutions: a list of lists of ComplexField-valued solutions. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x1,x2> = PolynomialRing(QQ,2) 

sage: test_sys = [x1^5-x1*x2^2-1, x2^5-x1*x2-1] 

sage: sol = phc.blackbox(test_sys, R2) # optional -- phc 

sage: len(sol.solutions()) # optional -- phc 

25 

""" 

try: 

return self.__solutions 

except AttributeError: 

pass 

solution_dicts = get_solution_dicts(self.output_file_contents, self.input_ring, get_failures = get_failures) 

self.__solution_dicts = solution_dicts 

solutions = [sol_dict.values() for sol_dict in solution_dicts] 

self.__solutions = solutions 

return solutions 

 

def variable_list(self): 

""" 

Returns the variables, as strings, in the order in which 

PHCpack has processed them. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x1,x2> = PolynomialRing(QQ,2) 

sage: test_sys = [x1^5-x1*x2^2-1, x2^5-x1*x2-1] 

sage: sol = phc.blackbox(test_sys, R2) # optional -- phc 

sage: sol.variable_list() # optional -- phc 

['x1', 'x2'] 

""" 

try: 

return self.__var_list 

except AttributeError: 

pass 

var_list = get_variable_list(self.output_file_contents) 

self.__var_list = var_list 

return var_list 

 

class PHC: 

""" 

A class to interface with PHCpack, for computing numerical 

homotopies and root counts. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import phc 

sage: R.<x,y> = PolynomialRing(CDF,2) 

sage: testsys = [x^2 + 1, x*y - 1] 

sage: phc.mixed_volume(testsys) # optional -- phc 

2 

sage: v = phc.blackbox(testsys, R) # optional -- phc 

sage: sols = v.solutions() # optional -- phc 

sage: sols.sort() # optional -- phc 

sage: sols # optional -- phc 

[[-1.00000000000000*I, 1.00000000000000*I], [1.00000000000000*I, -1.00000000000000*I]] 

sage: sol_dict = v.solution_dicts() # optional -- phc 

sage: x_sols_from_dict = [d[x] for d in sol_dict] # optional -- phc 

sage: x_sols_from_dict.sort(); x_sols_from_dict # optional -- phc 

[-1.00000000000000*I, 1.00000000000000*I] 

sage: residuals = [[test_equation.change_ring(CDF).subs(sol) for test_equation in testsys] for sol in v.solution_dicts()] # optional -- phc 

sage: residuals # optional -- phc 

[[0, 0], [0, 0]] 

""" 

 

def _output_from_command_list(self, command_list, polys, verbose = False): 

""" 

A pexpect interface to phcpack, given a command list for 

interactive dialogs. The input file is supplied from the 

polynomial list, output file is also supplied. This is 

only used as a building block for the interface. 

 

INPUT: 

 

- command_list -- a list of commands to phc 

- polys -- a polynomial system as a list of polynomials 

 

OUTPUT: 

 

- an output string from phc 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [(x-1)^2+(y-1)-1, x^2+y^2-1] # optional -- phc 

sage: a = phc._output_from_command_list(['phc -m','4','n','n','n'], start_sys) # optional -- phc 

""" 

# Get temporary file names (these will be in SAGE_HOME/.sage/tmp/pid) 

input_filename = tmp_filename() 

output_filename = tmp_filename() 

 

# Get the input polynomial text 

input = self._input_file(polys) 

if verbose: 

print("Writing the input file to %s" % input_filename) 

open(input_filename, 'w').write(input) 

 

if verbose: 

print("The following file will be the input polynomial file to phc.") 

print(input) 

 

# Create a phc process 

child_phc = pexpect.spawn(command_list[0]) 

# feed it the commands 

child_phc.sendline('y') 

child_phc.sendline(input_filename) 

child_phc.sendline(output_filename) 

for command_string in command_list[1:]: 

if verbose: 

print(command_string) 

child_phc.sendline(command_string) 

child_phc.expect('results') 

read_stuff = child_phc.read() 

if verbose: 

print(read_stuff) 

child_phc.close() 

if not os.path.exists(output_filename): 

raise RuntimeError("The output file does not exist; something went wrong running phc.") 

 

# Delete the input file 

os.unlink(input_filename) 

 

# Return the output filename 

return output_filename 

 

def _input_file(self, polys): 

""" 

This is used internally to implement the PHC interface. 

 

INPUT: 

 

- polys -- a list of polynomials in a Sage polynomial ring 

over a field that embeds into the complex 

numbers. 

 

OUTPUT: 

 

- a PHC input file (as a text string) that describes these - 

polynomials. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [(x-1)^2+(y-1)-1, x^2+y^2-1] 

sage: phc._input_file(start_sys) # optional -- phc 

'2\nx^2 - 2*x + y - 1;\nx^2 + y^2 - 1;\n' 

""" 

if not isinstance(polys, (list, tuple)): 

raise TypeError('polys must be a list or tuple') 

s = '%s\n'%len(polys) 

for f in polys: 

s += f._repr_() + ';\n' # note the semicolon *terminators* 

 

return s 

 

def _parse_path_file(self, input_filename, verbose = False): 

""" 

Takes a phpack output file containing path tracking information 

and parses it into a list of lists of dictionaries - i.e. a 

list of solutions paths, where each solution path is a list of 

dictionaries of variable and homotopy parameter values. 

 

INPUT: 

 

- input_filename -- file must have path-tracking information 

 

OUTPUT: 

 

- a list of lists of dictionaries, described above 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^5-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: start_save = sol.save_as_start() # optional -- phc 

sage: end_sys = [x^5-2,y^5-x^2] # optional -- phc 

sage: path_track_filename = phc._path_track_file(start_save, end_sys, R2, c_skew = .001) # optional -- phc 

sage: sol_paths = phc._parse_path_file(path_track_filename) # optional -- phc 

sage: len(sol_paths) # optional -- phc 

25 

""" 

 

if not os.path.exists(input_filename): 

raise RuntimeError("The file containing output from phc (" + input_filename + ") cannot be found") 

 

fh = open(input_filename) 

line_idx = 0 

begin = 0 

count = 0 

solutions_dicts = [] 

steps_dicts = [] 

 

# regular expressions for matching certain output types 

var_cnt_regex = re.compile('^ +([0-9]+)') 

output_regex = re.compile('^OUTPUT INFORMATION DURING') 

t_regex = re.compile('(^t +: +(-{0,1}[0-9]+\.[0-9]+E[-+][0-9]+) +(-{0,1}[0-9]+\.[0-9]+E[-+][0-9]+)$)', re.IGNORECASE) 

sols_regex = re.compile('(^ *(([a-z]|[0-9])+) +: +(-?[0-9]+\.[0-9]+E[-+][0-9]+) +(-?[0-9]+\.[0-9]+E[-+][0-9]+)$)', re.IGNORECASE) 

complete_regex= re.compile('^TIMING INFORMATION') 

 

breakfast = False 

a_line = fh.readline() 

end_test = '' 

while a_line: 

# processing.... 

a_line = a_line.replace("\n", '') 

if line_idx == 0: 

m = var_cnt_regex.match(a_line) 

if m: 

count = Integer(m.group(1)) 

if count > 0: 

m = output_regex.match(a_line) 

if m: 

begin = 1 

if begin: 

m = t_regex.match(a_line) 

if m: 

# put the t-values into a dict 

# m.group(2) contains the real val 

# m.group(3) contains the imaginary val 

# fh_w.write( "T=> G1(" + m.group(2) + '),G2(' + m.group(3) + ")\n") 

# read off two lines - this should be 'm' and 'the solution for t :' 

a_line = fh.readline() 

end_test = a_line # store this to check for end of solution 

a_line = fh.readline() 

t_val = CC(m.group(2), m.group(3)) 

temp_dict = {} 

temp_dict["t"] = t_val 

for i in range(0, count): 

a_line = fh.readline() 

m = sols_regex.match(a_line) 

if m: 

# m.group(2) contains our var name 

# m.group(4) contains our real val 

# m.group(5) contains our imaginary val 

temp_dict[m.group(2)] = CC(m.group(4),m.group(5)) 

steps_dicts.append(temp_dict) 

# check if its the end of a solution 

if end_test.find('Length of path') != -1: 

if verbose: 

print("recording sol") 

if steps_dicts != []: 

solutions_dicts.append(steps_dicts) 

steps_dicts = [] 

m = complete_regex.match(a_line) 

if m: 

breakfast = True 

if breakfast: 

break 

line_idx += 1 

a_line = fh.readline() 

fh.close() 

return solutions_dicts 

 

def _path_track_file(self, start_filename_or_string, polys, input_ring, c_skew = 0.001, verbose = False): 

""" 

Returns the filename which contains path tracking output. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^6-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: start_save = sol.save_as_start() # optional -- phc 

sage: end_sys = [x^7-2,y^5-x^2] # optional -- phc 

sage: path_track_filename = phc._path_track_file(start_save, end_sys, R2, c_skew = .001) # optional -- phc 

sage: sol_paths = phc._parse_path_file(path_track_filename) # optional -- phc 

sage: len(sol_paths) # optional -- phc 

30 

""" 

# Probably unnecessarily redundant from the start_from function 

if start_filename_or_string.find('THE SOLUTIONS') != -1: 

start_filename = tmp_filename() 

start_file = open(start_filename, 'w') 

start_file.write(start_filename_or_string) 

start_file.close() 

elif os.path.exists(start_filename_or_string): 

start_filename = start_filename_or_string 

else: 

raise RuntimeError("There is something wrong with your start string or filename") 

 

return self._output_from_command_list(['phc','0','0','A',start_filename, 'y','1','0','n','k','2','a','1',str(c_skew),'0','0','2'], polys, verbose = verbose) 

 

def path_track(self, start_sys, end_sys, input_ring, c_skew = .001, saved_start = None): 

""" 

This function computes homotopy paths between the solutions of start_sys and end_sys. 

 

INPUT: 

 

- start_sys -- a square polynomial system, given as a list of polynomials 

- end_sys -- same type as start_sys 

- input_ring -- for coercion of the variables into the desired ring. 

- c_skew -- optional. the imaginary part of homotopy multiplier; nonzero values 

are often necessary to avoid intermediate path collisions 

- saved_start -- optional. A phc output file. If not given, start system solutions 

are computed via the phc.blackbox function. 

 

OUTPUT: 

 

- a list of paths as dictionaries, with the keys variables and t-values on the path. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^6-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: start_save = sol.save_as_start() # optional -- phc 

sage: end_sys = [x^7-2,y^5-x^2] # optional -- phc 

sage: sol_paths = phc.path_track(start_sys, end_sys, R2, saved_start = start_save) # optional -- phc 

sage: len(sol_paths) # optional -- phc 

30 

""" 

 

if not saved_start: 

sol = phc.blackbox(start_sys, input_ring) 

saved_start = sol.save_as_start() 

path_track_filename = phc._path_track_file(saved_start, end_sys, input_ring = input_ring, c_skew = c_skew) 

sol_paths = phc._parse_path_file(path_track_filename) 

os.unlink(path_track_filename) 

return sol_paths 

 

def plot_paths_2d(self, start_sys, end_sys, input_ring, c_skew = .001, endpoints = True, saved_start = None, rand_colors = False): 

""" 

This returns a graphics object of solution paths in the complex plane. 

 

INPUT: 

 

- start_sys -- a square polynomial system, given as a list of polynomials 

- end_sys -- same type as start_sys 

- input_ring -- for coercion of the variables into the desired ring. 

- c_skew -- optional. the imaginary part of homotopy multiplier; nonzero values 

are often necessary to avoid intermediate path collisions 

- endpoints -- optional. Whether to draw in the ends of paths as points. 

- saved_start -- optional. A phc output file. If not given, start system solutions 

are computed via the phc.blackbox function. 

 

OUTPUT: 

 

- lines and points of solution paths 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: from sage.structure.sage_object import SageObject 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^5-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: start_save = sol.save_as_start() # optional -- phc 

sage: end_sys = [x^5-25,y^5-x^2] # optional -- phc 

sage: testing = phc.plot_paths_2d(start_sys, end_sys, R2) # optional -- phc 

sage: type(testing) # optional -- phc (normally use plot here) 

<class 'sage.plot.graphics.Graphics'> 

""" 

paths = phc.path_track(start_sys, end_sys, input_ring, c_skew = c_skew, saved_start = saved_start) 

path_lines = [] 

sol_pts = [] 

if rand_colors: 

r_color = {} 

for a_var in input_ring.gens(): 

var_name = str(a_var) 

r_color[var_name] = (random(),random(),random()) 

for a_sol in paths: 

for a_var in input_ring.gens(): 

var_name = str(a_var) 

temp_line = [] 

for data in a_sol: 

temp_line.append([data[var_name].real(), data[var_name].imag()]) 

if rand_colors: 

path_lines.append(line(temp_line, rgbcolor = r_color[var_name])) 

else: 

path_lines.append(line(temp_line)) 

if endpoints: 

sol_pts = [] 

for a_sol in paths: 

for a_var in input_ring.gens(): 

var_name = str(a_var) 

sol_pts.append(point([a_sol[0][var_name].real(), a_sol[0][var_name].imag()])) 

sol_pts.append(point([a_sol[-1][var_name].real(), a_sol[-1][var_name].imag()])) 

return sum(sol_pts) + sum(path_lines) 

else: 

return sum(path_lines) 

 

def mixed_volume(self, polys, verbose=False): 

""" 

Computes the mixed volume of the polynomial system given by the input polys. 

 

INPUT: 

 

- polys -- a list of multivariate polynomials (elements of a multivariate 

polynomial ring). 

- verbose -- print lots of verbose information about what this function does. 

 

OUTPUT: 

 

- The mixed volume. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y,z> = PolynomialRing(QQ,3) 

sage: test_sys = [(x+y+z)^2-1,x^2-x,y^2-1] 

sage: phc.mixed_volume(test_sys) # optional -- phc 

4 

""" 

output_filename = self._output_from_command_list(['phc -m','4','n','n','n'], polys, verbose = verbose) 

 

out = open(output_filename).read() 

# All done 

out_lines = out.split('\n') 

for a_line in out_lines: 

# the two conditions below are necessary because of changes in output format 

if a_line.find('The mixed volume equals :') == 0 or a_line.find('common mixed volume :') == 0: 

if verbose: 

print('found line: ' + a_line) 

mixed_vol = Integer(a_line.split(':')[1]) 

break 

 

try: 

return mixed_vol 

except NameError: 

raise RuntimeError("Mixed volume not found in output; something went wrong running phc.") 

 

 

def start_from(self, start_filename_or_string, polys, input_ring, path_track_file = None, verbose = False): 

""" 

This computes solutions starting from a phcpack solution file. 

 

INPUT: 

 

- start_filename_or_string -- the filename for a phcpack start system, 

or the contents of such a file as a string. Variable names must match 

the inputring variables. The value of the homotopy variable t should 

be 1, not 0. 

- polys -- a list of multivariate polynomials (elements of a multivariate 

polynomial ring). 

- input_ring: for coercion of the variables into the desired ring. 

- path_track_file: whether to save path-tracking information 

- verbose -- print lots of verbose information about what this function does. 

 

OUTPUT: 

 

- A solution in the form of a PHCObject. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^6-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: start_save = sol.save_as_start() # optional -- phc 

sage: end_sys = [x^7-2,y^5-x^2] # optional -- phc 

sage: sol = phc.start_from(start_save, end_sys, R2) # optional -- phc 

sage: len(sol.solutions()) # optional -- phc 

30 

""" 

input_filename = tmp_filename() 

output_filename = tmp_filename() 

 

if start_filename_or_string.find('THE SOLUTIONS') != -1: 

start_filename = tmp_filename() 

start_file = open(start_filename,'w') 

start_file.write(start_filename_or_string) 

start_file.close() 

elif os.path.exists(start_filename_or_string): 

start_filename = start_filename_or_string 

else: 

raise RuntimeError("There is something wrong with your start string or filename") 

 

# Get the input polynomial text 

input = self._input_file(polys) 

if verbose: 

print("Writing the input file to %s" % input_filename) 

open(input_filename, 'w').write(input) 

 

if verbose: 

print("The following file will be the input polynomial file to phc.") 

print(input) 

 

# Create a phc process 

child_phc = pexpect.spawn('phc') 

child_phc.sendline('y') 

child_phc.sendline(input_filename) 

child_phc.sendline(output_filename) 

child_phc.sendline('0') 

child_phc.sendline('0') 

child_phc.expect('Nonlinear Reduction') 

child_phc.sendline('A') 

child_phc.sendline(start_filename) 

child_phc.sendline('y') 

child_phc.sendline('1') 

child_phc.sendline('0') 

if verbose: 

phc_dialog = child_phc.read(size = 40) 

print(phc_dialog) 

child_phc.sendline('n') 

child_phc.sendline('0') 

if verbose: 

child_phc.expect('CURRENT CONTINUATION') 

phc_dialog = child_phc.read(size = 40) 

print(phc_dialog) 

child_phc.sendline('0') 

if path_track_file is None: 

child_phc.sendline('0') 

else: 

child_phc.sendline('2') 

child_phc.expect('results') 

dots = child_phc.read() 

if verbose: 

print("should be . : " + dots) 

 

#close down the process: 

child_phc.close() 

if not os.path.exists(output_filename): 

raise RuntimeError("The output file does not exist; something went wrong running phc.") 

 

# Read the output produced by PHC 

out = open(output_filename).read() 

 

# Delete the temporary files 

os.unlink(output_filename) 

os.unlink(input_filename) 

 

# All done 

return PHC_Object(out, input_ring) 

 

def blackbox(self, polys, input_ring, verbose = False): 

""" 

Returns as a string the result of running PHC with the given polynomials 

under blackbox mode (the '-b' option). 

 

INPUT: 

 

- polys -- a list of multivariate polynomials (elements of a multivariate 

polynomial ring). 

- input_ring -- for coercion of the variables into the desired ring. 

- verbose -- print lots of verbose information about what this function does. 

 

OUTPUT: 

 

- a PHC_Object object containing the phcpack output string. 

 

EXAMPLES:: 

 

sage: from sage.interfaces.phc import * 

sage: R2.<x,y> = PolynomialRing(QQ,2) 

sage: start_sys = [x^6-y^2,y^5-1] 

sage: sol = phc.blackbox(start_sys, R2) # optional -- phc 

sage: len(sol.solutions()) # optional -- phc 

30 

""" 

 

# Get three temporary file names (these will be in SAGE_HOME/.sage/tmp/pid) 

input_filename = tmp_filename() 

output_filename = input_filename + ".phc" 

log_filename = tmp_filename() 

 

# Get the input polynomial text 

input = self._input_file(polys) 

if verbose: 

print("Writing the input file to %s" % input_filename) 

open(input_filename, 'w').write(input) 

 

if verbose: 

print("The following file will be the input polynomial file to phc.") 

print(input) 

 

# Create the phc command line> 

cmd = 'phc -b %s %s'%(input_filename, output_filename) 

 

if verbose: 

print("The phc command line is:") 

print(cmd) 

 

# Do it -- make the system call. 

e = os.system(cmd) 

 

# Was there an error? 

if e: 

from sage.misc.sage_ostools import have_program 

if not have_program('phc'): 

print(os.system('which phc') + ' PHC needs to be installed and in your path') 

raise RuntimeError 

# todo -- why? etc. 

raise RuntimeError(open(log_filename).read() + "\nError running phc.") 

 

if not os.path.exists(output_filename): 

raise RuntimeError("The output file does not exist; something went wrong running phc.") 

 

# Read the output produced by PHC 

out = open(output_filename).read() 

 

# All done 

return PHC_Object(out, input_ring) 

 

 

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

 

# The unique phc interface instance. 

phc = PHC()