Ipopt Documentation  
 
Loading...
Searching...
No Matches
Ipopt::TNLPAdapter Class Reference

This class adapts the TNLP interface so it looks like an NLP interface. More...

#include <IpTNLPAdapter.hpp>

+ Inheritance diagram for Ipopt::TNLPAdapter:

Public Types

enum  FixedVariableTreatmentEnum { MAKE_PARAMETER = 0 , MAKE_PARAMETER_NODUAL , MAKE_CONSTRAINT , RELAX_BOUNDS }
 Enum for treatment of fixed variables option. More...
 
enum  DerivativeTestEnum { NO_TEST = 0 , FIRST_ORDER_TEST , SECOND_ORDER_TEST , ONLY_SECOND_ORDER_TEST }
 Enum for specifying which derivative test is to be performed. More...
 
enum  JacobianApproxEnum { JAC_EXACT = 0 , JAC_FINDIFF_VALUES }
 Enum for specifying technique for computing Jacobian. More...
 
enum  GradientApproxEnum { OBJGRAD_EXACT = 0 , OBJGRAD_FINDIFF_VALUES }
 Enum for specifying technique for computing objective Gradient. More...
 

Public Member Functions

virtual void GetQuasiNewtonApproximationSpaces (SmartPtr< VectorSpace > &approx_space, SmartPtr< Matrix > &P_approx)
 Method returning information on quasi-Newton approximation.
 
bool CheckDerivatives (DerivativeTestEnum deriv_test, Index deriv_test_start_index)
 Method for performing the derivative test.
 
SmartPtr< TNLPtnlp () const
 Accessor method for the underlying TNLP.
 
Constructors/Destructors
 TNLPAdapter (const SmartPtr< TNLP > tnlp, const SmartPtr< const Journalist > jnlst=NULL)
 Default constructor.
 
virtual ~TNLPAdapter ()
 Default destructor.
 
Exceptions
 DECLARE_STD_EXCEPTION (INVALID_TNLP)
 
 DECLARE_STD_EXCEPTION (ERROR_IN_TNLP_DERIVATIVE_TEST)
 
TNLPAdapter Initialization.
virtual bool ProcessOptions (const OptionsList &options, const std::string &prefix)
 Overload if you want the chance to process options or parameters that may be specific to the NLP.
 
virtual bool GetSpaces (SmartPtr< const VectorSpace > &x_space, SmartPtr< const VectorSpace > &c_space, SmartPtr< const VectorSpace > &d_space, SmartPtr< const VectorSpace > &x_l_space, SmartPtr< const MatrixSpace > &px_l_space, SmartPtr< const VectorSpace > &x_u_space, SmartPtr< const MatrixSpace > &px_u_space, SmartPtr< const VectorSpace > &d_l_space, SmartPtr< const MatrixSpace > &pd_l_space, SmartPtr< const VectorSpace > &d_u_space, SmartPtr< const MatrixSpace > &pd_u_space, SmartPtr< const MatrixSpace > &Jac_c_space, SmartPtr< const MatrixSpace > &Jac_d_space, SmartPtr< const SymMatrixSpace > &Hess_lagrangian_space)
 Method for creating the derived vector / matrix types.
 
virtual bool GetBoundsInformation (const Matrix &Px_L, Vector &x_L, const Matrix &Px_U, Vector &x_U, const Matrix &Pd_L, Vector &d_L, const Matrix &Pd_U, Vector &d_U)
 Method for obtaining the bounds information.
 
virtual bool GetStartingPoint (SmartPtr< Vector > x, bool need_x, SmartPtr< Vector > y_c, bool need_y_c, SmartPtr< Vector > y_d, bool need_y_d, SmartPtr< Vector > z_L, bool need_z_L, SmartPtr< Vector > z_U, bool need_z_U)
 Method for obtaining the starting point for all the iterates.
 
virtual bool GetWarmStartIterate (IteratesVector &warm_start_iterate)
 Method for obtaining an entire iterate as a warmstart point.
 
TNLPAdapter evaluation routines.
virtual bool Eval_f (const Vector &x, Number &f)
 
virtual bool Eval_grad_f (const Vector &x, Vector &g_f)
 
virtual bool Eval_c (const Vector &x, Vector &c)
 
virtual bool Eval_jac_c (const Vector &x, Matrix &jac_c)
 
virtual bool Eval_d (const Vector &x, Vector &d)
 
virtual bool Eval_jac_d (const Vector &x, Matrix &jac_d)
 
virtual bool Eval_h (const Vector &x, Number obj_factor, const Vector &yc, const Vector &yd, SymMatrix &h)
 
virtual void GetScalingParameters (const SmartPtr< const VectorSpace > x_space, const SmartPtr< const VectorSpace > c_space, const SmartPtr< const VectorSpace > d_space, Number &obj_scaling, SmartPtr< Vector > &x_scaling, SmartPtr< Vector > &c_scaling, SmartPtr< Vector > &d_scaling) const
 Routines to get the scaling parameters.
 
Solution Reporting Methods
virtual void FinalizeSolution (SolverReturn status, const Vector &x, const Vector &z_L, const Vector &z_U, const Vector &c, const Vector &d, const Vector &y_c, const Vector &y_d, Number obj_value, const IpoptData *ip_data, IpoptCalculatedQuantities *ip_cq)
 This method is called at the very end of the optimization.
 
virtual bool IntermediateCallBack (AlgorithmMode mode, Index iter, Number obj_value, Number inf_pr, Number inf_du, Number mu, Number d_norm, Number regularization_size, Number alpha_du, Number alpha_pr, Index ls_trials, const IpoptData *ip_data, IpoptCalculatedQuantities *ip_cq)
 This method is called once per iteration, after the iteration summary output has been printed.
 
Methods for translating data for IpoptNLP into the TNLP data.

These methods are used to obtain the current (or final) data for the TNLP formulation from the IpoptNLP structure.

void ResortX (const Vector &x, Number *x_orig, bool usefixedvals=true)
 Sort the primal variables, and add the fixed values in x_orig.
 
void ResortG (const Vector &c, const Vector &d, Number *g_orig, bool correctrhs=false)
 Sort constraint values.
 
void ResortBounds (const Vector &x_L, Number *x_L_orig, const Vector &x_U, Number *x_U_orig)
 Provides values for lower and upper bounds on variables for given Ipopt-internal vectors.
 
IPOPT_DEPRECATED void ResortBnds (const Vector &x_L, Number *x_L_orig, const Vector &x_U, Number *x_U_orig, bool clearorig=true)
 Provides values for lower and upper bounds on variables for given Ipopt-internal vectors.
 
bool ResortBoundMultipliers (const Vector &x, const Vector &y_c, const Vector &y_d, const Vector &z_L, Number *z_L_orig, const Vector &z_U, Number *z_U_orig)
 Provides values for dual multipliers on lower and upper bounds on variables for given Ipopt-internal vectors.
 
void GetFullDimensions (Index &n, Index &m) const
 Get number of variables and number of constraints in TNLP.
 
void GetFixedVariables (Index &n_x_fixed, Index *&x_fixed_map, FixedVariableTreatmentEnum &fixed_variable_treatment) const
 Get number and indices of fixed variables.
 
void GetPermutationMatrices (SmartPtr< const ExpansionMatrix > &P_x_full_x, SmartPtr< const ExpansionMatrix > &P_x_x_L, SmartPtr< const ExpansionMatrix > &P_x_x_U, SmartPtr< const ExpansionMatrix > &P_c_g, SmartPtr< const ExpansionMatrix > &P_d_g) const
 Get mappings between TNLP indices and Ipopt internal indices for variables and constraints.
 
const NumberGetC_Rhs () const
 Get right-hand-sides that are added into c(x)
 
- Public Member Functions inherited from Ipopt::NLP
 NLP ()
 Default constructor.
 
virtual ~NLP ()
 Default destructor.
 
 DECLARE_STD_EXCEPTION (USER_SCALING_NOT_IMPLEMENTED)
 Exceptions.
 
 DECLARE_STD_EXCEPTION (INVALID_NLP)
 
- Public Member Functions inherited from Ipopt::ReferencedObject
 ReferencedObject ()
 
virtual ~ReferencedObject ()
 
Index ReferenceCount () const
 
void AddRef (const Referencer *referencer) const
 
void ReleaseRef (const Referencer *referencer) const
 

Static Public Member Functions

static void RegisterOptions (SmartPtr< RegisteredOptions > roptions)
 

Private Member Functions

Default Compiler Generated Methods

(Hidden to avoid implicit creation/calling).

These methods are not implemented and we do not want the compiler to implement them for us, so we declare them private and do not define them. This ensures that they will not be implicitly created/called.

 TNLPAdapter (const TNLPAdapter &)
 Copy Constructor.
 
void operator= (const TNLPAdapter &)
 Default Assignment Operator.
 
Methods to update the values in the local copies of vectors
bool update_local_x (const Vector &x)
 
bool update_local_lambda (const Vector &y_c, const Vector &y_d)
 
Internal routines for evaluating g and jac_g.

Values stored since they are used in both c and d routines.

bool internal_eval_g (bool new_x)
 
bool internal_eval_jac_g (bool new_x)
 
Internal methods for dealing with finite difference approximation
void initialize_findiff_jac (const Index *iRow, const Index *jCol)
 Initialize sparsity structure for finite difference Jacobian.
 

Private Attributes

TNLP::IndexStyleEnum index_style_
 Numbering style of variables and constraints.
 
Algorithmic parameters
Number nlp_lower_bound_inf_
 Value for a lower bound that denotes -infinity.
 
Number nlp_upper_bound_inf_
 Value for a upper bound that denotes infinity.
 
FixedVariableTreatmentEnum fixed_variable_treatment_
 Flag indicating how fixed variables should be handled.
 
Number bound_relax_factor_
 Determines relaxation of fixing bound for RELAX_BOUNDS.
 
DerivativeTestEnum derivative_test_
 Maximal slack for one-sidedly bounded variables.
 
Number derivative_test_perturbation_
 Size of the perturbation for the derivative test.
 
Number derivative_test_tol_
 Relative threshold for marking deviation from finite difference test.
 
bool derivative_test_print_all_
 Flag indicating if all test values should be printed, or only those violating the threshold.
 
Index derivative_test_first_index_
 Index of first quantity to be checked.
 
bool warm_start_same_structure_
 Flag indicating whether the TNLP with identical structure has already been solved before.
 
HessianApproximationType hessian_approximation_
 Flag indicating what Hessian information is to be used.
 
Index num_linear_variables_
 Number of linear variables.
 
JacobianApproxEnum jacobian_approximation_
 Flag indicating how Jacobian is computed.
 
GradientApproxEnum gradient_approximation_
 Flag indicating how objective Gradient is computed.
 
Number findiff_perturbation_
 Size of the perturbation for the derivative approximation.
 
Number point_perturbation_radius_
 Maximal perturbation of the initial point.
 
bool dependency_detection_with_rhs_
 Flag indicating if rhs should be considered during dependency detection.
 
Number tol_
 Overall convergence tolerance.
 
Problem Size Data
Index n_full_x_
 full dimension of x (fixed + non-fixed)
 
Index n_full_g_
 full dimension of g (c + d)
 
Index nz_jac_c_
 non-zeros of the jacobian of c
 
Index nz_jac_c_no_extra_
 non-zeros of the jacobian of c without added constraints for fixed variables.
 
Index nz_jac_d_
 non-zeros of the jacobian of d
 
Index nz_full_jac_g_
 number of non-zeros in full-size Jacobian of g
 
Index nz_full_h_
 number of non-zeros in full-size Hessian
 
Index nz_h_
 number of non-zeros in the non-fixed-size Hessian
 
Index n_x_fixed_
 Number of fixed variables.
 
Local copy of spaces (for warm start)
SmartPtr< const VectorSpacex_space_
 
SmartPtr< const VectorSpacec_space_
 
SmartPtr< const VectorSpaced_space_
 
SmartPtr< const VectorSpacex_l_space_
 
SmartPtr< const MatrixSpacepx_l_space_
 
SmartPtr< const VectorSpacex_u_space_
 
SmartPtr< const MatrixSpacepx_u_space_
 
SmartPtr< const VectorSpaced_l_space_
 
SmartPtr< const MatrixSpacepd_l_space_
 
SmartPtr< const VectorSpaced_u_space_
 
SmartPtr< const MatrixSpacepd_u_space_
 
SmartPtr< const MatrixSpaceJac_c_space_
 
SmartPtr< const MatrixSpaceJac_d_space_
 
SmartPtr< const SymMatrixSpaceHess_lagrangian_space_
 
Local Copy of the Data
Numberfull_x_
 
Numberfull_lambda_
 copy of the full x vector (fixed & non-fixed)
 
Numberfull_g_
 copy of lambda (yc & yd)
 
Numberjac_g_
 copy of g (c & d)
 
Numberc_rhs_
 the values for the full jacobian of g
 
Tags for deciding when to update internal copies of vectors

the rhs values of c

TaggedObject::Tag x_tag_for_iterates_
 
TaggedObject::Tag y_c_tag_for_iterates_
 
TaggedObject::Tag y_d_tag_for_iterates_
 
TaggedObject::Tag x_tag_for_g_
 
TaggedObject::Tag x_tag_for_jac_g_
 
Internal Permutation Spaces and matrices
SmartPtr< ExpansionMatrixP_x_full_x_
 Expansion from fixed x (ipopt) to full x.
 
SmartPtr< ExpansionMatrixSpaceP_x_full_x_space_
 
SmartPtr< ExpansionMatrixP_x_x_L_
 Expansion from fixed x_L (ipopt) to full x.
 
SmartPtr< ExpansionMatrixSpaceP_x_x_L_space_
 
SmartPtr< ExpansionMatrixP_x_x_U_
 Expansion from fixed x_U (ipopt) to full x.
 
SmartPtr< ExpansionMatrixSpaceP_x_x_U_space_
 
SmartPtr< ExpansionMatrixSpaceP_c_g_space_
 Expansion from c only (ipopt) to full ampl c.
 
SmartPtr< ExpansionMatrixP_c_g_
 
SmartPtr< ExpansionMatrixSpaceP_d_g_space_
 Expansion from d only (ipopt) to full ampl d.
 
SmartPtr< ExpansionMatrixP_d_g_
 
Indexjac_idx_map_
 
Indexh_idx_map_
 
Indexx_fixed_map_
 Position of fixed variables.
 
std::vector< Indexjac_fixed_idx_map_
 Index mapping of Jacobian w.r.t.
 
std::vector< Indexjac_fixed_iRow_
 
std::vector< Indexjac_fixed_jCol_
 
Data for finite difference approximations of derivatives
Index findiff_jac_nnz_
 Number of unique nonzeros in constraint Jacobian.
 
Indexfindiff_jac_ia_
 Start position for nonzero indices in ja for each column of Jacobian.
 
Indexfindiff_jac_ja_
 Ordered by columns, for each column the row indices in Jacobian.
 
Indexfindiff_jac_postriplet_
 Position of entry in original triplet matrix.
 
Numberfindiff_x_l_
 Copy of the lower bounds.
 
Numberfindiff_x_u_
 Copy of the upper bounds.
 

Method implementing the detection of linearly dependent equality constraints

SmartPtr< TNLPtnlp_
 Pointer to the TNLP class (class specific to Number* vectors and triplet matrices)
 
SmartPtr< const Journalistjnlst_
 Journalist.
 
SmartPtr< TDependencyDetectordependency_detector_
 Object that can be used to detect linearly dependent rows in the equality constraint Jacobian.
 
bool DetermineDependentConstraints (Index n_x_var, const Index *x_not_fixed_map, const Number *x_l, const Number *x_u, const Number *g_l, const Number *g_u, Index n_c, const Index *c_map, std::list< Index > &c_deps)
 

Detailed Description

This class adapts the TNLP interface so it looks like an NLP interface.

This is an Adapter class (Design Patterns) that converts a TNLP to an NLP. This allows users to write to the "more convenient" TNLP interface.

Given a TNLP

\begin{eqnarray*} \mathrm{min} && f(x), \\ \mathrm{s.t.} && g_L \leq g(x) \leq g_U, &\qquad \lambda\\ && x_L \leq x \leq x_U, &\qquad z_L, z_U \end{eqnarray*}

let \(E = \{i : g_{L,i} = g_{U,i}\}\) and \(I = \{i : g_{L,i} \neq g_{U,i}\}\) be the indices of equality and inequality constraints, respectively. The dual variables for the constraints are \(\lambda\). The dual variables for the variable bounds are \(z_L\) and \(z_U\).

A TNLPAdapter represents the problem

\begin{eqnarray*} \mathrm{min} && f(x), \\ \mathrm{s.t.} && c(x) = 0, &\qquad y_c\\ && d_L \leq d(x) \leq d_U, &\qquad y_d \\ && x_L \leq x \leq x_U, &\qquad z_L, z_U \end{eqnarray*}

where \(c(x) = g_E(x) - g_{L,E}\), i.e., corresponding to equality constraints of TNLP, and \(d(x) = g_I(x)\), \(d_L = g_{L,I}\), \(d_U = g_{U,I}\), i.e., corresponding to inequality constraints of TNLP.

The dual variables for the constraints are \(y_c\) and \(y_d\). The dual variables for the bounds on slack and original variables are \(s_L\), \(s_U\), \(z_L\), \(z_U\).

Internally, Ipopt reformulates \(d_L \leq d(x) \leq d_U\) as

\begin{eqnarray*} && d(x) - s = 0 &\qquad y_d\\ && d_L \leq s \leq d_U &\qquad v_L, v_U \\ \end{eqnarray*}

If fixed variables are present ( \(x_{L,i} = x_{U,i}\)) in the TNLP and fixed_variable_treatment is set to make_parameter, then these variables are not made visible to Ipopt internally. If fixed_variable_treatment is set to make_constraint, then their bounds relaxed and equality constraints \(x_i - x_{L,i} = 0\) are added to the end of \(c(x) = 0\).

Definition at line 65 of file IpTNLPAdapter.hpp.

Member Enumeration Documentation

◆ FixedVariableTreatmentEnum

Enum for treatment of fixed variables option.

Enumerator
MAKE_PARAMETER 
MAKE_PARAMETER_NODUAL 
Since
3.14.0
MAKE_CONSTRAINT 
RELAX_BOUNDS 

Definition at line 240 of file IpTNLPAdapter.hpp.

◆ DerivativeTestEnum

Enum for specifying which derivative test is to be performed.

Enumerator
NO_TEST 
FIRST_ORDER_TEST 
SECOND_ORDER_TEST 
ONLY_SECOND_ORDER_TEST 

Definition at line 249 of file IpTNLPAdapter.hpp.

◆ JacobianApproxEnum

Enum for specifying technique for computing Jacobian.

Enumerator
JAC_EXACT 
JAC_FINDIFF_VALUES 

Definition at line 258 of file IpTNLPAdapter.hpp.

◆ GradientApproxEnum

Enum for specifying technique for computing objective Gradient.

Enumerator
OBJGRAD_EXACT 
OBJGRAD_FINDIFF_VALUES 

Definition at line 265 of file IpTNLPAdapter.hpp.

Constructor & Destructor Documentation

◆ TNLPAdapter() [1/2]

Ipopt::TNLPAdapter::TNLPAdapter ( const SmartPtr< TNLP tnlp,
const SmartPtr< const Journalist jnlst = NULL 
)

Default constructor.

◆ ~TNLPAdapter()

virtual Ipopt::TNLPAdapter::~TNLPAdapter ( )
virtual

Default destructor.

◆ TNLPAdapter() [2/2]

Ipopt::TNLPAdapter::TNLPAdapter ( const TNLPAdapter )
private

Copy Constructor.

Member Function Documentation

◆ DECLARE_STD_EXCEPTION() [1/2]

Ipopt::TNLPAdapter::DECLARE_STD_EXCEPTION ( INVALID_TNLP  )

◆ DECLARE_STD_EXCEPTION() [2/2]

Ipopt::TNLPAdapter::DECLARE_STD_EXCEPTION ( ERROR_IN_TNLP_DERIVATIVE_TEST  )

◆ ProcessOptions()

virtual bool Ipopt::TNLPAdapter::ProcessOptions ( const OptionsList ,
const std::string &   
)
virtual

Overload if you want the chance to process options or parameters that may be specific to the NLP.

Reimplemented from Ipopt::NLP.

◆ GetSpaces()

virtual bool Ipopt::TNLPAdapter::GetSpaces ( SmartPtr< const VectorSpace > &  x_space,
SmartPtr< const VectorSpace > &  c_space,
SmartPtr< const VectorSpace > &  d_space,
SmartPtr< const VectorSpace > &  x_l_space,
SmartPtr< const MatrixSpace > &  px_l_space,
SmartPtr< const VectorSpace > &  x_u_space,
SmartPtr< const MatrixSpace > &  px_u_space,
SmartPtr< const VectorSpace > &  d_l_space,
SmartPtr< const MatrixSpace > &  pd_l_space,
SmartPtr< const VectorSpace > &  d_u_space,
SmartPtr< const MatrixSpace > &  pd_u_space,
SmartPtr< const MatrixSpace > &  Jac_c_space,
SmartPtr< const MatrixSpace > &  Jac_d_space,
SmartPtr< const SymMatrixSpace > &  Hess_lagrangian_space 
)
virtual

Method for creating the derived vector / matrix types.

Do not delete these.

Implements Ipopt::NLP.

◆ GetBoundsInformation()

virtual bool Ipopt::TNLPAdapter::GetBoundsInformation ( const Matrix Px_L,
Vector x_L,
const Matrix Px_U,
Vector x_U,
const Matrix Pd_L,
Vector d_L,
const Matrix Pd_U,
Vector d_U 
)
virtual

Method for obtaining the bounds information.

Implements Ipopt::NLP.

◆ GetStartingPoint()

virtual bool Ipopt::TNLPAdapter::GetStartingPoint ( SmartPtr< Vector x,
bool  need_x,
SmartPtr< Vector y_c,
bool  need_y_c,
SmartPtr< Vector y_d,
bool  need_y_d,
SmartPtr< Vector z_L,
bool  need_z_L,
SmartPtr< Vector z_U,
bool  need_z_U 
)
virtual

Method for obtaining the starting point for all the iterates.

Implements Ipopt::NLP.

◆ GetWarmStartIterate()

virtual bool Ipopt::TNLPAdapter::GetWarmStartIterate ( IteratesVector warm_start_iterate)
virtual

Method for obtaining an entire iterate as a warmstart point.

The incoming IteratesVector has to be filled.

Reimplemented from Ipopt::NLP.

◆ Eval_f()

virtual bool Ipopt::TNLPAdapter::Eval_f ( const Vector x,
Number f 
)
virtual

Implements Ipopt::NLP.

◆ Eval_grad_f()

virtual bool Ipopt::TNLPAdapter::Eval_grad_f ( const Vector x,
Vector g_f 
)
virtual

Implements Ipopt::NLP.

◆ Eval_c()

virtual bool Ipopt::TNLPAdapter::Eval_c ( const Vector x,
Vector c 
)
virtual

Implements Ipopt::NLP.

◆ Eval_jac_c()

virtual bool Ipopt::TNLPAdapter::Eval_jac_c ( const Vector x,
Matrix jac_c 
)
virtual

Implements Ipopt::NLP.

◆ Eval_d()

virtual bool Ipopt::TNLPAdapter::Eval_d ( const Vector x,
Vector d 
)
virtual

Implements Ipopt::NLP.

◆ Eval_jac_d()

virtual bool Ipopt::TNLPAdapter::Eval_jac_d ( const Vector x,
Matrix jac_d 
)
virtual

Implements Ipopt::NLP.

◆ Eval_h()

virtual bool Ipopt::TNLPAdapter::Eval_h ( const Vector x,
Number  obj_factor,
const Vector yc,
const Vector yd,
SymMatrix h 
)
virtual

Implements Ipopt::NLP.

◆ GetScalingParameters()

virtual void Ipopt::TNLPAdapter::GetScalingParameters ( const SmartPtr< const VectorSpace ,
const SmartPtr< const VectorSpace ,
const SmartPtr< const VectorSpace ,
Number ,
SmartPtr< Vector > &  ,
SmartPtr< Vector > &  ,
SmartPtr< Vector > &   
) const
virtual

Routines to get the scaling parameters.

These do not need to be overloaded unless the options are set for user scaling.

Reimplemented from Ipopt::NLP.

◆ FinalizeSolution()

virtual void Ipopt::TNLPAdapter::FinalizeSolution ( SolverReturn  ,
const Vector ,
const Vector ,
const Vector ,
const Vector ,
const Vector ,
const Vector ,
const Vector ,
Number  ,
const IpoptData ,
IpoptCalculatedQuantities  
)
virtual

This method is called at the very end of the optimization.

It provides the final iterate to the user, so that it can be stored as the solution. The status flag indicates the outcome of the optimization, where SolverReturn is defined in IpAlgTypes.hpp.

Reimplemented from Ipopt::NLP.

◆ IntermediateCallBack()

virtual bool Ipopt::TNLPAdapter::IntermediateCallBack ( AlgorithmMode  ,
Index  ,
Number  ,
Number  ,
Number  ,
Number  ,
Number  ,
Number  ,
Number  ,
Number  ,
Index  ,
const IpoptData ,
IpoptCalculatedQuantities  
)
virtual

This method is called once per iteration, after the iteration summary output has been printed.

It provides the current information to the user to do with it anything she wants. It also allows the user to ask for a premature termination of the optimization by returning false, in which case Ipopt will terminate with a corresponding return status. The basic information provided in the argument list has the quantities values printed in the iteration summary line. If more information is required, a user can obtain it from the IpData and IpCalculatedQuantities objects. However, note that the provided quantities are all for the problem that Ipopt sees, i.e., the quantities might be scaled, fixed variables might be sorted out, etc. The status indicates things like whether the algorithm is in the restoration phase... In the restoration phase, the dual variables are probably not not changing.

Reimplemented from Ipopt::NLP.

◆ GetQuasiNewtonApproximationSpaces()

virtual void Ipopt::TNLPAdapter::GetQuasiNewtonApproximationSpaces ( SmartPtr< VectorSpace > &  approx_space,
SmartPtr< Matrix > &  P_approx 
)
virtual

Method returning information on quasi-Newton approximation.

Reimplemented from Ipopt::NLP.

◆ CheckDerivatives()

bool Ipopt::TNLPAdapter::CheckDerivatives ( DerivativeTestEnum  deriv_test,
Index  deriv_test_start_index 
)

Method for performing the derivative test.

◆ RegisterOptions()

static void Ipopt::TNLPAdapter::RegisterOptions ( SmartPtr< RegisteredOptions roptions)
static

◆ tnlp()

SmartPtr< TNLP > Ipopt::TNLPAdapter::tnlp ( ) const
inline

Accessor method for the underlying TNLP.

Definition at line 282 of file IpTNLPAdapter.hpp.

◆ ResortX()

void Ipopt::TNLPAdapter::ResortX ( const Vector x,
Number x_orig,
bool  usefixedvals = true 
)

Sort the primal variables, and add the fixed values in x_orig.

Parameters
xinternal values for primal variables x
x_origvector to fill with values from x
usefixedvalswhether to use stored variable fixings for fixed variables (true), or zero (false)
Since
3.14.0

◆ ResortG()

void Ipopt::TNLPAdapter::ResortG ( const Vector c,
const Vector d,
Number g_orig,
bool  correctrhs = false 
)

Sort constraint values.

To Ipopt, the equality constraints are presented with right hand side zero. Specifying correctrhs=true corrects for the original right hand side.

Parameters
cinternal activity for equality constraints
dinternal activity for inequality constraints
g_origvector to fill with values from c and d
correctrhswhether to add unscaled rhs-values for constraints that internally correspond to c(x)=0
Since
3.14.0

◆ ResortBounds()

void Ipopt::TNLPAdapter::ResortBounds ( const Vector x_L,
Number x_L_orig,
const Vector x_U,
Number x_U_orig 
)

Provides values for lower and upper bounds on variables for given Ipopt-internal vectors.

Similar to ResortX, but does so for two arrays and does not set any values for fixed variables.

Since
3.14.0
Parameters
x_Linternal values for lower bounds on x
x_L_origvector to fill with values from x_L
x_Uinternal values for upper bounds on x
x_U_origvector to fill with values from x_U

◆ ResortBnds()

IPOPT_DEPRECATED void Ipopt::TNLPAdapter::ResortBnds ( const Vector x_L,
Number x_L_orig,
const Vector x_U,
Number x_U_orig,
bool  clearorig = true 
)
inline

Provides values for lower and upper bounds on variables for given Ipopt-internal vectors.

Similar to ResortX, but does so for two arrays and does not set any values for fixed variables.

Deprecated:
Use ResortBounds() instead.
Parameters
x_Linternal values for lower bounds on x
x_L_origvector to fill with values from x_L
x_Uinternal values for upper bounds on x
x_U_origvector to fill with values from x_U
clearorigwhether to initialize complete x_L_orig and x_U_orig to 0.0 before setting values for non-fixed variables - ignored

Definition at line 331 of file IpTNLPAdapter.hpp.

◆ ResortBoundMultipliers()

bool Ipopt::TNLPAdapter::ResortBoundMultipliers ( const Vector x,
const Vector y_c,
const Vector y_d,
const Vector z_L,
Number z_L_orig,
const Vector z_U,
Number z_U_orig 
)

Provides values for dual multipliers on lower and upper bounds on variables for given Ipopt-internal vectors.

Similar to ResortBounds, but also provides dual values for fixed variables if fixed_variable_treatment is set to make_constraint or make_parameter.

Attention
If there are fixed variables and fixed_variable_treatment is make_parameter (the default), then the Gradient of f(x) and the Jacobian of g(x) may be reevaluated here (that's why the function needs x). Further, in this setting, this function only provides correct bound multipliers for fixed variables if x, y_c, and y_d correspond to the unscaled problem.
Returns
True, if bound multipliers could be assigned. False if there was an evaluation error when calculating bound multipliers for fixed variables.
Since
3.14.0
Parameters
xinternal values for primal variables x
y_cinternal values for equality constraint multipliers
y_dinternal values for inequality constraint multipliers
z_Linternal values for lower bound multipliers
z_L_origvector to fill with values from z_L
z_Uinternal values for upper bound multipliers
z_U_origvector to fill with values from z_U

◆ GetFullDimensions()

void Ipopt::TNLPAdapter::GetFullDimensions ( Index n,
Index m 
) const
inline

Get number of variables and number of constraints in TNLP.

Since
3.14.0
Parameters
nstorage for full dimension of x (fixed + non-fixed)
mstorage for full dimension of g (c + d)

Definition at line 368 of file IpTNLPAdapter.hpp.

◆ GetFixedVariables()

void Ipopt::TNLPAdapter::GetFixedVariables ( Index n_x_fixed,
Index *&  x_fixed_map,
FixedVariableTreatmentEnum fixed_variable_treatment 
) const
inline

Get number and indices of fixed variables.

Since
3.14.0
Parameters
n_x_fixedstorage for number of fixed variables in TNLP
x_fixed_mapstorage for pointer to array that holds indices of fixed variables (has length n_fixed_x, can be NULL if n_fixed_x=0)
fixed_variable_treatmenttreatment for fixed variables as used by TNLP

Definition at line 380 of file IpTNLPAdapter.hpp.

◆ GetPermutationMatrices()

void Ipopt::TNLPAdapter::GetPermutationMatrices ( SmartPtr< const ExpansionMatrix > &  P_x_full_x,
SmartPtr< const ExpansionMatrix > &  P_x_x_L,
SmartPtr< const ExpansionMatrix > &  P_x_x_U,
SmartPtr< const ExpansionMatrix > &  P_c_g,
SmartPtr< const ExpansionMatrix > &  P_d_g 
) const
inline

Get mappings between TNLP indices and Ipopt internal indices for variables and constraints.

See the various ResortXyz functions on usage.

Attention
P_x_full_x is set to NULL if there are no fixed variables or fixed_variable_treatment is not make_parameter
Since
3.14.0
Parameters
P_x_full_xmap from TNLP variable indices to Ipopt internal indices, filtered variables get mapped to -1
P_x_x_Lmap from indices on lower bounds on x to Ipopt internal indices for x
P_x_x_Umap from indices on upper bounds on x to Ipopt internal indices for x
P_c_gmap from indices on equality constraints (c(x)=0) into TNLP constraint indices (g_l <= g(x) <= g_u)
P_d_gmap from indices on inequality constraints (d(x)-s=0) into TNLP constraint indices (g_l <= g(x) <= g_u)

Definition at line 398 of file IpTNLPAdapter.hpp.

◆ GetC_Rhs()

const Number * Ipopt::TNLPAdapter::GetC_Rhs ( ) const
inline

Get right-hand-sides that are added into c(x)

Since
3.14.0

Definition at line 416 of file IpTNLPAdapter.hpp.

◆ operator=()

void Ipopt::TNLPAdapter::operator= ( const TNLPAdapter )
private

Default Assignment Operator.

◆ DetermineDependentConstraints()

bool Ipopt::TNLPAdapter::DetermineDependentConstraints ( Index  n_x_var,
const Index x_not_fixed_map,
const Number x_l,
const Number x_u,
const Number g_l,
const Number g_u,
Index  n_c,
const Index c_map,
std::list< Index > &  c_deps 
)
private

◆ update_local_x()

bool Ipopt::TNLPAdapter::update_local_x ( const Vector x)
private

◆ update_local_lambda()

bool Ipopt::TNLPAdapter::update_local_lambda ( const Vector y_c,
const Vector y_d 
)
private

◆ internal_eval_g()

bool Ipopt::TNLPAdapter::internal_eval_g ( bool  new_x)
private

◆ internal_eval_jac_g()

bool Ipopt::TNLPAdapter::internal_eval_jac_g ( bool  new_x)
private

◆ initialize_findiff_jac()

void Ipopt::TNLPAdapter::initialize_findiff_jac ( const Index iRow,
const Index jCol 
)
private

Initialize sparsity structure for finite difference Jacobian.

Member Data Documentation

◆ tnlp_

SmartPtr<TNLP> Ipopt::TNLPAdapter::tnlp_
private

Pointer to the TNLP class (class specific to Number* vectors and triplet matrices)

Definition at line 458 of file IpTNLPAdapter.hpp.

◆ jnlst_

SmartPtr<const Journalist> Ipopt::TNLPAdapter::jnlst_
private

Journalist.

Definition at line 461 of file IpTNLPAdapter.hpp.

◆ dependency_detector_

SmartPtr<TDependencyDetector> Ipopt::TNLPAdapter::dependency_detector_
private

Object that can be used to detect linearly dependent rows in the equality constraint Jacobian.

Definition at line 464 of file IpTNLPAdapter.hpp.

◆ nlp_lower_bound_inf_

Number Ipopt::TNLPAdapter::nlp_lower_bound_inf_
private

Value for a lower bound that denotes -infinity.

Definition at line 469 of file IpTNLPAdapter.hpp.

◆ nlp_upper_bound_inf_

Number Ipopt::TNLPAdapter::nlp_upper_bound_inf_
private

Value for a upper bound that denotes infinity.

Definition at line 471 of file IpTNLPAdapter.hpp.

◆ fixed_variable_treatment_

FixedVariableTreatmentEnum Ipopt::TNLPAdapter::fixed_variable_treatment_
private

Flag indicating how fixed variables should be handled.

Definition at line 473 of file IpTNLPAdapter.hpp.

◆ bound_relax_factor_

Number Ipopt::TNLPAdapter::bound_relax_factor_
private

Determines relaxation of fixing bound for RELAX_BOUNDS.

Definition at line 475 of file IpTNLPAdapter.hpp.

◆ derivative_test_

DerivativeTestEnum Ipopt::TNLPAdapter::derivative_test_
private

Maximal slack for one-sidedly bounded variables.

If a variable has only one bound, say a lower bound xL, then an upper bound xL + max_onesided_bound_slack_. If this value is zero, no upper bound is added. Whether and which derivative test should be performed at starting point

Definition at line 482 of file IpTNLPAdapter.hpp.

◆ derivative_test_perturbation_

Number Ipopt::TNLPAdapter::derivative_test_perturbation_
private

Size of the perturbation for the derivative test.

Definition at line 484 of file IpTNLPAdapter.hpp.

◆ derivative_test_tol_

Number Ipopt::TNLPAdapter::derivative_test_tol_
private

Relative threshold for marking deviation from finite difference test.

Definition at line 486 of file IpTNLPAdapter.hpp.

◆ derivative_test_print_all_

bool Ipopt::TNLPAdapter::derivative_test_print_all_
private

Flag indicating if all test values should be printed, or only those violating the threshold.

Definition at line 488 of file IpTNLPAdapter.hpp.

◆ derivative_test_first_index_

Index Ipopt::TNLPAdapter::derivative_test_first_index_
private

Index of first quantity to be checked.

Definition at line 490 of file IpTNLPAdapter.hpp.

◆ warm_start_same_structure_

bool Ipopt::TNLPAdapter::warm_start_same_structure_
private

Flag indicating whether the TNLP with identical structure has already been solved before.

Definition at line 492 of file IpTNLPAdapter.hpp.

◆ hessian_approximation_

HessianApproximationType Ipopt::TNLPAdapter::hessian_approximation_
private

Flag indicating what Hessian information is to be used.

Definition at line 494 of file IpTNLPAdapter.hpp.

◆ num_linear_variables_

Index Ipopt::TNLPAdapter::num_linear_variables_
private

Number of linear variables.

Definition at line 496 of file IpTNLPAdapter.hpp.

◆ jacobian_approximation_

JacobianApproxEnum Ipopt::TNLPAdapter::jacobian_approximation_
private

Flag indicating how Jacobian is computed.

Definition at line 498 of file IpTNLPAdapter.hpp.

◆ gradient_approximation_

GradientApproxEnum Ipopt::TNLPAdapter::gradient_approximation_
private

Flag indicating how objective Gradient is computed.

Definition at line 500 of file IpTNLPAdapter.hpp.

◆ findiff_perturbation_

Number Ipopt::TNLPAdapter::findiff_perturbation_
private

Size of the perturbation for the derivative approximation.

Definition at line 502 of file IpTNLPAdapter.hpp.

◆ point_perturbation_radius_

Number Ipopt::TNLPAdapter::point_perturbation_radius_
private

Maximal perturbation of the initial point.

Definition at line 504 of file IpTNLPAdapter.hpp.

◆ dependency_detection_with_rhs_

bool Ipopt::TNLPAdapter::dependency_detection_with_rhs_
private

Flag indicating if rhs should be considered during dependency detection.

Definition at line 506 of file IpTNLPAdapter.hpp.

◆ tol_

Number Ipopt::TNLPAdapter::tol_
private

Overall convergence tolerance.

Definition at line 509 of file IpTNLPAdapter.hpp.

◆ n_full_x_

Index Ipopt::TNLPAdapter::n_full_x_
private

full dimension of x (fixed + non-fixed)

Definition at line 515 of file IpTNLPAdapter.hpp.

◆ n_full_g_

Index Ipopt::TNLPAdapter::n_full_g_
private

full dimension of g (c + d)

Definition at line 517 of file IpTNLPAdapter.hpp.

◆ nz_jac_c_

Index Ipopt::TNLPAdapter::nz_jac_c_
private

non-zeros of the jacobian of c

Definition at line 519 of file IpTNLPAdapter.hpp.

◆ nz_jac_c_no_extra_

Index Ipopt::TNLPAdapter::nz_jac_c_no_extra_
private

non-zeros of the jacobian of c without added constraints for fixed variables.

Definition at line 521 of file IpTNLPAdapter.hpp.

◆ nz_jac_d_

Index Ipopt::TNLPAdapter::nz_jac_d_
private

non-zeros of the jacobian of d

Definition at line 523 of file IpTNLPAdapter.hpp.

◆ nz_full_jac_g_

Index Ipopt::TNLPAdapter::nz_full_jac_g_
private

number of non-zeros in full-size Jacobian of g

Definition at line 525 of file IpTNLPAdapter.hpp.

◆ nz_full_h_

Index Ipopt::TNLPAdapter::nz_full_h_
private

number of non-zeros in full-size Hessian

Definition at line 527 of file IpTNLPAdapter.hpp.

◆ nz_h_

Index Ipopt::TNLPAdapter::nz_h_
private

number of non-zeros in the non-fixed-size Hessian

Definition at line 529 of file IpTNLPAdapter.hpp.

◆ n_x_fixed_

Index Ipopt::TNLPAdapter::n_x_fixed_
private

Number of fixed variables.

Definition at line 531 of file IpTNLPAdapter.hpp.

◆ index_style_

TNLP::IndexStyleEnum Ipopt::TNLPAdapter::index_style_
private

Numbering style of variables and constraints.

Definition at line 535 of file IpTNLPAdapter.hpp.

◆ x_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::x_space_
private

Definition at line 539 of file IpTNLPAdapter.hpp.

◆ c_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::c_space_
private

Definition at line 540 of file IpTNLPAdapter.hpp.

◆ d_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::d_space_
private

Definition at line 541 of file IpTNLPAdapter.hpp.

◆ x_l_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::x_l_space_
private

Definition at line 542 of file IpTNLPAdapter.hpp.

◆ px_l_space_

SmartPtr<const MatrixSpace> Ipopt::TNLPAdapter::px_l_space_
private

Definition at line 543 of file IpTNLPAdapter.hpp.

◆ x_u_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::x_u_space_
private

Definition at line 544 of file IpTNLPAdapter.hpp.

◆ px_u_space_

SmartPtr<const MatrixSpace> Ipopt::TNLPAdapter::px_u_space_
private

Definition at line 545 of file IpTNLPAdapter.hpp.

◆ d_l_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::d_l_space_
private

Definition at line 546 of file IpTNLPAdapter.hpp.

◆ pd_l_space_

SmartPtr<const MatrixSpace> Ipopt::TNLPAdapter::pd_l_space_
private

Definition at line 547 of file IpTNLPAdapter.hpp.

◆ d_u_space_

SmartPtr<const VectorSpace> Ipopt::TNLPAdapter::d_u_space_
private

Definition at line 548 of file IpTNLPAdapter.hpp.

◆ pd_u_space_

SmartPtr<const MatrixSpace> Ipopt::TNLPAdapter::pd_u_space_
private

Definition at line 549 of file IpTNLPAdapter.hpp.

◆ Jac_c_space_

SmartPtr<const MatrixSpace> Ipopt::TNLPAdapter::Jac_c_space_
private

Definition at line 550 of file IpTNLPAdapter.hpp.

◆ Jac_d_space_

SmartPtr<const MatrixSpace> Ipopt::TNLPAdapter::Jac_d_space_
private

Definition at line 551 of file IpTNLPAdapter.hpp.

◆ Hess_lagrangian_space_

SmartPtr<const SymMatrixSpace> Ipopt::TNLPAdapter::Hess_lagrangian_space_
private

Definition at line 552 of file IpTNLPAdapter.hpp.

◆ full_x_

Number* Ipopt::TNLPAdapter::full_x_
private

Definition at line 557 of file IpTNLPAdapter.hpp.

◆ full_lambda_

Number* Ipopt::TNLPAdapter::full_lambda_
private

copy of the full x vector (fixed & non-fixed)

Definition at line 558 of file IpTNLPAdapter.hpp.

◆ full_g_

Number* Ipopt::TNLPAdapter::full_g_
private

copy of lambda (yc & yd)

Definition at line 559 of file IpTNLPAdapter.hpp.

◆ jac_g_

Number* Ipopt::TNLPAdapter::jac_g_
private

copy of g (c & d)

Definition at line 560 of file IpTNLPAdapter.hpp.

◆ c_rhs_

Number* Ipopt::TNLPAdapter::c_rhs_
private

the values for the full jacobian of g

Definition at line 561 of file IpTNLPAdapter.hpp.

◆ x_tag_for_iterates_

TaggedObject::Tag Ipopt::TNLPAdapter::x_tag_for_iterates_
private

Definition at line 566 of file IpTNLPAdapter.hpp.

◆ y_c_tag_for_iterates_

TaggedObject::Tag Ipopt::TNLPAdapter::y_c_tag_for_iterates_
private

Definition at line 567 of file IpTNLPAdapter.hpp.

◆ y_d_tag_for_iterates_

TaggedObject::Tag Ipopt::TNLPAdapter::y_d_tag_for_iterates_
private

Definition at line 568 of file IpTNLPAdapter.hpp.

◆ x_tag_for_g_

TaggedObject::Tag Ipopt::TNLPAdapter::x_tag_for_g_
private

Definition at line 569 of file IpTNLPAdapter.hpp.

◆ x_tag_for_jac_g_

TaggedObject::Tag Ipopt::TNLPAdapter::x_tag_for_jac_g_
private

Definition at line 570 of file IpTNLPAdapter.hpp.

◆ P_x_full_x_

SmartPtr<ExpansionMatrix> Ipopt::TNLPAdapter::P_x_full_x_
private

Expansion from fixed x (ipopt) to full x.

Definition at line 598 of file IpTNLPAdapter.hpp.

◆ P_x_full_x_space_

SmartPtr<ExpansionMatrixSpace> Ipopt::TNLPAdapter::P_x_full_x_space_
private

Definition at line 599 of file IpTNLPAdapter.hpp.

◆ P_x_x_L_

SmartPtr<ExpansionMatrix> Ipopt::TNLPAdapter::P_x_x_L_
private

Expansion from fixed x_L (ipopt) to full x.

Definition at line 602 of file IpTNLPAdapter.hpp.

◆ P_x_x_L_space_

SmartPtr<ExpansionMatrixSpace> Ipopt::TNLPAdapter::P_x_x_L_space_
private

Definition at line 603 of file IpTNLPAdapter.hpp.

◆ P_x_x_U_

SmartPtr<ExpansionMatrix> Ipopt::TNLPAdapter::P_x_x_U_
private

Expansion from fixed x_U (ipopt) to full x.

Definition at line 606 of file IpTNLPAdapter.hpp.

◆ P_x_x_U_space_

SmartPtr<ExpansionMatrixSpace> Ipopt::TNLPAdapter::P_x_x_U_space_
private

Definition at line 607 of file IpTNLPAdapter.hpp.

◆ P_c_g_space_

SmartPtr<ExpansionMatrixSpace> Ipopt::TNLPAdapter::P_c_g_space_
private

Expansion from c only (ipopt) to full ampl c.

Definition at line 610 of file IpTNLPAdapter.hpp.

◆ P_c_g_

SmartPtr<ExpansionMatrix> Ipopt::TNLPAdapter::P_c_g_
private

Definition at line 611 of file IpTNLPAdapter.hpp.

◆ P_d_g_space_

SmartPtr<ExpansionMatrixSpace> Ipopt::TNLPAdapter::P_d_g_space_
private

Expansion from d only (ipopt) to full ampl d.

Definition at line 614 of file IpTNLPAdapter.hpp.

◆ P_d_g_

SmartPtr<ExpansionMatrix> Ipopt::TNLPAdapter::P_d_g_
private

Definition at line 615 of file IpTNLPAdapter.hpp.

◆ jac_idx_map_

Index* Ipopt::TNLPAdapter::jac_idx_map_
private

Definition at line 617 of file IpTNLPAdapter.hpp.

◆ h_idx_map_

Index* Ipopt::TNLPAdapter::h_idx_map_
private

Definition at line 618 of file IpTNLPAdapter.hpp.

◆ x_fixed_map_

Index* Ipopt::TNLPAdapter::x_fixed_map_
private

Position of fixed variables.

This is required for a warm start

Definition at line 621 of file IpTNLPAdapter.hpp.

◆ jac_fixed_idx_map_

std::vector<Index> Ipopt::TNLPAdapter::jac_fixed_idx_map_
private

Index mapping of Jacobian w.r.t.

fixed variables.

Definition at line 624 of file IpTNLPAdapter.hpp.

◆ jac_fixed_iRow_

std::vector<Index> Ipopt::TNLPAdapter::jac_fixed_iRow_
private

Definition at line 625 of file IpTNLPAdapter.hpp.

◆ jac_fixed_jCol_

std::vector<Index> Ipopt::TNLPAdapter::jac_fixed_jCol_
private

Definition at line 626 of file IpTNLPAdapter.hpp.

◆ findiff_jac_nnz_

Index Ipopt::TNLPAdapter::findiff_jac_nnz_
private

Number of unique nonzeros in constraint Jacobian.

Definition at line 632 of file IpTNLPAdapter.hpp.

◆ findiff_jac_ia_

Index* Ipopt::TNLPAdapter::findiff_jac_ia_
private

Start position for nonzero indices in ja for each column of Jacobian.

Definition at line 634 of file IpTNLPAdapter.hpp.

◆ findiff_jac_ja_

Index* Ipopt::TNLPAdapter::findiff_jac_ja_
private

Ordered by columns, for each column the row indices in Jacobian.

Definition at line 636 of file IpTNLPAdapter.hpp.

◆ findiff_jac_postriplet_

Index* Ipopt::TNLPAdapter::findiff_jac_postriplet_
private

Position of entry in original triplet matrix.

Definition at line 638 of file IpTNLPAdapter.hpp.

◆ findiff_x_l_

Number* Ipopt::TNLPAdapter::findiff_x_l_
private

Copy of the lower bounds.

Definition at line 640 of file IpTNLPAdapter.hpp.

◆ findiff_x_u_

Number* Ipopt::TNLPAdapter::findiff_x_u_
private

Copy of the upper bounds.

Definition at line 642 of file IpTNLPAdapter.hpp.


The documentation for this class was generated from the following file: