/
customlm.py
1549 lines (1342 loc) · 77.8 KB
/
customlm.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
"""
Custom implementation of the Levenberg-Marquardt Algorithm
"""
#***************************************************************************************************
# Copyright 2015, 2019 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
# Under the terms of Contract DE-NA0003525 with NTESS, the U.S. Government retains certain rights
# in this software.
# Licensed under the Apache License, Version 2.0 (the "License"); you may not use this file except
# in compliance with the License. You may obtain a copy of the License at
# http://www.apache.org/licenses/LICENSE-2.0 or in the LICENSE file in the root pyGSTi directory.
#***************************************************************************************************
import signal as _signal
import time as _time
import numpy as _np
import scipy as _scipy
from pygsti.optimize import arraysinterface as _ari
from pygsti.optimize.customsolve import custom_solve as _custom_solve
from pygsti.baseobjs.verbosityprinter import VerbosityPrinter as _VerbosityPrinter
from pygsti.baseobjs.resourceallocation import ResourceAllocation as _ResourceAllocation
from pygsti.baseobjs.nicelyserializable import NicelySerializable as _NicelySerializable
# from scipy.optimize import OptimizeResult as _optResult
#Make sure SIGINT will generate a KeyboardInterrupt (even if we're launched in the background)
_signal.signal(_signal.SIGINT, _signal.default_int_handler)
#constants
_MACH_PRECISION = 1e-12
#MU_TOL1 = 1e10 # ??
#MU_TOL2 = 1e3 # ??
class OptimizerResult(object):
"""
The result from an optimization.
Parameters
----------
objective_func : ObjectiveFunction
The objective function that was optimized.
opt_x : numpy.ndarray
The optimal argument (x) value. Often a vector of parameters.
opt_f : numpy.ndarray
the optimal objective function (f) value. Often this is the least-squares
vector of objective function values.
opt_jtj : numpy.ndarray, optional
the optimial `dot(transpose(J),J)` value, where `J`
is the Jacobian matrix. This may be useful for computing
approximate error bars.
opt_unpenalized_f : numpy.ndarray, optional
the optimal objective function (f) value with any
penalty terms removed.
chi2_k_distributed_qty : float, optional
a value that is supposed to be chi2_k distributed.
optimizer_specific_qtys : dict, optional
a dictionary of additional optimization parameters.
"""
def __init__(self, objective_func, opt_x, opt_f=None, opt_jtj=None,
opt_unpenalized_f=None, chi2_k_distributed_qty=None,
optimizer_specific_qtys=None):
self.objective_func = objective_func
self.x = opt_x
self.f = opt_f
self.jtj = opt_jtj # jacobian.T * jacobian
self.f_no_penalties = opt_unpenalized_f
self.optimizer_specific_qtys = optimizer_specific_qtys
self.chi2_k_distributed_qty = chi2_k_distributed_qty
class Optimizer(_NicelySerializable):
"""
An optimizer. Optimizes an objective function.
"""
@classmethod
def cast(cls, obj):
"""
Cast `obj` to a :class:`Optimizer`.
If `obj` is already an `Optimizer` it is just returned,
otherwise this function tries to create a new object
using `obj` as a dictionary of constructor arguments.
Parameters
----------
obj : Optimizer or dict
The object to cast.
Returns
-------
Optimizer
"""
if isinstance(obj, cls):
return obj
else:
return cls(**obj) if obj else cls()
def __init__(self):
super().__init__()
class CustomLMOptimizer(Optimizer):
"""
A Levenberg-Marquardt optimizer customized for GST-like problems.
Parameters
----------
maxiter : int, optional
The maximum number of (outer) interations.
maxfev : int, optional
The maximum function evaluations.
tol : float or dict, optional
The tolerance, specified as a single float or as a dict
with keys `{'relx', 'relf', 'jac', 'maxdx'}`. A single
float sets the `'relf'` and `'jac'` elemments and leaves
the others at their default values.
fditer : int optional
Internally compute the Jacobian using a finite-difference method
for the first `fditer` iterations. This is useful when the initial
point lies at a special or singular point where the analytic Jacobian
is misleading.
first_fditer : int, optional
Number of finite-difference iterations applied to the first
stage of the optimization (only). Unused.
damping_mode : {'identity', 'JTJ', 'invJTJ', 'adaptive'}
How damping is applied. `'identity'` means that the damping parameter mu
multiplies the identity matrix. `'JTJ'` means that mu multiplies the
diagonal or singular values (depending on `scaling_mode`) of the JTJ
(Fischer information and approx. hessaian) matrix, whereas `'invJTJ'`
means mu multiplies the reciprocals of these values instead. The
`'adaptive'` mode adaptively chooses a damping strategy.
damping_basis : {'diagonal_values', 'singular_values'}
Whether the the diagonal or singular values of the JTJ matrix are used
during damping. If `'singular_values'` is selected, then a SVD of the
Jacobian (J) matrix is performed and damping is performed in the basis
of (right) singular vectors. If `'diagonal_values'` is selected, the
diagonal values of relevant matrices are used as a proxy for the the
singular values (saving the cost of performing a SVD).
damping_clip : tuple, optional
A 2-tuple giving upper and lower bounds for the values that mu multiplies.
If `damping_mode == "identity"` then this argument is ignored, as mu always
multiplies a 1.0 on the diagonal if the identity matrix. If None, then no
clipping is applied.
use_acceleration : bool, optional
Whether to include a geodesic acceleration term as suggested in
arXiv:1201.5885. This is supposed to increase the rate of
convergence with very little overhead. In practice we've seen
mixed results.
uphill_step_threshold : float, optional
Allows uphill steps when taking two consecutive steps in nearly
the same direction. The condition for accepting an uphill step
is that `(uphill_step_threshold-beta)*new_objective < old_objective`,
where `beta` is the cosine of the angle between successive steps.
If `uphill_step_threshold == 0` then no uphill steps are allowed,
otherwise it should take a value between 1.0 and 2.0, with 1.0 being
the most permissive to uphill steps.
init_munu : tuple, optional
If not None, a (mu, nu) tuple of 2 floats giving the initial values
for mu and nu.
oob_check_interval : int, optional
Every `oob_check_interval` outer iterations, the objective function
(`obj_fn`) is called with a second argument 'oob_check', set to True.
In this case, `obj_fn` can raise a ValueError exception to indicate
that it is Out Of Bounds. If `oob_check_interval` is 0 then this
check is never performed; if 1 then it is always performed.
oob_action : {"reject","stop"}
What to do when the objective function indicates (by raising a ValueError
as described above). `"reject"` means the step is rejected but the
optimization proceeds; `"stop"` means the optimization stops and returns
as converged at the last known-in-bounds point.
oob_check_mode : int, optional
An advanced option, expert use only. If 0 then the optimization is
halted as soon as an *attempt* is made to evaluate the function out of bounds.
If 1 then the optimization is halted only when a would-be *accepted* step
is out of bounds.
serial_solve_proc_threshold : int, optional
When there are fewer than this many processors, the optimizer will solve linear
systems serially, using SciPy on a single processor, rather than using a parallelized
Gaussian Elimination (with partial pivoting) algorithm coded in Python. Since SciPy's
implementation is more efficient, it's not worth using the parallel version until there
are many processors to spread the work among.
lsvec_mode : {'normal', 'percircuit'}
Whether the terms used in the least-squares optimization are the "elements" as computed
by the objective function's `.terms()` and `.lsvec()` methods (`'normal'` mode) or the
"per-circuit quantities" computed by the objective function's `.percircuit()` and
`.lsvec_percircuit()` methods (`'percircuit'` mode).
"""
def __init__(self, maxiter=100, maxfev=100, tol=1e-6, fditer=0, first_fditer=0, damping_mode="identity",
damping_basis="diagonal_values", damping_clip=None, use_acceleration=False,
uphill_step_threshold=0.0, init_munu="auto", oob_check_interval=0,
oob_action="reject", oob_check_mode=0, serial_solve_proc_threshold=100, lsvec_mode="normal"):
super().__init__()
if isinstance(tol, float): tol = {'relx': 1e-8, 'relf': tol, 'f': 1.0, 'jac': tol, 'maxdx': 1.0}
self.maxiter = maxiter
self.maxfev = maxfev
self.tol = tol
self.fditer = fditer
self.first_fditer = first_fditer
self.damping_mode = damping_mode
self.damping_basis = damping_basis
self.damping_clip = damping_clip
self.use_acceleration = use_acceleration
self.uphill_step_threshold = uphill_step_threshold
self.init_munu = init_munu
self.oob_check_interval = oob_check_interval
self.oob_action = oob_action
self.oob_check_mode = oob_check_mode
self.array_types = 3 * ('p',) + ('e', 'ep') # see custom_leastsq fn "-type"s -need to add 'jtj' type
self.called_objective_methods = ('lsvec', 'dlsvec') # the objective function methods we use (for mem estimate)
self.serial_solve_proc_threshold = serial_solve_proc_threshold
self.lsvec_mode = lsvec_mode
def _to_nice_serialization(self):
state = super()._to_nice_serialization()
state.update({
'maximum_iterations': self.maxiter,
'maximum_function_evaluations': self.maxfev,
'tolerance': self.tol,
'number_of_finite_difference_iterations': self.fditer,
'number_of_first_stage_finite_difference_iterations': self.first_fditer,
'damping_mode': self.damping_mode,
'damping_basis': self.damping_basis,
'damping_clip': self.damping_clip,
'use_acceleration': self.use_acceleration,
'uphill_step_threshold': self.uphill_step_threshold,
'initial_mu_and_nu': self.init_munu,
'out_of_bounds_check_interval': self.oob_check_interval,
'out_of_bounds_action': self.oob_action,
'out_of_bounds_check_mode': self.oob_check_mode,
'array_types': self.array_types,
'called_objective_function_methods': self.called_objective_methods,
'serial_solve_number_of_processors_threshold': self.serial_solve_proc_threshold,
'lsvec_mode': self.lsvec_mode
})
return state
@classmethod
def _from_nice_serialization(cls, state):
return cls(maxiter=state['maximum_iterations'],
maxfev=state['maximum_function_evaluations'],
tol=state['tolerance'],
fditer=state['number_of_finite_difference_iterations'],
first_fditer=state['number_of_first_stage_finite_difference_iterations'],
damping_mode=state['damping_mode'],
damping_basis=state['damping_basis'],
damping_clip=state['damping_clip'],
use_acceleration=state['use_acceleration'],
uphill_step_threshold=state['uphill_step_threshold'],
init_munu=state['initial_mu_and_nu'],
oob_check_interval=state['out_of_bounds_check_interval'],
oob_action=state['out_of_bounds_action'],
oob_check_mode=state['out_of_bounds_check_mode'],
serial_solve_proc_threshold=state['serial_solve_number_of_processors_threshold'],
lsvec_mode=state.get('lsvec_mode', 'normal'))
def run(self, objective, profiler, printer):
"""
Perform the optimization.
Parameters
----------
objective : ObjectiveFunction
The objective function to optimize.
profiler : Profiler
A profiler to track resource usage.
printer : VerbosityPrinter
printer to use for sending output to stdout.
"""
nExtra = objective.ex # number of additional "extra" elements
if self.lsvec_mode == 'normal':
objective_func = objective.lsvec
jacobian = objective.dlsvec
nEls = objective.layout.num_elements + nExtra # 'e' for array types
elif self.lsvec_mode == 'percircuit':
objective_func = objective.lsvec_percircuit
jacobian = objective.dlsvec_percircuit
nEls = objective.layout.num_circuits + nExtra # 'e' for array types
else:
raise ValueError("Invalid `lsvec_mode`: %s" % str(self.lsvec_mode))
x0 = objective.model.to_vector()
x_limits = objective.model.parameter_bounds
# x_limits should be a (num_params, 2)-shaped array, holding on each row the (min, max) values for the
# corresponding parameter (element of the "x" vector) or `None`. If `None`, then no limits are imposed.
# Check memory limit can handle what custom_leastsq will "allocate"
nP = len(x0) # 'p' for array types
objective.resource_alloc.check_can_allocate_memory(3 * nP + nEls + nEls * nP + nP * nP) # see array_types above
from ..layouts.distlayout import DistributableCOPALayout as _DL
ari = _ari.DistributedArraysInterface(objective.layout, self.lsvec_mode, nExtra) \
if isinstance(objective.layout, _DL) else _ari.UndistributedArraysInterface(nEls, nP)
opt_x, converged, msg, mu, nu, norm_f, f, opt_jtj = custom_leastsq(
objective_func, jacobian, x0,
max_iter=self.maxiter,
num_fd_iters=self.fditer,
f_norm2_tol=self.tol.get('f', 1.0),
jac_norm_tol=self.tol.get('jac', 1e-6),
rel_ftol=self.tol.get('relf', 1e-6),
rel_xtol=self.tol.get('relx', 1e-8),
max_dx_scale=self.tol.get('maxdx', 1.0),
damping_mode=self.damping_mode,
damping_basis=self.damping_basis,
damping_clip=self.damping_clip,
use_acceleration=self.use_acceleration,
uphill_step_threshold=self.uphill_step_threshold,
init_munu=self.init_munu,
oob_check_interval=self.oob_check_interval,
oob_action=self.oob_action,
oob_check_mode=self.oob_check_mode,
resource_alloc=objective.resource_alloc,
arrays_interface=ari,
serial_solve_proc_threshold=self.serial_solve_proc_threshold,
x_limits=x_limits,
verbosity=printer - 1, profiler=profiler)
printer.log("Least squares message = %s" % msg, 2)
assert(converged), "Failed to converge: %s" % msg
current_v = objective.model.to_vector()
if not _np.allclose(current_v, opt_x): # ensure the last model evaluation was at opt_x
objective_func(opt_x)
#objective.model.from_vector(opt_x) # performed within line above
#DEBUG CHECK SYNC between procs (especially for shared mem) - could REMOVE
# if objective.resource_alloc.comm is not None:
# comm = objective.resource_alloc.comm
# v_cmp = comm.bcast(objective.model.to_vector() if (comm.Get_rank() == 0) else None, root=0)
# v_matches_x = _np.allclose(objective.model.to_vector(), opt_x)
# same_as_root = _np.isclose(_np.linalg.norm(objective.model.to_vector() - v_cmp), 0.0)
# if not (v_matches_x and same_as_root):
# raise ValueError("Rank %d CUSTOMLM ERROR: END model vector-matches-x=%s and vector-is-same-as-root=%s"
# % (comm.rank, str(v_matches_x), str(same_as_root)))
# comm.barrier() # if we get past here, then *all* processors are OK
# if comm.rank == 0:
# print("OK - model vector == best_x and all vectors agree w/root proc's")
unpenalized_f = f[0:-objective.ex] if (objective.ex > 0) else f
unpenalized_normf = sum(unpenalized_f**2) # objective function without penalty factors
chi2k_qty = objective.chi2k_distributed_qty(norm_f)
return OptimizerResult(objective, opt_x, norm_f, opt_jtj, unpenalized_normf, chi2k_qty,
{'msg': msg, 'mu': mu, 'nu': nu, 'fvec': f})
#Scipy version...
# opt_x, _, _, msg, flag = \
# _spo.leastsq(objective_func, x0, xtol=tol['relx'], ftol=tol['relf'], gtol=tol['jac'],
# maxfev=maxfev * (len(x0) + 1), full_output=True, Dfun=jacobian) # pragma: no cover
# printer.log("Least squares message = %s; flag =%s" % (msg, flag), 2) # pragma: no cover
# opt_state = (msg,)
def custom_leastsq(obj_fn, jac_fn, x0, f_norm2_tol=1e-6, jac_norm_tol=1e-6,
rel_ftol=1e-6, rel_xtol=1e-6, max_iter=100, num_fd_iters=0,
max_dx_scale=1.0, damping_mode="identity", damping_basis="diagonal_values",
damping_clip=None, use_acceleration=False, uphill_step_threshold=0.0,
init_munu="auto", oob_check_interval=0, oob_action="reject", oob_check_mode=0,
resource_alloc=None, arrays_interface=None, serial_solve_proc_threshold=100,
x_limits=None, verbosity=0, profiler=None):
"""
An implementation of the Levenberg-Marquardt least-squares optimization algorithm customized for use within pyGSTi.
This general purpose routine mimic to a large extent the interface used by
`scipy.optimize.leastsq`, though it implements a newer (and more robust) version
of the algorithm.
Parameters
----------
obj_fn : function
The objective function. Must accept and return 1D numpy ndarrays of
length N and M respectively. Same form as scipy.optimize.leastsq.
jac_fn : function
The jacobian function (not optional!). Accepts a 1D array of length N
and returns an array of shape (M,N).
x0 : numpy.ndarray
Initial evaluation point.
f_norm2_tol : float, optional
Tolerace for `F^2` where `F = `norm( sum(obj_fn(x)**2) )` is the
least-squares residual. If `F**2 < f_norm2_tol`, then mark converged.
jac_norm_tol : float, optional
Tolerance for jacobian norm, namely if `infn(dot(J.T,f)) < jac_norm_tol`
then mark converged, where `infn` is the infinity-norm and
`f = obj_fn(x)`.
rel_ftol : float, optional
Tolerance on the relative reduction in `F^2`, that is, if
`d(F^2)/F^2 < rel_ftol` then mark converged.
rel_xtol : float, optional
Tolerance on the relative value of `|x|`, so that if
`d(|x|)/|x| < rel_xtol` then mark converged.
max_iter : int, optional
The maximum number of (outer) interations.
num_fd_iters : int optional
Internally compute the Jacobian using a finite-difference method
for the first `num_fd_iters` iterations. This is useful when `x0`
lies at a special or singular point where the analytic Jacobian is
misleading.
max_dx_scale : float, optional
If not None, impose a limit on the magnitude of the step, so that
`|dx|^2 < max_dx_scale^2 * len(dx)` (so elements of `dx` should be,
roughly, less than `max_dx_scale`).
damping_mode : {'identity', 'JTJ', 'invJTJ', 'adaptive'}
How damping is applied. `'identity'` means that the damping parameter mu
multiplies the identity matrix. `'JTJ'` means that mu multiplies the
diagonal or singular values (depending on `scaling_mode`) of the JTJ
(Fischer information and approx. hessaian) matrix, whereas `'invJTJ'`
means mu multiplies the reciprocals of these values instead. The
`'adaptive'` mode adaptively chooses a damping strategy.
damping_basis : {'diagonal_values', 'singular_values'}
Whether the the diagonal or singular values of the JTJ matrix are used
during damping. If `'singular_values'` is selected, then a SVD of the
Jacobian (J) matrix is performed and damping is performed in the basis
of (right) singular vectors. If `'diagonal_values'` is selected, the
diagonal values of relevant matrices are used as a proxy for the the
singular values (saving the cost of performing a SVD).
damping_clip : tuple, optional
A 2-tuple giving upper and lower bounds for the values that mu multiplies.
If `damping_mode == "identity"` then this argument is ignored, as mu always
multiplies a 1.0 on the diagonal if the identity matrix. If None, then no
clipping is applied.
use_acceleration : bool, optional
Whether to include a geodesic acceleration term as suggested in
arXiv:1201.5885. This is supposed to increase the rate of
convergence with very little overhead. In practice we've seen
mixed results.
uphill_step_threshold : float, optional
Allows uphill steps when taking two consecutive steps in nearly
the same direction. The condition for accepting an uphill step
is that `(uphill_step_threshold-beta)*new_objective < old_objective`,
where `beta` is the cosine of the angle between successive steps.
If `uphill_step_threshold == 0` then no uphill steps are allowed,
otherwise it should take a value between 1.0 and 2.0, with 1.0 being
the most permissive to uphill steps.
init_munu : tuple, optional
If not None, a (mu, nu) tuple of 2 floats giving the initial values
for mu and nu.
oob_check_interval : int, optional
Every `oob_check_interval` outer iterations, the objective function
(`obj_fn`) is called with a second argument 'oob_check', set to True.
In this case, `obj_fn` can raise a ValueError exception to indicate
that it is Out Of Bounds. If `oob_check_interval` is 0 then this
check is never performed; if 1 then it is always performed.
oob_action : {"reject","stop"}
What to do when the objective function indicates (by raising a ValueError
as described above). `"reject"` means the step is rejected but the
optimization proceeds; `"stop"` means the optimization stops and returns
as converged at the last known-in-bounds point.
oob_check_mode : int, optional
An advanced option, expert use only. If 0 then the optimization is
halted as soon as an *attempt* is made to evaluate the function out of bounds.
If 1 then the optimization is halted only when a would-be *accepted* step
is out of bounds.
resource_alloc : ResourceAllocation, optional
When not None, an resource allocation object used for distributing the computation
across multiple processors.
arrays_interface : ArraysInterface
An object that provides an interface for creating and manipulating data arrays.
serial_solve_proc_threshold : int optional
When there are fewer than this many processors, the optimizer will solve linear
systems serially, using SciPy on a single processor, rather than using a parallelized
Gaussian Elimination (with partial pivoting) algorithm coded in Python. Since SciPy's
implementation is more efficient, it's not worth using the parallel version until there
are many processors to spread the work among.
x_limits : numpy.ndarray, optional
A (num_params, 2)-shaped array, holding on each row the (min, max) values for the corresponding
parameter (element of the "x" vector). If `None`, then no limits are imposed.
verbosity : int, optional
Amount of detail to print to stdout.
profiler : Profiler, optional
A profiler object used for to track timing and memory usage.
Returns
-------
x : numpy.ndarray
The optimal solution.
converged : bool
Whether the solution converged.
msg : str
A message indicating why the solution converged (or didn't).
"""
resource_alloc = _ResourceAllocation.cast(resource_alloc)
comm = resource_alloc.comm
printer = _VerbosityPrinter.create_printer(verbosity, comm)
ari = arrays_interface # shorthand
# MEM from ..baseobjs.profiler import Profiler
# MEM debug_prof = Profiler(comm, True)
# MEM profiler = debug_prof
msg = ""
converged = False
global_x = x0.copy()
f = obj_fn(global_x) # 'E'-type array
norm_f = ari.norm2_f(f) # _np.linalg.norm(f)**2
half_max_nu = 2**62 # what should this be??
tau = 1e-3
alpha = 0.5 # for acceleration
nu = 2
mu = 1 # just a guess - initialized on 1st iter and only used if rejected
#Allocate potentially shared memory used in loop
JTJ = ari.allocate_jtj()
JTf = ari.allocate_jtf()
x = ari.allocate_jtf()
#x_for_jac = ari.allocate_x_for_jac()
if num_fd_iters > 0:
fdJac = ari.allocate_jac()
ari.allscatter_x(global_x, x)
if x_limits is not None:
x_lower_limits = ari.allocate_jtf()
x_upper_limits = ari.allocate_jtf()
ari.allscatter_x(x_limits[:, 0], x_lower_limits)
ari.allscatter_x(x_limits[:, 1], x_upper_limits)
if damping_basis == "singular_values":
Jac_V = ari.allocate_jtj()
if damping_mode == 'adaptive':
dx_lst = [ari.allocate_jtf(), ari.allocate_jtf(), ari.allocate_jtf()]
new_x_lst = [ari.allocate_jtf(), ari.allocate_jtf(), ari.allocate_jtf()]
global_new_x_lst = [global_x.copy() for i in range(3)]
else:
dx = ari.allocate_jtf()
new_x = ari.allocate_jtf()
global_new_x = global_x.copy()
if use_acceleration:
dx1 = ari.allocate_jtf()
dx2 = ari.allocate_jtf()
df2_x = ari.allocate_jtf()
JTdf2 = ari.allocate_jtf()
global_accel_x = global_x.copy()
# don't let any component change by more than ~max_dx_scale
if max_dx_scale:
max_norm_dx = (max_dx_scale**2) * len(global_x)
else: max_norm_dx = None
if not _np.isfinite(norm_f):
msg = "Infinite norm of objective function at initial point!"
if len(global_x) == 0: # a model with 0 parameters - nothing to optimize
msg = "No parameters to optimize"; converged = True
# DB: from ..tools import matrixtools as _mt
# DB: print("DB F0 (%s)=" % str(f.shape)); _mt.print_mx(f,prec=0,width=4)
#num_fd_iters = 1000000 # DEBUG: use finite difference iterations instead
# print("DEBUG: setting num_fd_iters == 0!"); num_fd_iters = 0 # DEBUG
last_accepted_dx = None
min_norm_f = 1e100 # sentinel
best_x = ari.allocate_jtf()
best_x[:] = x[:] # like x.copy() -the x-value corresponding to min_norm_f ('P'-type)
spow = 0.0 # for damping_mode == 'adaptive'
if damping_clip is not None:
def dclip(ar): return _np.clip(ar, damping_clip[0], damping_clip[1])
else:
def dclip(ar): return ar
if init_munu != "auto":
mu, nu = init_munu
best_x_state = (mu, nu, norm_f, f.copy(), spow, None) # need f.copy() b/c f is objfn mem
rawJTJ_scratch = None
jtj_buf = ari.allocate_jtj_shared_mem_buf()
try:
for k in range(max_iter): # outer loop
# assume global_x, x, f, fnorm hold valid values
if len(msg) > 0:
break # exit outer loop if an exit-message has been set
if norm_f < f_norm2_tol:
if oob_check_interval <= 1:
msg = "Sum of squares is at most %g" % f_norm2_tol
converged = True; break
else:
printer.log(("** Converged with out-of-bounds with check interval=%d, reverting to last "
"know in-bounds point and setting interval=1 **") % oob_check_interval, 2)
oob_check_interval = 1
x[:] = best_x[:]
mu, nu, norm_f, f[:], spow, _ = best_x_state
continue # can't make use of saved JTJ yet - recompute on nxt iter
#printer.log("--- Outer Iter %d: norm_f = %g, mu=%g" % (k,norm_f,mu))
if profiler: profiler.memory_check("custom_leastsq: begin outer iter *before de-alloc*")
Jac = None
if profiler: profiler.memory_check("custom_leastsq: begin outer iter")
# unnecessary b/c global_x is already valid: ari.allgather_x(x, global_x)
if k >= num_fd_iters:
Jac = jac_fn(global_x) # 'EP'-type, but doesn't actually allocate any more mem (!)
else:
# Note: x holds only number of "fine"-division params - need to use global_x, and
# Jac only holds a subset of the derivative and element columns and rows, respectively.
f_fixed = f.copy() # a static part of the distributed `f` resturned by obj_fn - MUST copy this.
pslice = ari.jac_param_slice(only_if_leader=True)
eps = 1e-7
#Don't do this: for ii, i in enumerate(range(pslice.start, pslice.stop)): (must keep procs in sync)
for i in range(len(global_x)):
x_plus_dx = global_x.copy()
x_plus_dx[i] += eps
fd = (obj_fn(x_plus_dx) - f_fixed) / eps
if pslice.start <= i < pslice.stop:
fdJac[:, i - pslice.start] = fd
#if comm is not None: comm.barrier() # overkill for shared memory leader host barrier
Jac = fdJac
#DEBUG: compare with analytic jacobian (need to uncomment num_fd_iters DEBUG line above too)
#Jac_analytic = jac_fn(x)
#if _np.linalg.norm(Jac_analytic-Jac) > 1e-6:
# print("JACDIFF = ",_np.linalg.norm(Jac_analytic-Jac)," per el=",
# _np.linalg.norm(Jac_analytic-Jac)/Jac.size," sz=",Jac.size)
# DB: from ..tools import matrixtools as _mt
# DB: print("DB JAC (%s)=" % str(Jac.shape)); _mt.print_mx(Jac,prec=0,width=4); assert(False)
if profiler: profiler.memory_check("custom_leastsq: after jacobian:"
+ "shape=%s, GB=%.2f" % (str(Jac.shape),
Jac.nbytes / (1024.0**3)))
Jnorm = _np.sqrt(ari.norm2_jac(Jac))
xnorm = _np.sqrt(ari.norm2_x(x))
printer.log("--- Outer Iter %d: norm_f = %g, mu=%g, |x|=%g, |J|=%g" % (k, norm_f, mu, xnorm, Jnorm))
#assert(_np.isfinite(Jac).all()), "Non-finite Jacobian!" # NaNs tracking
#assert(_np.isfinite(_np.linalg.norm(Jac))), "Finite Jacobian has inf norm!" # NaNs tracking
tm = _time.time()
#OLD MPI-enabled JTJ computation
##if my_mpidot_qtys is None:
## my_mpidot_qtys = _mpit.distribute_for_dot(Jac.T.shape, Jac.shape, resource_alloc)
#JTJ, JTJ_shm = _mpit.mpidot(Jac.T, Jac, my_mpidot_qtys[0], my_mpidot_qtys[1],
# my_mpidot_qtys[2], resource_alloc, JTJ, JTJ_shm) # _np.dot(Jac.T,Jac) 'PP'
ari.fill_jtj(Jac, JTJ, jtj_buf)
ari.fill_jtf(Jac, f, JTf) # 'P'-type
if profiler: profiler.add_time("custom_leastsq: dotprods", tm)
#assert(not _np.isnan(JTJ).any()), "NaN in JTJ!" # NaNs tracking
#assert(not _np.isinf(JTJ).any()), "inf in JTJ! norm Jac = %g" % _np.linalg.norm(Jac) # NaNs tracking
#assert(_np.isfinite(JTJ).all()), "Non-finite JTJ!" # NaNs tracking
#assert(_np.isfinite(JTf).all()), "Non-finite JTf!" # NaNs tracking
idiag = ari.jtj_diag_indices(JTJ)
norm_JTf = ari.infnorm_x(JTf)
norm_x = ari.norm2_x(x) # _np.linalg.norm(x)**2
undamped_JTJ_diag = JTJ[idiag].copy() # 'P'-type
#max_JTJ_diag = JTJ.diagonal().copy()
JTf *= -1.0; minus_JTf = JTf # use the same memory for -JTf below (shouldn't use JTf anymore)
#Maybe just have a minus_JTf variable?
# FUTURE TODO: keep tallying allocated memory, i.e. array_types (stopped here)
if damping_basis == "singular_values":
# Jac = U * s * Vh; J.T * J = conj(V) * s * U.T * U * s * Vh = conj(V) * s^2 * Vh
# Jac_U, Jac_s, Jac_Vh = _np.linalg.svd(Jac, full_matrices=False)
# Jac_V = _np.conjugate(Jac_Vh.T)
global_JTJ = ari.gather_jtj(JTJ)
if comm is None or comm.rank == 0:
global_Jac_s2, global_Jac_V = _np.linalg.eigh(global_JTJ)
ari.scatter_jtj(global_Jac_V, Jac_V)
comm.bcast(global_Jac_s2, root=0)
else:
ari.scatter_jtj(None, Jac_V)
global_Jac_s2 = comm.bcast(None, root=0)
#print("Rank %d: min s2 = %g" % (comm.rank, min(global_Jac_s2)))
#if min(global_Jac_s2) < -1e-4 and (comm is None or comm.rank == 0):
# print("WARNING: min Jac s^2 = %g (max = %g)" % (min(global_Jac_s2), max(global_Jac_s2)))
assert(min(global_Jac_s2) / abs(max(global_Jac_s2)) > -1e-6), "JTJ should be positive!"
global_Jac_s = _np.sqrt(_np.clip(global_Jac_s2, 1e-12, None)) # eigvals of JTJ must be >= 0
global_Jac_VT_mJTf = ari.global_svd_dot(Jac_V, minus_JTf) # = dot(Jac_V.T, minus_JTf)
#DEBUG
#num_large_svals = _np.count_nonzero(Jac_s > _np.max(Jac_s) / 1e2)
#Jac_Uproj = Jac_U[:,0:num_large_svals]
#JTJ_evals, JTJ_U = _np.linalg.eig(JTJ)
#printer.log("JTJ (dim=%d) eval min/max=%g, %g; %d large svals (of %d)" % (
# JTJ.shape[0], _np.min(_np.abs(JTJ_evals)), _np.max(_np.abs(JTJ_evals)),
# num_large_svals, len(Jac_s)))
if norm_JTf < jac_norm_tol:
if oob_check_interval <= 1:
msg = "norm(jacobian) is at most %g" % jac_norm_tol
converged = True; break
else:
printer.log(("** Converged with out-of-bounds with check interval=%d, reverting to last "
"know in-bounds point and setting interval=1 **") % oob_check_interval, 2)
oob_check_interval = 1
x[:] = best_x[:]
mu, nu, norm_f, f[:], spow, _ = best_x_state
continue # can't make use of saved JTJ yet - recompute on nxt iter
if k == 0:
if init_munu == "auto":
if damping_mode == 'identity':
mu = tau * ari.max_x(undamped_JTJ_diag) # initial damping element
#mu = min(mu, MU_TOL1)
else:
# initial multiplicative damping element
#mu = tau # initial damping element - but this seem to low, at least for termgap...
mu = min(1.0e5, ari.max_x(undamped_JTJ_diag) / norm_JTf) # Erik's heuristic
#tries to avoid making mu so large that dx is tiny and we declare victory prematurely
else:
mu, nu = init_munu
rawJTJ_scratch = JTJ.copy() # allocates the memory for a copy of JTJ so only update mem elsewhere
best_x_state = mu, nu, norm_f, f.copy(), spow, rawJTJ_scratch # update mu,nu,JTJ of initial best state
else:
#on all other iterations, update JTJ of best_x_state if best_x == x, i.e. if we've just evaluated
# a previously accepted step that was deemed the best we've seen so far
if _np.allclose(x, best_x):
rawJTJ_scratch[:, :] = JTJ[:, :] # use pre-allocated memory
rawJTJ_scratch[idiag] = undamped_JTJ_diag # no damping; the "raw" JTJ
best_x_state = best_x_state[0:5] + (rawJTJ_scratch,) # update mu,nu,JTJ of initial "best state"
#determing increment using adaptive damping
while True: # inner loop
if profiler: profiler.memory_check("custom_leastsq: begin inner iter")
#print("DB: Pre-damping JTJ diag = [",_np.min(_np.abs(JTJ[idiag])),_np.max(_np.abs(JTJ[idiag])),"]")
if damping_mode == 'identity':
assert(damping_clip is None), "damping_clip cannot be used with damping_mode == 'identity'"
if damping_basis == "singular_values":
reg_Jac_s = global_Jac_s + mu
#Notes:
#Previously we computed inv_JTJ here and below computed dx:
#inv_JTJ = _np.dot(Jac_V, _np.dot(_np.diag(1 / reg_Jac_s**2), Jac_V.T))
# dx = _np.dot(Jac_V, _np.diag(1 / reg_Jac_s**2), global_Jac_VT_mJTf
#But now we just compute reg_Jac_s here, and so the rest below.
else:
# ok if assume fine-param-proc.size == 1 (otherwise need to sync setting local JTJ)
JTJ[idiag] = undamped_JTJ_diag + mu # augment normal equations
elif damping_mode == 'JTJ':
if damping_basis == "singular_values":
reg_Jac_s = global_Jac_s + mu * dclip(global_Jac_s)
else:
add_to_diag = mu * dclip(undamped_JTJ_diag)
JTJ[idiag] = undamped_JTJ_diag + add_to_diag # ok if assume fine-param-proc.size == 1
elif damping_mode == 'invJTJ':
if damping_basis == "singular_values":
reg_Jac_s = global_Jac_s + mu * dclip(1.0 / global_Jac_s)
else:
add_to_diag = mu * dclip(1.0 / undamped_JTJ_diag)
JTJ[idiag] = undamped_JTJ_diag + add_to_diag # ok if assume fine-param-proc.size == 1
elif damping_mode == 'adaptive':
if damping_basis == "singular_values":
reg_Jac_s_lst = [global_Jac_s + mu * dclip(global_Jac_s**(spow + 0.1)),
global_Jac_s + mu * dclip(global_Jac_s**spow),
global_Jac_s + mu * dclip(global_Jac_s**(spow - 0.1))]
else:
add_to_diag_lst = [mu * dclip(undamped_JTJ_diag**(spow + 0.1)),
mu * dclip(undamped_JTJ_diag**spow),
mu * dclip(undamped_JTJ_diag**(spow - 0.1))]
else:
raise ValueError("Invalid damping mode: %s" % damping_mode)
#assert(_np.isfinite(JTJ).all()), "Non-finite JTJ (inner)!" # NaNs tracking
#assert(_np.isfinite(JTf).all()), "Non-finite JTf (inner)!" # NaNs tracking
try:
if profiler: profiler.memory_check("custom_leastsq: before linsolve")
tm = _time.time()
success = True
if damping_basis == 'diagonal_values':
if damping_mode == 'adaptive':
for ii, add_to_diag in enumerate(add_to_diag_lst):
JTJ[idiag] = undamped_JTJ_diag + add_to_diag # ok if assume fine-param-proc.size == 1
#dx_lst.append(_scipy.linalg.solve(JTJ, -JTf, sym_pos=True))
#dx_lst.append(custom_solve(JTJ, -JTf, resource_alloc))
_custom_solve(JTJ, minus_JTf, dx_lst[ii], ari, resource_alloc,
serial_solve_proc_threshold)
else:
#dx = _scipy.linalg.solve(JTJ, -JTf, sym_pos=True)
_custom_solve(JTJ, minus_JTf, dx, ari, resource_alloc, serial_solve_proc_threshold)
elif damping_basis == 'singular_values':
#Note: above solves JTJ*x = -JTf => x = inv_JTJ * (-JTf)
# but: J = U*s*Vh => JTJ = (VhT*s*UT)(U*s*Vh) = VhT*s^2*Vh, and inv_Vh = V b/c V is unitary
# so inv_JTJ = inv_Vh * 1/s^2 * inv_VhT = V * 1/s^2 * VT = (N,K)*(K,K)*(K,N) if use psuedoinv
if damping_mode == 'adaptive':
#dx_lst = [_np.dot(ijtj, minus_JTf) for ijtj in inv_JTJ_lst] # special case
for ii, s in enumerate(reg_Jac_s_lst):
ari.fill_dx_svd(Jac_V, (1 / s**2) * global_Jac_VT_mJTf, dx_lst[ii])
else:
# dx = _np.dot(inv_JTJ, minus_JTf)
ari.fill_dx_svd(Jac_V, (1 / reg_Jac_s**2) * global_Jac_VT_mJTf, dx)
else:
raise ValueError("Invalid damping_basis = '%s'" % damping_basis)
if profiler: profiler.add_time("custom_leastsq: linsolve", tm)
#except _np.linalg.LinAlgError:
except _scipy.linalg.LinAlgError: # DIST TODO - a different kind of exception caught?
success = False
if success and use_acceleration: # Find acceleration term:
assert(damping_mode != 'adaptive'), "Cannot use acceleration in adaptive mode (yet)"
assert(damping_basis != 'singular_values'), "Cannot use acceleration w/singular-value basis (yet)"
df2_eps = 1.0
try:
#df2 = (obj_fn(x + df2_dx) + obj_fn(x - df2_dx) - 2 * f) / \
# df2_eps**2 # 2nd deriv of f along dx direction
# Above line expanded to reuse shared memory
df2 = -2 * f
df2_x[:] = x + df2_eps * dx
ari.allgather_x(df2_x, global_accel_x)
df2 += obj_fn(global_accel_x)
df2_x[:] = x - df2_eps * dx
ari.allgather_x(df2_x, global_accel_x)
df2 += obj_fn(global_accel_x)
df2 /= df2_eps**2
f[:] = df2; df2 = f # use `f` as an appropriate shared-mem object for fill_jtf below
ari.fill_jtf(Jac, df2, JTdf2)
JTdf2 *= -0.5 # keep using JTdf2 memory in solve call below
#dx2 = _scipy.linalg.solve(JTJ, -0.5 * JTdf2, sym_pos=True) # Note: JTJ not init w/'adaptive'
_custom_solve(JTJ, JTdf2, dx2, ari, resource_alloc, serial_solve_proc_threshold)
dx1[:] = dx[:]
dx += dx2 # add acceleration term to dx
except _scipy.linalg.LinAlgError:
print("WARNING - linear solve failed for acceleration term!")
# but ok to continue - just stick with first order term
except ValueError:
print("WARNING - value error during computation of acceleration term!")
reject_msg = ""
if profiler: profiler.memory_check("custom_leastsq: after linsolve")
if success: # linear solve succeeded
#dx = _hack_dx(obj_fn, x, dx, Jac, JTJ, JTf, f, norm_f)
if damping_mode != 'adaptive':
new_x[:] = x + dx
norm_dx = ari.norm2_x(dx) # _np.linalg.norm(dx)**2
#ensure dx isn't too large - don't let any component change by more than ~max_dx_scale
if max_norm_dx and norm_dx > max_norm_dx:
dx *= _np.sqrt(max_norm_dx / norm_dx)
new_x[:] = x + dx
norm_dx = ari.norm2_x(dx) # _np.linalg.norm(dx)**2
#apply x limits (bounds)
if x_limits is not None:
# Approach 1: project x into valid space by simply clipping out-of-bounds values
for i, (x_el, lower, upper) in enumerate(zip(x, x_lower_limits, x_upper_limits)):
if new_x[i] < lower:
new_x[i] = lower
dx[i] = lower - x_el
elif new_x[i] > upper:
new_x[i] = upper
dx[i] = upper - x_el
norm_dx = ari.norm2_x(dx) # _np.linalg.norm(dx)**2
# Approach 2: by scaling back dx (seems less good, but here in case we want it later)
# # minimally reduce dx s.t. new_x = x + dx so that x_lower_limits <= x+dx <= x_upper_limits
# # x_lower_limits - x <= dx <= x_upper_limits - x. Note: use potentially updated dx from
# # max_norm_dx block above. For 0 <= scale <= 1,
# # 1) require x + scale*dx - x_upper_limits <= 0 => scale <= (x_upper_limits - x) / dx
# # [Note: above assumes dx > 0 b/c if not it moves x away from bound and scale < 0]
# # so if scale >= 0, then scale = min((x_upper_limits - x) / dx, 1.0)
# scale = None
# new_x[:] = (x_upper_limits - x) / dx
# new_x_min = ari.min_x(new_x)
# if 0 <= new_x_min < 1.0:
# scale = new_x_min
#
# # 2) require x + scale*dx - x_lower_limits <= 0 => scale <= (x - x_lower_limits) / (-dx)
# new_x[:] = (x_lower_limits - x) / dx
# new_x_min = ari.min_x(new_x)
# if 0 <= new_x_min < 1.0:
# scale = new_x_min if (scale is None) else min(new_x_min, scale)
#
# if scale is not None:
# dx *= scale
# new_x[:] = x + dx
# norm_dx = ari.norm2_x(dx) # _np.linalg.norm(dx)**2
else:
for dx, new_x in zip(dx_lst, new_x_lst):
new_x[:] = x + dx
norm_dx_lst = [ari.norm2_x(dx) for dx in dx_lst]
#ensure dx isn't too large - don't let any component change by more than ~max_dx_scale
if max_norm_dx:
for i, norm_dx in enumerate(norm_dx_lst):
if norm_dx > max_norm_dx:
dx_lst[i] *= _np.sqrt(max_norm_dx / norm_dx)
new_x_lst[i][:] = x + dx_lst[i]
norm_dx_lst[i] = ari.norm2_x(dx_lst[i])
#apply x limits (bounds)
if x_limits is not None:
for i, (dx, new_x) in enumerate(zip(dx_lst, new_x_lst)):
# Do same thing as above for each possible dx in dx_lst
# Approach 1:
for ii, (x_el, lower, upper) in enumerate(zip(x, x_lower_limits, x_upper_limits)):
if new_x[ii] < lower:
new_x[ii] = lower
dx[ii] = lower - x_el
elif new_x[ii] > upper:
new_x[ii] = upper
dx[ii] = upper - x_el
norm_dx_lst[i] = ari.norm2_x(dx) # _np.linalg.norm(dx)**2
# Approach 2:
# scale = None
# new_x[:] = (x_upper_limits - x) / dx
# new_x_min = ari.min_x(new_x)
# if 0 <= new_x_min < 1.0:
# scale = new_x_min
#
# new_x[:] = (x_lower_limits - x) / dx
# new_x_min = ari.min_x(new_x)
# if 0 <= new_x_min < 1.0:
# scale = new_x_min if (scale is None) else min(new_x_min, scale)
#
# if scale is not None:
# dx *= scale
# new_x[:] = x + dx
# norm_dx_lst[i] = ari.norm2_x(dx)
norm_dx = norm_dx_lst[1] # just use center value for printing & checks below
printer.log(" - Inner Loop: mu=%g, norm_dx=%g" % (mu, norm_dx), 2)
#MEM if profiler: profiler.memory_check("custom_leastsq: mid inner loop")
#print("DB: new_x = ", new_x)
if norm_dx < (rel_xtol**2) * norm_x: # and mu < MU_TOL2:
if oob_check_interval <= 1:
msg = "Relative change, |dx|/|x|, is at most %g" % rel_xtol
converged = True; break
else:
printer.log(("** Converged with out-of-bounds with check interval=%d, reverting to last "
"know in-bounds point and setting interval=1 **") % oob_check_interval, 2)
oob_check_interval = 1
x[:] = best_x[:]