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

import numpy as np 

import copy 

 

try: 

from functools import lru_cache # For Python 3 only 

except ImportError: # Python 2: 

import time 

import functools 

import collections 

 

# Note to avoid using external packages such as functools32 we use this code 

# only using the standard library 

def lru_cache(maxsize=255, timeout=None): 

""" 

Thanks to ilialuk @ https://stackoverflow.com/users/2121105/ilialuk for 

this code snippet. Modifications by S. Endres 

""" 

 

class LruCacheClass(object): 

def __init__(self, input_func, max_size, timeout): 

self._input_func = input_func 

self._max_size = max_size 

self._timeout = timeout 

 

# This will store the cache for this function, 

# format - {caller1 : [OrderedDict1, last_refresh_time1], 

# caller2 : [OrderedDict2, last_refresh_time2]}. 

# In case of an instance method - the caller is the instance, 

# in case called from a regular function - the caller is None. 

self._caches_dict = {} 

 

def cache_clear(self, caller=None): 

# Remove the cache for the caller, only if exists: 

if caller in self._caches_dict: 

del self._caches_dict[caller] 

self._caches_dict[caller] = [collections.OrderedDict(), 

time.time()] 

 

def __get__(self, obj, objtype): 

""" Called for instance methods """ 

return_func = functools.partial(self._cache_wrapper, obj) 

return_func.cache_clear = functools.partial(self.cache_clear, 

obj) 

# Return the wrapped function and wraps it to maintain the 

# docstring and the name of the original function: 

return functools.wraps(self._input_func)(return_func) 

 

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

""" Called for regular functions """ 

return self._cache_wrapper(None, *args, **kwargs) 

 

# Set the cache_clear function in the __call__ operator: 

__call__.cache_clear = cache_clear 

 

def _cache_wrapper(self, caller, *args, **kwargs): 

# Create a unique key including the types (in order to 

# differentiate between 1 and '1'): 

kwargs_key = "".join(map( 

lambda x: str(x) + str(type(kwargs[x])) + str(kwargs[x]), 

sorted(kwargs))) 

key = "".join( 

map(lambda x: str(type(x)) + str(x), args)) + kwargs_key 

 

# Check if caller exists, if not create one: 

if caller not in self._caches_dict: 

self._caches_dict[caller] = [collections.OrderedDict(), 

time.time()] 

else: 

# Validate in case the refresh time has passed: 

if self._timeout is not None: 

if (time.time() - self._caches_dict[caller][1] 

> self._timeout): 

self.cache_clear(caller) 

 

# Check if the key exists, if so - return it: 

cur_caller_cache_dict = self._caches_dict[caller][0] 

if key in cur_caller_cache_dict: 

return cur_caller_cache_dict[key] 

 

# Validate we didn't exceed the max_size: 

if len(cur_caller_cache_dict) >= self._max_size: 

# Delete the first item in the dict: 

try: 

cur_caller_cache_dict.popitem(False) 

except KeyError: 

pass 

# Call the function and store the data in the cache (call it 

# with the caller in case it's an instance function) 

if caller is not None: 

args = (caller,) + args 

cur_caller_cache_dict[key] = self._input_func(*args, **kwargs) 

 

return cur_caller_cache_dict[key] 

 

# Return the decorator wrapping the class (also wraps the instance to 

# maintain the docstring and the name of the original function): 

return (lambda input_func: functools.wraps(input_func)( 

LruCacheClass(input_func, maxsize, timeout))) 

 

 

class Complex: 

def __init__(self, dim, func, func_args=(), symmetry=False, bounds=None, 

g_cons=None, g_args=()): 

self.dim = dim 

self.bounds = bounds 

self.symmetry = symmetry # TODO: Define the functions to be used 

# here in init to avoid if checks 

self.gen = 0 

self.perm_cycle = 0 

 

# Every cell is stored in a list of its generation, 

# ex. the initial cell is stored in self.H[0] 

# 1st get new cells are stored in self.H[1] etc. 

# When a cell is subgenerated it is removed from this list 

 

self.H = [] # Storage structure of cells 

# Cache of all vertices 

self.V = VertexCache(func, func_args, bounds, g_cons, g_args) 

 

# Generate n-cube here: 

self.n_cube(dim, symmetry=symmetry) 

 

# TODO: Assign functions to a the complex instead 

if symmetry: 

self.generation_cycle = 1 

# self.centroid = self.C0()[-1].x 

# self.C0.centroid = self.centroid 

else: 

self.add_centroid() 

 

self.H.append([]) 

self.H[0].append(self.C0) 

self.hgr = self.C0.homology_group_rank() 

self.hgrd = 0 # Complex group rank differential 

# self.hgr = self.C0.hg_n 

 

# Build initial graph 

self.graph_map() 

 

self.performance = [] 

self.performance.append(0) 

self.performance.append(0) 

 

def __call__(self): 

return self.H 

 

def n_cube(self, dim, symmetry=False, printout=False): 

""" 

Generate the simplicial triangulation of the n dimensional hypercube 

containing 2**n vertices 

""" 

origin = list(np.zeros(dim, dtype=int)) 

self.origin = origin 

supremum = list(np.ones(dim, dtype=int)) 

self.supremum = supremum 

 

# tuple versions for indexing 

origintuple = tuple(origin) 

supremumtuple = tuple(supremum) 

 

x_parents = [origintuple] 

 

if symmetry: 

self.C0 = Simplex(0, 0, 0, self.dim) # Initial cell object 

self.C0.add_vertex(self.V[origintuple]) 

 

i_s = 0 

self.perm_symmetry(i_s, x_parents, origin) 

self.C0.add_vertex(self.V[supremumtuple]) 

else: 

self.C0 = Cell(0, 0, origin, supremum) # Initial cell object 

self.C0.add_vertex(self.V[origintuple]) 

self.C0.add_vertex(self.V[supremumtuple]) 

 

i_parents = [] 

self.perm(i_parents, x_parents, origin) 

 

if printout: 

print("Initial hyper cube:") 

for v in self.C0(): 

v.print_out() 

 

def perm(self, i_parents, x_parents, xi): 

# TODO: Cut out of for if outside linear constraint cutting planes 

xi_t = tuple(xi) 

 

# Construct required iterator 

iter_range = [x for x in range(self.dim) if x not in i_parents] 

 

for i in iter_range: 

i2_parents = copy.copy(i_parents) 

i2_parents.append(i) 

xi2 = copy.copy(xi) 

xi2[i] = 1 

# Make new vertex list a hashable tuple 

xi2_t = tuple(xi2) 

# Append to cell 

self.C0.add_vertex(self.V[xi2_t]) 

# Connect neighbours and vice versa 

# Parent point 

self.V[xi2_t].connect(self.V[xi_t]) 

 

# Connect all family of simplices in parent containers 

for x_ip in x_parents: 

self.V[xi2_t].connect(self.V[x_ip]) 

 

x_parents2 = copy.copy(x_parents) 

x_parents2.append(xi_t) 

 

# Permutate 

self.perm(i2_parents, x_parents2, xi2) 

 

def perm_symmetry(self, i_s, x_parents, xi): 

# TODO: Cut out of for if outside linear constraint cutting planes 

xi_t = tuple(xi) 

xi2 = copy.copy(xi) 

xi2[i_s] = 1 

# Make new vertex list a hashable tuple 

xi2_t = tuple(xi2) 

# Append to cell 

self.C0.add_vertex(self.V[xi2_t]) 

# Connect neighbours and vice versa 

# Parent point 

self.V[xi2_t].connect(self.V[xi_t]) 

 

# Connect all family of simplices in parent containers 

for x_ip in x_parents: 

self.V[xi2_t].connect(self.V[x_ip]) 

 

x_parents2 = copy.copy(x_parents) 

x_parents2.append(xi_t) 

 

i_s += 1 

if i_s == self.dim: 

return 

# Permutate 

self.perm_symmetry(i_s, x_parents2, xi2) 

 

def add_centroid(self): 

"""Split the central edge between the origin and supremum of 

a cell and add the new vertex to the complex""" 

self.centroid = list( 

(np.array(self.origin) + np.array(self.supremum)) / 2.0) 

self.C0.add_vertex(self.V[tuple(self.centroid)]) 

self.C0.centroid = self.centroid 

 

# Disconnect origin and supremum 

self.V[tuple(self.origin)].disconnect(self.V[tuple(self.supremum)]) 

 

# Connect centroid to all other vertices 

for v in self.C0(): 

self.V[tuple(self.centroid)].connect(self.V[tuple(v.x)]) 

 

self.centroid_added = True 

return 

 

# Construct incidence array: 

def incidence(self): 

if self.centroid_added: 

self.structure = np.zeros([2 ** self.dim + 1, 2 ** self.dim + 1], 

dtype=int) 

else: 

self.structure = np.zeros([2 ** self.dim, 2 ** self.dim], 

dtype=int) 

 

for v in self.HC.C0(): 

for v2 in v.nn: 

self.structure[v.index, v2.index] = 1 

 

return 

 

# A more sparse incidence generator: 

def graph_map(self): 

""" Make a list of size 2**n + 1 where an entry is a vertex 

incidence, each list element contains a list of indexes 

corresponding to that entries neighbours""" 

 

self.graph = [[v2.index for v2 in v.nn] for v in self.C0()] 

 

# Graph structure method: 

# 0. Capture the indices of the initial cell. 

# 1. Generate new origin and supremum scalars based on current generation 

# 2. Generate a new set of vertices corresponding to a new 

# "origin" and "supremum" 

# 3. Connected based on the indices of the previous graph structure 

# 4. Disconnect the edges in the original cell 

 

def sub_generate_cell(self, C_i, gen): 

"""Subgenerate a cell `C_i` of generation `gen` and 

homology group rank `hgr`.""" 

origin_new = tuple(C_i.centroid) 

centroid_index = len(C_i()) - 1 

 

# If not gen append 

try: 

self.H[gen] 

except IndexError: 

self.H.append([]) 

 

# Generate subcubes using every extreme vertex in C_i as a supremum 

# and the centroid of C_i as the origin 

H_new = [] # list storing all the new cubes split from C_i 

for i, v in enumerate(C_i()[:-1]): 

supremum = tuple(v.x) 

H_new.append( 

self.construct_hypercube(origin_new, supremum, gen, C_i.hg_n)) 

 

for i, connections in enumerate(self.graph): 

# Present vertex V_new[i]; connect to all connections: 

if i == centroid_index: # Break out of centroid 

break 

 

for j in connections: 

C_i()[i].disconnect(C_i()[j]) 

 

# Destroy the old cell 

if C_i is not self.C0: # Garbage collector does this anyway; not needed 

del C_i 

 

# TODO: Recalculate all the homology group ranks of each cell 

return H_new 

 

def split_generation(self): 

""" 

Run sub_generate_cell for every cell in the current complex self.gen 

""" 

no_splits = False # USED IN SHGO 

try: 

for c in self.H[self.gen]: 

if self.symmetry: 

# self.sub_generate_cell_symmetry(c, self.gen + 1) 

self.split_simplex_symmetry(c, self.gen + 1) 

else: 

self.sub_generate_cell(c, self.gen + 1) 

except IndexError: 

no_splits = True # USED IN SHGO 

 

self.gen += 1 

return no_splits # USED IN SHGO 

 

# @lru_cache(maxsize=None) 

def construct_hypercube(self, origin, supremum, gen, hgr, 

printout=False): 

""" 

Build a hypercube with triangulations symmetric to C0. 

 

Parameters 

---------- 

origin : vec 

supremum : vec (tuple) 

gen : generation 

hgr : parent homology group rank 

""" 

 

# Initiate new cell 

C_new = Cell(gen, hgr, origin, supremum) 

C_new.centroid = tuple( 

(np.array(origin) + np.array(supremum)) / 2.0) 

 

# Build new indexed vertex list 

V_new = [] 

 

# Cached calculation 

for i, v in enumerate(self.C0()[:-1]): 

t1 = self.generate_sub_cell_t1(origin, v.x) 

t2 = self.generate_sub_cell_t2(supremum, v.x) 

 

vec = t1 + t2 

 

vec = tuple(vec) 

C_new.add_vertex(self.V[vec]) 

V_new.append(vec) 

 

# Add new centroid 

C_new.add_vertex(self.V[C_new.centroid]) 

V_new.append(C_new.centroid) 

 

# Connect new vertices #TODO: Thread into other loop; no need for V_new 

for i, connections in enumerate(self.graph): 

# Present vertex V_new[i]; connect to all connections: 

for j in connections: 

self.V[V_new[i]].connect(self.V[V_new[j]]) 

 

if printout: 

print("A sub hyper cube with:") 

print("origin: {}".format(origin)) 

print("supremum: {}".format(supremum)) 

for v in C_new(): 

v.print_out() 

 

# Append the new cell to the to complex 

self.H[gen].append(C_new) 

 

return C_new 

 

def split_simplex_symmetry(self, S, gen): 

""" 

Split a hypersimplex S into two sub simplices by building a hyperplane 

which connects to a new vertex on an edge (the longest edge in 

dim = {2, 3}) and every other vertex in the simplex that is not 

connected to the edge being split. 

 

This function utilizes the knowledge that the problem is specified 

with symmetric constraints 

 

The longest edge is tracked by an ordering of the 

vertices in every simplices, the edge between first and second 

vertex is the longest edge to be split in the next iteration. 

""" 

# If not gen append 

try: 

self.H[gen] 

except IndexError: 

self.H.append([]) 

 

# Find new vertex. 

# V_new_x = tuple((np.array(C()[0].x) + np.array(C()[1].x)) / 2.0) 

s = S() 

firstx = s[0].x 

lastx = s[-1].x 

V_new = self.V[tuple((np.array(firstx) + np.array(lastx)) / 2.0)] 

 

# Disconnect old longest edge 

self.V[firstx].disconnect(self.V[lastx]) 

 

# Connect new vertices to all other vertices 

for v in s[:]: 

v.connect(self.V[V_new.x]) 

 

# New "lower" simplex 

S_new_l = Simplex(gen, S.hg_n, self.generation_cycle, 

self.dim) 

S_new_l.add_vertex(s[0]) 

S_new_l.add_vertex(V_new) # Add new vertex 

for v in s[1:-1]: # Add all other vertices 

S_new_l.add_vertex(v) 

 

# New "upper" simplex 

S_new_u = Simplex(gen, S.hg_n, S.generation_cycle, self.dim) 

 

# First vertex on new long edge 

S_new_u.add_vertex(s[S_new_u.generation_cycle + 1]) 

 

for v in s[1:-1]: # Remaining vertices 

S_new_u.add_vertex(v) 

 

for k, v in enumerate(s[1:-1]): # iterate through inner vertices 

if k == S.generation_cycle: 

S_new_u.add_vertex(V_new) 

else: 

S_new_u.add_vertex(v) 

 

S_new_u.add_vertex(s[-1]) # Second vertex on new long edge 

 

self.H[gen].append(S_new_l) 

self.H[gen].append(S_new_u) 

 

return 

 

@lru_cache(maxsize=None) 

def generate_sub_cell_2(self, origin, supremum, v_x_t): # No hits 

""" 

Use the origin and supremum vectors to find a new cell in that 

subspace direction 

 

NOTE: NOT CURRENTLY IN USE! 

 

Parameters 

---------- 

origin : tuple vector (hashable) 

supremum : tuple vector (hashable) 

 

Returns 

------- 

 

""" 

t1 = self.generate_sub_cell_t1(origin, v_x_t) 

t2 = self.generate_sub_cell_t2(supremum, v_x_t) 

vec = t1 + t2 

return tuple(vec) 

 

@lru_cache(maxsize=None) 

def generate_sub_cell_t1(self, origin, v_x): 

# TODO: Calc these arrays outside 

v_o = np.array(origin) 

return v_o - v_o * np.array(v_x) 

 

@lru_cache(maxsize=None) 

def generate_sub_cell_t2(self, supremum, v_x): 

v_s = np.array(supremum) 

return v_s * np.array(v_x) 

 

# Plots 

def plot_complex(self): 

""" 

Here C is the LIST of simplexes S in the 

2 or 3 dimensional complex 

 

To plot a single simplex S in a set C, use ex. [C[0]] 

""" 

from matplotlib import pyplot 

if self.dim == 2: 

pyplot.figure() 

for C in self.H: 

for c in C: 

for v in c(): 

if self.bounds is None: 

x_a = np.array(v.x, dtype=float) 

else: 

x_a = np.array(v.x, dtype=float) 

for i in range(len(self.bounds)): 

x_a[i] = (x_a[i] * (self.bounds[i][1] 

- self.bounds[i][0]) 

+ self.bounds[i][0]) 

 

# logging.info('v.x_a = {}'.format(x_a)) 

 

pyplot.plot([x_a[0]], [x_a[1]], 'o') 

 

xlines = [] 

ylines = [] 

for vn in v.nn: 

if self.bounds is None: 

xn_a = np.array(vn.x, dtype=float) 

else: 

xn_a = np.array(vn.x, dtype=float) 

for i in range(len(self.bounds)): 

xn_a[i] = (xn_a[i] * (self.bounds[i][1] 

- self.bounds[i][0]) 

+ self.bounds[i][0]) 

 

# logging.info('vn.x = {}'.format(vn.x)) 

 

xlines.append(xn_a[0]) 

ylines.append(xn_a[1]) 

xlines.append(x_a[0]) 

ylines.append(x_a[1]) 

 

pyplot.plot(xlines, ylines) 

 

if self.bounds is None: 

pyplot.ylim([-1e-2, 1 + 1e-2]) 

pyplot.xlim([-1e-2, 1 + 1e-2]) 

else: 

pyplot.ylim( 

[self.bounds[1][0] - 1e-2, self.bounds[1][1] + 1e-2]) 

pyplot.xlim( 

[self.bounds[0][0] - 1e-2, self.bounds[0][1] + 1e-2]) 

 

pyplot.show() 

 

elif self.dim == 3: 

fig = pyplot.figure() 

ax = fig.add_subplot(111, projection='3d') 

 

for C in self.H: 

for c in C: 

for v in c(): 

x = [] 

y = [] 

z = [] 

# logging.info('v.x = {}'.format(v.x)) 

x.append(v.x[0]) 

y.append(v.x[1]) 

z.append(v.x[2]) 

for vn in v.nn: 

x.append(vn.x[0]) 

y.append(vn.x[1]) 

z.append(vn.x[2]) 

x.append(v.x[0]) 

y.append(v.x[1]) 

z.append(v.x[2]) 

# logging.info('vn.x = {}'.format(vn.x)) 

 

ax.plot(x, y, z, label='simplex') 

 

pyplot.show() 

else: 

print("dimension higher than 3 or wrong complex format") 

return 

 

 

class VertexGroup(object): 

def __init__(self, p_gen, p_hgr): 

self.p_gen = p_gen # parent generation 

self.p_hgr = p_hgr # parent homology group rank 

self.hg_n = None 

self.hg_d = None 

 

# Maybe add parent homology group rank total history 

# This is the sum off all previously split cells 

# cumulatively throughout its entire history 

self.C = [] 

 

def __call__(self): 

return self.C 

 

def add_vertex(self, V): 

if V not in self.C: 

self.C.append(V) 

 

def homology_group_rank(self): 

""" 

Returns the homology group order of the current cell 

""" 

if self.hg_n is None: 

self.hg_n = sum(1 for v in self.C if v.minimiser()) 

 

return self.hg_n 

 

def homology_group_differential(self): 

""" 

Returns the difference between the current homology group of the 

cell and it's parent group 

""" 

if self.hg_d is None: 

self.hgd = self.hg_n - self.p_hgr 

 

return self.hgd 

 

def polytopial_sperner_lemma(self): 

""" 

Returns the number of stationary points theoretically contained in the 

cell based information currently known about the cell 

""" 

pass 

 

def print_out(self): 

""" 

Print the current cell to console 

""" 

for v in self(): 

v.print_out() 

 

 

class Cell(VertexGroup): 

""" 

Contains a cell that is symmetric to the initial hypercube triangulation 

""" 

 

def __init__(self, p_gen, p_hgr, origin, supremum): 

super(Cell, self).__init__(p_gen, p_hgr) 

 

self.origin = origin 

self.supremum = supremum 

self.centroid = None # (Not always used) 

# TODO: self.bounds 

 

 

class Simplex(VertexGroup): 

""" 

Contains a simplex that is symmetric to the initial symmetry constrained 

hypersimplex triangulation 

""" 

 

def __init__(self, p_gen, p_hgr, generation_cycle, dim): 

super(Simplex, self).__init__(p_gen, p_hgr) 

 

self.generation_cycle = (generation_cycle + 1) % (dim - 1) 

 

 

class Vertex: 

def __init__(self, x, bounds=None, func=None, func_args=(), g_cons=None, 

g_cons_args=(), nn=None, index=None): 

self.x = x 

self.order = sum(x) 

x_a = np.array(x, dtype=float) 

if bounds is not None: 

for i, (lb, ub) in enumerate(bounds): 

x_a[i] = x_a[i] * (ub - lb) + lb 

 

# TODO: Make saving the array structure optional 

self.x_a = x_a 

 

# Note Vertex is only initiated once for all x so only 

# evaluated once 

if func is not None: 

self.feasible = True 

if g_cons is not None: 

for g, args in zip(g_cons, g_cons_args): 

if g(self.x_a, *args) < 0.0: 

self.f = np.inf 

self.feasible = False 

break 

if self.feasible: 

self.f = func(x_a, *func_args) 

 

if nn is not None: 

self.nn = nn 

else: 

self.nn = set() 

 

self.fval = None 

self.check_min = True 

 

# Index: 

if index is not None: 

self.index = index 

 

def __hash__(self): 

return hash(self.x) 

 

def connect(self, v): 

if v is not self and v not in self.nn: 

self.nn.add(v) 

v.nn.add(self) 

 

if self.minimiser(): 

v._min = False 

v.check_min = False 

 

# TEMPORARY 

self.check_min = True 

v.check_min = True 

 

def disconnect(self, v): 

if v in self.nn: 

self.nn.remove(v) 

v.nn.remove(self) 

self.check_min = True 

v.check_min = True 

 

def minimiser(self): 

"""Check whether this vertex is strictly less than all its neighbours""" 

if self.check_min: 

self._min = all(self.f < v.f for v in self.nn) 

self.check_min = False 

 

return self._min 

 

def print_out(self): 

print("Vertex: {}".format(self.x)) 

constr = 'Connections: ' 

for vc in self.nn: 

constr += '{} '.format(vc.x) 

 

print(constr) 

print('Order = {}'.format(self.order)) 

 

 

class VertexCache: 

def __init__(self, func, func_args=(), bounds=None, g_cons=None, 

g_cons_args=(), indexed=True): 

 

self.cache = {} 

self.func = func 

self.g_cons = g_cons 

self.g_cons_args = g_cons_args 

self.func_args = func_args 

self.bounds = bounds 

self.nfev = 0 

self.size = 0 

 

if indexed: 

self.index = -1 

 

def __getitem__(self, x, indexed=True): 

try: 

return self.cache[x] 

except KeyError: 

if indexed: 

self.index += 1 

xval = Vertex(x, bounds=self.bounds, 

func=self.func, func_args=self.func_args, 

g_cons=self.g_cons, 

g_cons_args=self.g_cons_args, 

index=self.index) 

else: 

xval = Vertex(x, bounds=self.bounds, 

func=self.func, func_args=self.func_args, 

g_cons=self.g_cons, 

g_cons_args=self.g_cons_args) 

 

# logging.info("New generated vertex at x = {}".format(x)) 

# NOTE: Surprisingly high performance increase if logging is commented out 

self.cache[x] = xval 

 

# TODO: Check 

if self.func is not None: 

if self.g_cons is not None: 

if xval.feasible: 

self.nfev += 1 

self.size += 1 

else: 

self.size += 1 

else: 

self.nfev += 1 

self.size += 1 

 

return self.cache[x]