OPTPP::OptConstrNewtonLike Class Reference

#include <OptConstrNewtonLike.h>

Inheritance diagram for OPTPP::OptConstrNewtonLike:

OPTPP::OptimizeClass OPTPP::OptConstrNewton1Deriv OPTPP::OptConstrNewton2Deriv OPTPP::OptNIPSLike OPTPP::OptConstrFDNewton OPTPP::OptConstrQNewton OPTPP::OptConstrNewton OPTPP::OptDHNIPS OPTPP::OptFDNIPS OPTPP::OptNIPS OPTPP::OptQNIPS List of all members.

Public Member Functions

 OptConstrNewtonLike ()
 OptConstrNewtonLike (int n)
 OptConstrNewtonLike (int n, UPDATEFCN u)
 OptConstrNewtonLike (int n, TOLS t)
virtual ~OptConstrNewtonLike ()
virtual void acceptStep (int k, int step_type)
 Accept this step and update the nonlinear model.
virtual NEWMAT::ColumnVector computeSearch (NEWMAT::SymmetricMatrix &H)
 Solve for the Newton direction.
virtual void updateModel (int k, int ndim, NEWMAT::ColumnVector x)
virtual void reset ()
virtual int checkConvg ()
 Check to see if algorithm satisfies the convergence criterion.
virtual int checkDeriv ()
 Compare the analytic gradient with the finite-difference approximation.
virtual int computeStep (NEWMAT::ColumnVector sk)
 If an acceptable step not found, returns an error code = -1.
virtual real computeMaxStep (NEWMAT::ColumnVector sk)
virtual void initOpt ()
 Initialize algorithmic parameters.
virtual void initHessian ()
 Compute the Hessian or its approximation at the initial point.
virtual double initTrustRegionSize () const
 Initialize the radius of the trust-region.
virtual void optimize ()
 Invoke a constrained Newton's method.
virtual void readOptInput ()
 Read user-specified input options from a file.
int checkAnalyticFDGrad ()
 Check finite difference approximation with analytic derivatives.
int getFevals () const
int getGevals () const
real getTRSize () const
void setTRSize (real delta)
 Set radius of trust-region.
real getGradMult () const
void setGradMult (real tau)
 Set gradient multiplier which is used to compute radius of trust-region.
int getSearchSize () const
void setSearchSize (int sss)
 Set number of points in search scheme for trustpds strategy.
bool getWarmStart () const
void UseWarmStart (NEWMAT::SymmetricMatrix &H)
void setSearchStrategy (SearchStrategy search)
SearchStrategy getSearchStrategy () const
 Set globalization strategy for optimization algorithm.
void setDerivOption (DerivOption d)
DerivOption getDerivOption () const
 Set type of finite difference approximation (forward, backward, or central).
NEWMAT::ColumnVector getEqualityMultiplier () const
void setEqualityMultiplier (const NEWMAT::ColumnVector &ymult)
 Store the current Lagrangian multipliers corresponding to equality constraints.
NEWMAT::ColumnVector getInequalityMultiplier () const
void setInequalityMultiplier (const NEWMAT::ColumnVector &zmult)
 Store the current Lagrangian multipliers corresponding to inequality constraints.
NEWMAT::ColumnVector getSlacks () const
void setSlacks (const NEWMAT::ColumnVector &slackVar)
 Store the current slack vector.
MeritFcn getMeritFcn () const
virtual void setMeritFcn (MeritFcn option)
 Specify the merit function to used in step acceptance test.
NEWMAT::ColumnVector getGradL () const
virtual void setGradL (NEWMAT::ColumnVector gradl_value)
 Store the gradient of the Lagrangian at the current iteration.
NEWMAT::ColumnVector getGradLPrev () const
virtual void setGradLPrev (NEWMAT::ColumnVector gradl_value)
 Store the gradient of the Lagrangian at the previous iteration.
NEWMAT::ColumnVector getConstraintResidual () const
virtual void setConstraintResidual (const NEWMAT::ColumnVector &constraint_value)
 Store the residuals of the constraints at the current iteration.
NEWMAT::Matrix getConstraintGradient () const
virtual void setConstraintGradient (const NEWMAT::Matrix &constraint_grad)
 Store the current gradients of the constraints.
real getCost () const
void setCost (real value)
 Store current value of merit function.
bool getFeasibilityRecoveryFlag () const
void setFeasibilityRecoveryFlag (bool flag)
 Set switch to turn on feasibility recovery method.
int getMaxFeasIter () const
void setMaxFeasIter (int k)
 Set maximum number of iterations for feasibility recovery method of trust-region.
NEWMAT::SymmetricMatrix getHessian () const
void setHessian (NEWMAT::SymmetricMatrix &H)
 Store current Hessian matrix.
virtual NEWMAT::SymmetricMatrix updateH (NEWMAT::SymmetricMatrix &H, int k)=0
 Compute the Hessian of the Lagrangrian or its approximation at iteration k.
NEWMAT::ColumnVector computeFFK1Ind (const NEWMAT::ColumnVector &xc)
 Compute the active set according to Facchinei, Fischer, and Kanzow indicator.
NEWMAT::ColumnVector computeFFK2Ind (const NEWMAT::ColumnVector &xc)
 Compute the active set according to Facchinei, Fischer, and Kanzow indicator.
NEWMAT::ColumnVector computeTapiaInd (const NEWMAT::ColumnVector &step)
 Compute the Tapia indicators for the constrained problem.
void printStatus (char *)
 Print status of the constrained Newton's method.
void printMultipliers (char *)
 Output the Lagrangian multipliers to the screen.
void fPrintMultipliers (ostream *nlpout, char *)
 Print the Lagrangian multipliers to a file.
void fPrintSecSuff (ostream *nlpout, NEWMAT::ColumnVector &info)
 Print second order sufficiency information to a file.

Protected Member Functions

virtual NLP1nlprob () const =0
 returns an NLP1 pointer
void defaultAcceptStep (int, int)
NEWMAT::ColumnVector defaultComputeSearch (NEWMAT::SymmetricMatrix &)

Protected Attributes

int me
 Number of equality constraints.
int mi
 Number of inequality constraints.
int grad_evals
 Number of gradient evaluations.
NEWMAT::ColumnVector gprev
 Gradient at previous iteration.
NEWMAT::ColumnVector z
 Lagrange mult. corr. to ineq constraints.
NEWMAT::ColumnVector y
 Lagrange mult. corr. to eq constraints.
NEWMAT::ColumnVector s
 Slack variables corr. to ineq. constraints.
NEWMAT::ColumnVector constrType
 Vector which contains the type of constraint.
NEWMAT::ColumnVector constraintResidual
 Constraint residual at xc.
NEWMAT::ColumnVector gradl
 Current gradient of the lagrangian.
NEWMAT::ColumnVector gradlprev
 Previous gradient of the lagrangian.
NEWMAT::Matrix constraintGradient
 Current constraint gradient.
NEWMAT::Matrix constraintGradientPrev
 Previous constraint gradient.
NEWMAT::SymmetricMatrix Hessian
 Current Hessian.
NEWMAT::SymmetricMatrix hessl
 Current Hessian of the lagrangian.
SearchStrategy strategy
 User-specified globalization strategy.
DerivOption finitediff
 User-specified derivative option.
MeritFcn mfcn
 Merit function option.
real TR_size
 Size of the trust region radius.
real gradMult
 Gradient multiplier to compute TR_size.
int searchSize
 Search pattern size for TRPDS.
real cost
 Value of the merit function.
bool WarmStart
bool feas_flag
 Switch to turn on the feasibility recovery method.
int max_feas_iter
 Maximize number of iterations to achieve feasibility.

Friends

int trustregion (NLP1 *, ostream *, NEWMAT::SymmetricMatrix &, NEWMAT::ColumnVector &, NEWMAT::ColumnVector &, real &, real &, real stpmax, real stpmin)
int trustpds (NLP1 *, ostream *, NEWMAT::SymmetricMatrix &, NEWMAT::ColumnVector &, NEWMAT::ColumnVector &, real &, real &, real stpmax, real stpmin, int)

Detailed Description

Constrained Newton abstract data classes OptConstrNewtonLike OptConstrNewton1Deriv OptConstrNewton2Deriv

OptConstrNewtonLike provides common data and functionality for OptConstrFDNewton, OptConstrQNewton, OptConstrNewton, OptFDNIPS, OptQNIPS, and OptNIPS methods.

Author:
J.C. Meza, Lawrence Berkeley National Laboratory
Note:
Modified by P.J. Williams
Date:
03/27/2007


Constructor & Destructor Documentation

OPTPP::OptConstrNewtonLike::OptConstrNewtonLike (  )  [inline]

Default Constructor

See also:
OptConstrNewtonLike(int n)

OptConstrNewtonLike(int n, UPDATEFCN u)

OptConstrNewtonLike(int n, TOLS t)

OPTPP::OptConstrNewtonLike::OptConstrNewtonLike ( int  n  )  [inline]

Parameters:
n an integer argument
See also:
OptConstrNewtonLike(int n, UPDATEFCN u)

OptConstrNewtonLike(int n, TOLS t)

OPTPP::OptConstrNewtonLike::OptConstrNewtonLike ( int  n,
UPDATEFCN  u 
) [inline]

Parameters:
n an integer argument
u a function pointer.
See also:
OptConstrNewtonLike(int n)

OptConstrNewtonLike(int n, TOLS t)

OPTPP::OptConstrNewtonLike::OptConstrNewtonLike ( int  n,
TOLS  t 
) [inline]

Parameters:
n an integer argument
t tolerance class reference.
See also:
OptConstrNewtonLike(int n)

OptConstrNewtonLike(int n, UPDATEFCN u)

virtual OPTPP::OptConstrNewtonLike::~OptConstrNewtonLike (  )  [inline, virtual]

Destructor


Member Function Documentation

virtual void OPTPP::OptConstrNewtonLike::acceptStep ( int  k,
int  step_type 
) [inline, virtual]

Accept this step and update the nonlinear model.

Implements OPTPP::OptimizeClass.

int OPTPP::OptConstrNewtonLike::checkAnalyticFDGrad (  ) 

Check finite difference approximation with analytic derivatives.

int OPTPP::OptConstrNewtonLike::checkConvg (  )  [virtual]

Check to see if algorithm satisfies the convergence criterion.

Implements OPTPP::OptimizeClass.

Reimplemented in OPTPP::OptNIPSLike.

int OPTPP::OptConstrNewtonLike::checkDeriv (  )  [virtual]

Compare the analytic gradient with the finite-difference approximation.

Reimplemented in OPTPP::OptConstrFDNewton, OPTPP::OptConstrNewton, OPTPP::OptConstrQNewton, OPTPP::OptFDNIPS, OPTPP::OptNIPSLike, and OPTPP::OptQNIPS.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::computeFFK1Ind ( const NEWMAT::ColumnVector &  xc  ) 

Compute the active set according to Facchinei, Fischer, and Kanzow indicator.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::computeFFK2Ind ( const NEWMAT::ColumnVector &  xc  ) 

Compute the active set according to Facchinei, Fischer, and Kanzow indicator.

virtual real OPTPP::OptConstrNewtonLike::computeMaxStep ( NEWMAT::ColumnVector  sk  )  [virtual]

virtual NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::computeSearch ( NEWMAT::SymmetricMatrix &  H  )  [inline, virtual]

Solve for the Newton direction.

Implements OPTPP::OptimizeClass.

virtual int OPTPP::OptConstrNewtonLike::computeStep ( NEWMAT::ColumnVector  sk  )  [virtual]

If an acceptable step not found, returns an error code = -1.

Reimplemented in OPTPP::OptNIPSLike.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::computeTapiaInd ( const NEWMAT::ColumnVector &  step  ) 

Compute the Tapia indicators for the constrained problem.

void OPTPP::OptConstrNewtonLike::defaultAcceptStep ( int  iter,
int  step_type 
) [protected, virtual]

Provide default implementation of AcceptStep

Reimplemented from OPTPP::OptimizeClass.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::defaultComputeSearch ( NEWMAT::SymmetricMatrix &   )  [protected]

void OPTPP::OptConstrNewtonLike::fPrintMultipliers ( ostream *  nlpout,
char *   
)

Print the Lagrangian multipliers to a file.

void OPTPP::OptConstrNewtonLike::fPrintSecSuff ( ostream *  nlpout,
NEWMAT::ColumnVector &  info 
)

Print second order sufficiency information to a file.

NEWMAT::Matrix OPTPP::OptConstrNewtonLike::getConstraintGradient (  )  const [inline]

Returns:
Gradient of the constraints at the current iteration

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::getConstraintResidual (  )  const [inline]

Returns:
Residuals of the constraints at the current iteration

real OPTPP::OptConstrNewtonLike::getCost (  )  const [inline]

Returns:
Value of merit function

DerivOption OPTPP::OptConstrNewtonLike::getDerivOption (  )  const [inline]

Set type of finite difference approximation (forward, backward, or central).

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::getEqualityMultiplier (  )  const [inline]

Returns:
Lagrangian multipliers corresponding to equality constraints

bool OPTPP::OptConstrNewtonLike::getFeasibilityRecoveryFlag (  )  const [inline]

Returns:
Switch to turn on feasibility recovery method

int OPTPP::OptConstrNewtonLike::getFevals (  )  const [inline]

Returns:
The number of function evaluations

int OPTPP::OptConstrNewtonLike::getGevals (  )  const [inline]

Returns:
The number of gradient evaluations

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::getGradL (  )  const [inline]

Returns:
Gradient of the Lagrangian at the current iteration

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::getGradLPrev (  )  const [inline]

Returns:
Gradient of the Lagrangian at the previous iteration

real OPTPP::OptConstrNewtonLike::getGradMult (  )  const [inline]

Returns:
Gradient multiplier to compute radius of trust-region

NEWMAT::SymmetricMatrix OPTPP::OptConstrNewtonLike::getHessian (  )  const [inline]

Returns:
Hessian matrix

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::getInequalityMultiplier (  )  const [inline]

Returns:
Lagrangian multipliers corresponding to inequality constraints

int OPTPP::OptConstrNewtonLike::getMaxFeasIter (  )  const [inline]

Returns:
Maximum number of iterations for feasibility recovery method

MeritFcn OPTPP::OptConstrNewtonLike::getMeritFcn (  )  const [inline]

Returns:
Merit function

int OPTPP::OptConstrNewtonLike::getSearchSize (  )  const [inline]

Returns:
Number of points in search scheme (only relevant when trustpds search strategy is selected)

SearchStrategy OPTPP::OptConstrNewtonLike::getSearchStrategy (  )  const [inline]

Set globalization strategy for optimization algorithm.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::getSlacks (  )  const [inline]

Returns:
Slack variables associated with inequality constraints

real OPTPP::OptConstrNewtonLike::getTRSize (  )  const [inline]

Returns:
Trust region radius

bool OPTPP::OptConstrNewtonLike::getWarmStart (  )  const [inline]

void OPTPP::OptConstrNewtonLike::initHessian (  )  [virtual]

Compute the Hessian or its approximation at the initial point.

Reimplemented in OPTPP::OptConstrNewton, OPTPP::OptDHNIPS, OPTPP::OptNIPS, and OPTPP::OptNIPSLike.

void OPTPP::OptConstrNewtonLike::initOpt (  )  [virtual]

Initialize algorithmic parameters.

Reimplemented in OPTPP::OptNIPSLike.

double OPTPP::OptConstrNewtonLike::initTrustRegionSize (  )  const [virtual]

Initialize the radius of the trust-region.

virtual NLP1* OPTPP::OptConstrNewtonLike::nlprob (  )  const [protected, pure virtual]

returns an NLP1 pointer

Implemented in OPTPP::OptConstrNewton1Deriv, OPTPP::OptConstrNewton2Deriv, OPTPP::OptDHNIPS, OPTPP::OptFDNIPS, OPTPP::OptNIPS, OPTPP::OptNIPSLike, and OPTPP::OptQNIPS.

void OPTPP::OptConstrNewtonLike::optimize (  )  [virtual]

Invoke a constrained Newton's method.

Implements OPTPP::OptimizeClass.

Reimplemented in OPTPP::OptNIPSLike.

void OPTPP::OptConstrNewtonLike::printMultipliers ( char *   ) 

Output the Lagrangian multipliers to the screen.

void OPTPP::OptConstrNewtonLike::printStatus ( char *   )  [virtual]

Print status of the constrained Newton's method.

Implements OPTPP::OptimizeClass.

Reimplemented in OPTPP::OptConstrNewton, OPTPP::OptDHNIPS, OPTPP::OptNIPS, and OPTPP::OptNIPSLike.

void OPTPP::OptConstrNewtonLike::readOptInput (  )  [virtual]

Read user-specified input options from a file.

Implements OPTPP::OptimizeClass.

Reimplemented in OPTPP::OptNIPSLike.

void OPTPP::OptConstrNewtonLike::reset (  )  [virtual]

Implements OPTPP::OptimizeClass.

Reimplemented in OPTPP::OptDHNIPS, and OPTPP::OptNIPSLike.

virtual void OPTPP::OptConstrNewtonLike::setConstraintGradient ( const NEWMAT::Matrix &  constraint_grad  )  [inline, virtual]

Store the current gradients of the constraints.

virtual void OPTPP::OptConstrNewtonLike::setConstraintResidual ( const NEWMAT::ColumnVector &  constraint_value  )  [inline, virtual]

Store the residuals of the constraints at the current iteration.

void OPTPP::OptConstrNewtonLike::setCost ( real  value  )  [inline]

Store current value of merit function.

void OPTPP::OptConstrNewtonLike::setDerivOption ( DerivOption  d  )  [inline]

Returns:
Type of finite difference approximation

void OPTPP::OptConstrNewtonLike::setEqualityMultiplier ( const NEWMAT::ColumnVector &  ymult  )  [inline]

Store the current Lagrangian multipliers corresponding to equality constraints.

void OPTPP::OptConstrNewtonLike::setFeasibilityRecoveryFlag ( bool  flag  )  [inline]

Set switch to turn on feasibility recovery method.

virtual void OPTPP::OptConstrNewtonLike::setGradL ( NEWMAT::ColumnVector  gradl_value  )  [inline, virtual]

Store the gradient of the Lagrangian at the current iteration.

virtual void OPTPP::OptConstrNewtonLike::setGradLPrev ( NEWMAT::ColumnVector  gradl_value  )  [inline, virtual]

Store the gradient of the Lagrangian at the previous iteration.

void OPTPP::OptConstrNewtonLike::setGradMult ( real  tau  )  [inline]

Set gradient multiplier which is used to compute radius of trust-region.

void OPTPP::OptConstrNewtonLike::setHessian ( NEWMAT::SymmetricMatrix &  H  )  [inline]

Store current Hessian matrix.

void OPTPP::OptConstrNewtonLike::setInequalityMultiplier ( const NEWMAT::ColumnVector &  zmult  )  [inline]

Store the current Lagrangian multipliers corresponding to inequality constraints.

void OPTPP::OptConstrNewtonLike::setMaxFeasIter ( int  k  )  [inline]

Set maximum number of iterations for feasibility recovery method of trust-region.

virtual void OPTPP::OptConstrNewtonLike::setMeritFcn ( MeritFcn  option  )  [inline, virtual]

Specify the merit function to used in step acceptance test.

Reimplemented in OPTPP::OptNIPSLike.

void OPTPP::OptConstrNewtonLike::setSearchSize ( int  sss  )  [inline]

Set number of points in search scheme for trustpds strategy.

void OPTPP::OptConstrNewtonLike::setSearchStrategy ( SearchStrategy  search  )  [inline]

Returns:
Globalization strategy for optimization algorithm

void OPTPP::OptConstrNewtonLike::setSlacks ( const NEWMAT::ColumnVector &  slackVar  )  [inline]

Store the current slack vector.

void OPTPP::OptConstrNewtonLike::setTRSize ( real  delta  )  [inline]

Set radius of trust-region.

virtual NEWMAT::SymmetricMatrix OPTPP::OptConstrNewtonLike::updateH ( NEWMAT::SymmetricMatrix &  H,
int  k 
) [pure virtual]

Compute the Hessian of the Lagrangrian or its approximation at iteration k.

Implemented in OPTPP::OptConstrFDNewton, OPTPP::OptConstrNewton, OPTPP::OptConstrQNewton, OPTPP::OptDHNIPS, OPTPP::OptFDNIPS, OPTPP::OptNIPS, OPTPP::OptNIPSLike, and OPTPP::OptQNIPS.

virtual void OPTPP::OptConstrNewtonLike::updateModel ( int  k,
int  ndim,
NEWMAT::ColumnVector  x 
) [inline, virtual]

Implements OPTPP::OptimizeClass.

void OPTPP::OptConstrNewtonLike::UseWarmStart ( NEWMAT::SymmetricMatrix &  H  )  [inline]


Friends And Related Function Documentation

int trustpds ( NLP1 ,
ostream *  ,
NEWMAT::SymmetricMatrix &  ,
NEWMAT::ColumnVector &  ,
NEWMAT::ColumnVector &  ,
real ,
real ,
real  stpmax,
real  stpmin,
int   
) [friend]

int trustregion ( NLP1 ,
ostream *  ,
NEWMAT::SymmetricMatrix &  ,
NEWMAT::ColumnVector &  ,
NEWMAT::ColumnVector &  ,
real ,
real ,
real  stpmax,
real  stpmin 
) [friend]


Member Data Documentation

NEWMAT::Matrix OPTPP::OptConstrNewtonLike::constraintGradient [protected]

Current constraint gradient.

NEWMAT::Matrix OPTPP::OptConstrNewtonLike::constraintGradientPrev [protected]

Previous constraint gradient.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::constraintResidual [protected]

Constraint residual at xc.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::constrType [protected]

Vector which contains the type of constraint.

real OPTPP::OptConstrNewtonLike::cost [protected]

Value of the merit function.

bool OPTPP::OptConstrNewtonLike::feas_flag [protected]

Switch to turn on the feasibility recovery method.

DerivOption OPTPP::OptConstrNewtonLike::finitediff [protected]

User-specified derivative option.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::gprev [protected]

Gradient at previous iteration.

int OPTPP::OptConstrNewtonLike::grad_evals [protected]

Number of gradient evaluations.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::gradl [protected]

Current gradient of the lagrangian.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::gradlprev [protected]

Previous gradient of the lagrangian.

real OPTPP::OptConstrNewtonLike::gradMult [protected]

Gradient multiplier to compute TR_size.

NEWMAT::SymmetricMatrix OPTPP::OptConstrNewtonLike::Hessian [protected]

Current Hessian.

NEWMAT::SymmetricMatrix OPTPP::OptConstrNewtonLike::hessl [protected]

Current Hessian of the lagrangian.

int OPTPP::OptConstrNewtonLike::max_feas_iter [protected]

Maximize number of iterations to achieve feasibility.

int OPTPP::OptConstrNewtonLike::me [protected]

Number of equality constraints.

MeritFcn OPTPP::OptConstrNewtonLike::mfcn [protected]

Merit function option.

int OPTPP::OptConstrNewtonLike::mi [protected]

Number of inequality constraints.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::s [protected]

Slack variables corr. to ineq. constraints.

int OPTPP::OptConstrNewtonLike::searchSize [protected]

Search pattern size for TRPDS.

SearchStrategy OPTPP::OptConstrNewtonLike::strategy [protected]

User-specified globalization strategy.

real OPTPP::OptConstrNewtonLike::TR_size [protected]

Size of the trust region radius.

bool OPTPP::OptConstrNewtonLike::WarmStart [protected]

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::y [protected]

Lagrange mult. corr. to eq constraints.

NEWMAT::ColumnVector OPTPP::OptConstrNewtonLike::z [protected]

Lagrange mult. corr. to ineq constraints.


The documentation for this class was generated from the following files:
Bug Reports    OPT++ Developers    Copyright Information    GNU Lesser General Public License
Documentation, generated by , last revised August 30, 2006.