CompositePrimalVector

class kona.linalg.vectors.composite.CompositePrimalVector(primal_vec, dual_ineq)[source]

Bases: kona.linalg.vectors.composite.CompositeVector

A composite vector representing a combined design and slack vectors..

Parameters:
  • _memory (KonaMemory) – All-knowing Kona memory manager.
  • design (DesignVector) – Design component of the composite vector.
  • slack (DualVectorINEQ) – Slack components of the composite vector.
convert_to_dual(dual_vector)[source]
equals_init_design()[source]
equals_lagrangian_total_gradient(at_primal, at_state, at_dual, at_adjoint, obj_scale=1.0, cnstr_scale=1.0)[source]

Computes the total primal derivative of the Lagrangian.

In this case, the primal derivative includes the slack derivative.

\[\begin{split}\nabla_{primal} \mathcal{L} = \begin{bmatrix} \nabla_x f(x, u) + \nabla_x c_{eq}(x, u)^T \lambda_{eq} + \nabla_x c_{inq}(x, u)^T \lambda_{ineq} \\ \muS^{-1}e - \lambda_{ineq} \end{bmatrix}\end{split}\]
Parameters:
  • at_primal (CompositePrimalVector) – The design/slack vector at which the derivative is computed.
  • at_state (StateVector) – State variables at which the derivative is computed.
  • at_dual (DualVector) – Lagrange multipliers at which the derivative is computed.
  • at_adjoint (StateVector) – Pre-computed adjoint variables for the Lagrangian.
  • obj_scale (float, optional) – Scaling for the objective function.
  • cnstr_scale (float, optional) – Scaling for the constraints.
init_slack = 1.0
restrict_to_design()[source]
restrict_to_target()[source]

CompositeDualVector

class kona.linalg.vectors.composite.CompositeDualVector(dual_eq, dual_ineq)[source]

Bases: kona.linalg.vectors.composite.CompositeVector

A composite vector representing a combined equality and inequality constraints.

Parameters:
convert_to_design(primal_vector)[source]
equals_constraints(at_primal, at_state, scale=1.0)[source]

Evaluate equality and inequality constraints in-place.

Parameters:
restrict_to_idf()[source]
restrict_to_regular()[source]

PrimalDualVector

class kona.linalg.vectors.composite.PrimalDualVector(primal_vec, eq_vec=None, ineq_vec=None)[source]

Bases: kona.linalg.vectors.composite.CompositeVector

A composite vector made up of primal, dual equality, and dual inequality vectors. :param _memory: All-knowing Kona memory manager. :type _memory: KonaMemory :param primal: Primal component of the composite vector. :type primal: DesignVector :param eq: Dual component corresponding to the equality constraints. :type eq: DualVectorEQ :param ineq: Dual component corresponding to the inequality constraints. :type ineq: DualVectorINEQ

equals_KKT_conditions(x, state, adjoint, obj_scale=1.0, cnstr_scale=1.0)[source]

Calculates the total derivative of the Lagrangian \(\mathcal{L}(x, u) = f(x, u)+ \lambda_{h}^T h(x, u) + \lambda_{g}^T g(x, u)\) with respect to \(\begin{pmatrix}x && \lambda_{h} && \lambda_{g}\end{pmatrix}^T\), where \(h\) denotes the equality constraints (if any) and \(g\) denotes the inequality constraints (if any). Note that these (total) derivatives do not represent the complete set of first-order optimality conditions in the case of inequality constraints. :param x: Evaluate derivatives at this primal-dual point. :type x: PrimalDualVector :param state: Evaluate derivatives at this state point. :type state: StateVector :param adjoint: Evaluate derivatives using this adjoint vector. :type adjoint: StateVector :param obj_scale: Scaling for the objective function. :type obj_scale: float, optional :param cnstr_scale: Scaling for the constraints. :type cnstr_scale: float, optional

equals_homotopy_residual(dLdx, x, init, mu=1.0)[source]

Using dLdx=:math:begin{pmatrix} nabla_x L && h && g end{pmatrix}, which can be obtained from the method equals_KKT_conditions, as well as the initial values init=:math:begin{pmatrix} x_0 && h(x_0,u_0) && g(x_0,u_0) end{pmatrix} and the current point x=:math:begin{pmatrix} x && lambda_h && lambda_g end{pmatrix}, this method computes the following nonlinear vector function: .. math:: r(x,lambda_h,lambda_g;mu) = begin{bmatrix} muleft[nabla_x f(x, u) - nabla_x h(x, u)^T lambda_{h} - nabla_x g(x, u)^T lambda_{g}right] + (1 - mu)(x - x_0) \ -mu h(x,u) - (1-mu)lambda_h \ -|g(x,u) - (1-mu)*g_0 - lambda_g|^3 + (g(x,u) - (1-mu)g_0)^3 + lambda_g^3 - (1-mu)hat{g} end{bmatrix} where \(h(x,u)\) are the equality constraints, and \(g(x,u)\) are the inequality constraints. The vectors \(\lambda_h\) and \(\lambda_g\) are the associated Lagrange multipliers. When mu=1.0, we recover a set of nonlinear algebraic equations equivalent to the first-order optimality conditions. :param dLdx: The total derivative of the Lagranginan with respect to the primal and dual variables. :type dLdx: PrimalDualVector :param x: The current solution vector value corresponding to dLdx. :type x: PrimalDualVector :param init: The initial primal variable, as well as the initial constraint values. :type init: PrimalDualVector :param mu: Homotopy parameter; must be between 0 and 1. :type mu: float

equals_init_guess()[source]

Sets the primal-dual vector to the initial guess, using the initial design.

equals_predictor_rhs(dLdx, x, init, mu=1.0)[source]

Using dLdx=:math:begin{pmatrix} nabla_x L && h && g end{pmatrix}, which can be obtained from the method equals_KKT_conditions, as well as the initial values init=:math:begin{pmatrix} x_0 && h(x_0,u_0) && g(x_0,u_0) end{pmatrix} and the current point x=:math:begin{pmatrix} x && lambda_h && lambda_g end{pmatrix}, this method computes the right-hand-side for the homotopy-predictor step, that is .. math:: partial r/partial mu = begin{bmatrix} left[nabla_x f(x, u) - nabla_x h(x, u)^T lambda_{h} - nabla_x g(x, u)^T lambda_{g}right] - (x - x_0) \ -h(x,u) + lambda_h \ -3*(g_0)|g(x,u) - (1-mu)*g_0 - lambda_g|^2 + 3*g_0*(g(x,u) - (1-mu)g_0)^2 + hat{g} end{bmatrix} where \(h(x,u)\) are the equality constraints, and \(g(x,u)\) are the inequality constraints. The vectors \(\lambda_h\) and \(\lambda_g\) are the associated Lagrange multipliers. :param dLdx: The total derivative of the Lagranginan with respect to the primal and dual variables. :type dLdx: PrimalDualVector :param x: The current solution vector value corresponding to dLdx. :type x: PrimalDualVector :param init: The initial primal variable, as well as the initial constraint values. :type init: PrimalDualVector :param mu: Homotopy parameter; must be between 0 and 1. :type mu: float

get_base_data(A)[source]

Inserts the PrimalDualVector’s underlying data into the given array :param A: Array into which data is inserted. :type A: numpy array

get_dual()[source]

Returns 1) eq or ineq if only one is None, 2) a CompositeDualVector if neither is none or 3) None if both are none

get_num_var()[source]

Returns the total number of variables in the PrimalDualVector

get_optimality_and_feasiblity()[source]

Returns the norm of the primal (opt) the dual parts of the vector (feas). If the dual parts of the vector are both None, then feas is returned as zero.

init_dual = 0.0
set_base_data(A)[source]

Copies the given array into the PrimalDualVector’s underlying data :param A: Array that is copied into the PrimalDualVector. :type A: numpy array

ReducedKKTVector

class kona.linalg.vectors.composite.ReducedKKTVector(primal_vec, dual_vec)[source]

Bases: kona.linalg.vectors.composite.CompositeVector

A composite vector representing a combined primal and dual vectors.

Parameters:
  • _memory (KonaMemory) – All-knowing Kona memory manager.
  • _primal (DesignVector or CompositePrimalVector) – Primal component of the composite vector.
  • _dual (DualVector) – Dual components of the composite vector.
equals_KKT_conditions(x, state, adjoint, barrier=None, obj_scale=1.0, cnstr_scale=1.0)[source]

Calculates the total derivative of the Lagrangian \(\mathcal{L}(x, u) = f(x, u)+ \lambda_{eq}^T c_{eq}(x, u) + \lambda_{ineq}^T (c_{ineq}(x, u) - s)\) with respect to \(\begin{pmatrix}x && s && \lambda_{eq} && \lambda_{ineq}\end{pmatrix}^T\). This total derivative represents the Karush-Kuhn-Tucker (KKT) convergence conditions for the optimization problem defined by \(\mathcal{L}(x, s, \lambda_{eq}, \lambda_{ineq})\) where the stat variables \(u(x)\) are treated as implicit functions of the design.

The full expression of the KKT conditions are:

\[\begin{split}\nabla \mathcal{L} = \begin{bmatrix} \nabla_x f(x, u) + \nabla_x c_{eq}(x, u)^T \lambda_{eq} + \nabla_x c_{inq}(x, u)^T \lambda_{ineq} \\ \muS^{-1}e - \lambda_{ineq} \\ c_{eq}(x, u) \\ c_{ineq}(x, u) - s \end{bmatrix}\end{split}\]
Parameters:
  • x (ReducedKKTVector) – Evaluate KKT conditions at this primal-dual point.
  • state (StateVector) – Evaluate KKT conditions at this state point.
  • adjoint (StateVector) – Evaluate KKT conditions using this adjoint vector.
  • barrier (float, optional) – Log barrier coefficient for slack variable non-negativity.
  • obj_scale (float, optional) – Scaling for the objective function.
  • cnstr_scale (float, optional) – Scaling for the constraints.
equals_init_guess()[source]

Sets the KKT vector to the initial guess, using the initial design.

init_dual = 0.0