MinimizationReporter

class openmm.openmm.MinimizationReporter(*args)

A MinimizationReporter can be passed to LocalEnergyMinimizer::minimize() to provide periodic information on the progress of minimization, and to give you the chance to stop minimization early. Define a subclass that overrides report() and implement it to take whatever action you want.

To correctly interpret the information passed to the reporter, you need to know a bit about how the minimizer works. The L-BFGS algorithm used by the minimizer does not support constraints. The minimizer therefore replaces all constraints with harmonic restraints, then performs unconstrained minimization of a combined objective function that is the sum of the system’s potential energy and the restraint energy. Once minimization completes, it checks whether all constraints are satisfied to an acceptable tolerance. It not, it increases the strength of the harmonic restraints and performs additional minimization. If the error in constrained distances is especially large, it may choose to throw out all work that has been done so far and start over with stronger restraints. This has several important consequences.

  • The objective function being minimized not actually the same as the potential energy.

  • The objective function and the potential energy can both increase between iterations.

  • The total number of iterations performed could be larger than the number specified by the maxIterations argument, if that many iterations leaves unacceptable constraint errors.

  • All work is provisional. It is possible for the minimizer to throw it out and start over.

__init__(self)MinimizationReporter
__init__(self, other)MinimizationReporter

A MinimizationReporter can be passed to LocalEnergyMinimizer::minimize() to provide periodic information on the progress of minimization, and to give you the chance to stop minimization early. Define a subclass that overrides report() and implement it to take whatever action you want.

To correctly interpret the information passed to the reporter, you need to know a bit about how the minimizer works. The L-BFGS algorithm used by the minimizer does not support constraints. The minimizer therefore replaces all constraints with harmonic restraints, then performs unconstrained minimization of a combined objective function that is the sum of the system’s potential energy and the restraint energy. Once minimization completes, it checks whether all constraints are satisfied to an acceptable tolerance. It not, it increases the strength of the harmonic restraints and performs additional minimization. If the error in constrained distances is especially large, it may choose to throw out all work that has been done so far and start over with stronger restraints. This has several important consequences.

  • The objective function being minimized not actually the same as the potential energy.

  • The objective function and the potential energy can both increase between iterations.

  • The total number of iterations performed could be larger than the number specified by the maxIterations argument, if that many iterations leaves unacceptable constraint errors.

  • All work is provisional. It is possible for the minimizer to throw it out and start over.

Methods

__init__(-> MinimizationReporter)

A MinimizationReporter can be passed to LocalEnergyMinimizer::minimize() to provide periodic information on the progress of minimization, and to give you the chance to stop minimization early.

report(self, iteration, x, grad)

This is called after each iteration to provide information about the current status of minimization.

Attributes

thisown

The membership flag

property thisown

The membership flag

report(self, iteration, x, grad)bool

This is called after each iteration to provide information about the current status of minimization. It receives the current particle coordinates, the gradient of the objective function with respect to them, and a set of useful statistics. In particular, args contains these values:

“system energy”: the current potential energy of the system

“restraint energy”: the energy of the harmonic restraints

“restraint strength”: the force constant of the restraints (in kJ/mol/nm^2)

“max constraint error”: the maximum relative error in the length of any constraint

If this function returns true, it will cause the L-BFGS optimizer to immediately exit. If all constrained distances are sufficiently close to their target values, minimize() will return. If any constraint error is unacceptably large, it will instead cause the minimizer to immediately increase the strength of the harmonic restraints and perform additional optimization.

Parameters
  • iteration (int) – the index of the current iteration. This refers to the current call to the L-BFGS optimizer. Each time the minimizer increases the restraint strength, the iteration index is reset to 0.

  • x (vector< double >) – the current particle positions in flattened order: the three coordinates of the first particle, then the three coordinates of the second particle, etc.

  • grad (vector< double >) – the current gradient of the objective function (potential energy plus restraint energy) with respect to the particle coordinates, in flattened order

  • args (map< std::string, double >) – additional statistics described above about the current state of minimization

Returns

whether to immediately stop minimization

Return type

bool