-
Notifications
You must be signed in to change notification settings - Fork 24
/
persistent_aposmm.py
1194 lines (913 loc) · 48.6 KB
/
persistent_aposmm.py
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
"""
This module contains methods used our implementation of the Asynchronously
Parallel Optimization Solver for finding Multiple Minima (APOSMM) method
described in detail in the paper
`https://doi.org/10.1007/s12532-017-0131-4 <https://doi.org/10.1007/s12532-017-0131-4>`_
This implementation of APOSMM was developed by Kaushik Kulkarni and Jeffrey
Larson in the summer of 2019.
"""
__all__ = ['initialize_APOSMM', 'decide_where_to_start_localopt', 'update_history_dist']
import numpy as np
from scipy.spatial.distance import cdist
from scipy import optimize as sp_opt
from petsc4py import PETSc
from mpi4py import MPI
from math import log, gamma, pi, sqrt
import nlopt
import dfols
from libensemble.message_numbers import STOP_TAG, PERSIS_STOP, FINISHED_PERSISTENT_GEN_TAG
from libensemble.tools.gen_support import send_mgr_worker_msg
from libensemble.tools.gen_support import get_mgr_worker_msg
from multiprocessing import Event, Process, Queue
class APOSMMException(Exception):
"Raised for any exception in APOSMM"
class ConvergedMsg(object):
"""
Message communicated when a local optimization is converged.
"""
def __init__(self, x, opt_flag):
self.x = x
self.opt_flag = opt_flag
class ErrorMsg(object):
"""
Message communicated when a local optimization has an exception.
"""
def __init__(self, x):
self.x = x
def aposmm(H, persis_info, gen_specs, libE_info):
"""
APOSMM coordinates multiple local optimization runs, starting from points
which do not have a better point nearby (within a distance ``r_k``). This
generation function produces/requires the following fields in ``H``:
- ``'x' [n floats]``: Parameters being optimized over
- ``'x_on_cube' [n floats]``: Parameters scaled to the unit cube
- ``'f' [float]``: Objective function being minimized
- ``'local_pt' [bool]``: True if point from a local optimization run
- ``'dist_to_unit_bounds' [float]``: Distance to domain boundary
- ``'dist_to_better_l' [float]``: Dist to closest better local opt point
- ``'dist_to_better_s' [float]``: Dist to closest better sample point
- ``'ind_of_better_l' [int]``: Index of point ``'dist_to_better_l``' away
- ``'ind_of_better_s' [int]``: Index of point ``'dist_to_better_s``' away
- ``'started_run' [bool]``: True if point has started a local opt run
- ``'num_active_runs' [int]``: Number of active local runs point is in
- ``'local_min' [float]``: True if point has been ruled a local minima
- ``'sim_id' [int]``: Row number of entry in history
and optionally
- ``'fvec' [m floats]``: All objective components (if calculated together)
- ``'obj_component' [int]``: Index corresponding to value in ``'f_i``'
- ``'pt_id' [int]``: Identify the point (useful when evaluating different
objective components for a given ``'x'``)
When using libEnsemble to do individual objective component evaluations,
APOSMM will return ``gen_specs['user']['components']`` copies of each point, but
the component=0 entry of each point will only be considered when
- deciding where to start a run,
- best nearby point,
- storing the order of the points is the run
- storing the combined objective function value
- etc
Necessary quantities in ``gen_specs['user']`` are:
- ``'lb' [n floats]``: Lower bound on search domain
- ``'ub' [n floats]``: Upper bound on search domain
- ``'initial_sample_size' [int]``: Number of uniformly sampled points
must be returned (non-nan value) before a local opt run is started
- ``'localopt_method' [str]``: Name of an NLopt, PETSc/TAO, or SciPy method
(see 'advance_local_run' below for supported methods)
Optional ``gen_specs['user']`` entries are:
- ``'sample_points' [numpy array]``: Points to be sampled (original domain)
- ``'combine_component_func' [func]``: Function to combine obj components
- ``'components' [int]``: Number of objective components
- ``'dist_to_bound_multiple' [float in (0,1]]``: What fraction of the
distance to the nearest boundary should the initial step size be in
localopt runs
- ``'lhs_divisions' [int]``: Number of Latin hypercube sampling partitions
(0 or 1 results in uniform sampling)
- ``'mu' [float]``: Distance from the boundary that all localopt starting
points must satisfy
- ``'nu' [float]``: Distance from identified minima that all starting
points must satisfy
- ``'rk_const' [float]``: Multiplier in front of the r_k value
- ``'max_active_runs' [int]``: Bound on number of runs APOSMM is advancing
And ``gen_specs['user']`` convergence tolerances for NLopt, PETSc/TAO, SciPy
- ``'fatol' [float]``:
- ``'ftol_abs' [float]``:
- ``'ftol_rel' [float]``:
- ``'gatol' [float]``:
- ``'grtol' [float]``:
- ``'xtol_abs' [float]``:
- ``'xtol_rel' [float]``:
- ``'tol' [float]``:
As a default, APOSMM starts a local optimization runs from a point that:
- is not in an active local optimization run,
- is more than ``mu`` from the boundary (in the unit-cube domain),
- is more than ``nu`` from identified minima (in the unit-cube domain),
- does not have a better point within a distance ``r_k`` of it.
If the above results in more than ``'max_active_runs'`` being advanced, the
best point in each run is determined and the dist_to_better is computed
(with inf being the value for the best run). Then those
``'max_active_runs'`` runs with largest dist_to_better are advanced
(breaking ties arbitrarily).
:Note:
``gen_specs['user']['combine_component_func']`` must be defined when there are
multiple objective components.
:Note:
APOSMM critically uses ``persis_info`` to store information about
active runs, order of points in each run, etc. The allocation function
must ensure it's always given.
.. seealso::
`test_branin_aposmm_nlopt_and_then_scipy.py <https://github.com/Libensemble/libensemble/blob/develop/libensemble/tests/regression_tests/test_branin_aposmm_nlopt_and_then_scipy.py>`_
for basic APOSMM usage.
.. seealso::
`test_chwirut_aposmm_one_residual_at_a_time.py <https://github.com/Libensemble/libensemble/blob/develop/libensemble/tests/regression_tests/test_chwirut_aposmm_one_residual_at_a_time.py>`_
for an example of APOSMM coordinating multiple local optimization runs
for an objective with more than one component.
"""
"""
Description of intermediate variables in aposmm_logic:
n: domain dimension
n_s: the number of complete evaluations of sampled points
updated_inds: indices of H that have been updated (and so all their
information must be sent back to libE manager to update)
O: new points to be sent back to the history
When re-running a local opt method to get the next point:
advance_local_run.x_new: stores the first new point requested by
a local optimization method
advance_local_run.pt_in_run: counts function evaluations to know
when a new point is given
starting_inds: indices where a runs should be started.
active_runs: indices of active local optimization runs
sorted_run_inds: indices of the considered run (in the order they were
requested by the localopt method)
x_opt: the reported minimum from a localopt run (disregarded
unless exit_code isn't 0)
exit_code: 0 if a new localopt point has been found, otherwise it's
the NLopt/TAO/SciPy code
samples_needed: Number of additional uniformly drawn samples needed
Description of persistent variables used to maintain the state of APOSMM
persis_info['total_runs']: Running count of started/completed localopt runs
persis_info['run_order']: Sequence of indices of points in unfinished runs
persis_info['old_runs']: Sequence of indices of points in finished runs
"""
user_specs = gen_specs['user']
n, n_s, rk_const, ld, mu, nu, comm, local_H = initialize_APOSMM(H, user_specs, libE_info)
local_opters, sim_id_to_child_inds, run_order, run_pts, total_runs, fields_to_pass = initialize_children(user_specs)
if user_specs['initial_sample_size'] != 0:
# Send our initial sample. We don't need to check that n_s is large enough:
# the alloc_func only returns when the initial sample has function values.
persis_info = add_k_sample_points_to_local_H(user_specs['initial_sample_size'], user_specs,
persis_info, n, comm, local_H,
sim_id_to_child_inds)
send_mgr_worker_msg(comm, local_H[-user_specs['initial_sample_size']:][[i[0] for i in gen_specs['out']]])
something_sent = True
else:
something_sent = False
tag = None
first_pass = True
while 1:
new_opt_inds_to_send_mgr = []
new_inds_to_send_mgr = []
if something_sent:
tag, Work, calc_in = get_mgr_worker_msg(comm)
if tag in [STOP_TAG, PERSIS_STOP]:
clean_up_and_stop(local_H, local_opters, run_order)
persis_info['run_order'] = run_order
break
n_s, n_r = update_local_H_after_receiving(local_H, n, n_s, user_specs, Work, calc_in, fields_to_pass)
for row in calc_in:
if sim_id_to_child_inds.get(row['sim_id']):
# Point came from a child local opt run
for child_idx in sim_id_to_child_inds[row['sim_id']]:
x_new = local_opters[child_idx].iterate(row[fields_to_pass])
if isinstance(x_new, ConvergedMsg):
x_opt = x_new.x
opt_flag = x_new.opt_flag
opt_ind = update_history_optimal(x_opt, opt_flag, local_H, run_order[child_idx])
new_opt_inds_to_send_mgr.append(opt_ind)
local_opters.pop(child_idx)
else:
add_to_local_H(local_H, x_new, user_specs, local_flag=1, on_cube=True)
new_inds_to_send_mgr.append(len(local_H)-1)
run_order[child_idx].append(local_H[-1]['sim_id'])
run_pts[child_idx].append(x_new)
if local_H[-1]['sim_id'] in sim_id_to_child_inds:
sim_id_to_child_inds[local_H[-1]['sim_id']] += (child_idx, )
else:
sim_id_to_child_inds[local_H[-1]['sim_id']] = (child_idx, )
starting_inds = decide_where_to_start_localopt(local_H, n, n_s, rk_const, ld, mu, nu)
for ind in starting_inds:
if len([p for p in local_opters.values() if p.is_running]) < user_specs.get('max_active_runs', np.inf):
local_H['started_run'][ind] = 1
# Initialize a local opt run
local_opter = LocalOptInterfacer(user_specs, local_H[ind]['x_on_cube'],
local_H[ind]['f'] if 'f' in fields_to_pass else local_H[ind]['fvec'],
local_H[ind]['grad'] if 'grad' in fields_to_pass else None)
local_opters[total_runs] = local_opter
x_new = local_opter.iterate(local_H[ind][fields_to_pass]) # Assuming the second point can't be ruled optimal
add_to_local_H(local_H, x_new, user_specs, local_flag=1, on_cube=True)
new_inds_to_send_mgr.append(len(local_H)-1)
run_order[total_runs] = [ind, local_H[-1]['sim_id']]
run_pts[total_runs] = [local_H['x_on_cube'], x_new]
if local_H[-1]['sim_id'] in sim_id_to_child_inds:
sim_id_to_child_inds[local_H[-1]['sim_id']] += (total_runs, )
else:
sim_id_to_child_inds[local_H[-1]['sim_id']] = (total_runs, )
total_runs += 1
if first_pass:
num_samples_needed = persis_info['nworkers'] - 1 - len(new_inds_to_send_mgr)
first_pass = False
else:
num_samples_needed = n_r-len(new_inds_to_send_mgr)
if num_samples_needed > 0:
persis_info = add_k_sample_points_to_local_H(num_samples_needed, user_specs, persis_info, n, comm, local_H, sim_id_to_child_inds)
new_inds_to_send_mgr = new_inds_to_send_mgr + list(range(len(local_H)-num_samples_needed, len(local_H)))
send_mgr_worker_msg(comm, local_H[new_inds_to_send_mgr + new_opt_inds_to_send_mgr][[i[0] for i in gen_specs['out']]])
something_sent = True
return [], persis_info, FINISHED_PERSISTENT_GEN_TAG
class LocalOptInterfacer(object):
def __init__(self, user_specs, x0, f0, grad0=None):
"""
:param x0: A numpy array of the initial guess solution. This guess
should be scaled to a unit cube.
:param f0: A numpy array of the initial function value.
.. warning:: In order to have correct functioning of the local
optimization child processes. ~self.iterate~ should be called
immediately after creating the class.
"""
self.parent_can_read = Event()
self.comm_queue = Queue()
self.child_can_read = Event()
self.x0 = x0.copy()
self.f0 = f0.copy()
if grad0 is not None:
self.grad0 = grad0.copy()
else:
self.grad0 = None
# Setting the local optimization method
if user_specs['localopt_method'] in ['LN_SBPLX', 'LN_BOBYQA', 'LN_COBYLA', 'LN_NEWUOA', 'LN_NELDERMEAD', 'LD_MMA']:
run_local_opt = run_local_nlopt
elif user_specs['localopt_method'] in ['pounders', 'blmvm', 'nm']:
run_local_opt = run_local_tao
elif user_specs['localopt_method'] in ['scipy_Nelder-Mead', 'scipy_COBYLA']:
run_local_opt = run_local_scipy_opt
elif user_specs['localopt_method'] in ['dfols']:
run_local_opt = run_local_dfols
elif user_specs['localopt_method'] in ['external_localopt']:
run_local_opt = run_external_localopt
self.parent_can_read.clear()
self.process = Process(target=opt_runner, args=(run_local_opt, user_specs,
self.comm_queue, x0, f0, self.child_can_read,
self.parent_can_read))
self.process.start()
self.is_running = True
self.parent_can_read.wait()
x_new = self.comm_queue.get()
if isinstance(x_new, ErrorMsg):
raise APOSMMException(x_new.x)
assert np.allclose(x_new, x0, rtol=1e-15, atol=1e-15), "The first point requested by this run does not match the starting point. Exiting"
def iterate(self, data):
"""
Returns an instance of either :class:`numpy.ndarray` corresponding to the next
iterative guess or :class:`ConvergedMsg` when the solver is converged.
:param f: A numpy array of the function evaluation.
:param grad: A numpy array of the function's gradient.
"""
self.parent_can_read.clear()
if 'grad' in data.dtype.names:
self.comm_queue.put((data['x_on_cube'], data['f'], data['grad']))
elif 'fvec' in data.dtype.names:
self.comm_queue.put((data['x_on_cube'], data['fvec']))
else:
self.comm_queue.put((data['x_on_cube'], data['f'], ))
self.child_can_read.set()
self.parent_can_read.wait()
x_new = self.comm_queue.get()
if isinstance(x_new, ErrorMsg):
raise APOSMMException(x_new.x)
elif isinstance(x_new, ConvergedMsg):
self.process.join()
self.comm_queue.close()
self.comm_queue.join_thread()
self.is_running = False
else:
x_new = np.atleast_2d(x_new)
return x_new
def destroy(self, previous_x):
while not isinstance(previous_x, ConvergedMsg):
self.parent_can_read.clear()
if self.grad0 is None:
self.comm_queue.put((previous_x, 0*np.ones_like(self.f0),))
else:
self.comm_queue.put((previous_x, 0*np.ones_like(self.f0), np.zeros_like(self.grad0)))
self.child_can_read.set()
self.parent_can_read.wait()
previous_x = self.comm_queue.get()
assert isinstance(previous_x, ConvergedMsg)
self.process.join()
self.comm_queue.close()
self.comm_queue.join_thread()
self.is_running = False
def run_local_nlopt(user_specs, comm_queue, x0, f0, child_can_read, parent_can_read):
# print('[Child]: Started local opt at {}.'.format(x0), flush=True)
n = len(user_specs['ub'])
opt = nlopt.opt(getattr(nlopt, user_specs['localopt_method']), n)
lb = np.zeros(n)
ub = np.ones(n)
if not user_specs.get('periodic'):
opt.set_lower_bounds(lb)
opt.set_upper_bounds(ub)
# Care must be taken here because a too-large initial step causes nlopt to move the starting point!
dist_to_bound = min(min(ub-x0), min(x0-lb))
assert dist_to_bound > np.finfo(np.float32).eps, "The distance to the boundary is too small for NLopt to handle"
if 'dist_to_bound_multiple' in user_specs:
opt.set_initial_step(dist_to_bound*user_specs['dist_to_bound_multiple'])
else:
opt.set_initial_step(dist_to_bound)
run_max_eval = user_specs.get('run_max_eval', 1000*n)
opt.set_maxeval(run_max_eval)
opt.set_min_objective(lambda x, grad: nlopt_callback_fun(x, grad,
comm_queue, child_can_read, parent_can_read,
user_specs))
if 'xtol_rel' in user_specs:
opt.set_xtol_rel(user_specs['xtol_rel'])
if 'ftol_rel' in user_specs:
opt.set_ftol_rel(user_specs['ftol_rel'])
if 'xtol_abs' in user_specs:
opt.set_xtol_abs(user_specs['xtol_abs'])
if 'ftol_abs' in user_specs:
opt.set_ftol_abs(user_specs['ftol_abs'])
# FIXME: Do we need to do something of the final 'x_opt'?
# print('[Child]: Started my optimization', flush=True)
x_opt = opt.optimize(x0)
return_val = opt.last_optimize_result()
if return_val >= 1 and return_val <= 4:
# These return values correspond to an optimium being identified
# https://nlopt.readthedocs.io/en/latest/NLopt_Reference/#return-values
opt_flag = 1
elif return_val >= 5:
print("The run started from " + str(x0) + " reached it maximum number "
"of function evaluations: " + str(run_max_eval) + ". No point from "
"this run will be ruled as a minimum! APOSMM may start a new run "
"from some point in this run.")
opt_flag = 0
else:
print("NLopt returned with a negative return value, which indicates an error")
opt_flag = 0
if user_specs.get('periodic'):
x_opt = x_opt % 1 # Shift x_opt to be in the correct location in the unit cube (not the domain user_specs['lb'] - user_specs['ub'])
finish_queue(x_opt, opt_flag, comm_queue, parent_can_read, user_specs)
def run_local_scipy_opt(user_specs, comm_queue, x0, f0, child_can_read, parent_can_read):
# Construct the bounds in the form of constraints
cons = []
for factor in range(len(x0)):
lo = {'type': 'ineq',
'fun': lambda x, lb=user_specs['lb'][factor], i=factor: x[i]-lb}
up = {'type': 'ineq',
'fun': lambda x, ub=user_specs['ub'][factor], i=factor: ub-x[i]}
cons.append(lo)
cons.append(up)
method = user_specs['localopt_method'][6:]
# print('[Child]: Started my optimization', flush=True)
res = sp_opt.minimize(lambda x: scipy_dfols_callback_fun(x, comm_queue,
child_can_read, parent_can_read, user_specs), x0,
# constraints=cons,
method=method, **user_specs.get('scipy_kwargs', {}))
if res['status'] == 0:
opt_flag = 1
else:
print("The SciPy localopt run started from " + str(x0) + " stopped "
" without finding a local min. The status of the run is " + str(res['status']) +
". No point from this run will be ruled as a minimum! APOSMM may "
"start a new run from some point in this run.")
opt_flag = 0
if user_specs.get('periodic'):
x_opt = res['x'] % 1
else:
x_opt = res['x']
finish_queue(x_opt, opt_flag, comm_queue, parent_can_read, user_specs)
def run_external_localopt(user_specs, comm_queue, x0, f0, child_can_read, parent_can_read):
import subprocess
import os
from uuid import uuid4
run_id = uuid4().hex
x_file = 'x_' + run_id + '.txt'
y_file = 'y_' + run_id + '.txt'
x_done_file = 'x_done_' + run_id + '.txt'
y_done_file = 'y_done_' + run_id + '.txt'
opt_file = 'opt_' + run_id + '.txt'
# cmd = ["matlab", "-nodisplay", "-nodesktop", "-nojvm", "-nosplash", "-r",
cmd = ["octave", "--no-window-system", "--eval",
"x0=[" + " ".join(["{:18.18f}".format(x) for x in x0]) + "];"
"opt_file='" + opt_file + "';"
"x_file='" + x_file + "';"
"y_file='" + y_file + "';"
"x_done_file='" + x_done_file + "';"
"y_done_file='" + y_done_file + "';"
"call_matlab_octave_script"]
p = subprocess.Popen(cmd, shell=False, stdout=subprocess.DEVNULL)
while p.poll() is None: # Process still going
if os.path.isfile(x_done_file): # x file exists
x = np.loadtxt(x_file)
os.remove(x_done_file)
x_recv, f_recv = put_set_wait_get(x, comm_queue, parent_can_read, child_can_read, user_specs)
np.savetxt(y_file, [f_recv])
open(y_done_file, 'w').close()
x_opt = np.loadtxt(opt_file)
opt_flag = np.loadtxt(opt_file + "_flag")
for f in [x_file, y_file, opt_file]:
os.remove(f)
finish_queue(x_opt, opt_flag, comm_queue, parent_can_read, user_specs)
def run_local_dfols(user_specs, comm_queue, x0, f0, child_can_read, parent_can_read):
# Define bound constraints (lower <= x <= upper)
lb = np.zeros(len(x0))
ub = np.ones(len(x0))
# Set random seed (for reproducibility)
np.random.seed(0)
# Care must be taken here because a too-large initial step causes DFO-LS to move the starting point!
dist_to_bound = min(min(ub-x0), min(x0-lb))
assert dist_to_bound > np.finfo(np.float32).eps, "The distance to the boundary is too small"
assert 'bounds' not in user_specs.get('dfols_kwargs', {}), "APOSMM must set the bounds for DFO-LS"
assert 'rhobeg' not in user_specs.get('dfols_kwargs', {}), "APOSMM must set rhobeg for DFO-LS"
assert 'x0' not in user_specs.get('dfols_kwargs', {}), "APOSMM must set x0 for DFO-LS"
# Call DFO-LS
soln = dfols.solve(lambda x: scipy_dfols_callback_fun(x, comm_queue, child_can_read, parent_can_read, user_specs),
x0, bounds=(lb, ub), rhobeg=0.5*dist_to_bound, **user_specs.get('dfols_kwargs', {}))
x_opt = soln.x
if soln.flag == soln.EXIT_SUCCESS:
opt_flag = 1
else:
print("The DFO-LS run started from " + str(x0) + " stopped with an exit "
"flag of " + str(soln.flag) + ". No point from this run will be "
"ruled as a minimum! APOSMM may start a new run from some point "
"in this run.")
opt_flag = 0
finish_queue(x_opt, opt_flag, comm_queue, parent_can_read, user_specs)
def run_local_tao(user_specs, comm_queue, x0, f0, child_can_read, parent_can_read):
assert isinstance(x0, np.ndarray)
tao_comm = MPI.COMM_SELF
n, = x0.shape
if f0.shape == ():
m = 1
else:
m, = f0.shape
# Create starting point, bounds, and tao object
x = PETSc.Vec().create(tao_comm)
x.setSizes(n)
x.setFromOptions()
x.array = x0
lb = x.duplicate()
ub = x.duplicate()
lb.array = 0*np.ones(n)
ub.array = 1*np.ones(n)
tao = PETSc.TAO().create(tao_comm)
tao.setType(user_specs['localopt_method'])
if user_specs['localopt_method'] == 'pounders':
f = PETSc.Vec().create(tao_comm)
f.setSizes(m)
f.setFromOptions()
if hasattr(tao, 'setResidual'):
tao.setResidual(lambda tao, x, f: tao_callback_fun_pounders(tao, x, f, comm_queue, child_can_read, parent_can_read, user_specs), f)
else:
tao.setSeparableObjective(lambda tao, x, f: tao_callback_fun_pounders(tao, x, f, comm_queue, child_can_read, parent_can_read, user_specs), f)
delta_0 = user_specs['dist_to_bound_multiple']*np.min([np.min(ub.array-x.array), np.min(x.array-lb.array)])
PETSc.Options().setValue('-tao_pounders_delta', str(delta_0))
elif user_specs['localopt_method'] == 'nm':
tao.setObjective(lambda tao, x: tao_callback_fun_nm(tao, x, comm_queue, child_can_read, parent_can_read, user_specs))
elif user_specs['localopt_method'] == 'blmvm':
g = PETSc.Vec().create(tao_comm)
g.setSizes(n)
g.setFromOptions()
tao.setObjectiveGradient(lambda tao, x, g: tao_callback_fun_grad(tao, x, g, comm_queue, child_can_read, parent_can_read, user_specs))
# Set everything for tao before solving
# FIXME: Hard-coding 100 as the max funcs as couldn't find any other
# sensible value.
PETSc.Options().setValue('-tao_max_funcs', str(user_specs.get('run_max_eval', 1000*n)))
tao.setFromOptions()
tao.setVariableBounds((lb, ub))
tao.setTolerances(grtol=user_specs.get('grtol', 1e-8), gatol=user_specs.get('gatol', 1e-8))
tao.setInitial(x)
# print('[Child]: Started my optimization', flush=True)
tao.solve(x)
x_opt = tao.getSolution().getArray()
exit_code = tao.getConvergedReason()
if exit_code > 0:
opt_flag = 1
else:
# https://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Tao/TaoGetConvergedReason.html
print("The run started from " + str(x0) + " exited with a nonpositive reason. No point from "
"this run will be ruled as a minimum! APOSMM may start a new run from some point in this run.")
opt_flag = 0
if user_specs['localopt_method'] == 'pounders':
f.destroy()
elif user_specs['localopt_method'] == 'blmvm':
g.destroy()
lb.destroy()
ub.destroy()
x.destroy()
tao.destroy()
finish_queue(x_opt, opt_flag, comm_queue, parent_can_read, user_specs)
def opt_runner(run_local_opt, user_specs, comm_queue, x0, f0, child_can_read, parent_can_read):
try:
run_local_opt(user_specs, comm_queue, x0, f0, child_can_read, parent_can_read)
except Exception as e:
comm_queue.put(ErrorMsg(e))
parent_can_read.set()
# Callback functions and routines
def nlopt_callback_fun(x, grad, comm_queue, child_can_read, parent_can_read, user_specs):
if user_specs['localopt_method'] in ['LD_MMA']:
x_recv, f_recv, grad_recv = put_set_wait_get(x, comm_queue, parent_can_read, child_can_read, user_specs)
grad[:] = grad_recv
else:
assert user_specs['localopt_method'] in ['LN_SBPLX', 'LN_BOBYQA', 'LN_NEWUOA',
'LN_COBYLA', 'LN_NELDERMEAD', 'LD_MMA']
x_recv, f_recv = put_set_wait_get(x, comm_queue, parent_can_read, child_can_read, user_specs)
return f_recv
def scipy_dfols_callback_fun(x, comm_queue, child_can_read, parent_can_read, user_specs):
x_recv, f_x_recv, = put_set_wait_get(x, comm_queue, parent_can_read, child_can_read, user_specs)
return f_x_recv
def tao_callback_fun_nm(tao, x, comm_queue, child_can_read, parent_can_read, user_specs):
x_recv, f_recv, = put_set_wait_get(x.array_r, comm_queue, parent_can_read, child_can_read, user_specs)
return f_recv
def tao_callback_fun_pounders(tao, x, f, comm_queue, child_can_read, parent_can_read, user_specs):
x_recv, f_recv, = put_set_wait_get(x.array_r, comm_queue, parent_can_read, child_can_read, user_specs)
f.array[:] = f_recv
return f
def tao_callback_fun_grad(tao, x, g, comm_queue, child_can_read, parent_can_read, user_specs):
x_recv, f_recv, grad_recv = put_set_wait_get(x.array_r, comm_queue, parent_can_read, child_can_read, user_specs)
g.array[:] = grad_recv
return f_recv
def finish_queue(x_opt, opt_flag, comm_queue, parent_can_read, user_specs):
if user_specs.get('print') and opt_flag:
print('Local optimum on the [0,1]^n domain', x_opt, flush=True)
comm_queue.put(ConvergedMsg(x_opt, opt_flag))
parent_can_read.set()
def put_set_wait_get(x, comm_queue, parent_can_read, child_can_read, user_specs):
"""This routine is used by children process callback functions. It:
- puts x into a comm_queue,
- tells the parent it can read,
- tells the child to wait
- receives the values put in the comm_queue by the parent
- checks that the first value received matches x
- removes the wait on the child
- returns values"""
comm_queue.put(x)
# print('[Child]: I just put x_on_cube =', x.array, flush=True)
# print('[Child]: Parent should no longer wait.', flush=True)
parent_can_read.set()
# print('[Child]: I have started waiting', flush=True)
child_can_read.wait()
# print('[Child]: Wohooo.. I am free folks', flush=True)
values = comm_queue.get()
child_can_read.clear()
if user_specs.get('periodic'):
assert np.allclose(x % 1, values[0] % 1, rtol=1e-15, atol=1e-15), "The point I gave is not the point I got back!"
else:
assert np.allclose(x, values[0], rtol=1e-15, atol=1e-15), "The point I gave is not the point I got back!"
return values
def update_local_H_after_receiving(local_H, n, n_s, user_specs, Work, calc_in, fields_to_pass):
for name in ['f', 'x_on_cube', 'grad', 'fvec']:
if name in fields_to_pass:
assert name in calc_in.dtype.names, name + " must be returned to persistent_aposmm for localopt_method" + user_specs['localopt_method']
for name in calc_in.dtype.names:
local_H[name][Work['libE_info']['H_rows']] = calc_in[name]
local_H['returned'][Work['libE_info']['H_rows']] = True
n_s += np.sum(~local_H[Work['libE_info']['H_rows']]['local_pt'])
n_r = len(Work['libE_info']['H_rows'])
# dist -> distance
update_history_dist(local_H, n)
return n_s, n_r
def add_to_local_H(local_H, pts, user_specs, local_flag=0, sorted_run_inds=[], run=[], on_cube=True):
"""
Adds points to O, the numpy structured array to be sent back to the manager
"""
assert not local_flag or len(pts) == 1, "Can't > 1 local points"
len_local_H = len(local_H)
ub = user_specs['ub']
lb = user_specs['lb']
num_pts = len(pts)
local_H.resize(len(local_H)+num_pts, refcheck=False) # Adds num_pts rows of zeros to O
if on_cube:
local_H['x_on_cube'][-num_pts:] = pts
local_H['x'][-num_pts:] = pts*(ub-lb)+lb
else:
local_H['x_on_cube'][-num_pts:] = (pts-lb)/(ub-lb)
local_H['x'][-num_pts:] = pts
if user_specs.get('periodic'):
local_H['x_on_cube'][-num_pts:] = local_H['x_on_cube'][-num_pts:] % 1
local_H['sim_id'][-num_pts:] = np.arange(len_local_H, len_local_H+num_pts)
local_H['local_pt'][-num_pts:] = local_flag
initialize_dists_and_inds(local_H, num_pts)
if local_flag:
local_H['num_active_runs'][-num_pts] += 1
def initialize_dists_and_inds(local_H, num_pts):
local_H['dist_to_unit_bounds'][-num_pts:] = np.inf
local_H['dist_to_better_l'][-num_pts:] = np.inf
local_H['dist_to_better_s'][-num_pts:] = np.inf
local_H['ind_of_better_l'][-num_pts:] = -1
local_H['ind_of_better_s'][-num_pts:] = -1
def update_history_dist(H, n):
"""
Updates distances/indices after new points that have been evaluated.
.. seealso::
`start_persistent_local_opt_gens.py <https://github.com/Libensemble/libensemble/blob/develop/libensemble/alloc_funcs/start_persistent_local_opt_gens.py>`_
"""
new_inds = np.where(~H['known_to_aposmm'])[0]
p = np.logical_and.reduce((H['returned'], ~np.isnan(H['f'])))
for new_ind in new_inds:
# Loop over new returned points and update their distances
if p[new_ind]:
H['known_to_aposmm'][new_ind] = True
# Compute distance to boundary
H['dist_to_unit_bounds'][new_ind] = min(min(np.ones(n)-H['x_on_cube'][new_ind]), min(H['x_on_cube'][new_ind]-np.zeros(n)))
dist_to_all = cdist(H['x_on_cube'][[new_ind]], H['x_on_cube'][p], 'euclidean').flatten()
new_better_than = H['f'][new_ind] < H['f'][p]
# Update any other points if new_ind is closer and better
if H['local_pt'][new_ind]:
inds_of_p = np.logical_and(dist_to_all < H['dist_to_better_l'][p], new_better_than)
updates = np.where(p)[0][inds_of_p]
H['dist_to_better_l'][updates] = dist_to_all[inds_of_p]
H['ind_of_better_l'][updates] = new_ind
else:
inds_of_p = np.logical_and(dist_to_all < H['dist_to_better_s'][p], new_better_than)
updates = np.where(p)[0][inds_of_p]
H['dist_to_better_s'][updates] = dist_to_all[inds_of_p]
H['ind_of_better_s'][updates] = new_ind
# Since we allow equality when deciding better_than_new_l and
# better_than_new_s, we have to prevent new_ind from being its own
# better point.
better_than_new_l = np.logical_and.reduce((~new_better_than, H['local_pt'][p], H['sim_id'][p] != new_ind))
better_than_new_s = np.logical_and.reduce((~new_better_than, ~H['local_pt'][p], H['sim_id'][p] != new_ind))
# Who is closest to ind and better
if np.any(better_than_new_l):
ind = dist_to_all[better_than_new_l].argmin()
H['ind_of_better_l'][new_ind] = H['sim_id'][p][np.nonzero(better_than_new_l)[0][ind]]
H['dist_to_better_l'][new_ind] = dist_to_all[better_than_new_l][ind]
if np.any(better_than_new_s):
ind = dist_to_all[better_than_new_s].argmin()
H['ind_of_better_s'][new_ind] = H['sim_id'][p][np.nonzero(better_than_new_s)[0][ind]]
H['dist_to_better_s'][new_ind] = dist_to_all[better_than_new_s][ind]
# if not ignore_L8:
# r_k = calc_rk(len(H['x_on_cube'][0]), n_s, rk_const, lhs_divisions)
# H['worse_within_rk'][new_ind][p] = np.logical_and.reduce((H['f'][new_ind] <= H['f'][p], dist_to_all <= r_k))
# # Add trues if new point is 'worse_within_rk'
# inds_to_change = np.logical_and.reduce((H['dist_to_all'][p,new_ind] <= r_k, H['f'][new_ind] >= H['f'][p], H['sim_id'][p] != new_ind))
# H['worse_within_rk'][inds_to_change,new_ind] = True
# if not H['local_pt'][new_ind]:
# H['worse_within_rk'][H['dist_to_all'] > r_k] = False
if np.any(~H['local_pt']) and not np.any(np.isinf(H['dist_to_better_s'][~H['local_pt']])):
# Our best sample point was not identified because the min was not unique.
min_inds = H['f'][~H['local_pt']] == np.min(H['f'][~H['local_pt']])
assert len(min_inds) >= 2, "Check this"
# Take the first point with this value to be the best sample point
best_samp = H['sim_id'][~H['local_pt']][min_inds][0]
H['dist_to_better_s'][best_samp] = np.inf
H['ind_of_better_s'][best_samp] = -1
# if np.any(H['local_pt']) and not np.any(np.isinf(H['dist_to_better_l'][H['local_pt']])):
# # Our best sample point was not identified because the min was not unique.
# min_inds = H['f'][H['local_pt']] == np.min(H['f'][H['local_pt']])
# assert len(min_inds) >= 2, "Check this"
# # Take the first point with this value to be the best sample point
# best_local = H['sim_id'][H['local_pt']][min_inds][0]
# H['dist_to_better_l'][best_local] = np.inf
# H['ind_of_better_l'][best_local] = -1
def update_history_optimal(x_opt, opt_flag, H, run_inds):
"""
Updated the history after any point has been declared a local minimum
"""
# opt_ind = np.where(np.logical_and(np.equal(x_opt,H['x_on_cube']).all(1),~np.isinf(H['f'])))[0] # This fails on some problems. x_opt is 1e-16 away from the point that was given and opt_ind is empty
run_inds = np.unique(run_inds)
dists = np.linalg.norm(H['x_on_cube'][run_inds]-x_opt, axis=1)
ind = np.argmin(dists)
opt_ind = run_inds[ind]
tol_x1 = 1e-15
# Instead of failing, we accept x_opt that is slightly different from its value in H
# assert dists[ind] <= tol_x1, "Closest point to x_opt not within {}?".format(tol_x1)
if dists[ind] > tol_x1:
print("Dist from reported x_opt to closest evaluated point is: " + str(dists[ind]) + "\n" +
"Check that the local optimizer is working correctly\n", x_opt, run_inds, flush=True)
tol_x2 = 1e-8
failsafe = np.logical_and(H['f'][run_inds] < H['f'][opt_ind], dists < tol_x2)
if opt_flag:
if np.any(failsafe):
print("This run has {} point(s) with smaller 'f' value within {} of "
"the point ruled to be the run minimum. \nMarking all as being "
"a 'local_min' to prevent APOSMM from starting another run "
"immediately from these points.".format(sum(failsafe), tol_x2))
print("Sim_ids to be marked optimal: ", opt_ind, run_inds[failsafe])
print("Check that the local optimizer is working correctly", flush=True)
H['local_min'][run_inds[failsafe]] = 1
H['local_min'][opt_ind] = 1
H['num_active_runs'][run_inds] -= 1
return opt_ind
def decide_where_to_start_localopt(H, n, n_s, rk_const, ld=0, mu=0, nu=0):
"""
Finds points in the history that satisfy the conditions (S1-S5 and L1-L8) in
Table 1 of the `APOSMM paper <https://doi.org/10.1007/s12532-017-0131-4>`_
This method first identifies sample points satisfying S2-S5, and then
identifies all localopt points that satisfy L1-L7.
We then start from any sample point also satisfying S1.
We do not check condition L8 currently.
We don't consider points in the history that have not returned from
computation, or that have a ``nan`` value. Also, note that ``mu`` and ``nu``
implicitly depend on the scaling that is happening with the domain. That
is, adjusting the initial domain can make a run start (or not start) at
a point that didn't (or did) previously.
Parameters
----------
H: numpy structured array
History array storing rows for each point.
r_k_const: float
Radius for deciding when to start runs
lhs_divisions: integer
Number of Latin hypercube sampling divisions (0 or 1 means uniform
random sampling over the domain)
mu: nonnegative float
Distance from the boundary that all starting points must satisfy
nu: nonnegative float
Distance from identified minima that all starting points must satisfy
gamma_quantile: float in (0,1]
Only sample points whose function values are in the lower
gamma_quantile can start localopt runs
Returns
----------
start_inds: list
Indices where a local opt run should be started
.. seealso::
`start_persistent_local_opt_gens.py <https://github.com/Libensemble/libensemble/blob/develop/libensemble/alloc_funcs/start_persistent_local_opt_gens.py>`_
"""
r_k = calc_rk(n, n_s, rk_const, ld)
if nu > 0:
test_2_through_5 = np.logical_and.reduce((
H['returned'] == 1, # have a returned function value
H['dist_to_better_s'] >
r_k, # no better sample point within r_k (L2)
~H['started_run'], # have not started a run (L3)
H['dist_to_unit_bounds'] >=
mu, # have all components at least mu away from bounds (L4)
np.all(
cdist(H['x_on_cube'], H['x_on_cube'][H['local_min']]) >= nu,
axis=1) # distance nu away from known local mins (L5)
))
else:
test_2_through_5 = np.logical_and.reduce((
H['returned'] == 1, # have a returned function value
H['dist_to_better_s'] >
r_k, # no better sample point within r_k (L2)
~H['started_run'], # have not started a run (L3)
H['dist_to_unit_bounds'] >=
mu, # have all components at least mu away from bounds (L4)
)) # (L5) is always true when nu = 0
# assert gamma_quantile == 1, "This is not supported yet. What is the best way to decide this when there are NaNs present in H['f']?"
# if gamma_quantile < 1:
# cut_off_value = np.sort(H['f'][~H['local_pt']])[np.floor(gamma_quantile*(sum(~H['local_pt'])-1)).astype(int)]
# else:
# cut_off_value = np.inf
# Find the indices of points that...
sample_seeds = np.logical_and.reduce((
~H['local_pt'], # are not localopt points