2.8. Optimization

2.8.1. Site constraints

2.8.2. Minimum distance constraints

class opentidalfarm.optimisation_helpers.ConvexPolygonSiteConstraint(farm, vertices)[source]

Bases: InequalityConstraint

Generates the inequality constraints for generic polygon constraints. The parameter polygon must be a list of point coordinates that describes the site edges in anti-clockwise order.

function(m)[source]
jacobian(m)[source]
length()[source]
output_workspace()[source]
class opentidalfarm.optimisation_helpers.DomainRestrictionConstraints(config, feasible_area, attraction_center)[source]

Bases: InequalityConstraint

function(m)[source]
jacobian(m)[source]
length()[source]
class opentidalfarm.optimisation_helpers.MinimumDistanceConstraints(turbine_positions, minimum_distance, controls)[source]

Bases: InequalityConstraint

This class implements minimum distance constraints between turbines.

Note

This class subclasses dolfin_adjoint.InequalityConstraint. The following method names must not change:

  • length(self)

  • function(self, m)

  • jacobian(self, m)

function(m)[source]

Return an object which must be zero for the point to be feasible.

Parameters:

m – The serialized paramaterisation of the turbines.

Tpye m:

numpy.ndarray.

Returns:

numpy.ndarray – each entry must be zero for a poinst to be feasible.

jacobian(m)[source]

Returns the gradient of the constraint function.

Return a list of vector-like objects representing the gradient of the constraint function with respect to the parameter m.

Parameters:

m – The serialized paramaterisation of the turbines.

Tpye m:

numpy.ndarray.

Returns:

numpy.ndarray – the gradient of the constraint function with respect to each input parameter m.

length()[source]

Returns the number of constraints len(function(m)).

class opentidalfarm.optimisation_helpers.MinimumDistanceConstraintsLargeArrays(turbine_positions, minimum_distance, controls)[source]

Bases: InequalityConstraint

This class implements minimum distance constraints between turbines which is well suited for large arrays.

It enforces the minimum distance with only one enquality constraint of the form:

\[\sum_{i,j=0, i>j}^{N-1} \min(||p_i - p_j||^2 - D^2, 0) \ge 0\]

where N is the number of turbines, \(p_i\) is the position of the i’th turbine and D is the minimum distance between two turbines.

Note

This class subclasses dolfin_adjoint.InequalityConstraint. The following method names must not change:

  • length(self)

  • function(self, m)

  • jacobian(self, m)

function(m)[source]

Return an object which must be >0 for the point to be feasible.

Parameters:

m – The serialized paramaterisation of the turbines.

Tpye m:

numpy.ndarray.

Returns:

numpy.ndarray – each entry must be zero for a poinst to be feasible.

jacobian(m)[source]

Returns the gradient of the constraint function.

Return a list of vector-like objects representing the gradient of the constraint function with respect to the parameter m.

Parameters:

m – The serialized paramaterisation of the turbines.

Tpye m:

numpy.ndarray.

Returns:

numpy.ndarray – the gradient of the constraint function with respect to each input parameter m.

length()[source]

Returns the number of constraints len(function(m)).

opentidalfarm.optimisation_helpers.friction_constraints(config, lb=0.0, ub=None)[source]

This function returns the constraints to ensure that the turbine friction controls remain reasonable.

opentidalfarm.optimisation_helpers.get_distance_function(config, domains)[source]
opentidalfarm.optimisation_helpers.get_domain_constraints(config, feasible_area, attraction_center)[source]
opentidalfarm.optimisation_helpers.position_constraints(config)[source]

This function returns the constraints to ensure that the turbine positions remain inside the domain.