BrentSolver#
BrentSolver is an implementation of Scipy’s Brentq, which is based on the Wijngaarden-Dekker-Brent method. It can solve for a single state bracketed between a lower and upper bound using the bisection method, where the value at the two bounding points must have opposite signs. Convergence is guaranteed for generally well-behaved problems on the interval. Derivatives are not required for this solver, though it is limited to solving a single state. Brent can be nested with other solvers, including other Brent solvers.
BrentSolver Options#
| Option | Default | Acceptable Values | Acceptable Types | Description |
|---|---|---|---|---|
| atol | 1e-10 | N/A | N/A | absolute error tolerance |
| debug_print | False | [True, False] | ['bool'] | If True, the values of input and output variables at the start of iteration are printed and written to a file after a failure to converge or when encountering aninvalid value in the residual. |
| err_on_non_converge | False | [True, False] | ['bool'] | When True, AnalysisError will be raised if we don't converge. |
| iprint | 1 | N/A | ['int'] | whether to print output |
| lower_bound | 0.0 | N/A | N/A | Lower bound for the root search |
| lower_bound_target | N/A | N/A | N/A | Openmdao path to the lower bound. When specified, this takes precedence over the value specified in lower_bound. |
| maxiter | 10 | N/A | ['int'] | maximum number of iterations |
| restart_from_successful | False | [True, False] | ['bool'] | If True, the states are cached after a successful solve and used to restart the solver in the case of a failed solve. |
| rtol | 1e-10 | N/A | N/A | relative error tolerance |
| state_target | N/A | N/A | ['str'] | Name of the implicit state to be solved |
| upper_bound | 100.0 | N/A | N/A | Upper bound for the root search |
| upper_bound_target | N/A | N/A | N/A | Openmdao path to the upper bound. When specified, this takes precedence over the value specified in upper_bound. |
BrentSolver Constructor#
The call signature for the BrentSolver constructor is:
- BrentSolver.__init__(**kwargs)[source]
Initialize all attributes.
BrentSolver Example#
The following simple example shows the use of the Brent solver to find the output of an ImplicitComponent that models the equation x = a*x**n + b*x - c.
import openmdao.api as om
class CompTest(om.ImplicitComponent):
def setup(self):
self.add_input('a', val=1.)
self.add_input('b', val=1.)
self.add_input('c', val=10.)
self.add_input('n', val=77.0/27.0)
self.add_output('x', val=2., lower=0, upper=100)
def apply_nonlinear(self, inputs, outputs, residuals, discrete_inputs=None, discrete_outputs=None):
a = inputs['a']
b = inputs['b']
c = inputs['c']
n = inputs['n']
x = outputs['x']
# Can't take fractional power of negative number
if x >= 0.0:
fact = x ** n
else:
fact = - (-x) ** n
residuals['x'] = a * fact + b * x - c
prob = om.Problem()
model = prob.model
model.add_subsystem('comp', CompTest(), promotes=['*'])
model.nonlinear_solver = om.BrentSolver(
state_target='x',
lower_bound=0.0,
upper_bound=80.0,
maxiter=100,
atol=1e-8,
rtol=1e-8,
)
prob.setup()
prob.set_solver_print(2)
prob.run_model()
NL: BRENT 0 ; 0.780735554 1
NL: BRENT 1 ; 10 12.8084342
NL: BRENT 2 ; 267574.261 342720.732
NL: BRENT 3 ; 9.99701022 12.8046048
NL: BRENT 4 ; 710.851762 910.489806
NL: BRENT 5 ; 9.85458137 12.6221757
NL: BRENT 6 ; 97.5661059 124.966905
NL: BRENT 7 ; 9.18003596 11.7581887
NL: BRENT 8 ; 12.3023868 15.7574312
NL: BRENT 9 ; 4.95794295 6.35034862
NL: BRENT 10 ; 3.79808278 4.86474935
NL: BRENT 11 ; 0.739762487 0.947519917
NL: BRENT 12 ; 0.0830067028 0.106318589
NL: BRENT 13 ; 0.000386953714 0.00049562712
NL: BRENT 14 ; 1.10893672e-06 1.42037431e-06
NL: BRENT 15 ; 1.47455381e-11 1.88867255e-11
NL: BRENT 16 ; 1.83156308e-07 2.34594553e-07
NL: BRENT 17 ; 1.47455381e-11 1.88867255e-11
NL: BRENT Converged
print(prob.get_val('x'))
[2.06720359]
The convergence history clearly shows the bisection and bracketing process, where the model is evaluated at closer points to the eventual solution.
BrentSolver Option Example: Upper and Lower Bounds from the Model#
The BrentSolver also allows the lower and upper bounds to be sourced from the model output, which is useful in cases where the interval of interest can be determined ahead of time by an upstream calculation. Note that the bounds are only queried at the start of the Brent iteration.
import openmdao.api as om
prob = om.Problem()
model = prob.model
model.add_subsystem('lower', om.ExecComp('low = 2*a'), promotes=['*'])
model.add_subsystem('upper', om.ExecComp('high = 2*b'), promotes=['*'])
model.add_subsystem('comp', CompTest(), promotes=['*'])
model.nonlinear_solver = om.BrentSolver(
state_target='x',
maxiter=100,
atol=1e-8,
rtol=1e-8,
lower_bound_target='low',
upper_bound_target='high',
)
prob.setup()
prob.set_solver_print(2)
prob.setup()
prob.set_val('a', -5.0)
prob.set_val('b', 55.0)
prob.run_model()
NL: BRENT 0 ; 63.9036778 1
NL: BRENT 1 ; 2994.85472 46.8651386
NL: BRENT 2 ; 3310761.75 51808.6261
NL: BRENT 3 ; 2891.97277 45.2551851
NL: BRENT 4 ; 819.781457 12.8283924
NL: BRENT 5 ; 379685.178 5941.52311
NL: BRENT 6 ; 764.522075 11.9636631
NL: BRENT 7 ; 203.336544 3.18192241
NL: BRENT 8 ; 38387.5788 600.710008
NL: BRENT 9 ; 171.104403 2.67753608
NL: BRENT 10 ; 36.6139392 0.572955117
NL: BRENT 11 ; 2669.37076 41.7717861
NL: BRENT 12 ; 13.1605549 0.205943623
NL: BRENT 13 ; 1.31446198 0.0205694262
NL: BRENT 14 ; 0.0566167805 0.000885970612
NL: BRENT 15 ; 0.000262898643 4.11398299e-06
NL: BRENT 16 ; 5.30071702e-08 8.29485439e-10
NL: BRENT 17 ; 2.54417273e-06 3.98126183e-08
NL: BRENT 18 ; 5.30071702e-08 8.29485439e-10
NL: BRENT Converged
print(prob.get_val('x'))
[-3.74515373]