NLPInterfacePack: C++ Interfaces and Implementation for Non-Linear Programs  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | List of all members
NLPInterfacePack::NLP Class Referenceabstract

NLP interface class {abstract}. More...

#include <NLPInterfacePack_NLP.hpp>

Inheritance diagram for NLPInterfacePack::NLP:
Inheritance graph
[legend]

Classes

class  IncompatibleType
 Thrown if an incompatible object is used. More...
 
class  InvalidInitialization
 Thrown from initialize() if some logical error occured. More...
 
class  NoBounds
 Thrown some bounds do not existe. More...
 
class  UnInitialized
 Thrown if any member functions are called before initialize() has been called. More...
 
struct  ZeroOrderInfo
 Struct for objective and constriants (pointer). More...
 

Public Types

typedef Teuchos::RCP< const
VectorSpace > 
vec_space_ptr_t
 
typedef Teuchos::RCP< const
OptionsFromStreamPack::OptionsFromStream
options_ptr_t
 

Public Member Functions

const ZeroOrderInfo zero_order_info () const
 Return pointer to set quantities. More...
 
const ZeroOrderInfo zero_order_info_breve () const
 Return pointer to set hat quantities. More...
 

Static Public Member Functions

static value_type infinite_bound ()
 Value for an infinite bound. More...
 

Protected Member Functions

template<class T >
void assert_ref_set (T *p, std::string info) const
 Assert referece has been set for a quanity. More...
 

Constructors, Destructor

 NLP ()
 Initialize to no reference set to calculation quanities. More...
 
virtual ~NLP ()
 Destructor that cleans all the memory it owns. More...
 

NLP initialization

virtual void force_xinit_in_bounds (bool force_xinit_in_bounds)=0
 Set if the initial point must be within the bounds. More...
 
virtual bool force_xinit_in_bounds () const =0
 Returns if the initial point must be within the bounds. More...
 
virtual void set_options (const options_ptr_t &options)
 Set the options that this NLP may be interested in. More...
 
virtual const options_ptr_tget_options () const
 Get the OptionsFromStream object being used to extract the options from. More...
 
virtual void initialize (bool test_setup=false)
 Initialize the NLP before it is used. More...
 
virtual bool is_initialized () const =0
 Return if this is initialized. More...
 

Dimensionality.

virtual size_type n () const
 Return the number of variables. More...
 
virtual size_type m () const
 Return the number of general equality constraints. More...
 

Vector space objects

virtual vec_space_ptr_t space_x () const =0
 Vector space object for unknown variables x (dimension n). More...
 
virtual vec_space_ptr_t space_c () const =0
 Vector space object for general equality constraints c(x) (dimension m). More...
 

Bounds on the unknown variables x.

virtual size_type num_bounded_x () const =0
 Returns the number of variables in x(i) for which xl(i)> -infinite_bound() or xu(i) < +infinite_bound(). More...
 
virtual const Vectorxl () const =0
 Returns the lower bounds on the variables x. More...
 
virtual const Vectorxu () const =0
 Returns a reference to the vector of upper bounds on the variables x. More...
 
virtual value_type max_var_bounds_viol () const =0
 Set the maximum absolute value for which the variable bounds may be violated by when computing function and gradient values. More...
 

Initial guess of NLP solution

virtual const Vectorxinit () const =0
 Returns a reference to the vector of the initial guess for the solution x. More...
 
virtual void get_init_lagrange_mult (VectorMutable *lambda, VectorMutable *nu) const
 Get the initial value of the Lagrange multipliers lambda. More...
 

Set and access storage for the objective function value f(x).

virtual void set_f (value_type *f)
 Set a pointer to an value to be updated when this->calc_f() is called. More...
 
virtual value_type * get_f ()
 Return pointer passed to this->set_f(). More...
 
virtual value_type & f ()
 Returns non-const *this->get_f(). More...
 
virtual const value_type & f () const
 Returns const *this->get_f(). More...
 

Set and access storage for the residual of the general equality constriants c(x).

virtual void set_c (VectorMutable *c)
 Set a pointer to a vector to be updated when this->calc_c() is called. More...
 
virtual VectorMutableget_c ()
 Return pointer passed to this->set_c(). More...
 
virtual VectorMutablec ()
 Returns non-const *this->get_c(). More...
 
virtual const Vectorc () const
 Returns const *this->get_c(). More...
 

Unset calculation quantities

virtual void unset_quantities ()
 Call to unset all storage quantities (both in this class and all subclasses). More...
 

Calculation members

virtual void scale_f (value_type scale_f)=0
 Set the scaling of the objective function. More...
 
virtual value_type scale_f () const =0
 Return the scaling being used for the objective function. More...
 
virtual void calc_f (const Vector &x, bool newx=true) const
 Update the value for the objective f at the point x and put it in the stored reference. More...
 
virtual void calc_c (const Vector &x, bool newx=true) const
 Update the constraint residual vector for c at the point x and put it in the stored reference. More...
 

Report final solution

virtual void report_final_solution (const Vector &x, const Vector *lambda, const Vector *nu, bool is_optimal)
 Used by the solver to report the final solution and multipliers. More...
 

Objective and constraint function evaluation counts.

virtual size_type num_f_evals () const
 Gives the number of object function f(x) evaluations called by the solver since initialize() was called. More...
 
virtual size_type num_c_evals () const
 Gives the number of constraint function c(x) evaluations called by the solver since initialize() was called. Throws exception if this->m() == 0. More...
 

General inequalities and slack variables

virtual size_type ns () const
 Return the number of slack variables (i.e. number of general inequalities). More...
 
virtual vec_space_ptr_t space_c_breve () const
 Vector space object for the original equalities c_breve(x_breve) More...
 
virtual vec_space_ptr_t space_h_breve () const
 Vector space object for the original inequalities h_breve(x_breve) More...
 
virtual const Vectorhl_breve () const
 Returns a reference to the vector of lower bounds on the general inequality constraints h_breve(x_breve). More...
 
virtual const Vectorhu_breve () const
 Returns a reference to the vector of upper bounds on the general inequality constraints h_breve(x_breve). More...
 
virtual void set_c_breve (VectorMutable *c_breve)
 Set a pointer to a vector to be updated when this->calc_c_breve() is called. More...
 
virtual VectorMutableget_c_breve ()
 Return pointer passed to this->set_c_breve(). More...
 
virtual VectorMutablec_breve ()
 Returns non-const *this->get_c_breve(). More...
 
virtual const Vectorc_breve () const
 Returns const *this->get_c_breve(). More...
 
virtual void set_h_breve (VectorMutable *h_breve)
 Set a pointer to a vector to be updated when this->calc_h_breve() is called. More...
 
virtual VectorMutableget_h_breve ()
 Return pointer passed to this->set_h_breve(). More...
 
virtual VectorMutableh_breve ()
 Returns non-const *this->get_h_breve(). More...
 
virtual const Vectorh_breve () const
 Returns const *this->get_h_breve(). More...
 
virtual const Permutation & P_var () const
 Return the permutation object for the variables. More...
 
virtual const Permutation & P_equ () const
 Return the permutation object for the constraints. More...
 
virtual void calc_c_breve (const Vector &x, bool newx=true) const
 Update the constraint residual vector for c_breve at the point x and put it in the stored reference. More...
 
virtual void calc_h_breve (const Vector &x, bool newx=true) const
 Update the constraint residual vector for h_breve at the point x and put it in the stored reference. More...
 

Protected methods to be overridden by subclasses

virtual void imp_calc_f (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info) const =0
 Overridden to compute f(x) (and perhaps other quantities if set). More...
 
virtual void imp_calc_c (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info) const =0
 Overridden to compute c(x) and perhaps f(x) and/or h(x) (if multiple calculaiton = true). More...
 
virtual void imp_calc_c_breve (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info_breve) const
 Overridden to compute c_breve(x_breve) and perhaps f(x) and/or h_breve(x_breve) More...
 
virtual void imp_calc_h_breve (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info_breve) const
 Overridden to compute h_breve(x_breve) and perhaps f(x) and/or c_breve(x_breve). More...
 

Detailed Description

NLP interface class {abstract}.

Overview:

This class represents an abstract interface to a general nonlinear programming problem of the form (in mathematical and ASCII notation):

\[ \begin{array}{lcl} \mbox{min} & & f(x) \\ \mbox{s.t.} & & c(x) = 0 \\ & & x^L \leq x \leq x^U \end{array} \]

where:

    min    f(x)
    s.t.   c(x) = 0
           xl <= x <= xu
 where:
           x    <: space_x
           f(x) <: space_x -> R^1
           c(x) <: space_x -> space_c 
           space_x <: R^n
           space_c <: R^n -> R^m

The NLP is defined in terms of vector spaces for the unknowns x (space_x), the equality constraints c (space_c)and nonlinear operator functions f(x) and c(x). In the above form, none of the variables are fixed between bounds (strictly xl < xu). It is allowed however for m == 0 for the elimination of general constriants. It is also allowed for n == m in which case this represents a fully determined system of nonlinear equaltions. In any case, an objective function is always included in the formutation and will impact solution algorithms.

Special types of NLPs are identified as:

  1. Fully general NLP :
    • ( xl != -Inf || xu != +Inf ) && ( m > 0 )
  2. General equality only constrained NLP :
    • ( xl == -Inf && xu == +Inf ) && ( m > 0 )
  3. Bound constrained NLP :
    • ( xl != -Inf || xu != +Inf ) && ( m == 0 )
  4. Unconstrained NLP :
    • ( xl == -Inf && xu == +Inf ) && ( m == 0 )
  5. Nonlinear Equations (NLE) :
    • n == m

If n==m but some of the equations in c(x) are dependent (but consistent) then the problem is actually an NLP and not an NLE and it is possible that some of the variable bounds my be active at the solution but in general they can not be. An optimization algorithm may refuse to solve some of the above problems but this interface allows all of these different types of mathematical programming problems to be represented using this interface.

The Lagrangian for this problem is defined by:

 L = f(x) + lambda' * c(x) + nul' * ( xl - x ) + nuu' * ( x - xu )

The optimality conditions are given by:

 del(L,x) = del(f,x) + del(c,x) * lambda + nu = 0
 c(x) = 0
 nuu(i) * ( x(i) - xu(i) ) = 0,      for i = 1...n
 nuu(i) * ( x(i) - xu(i) ) = 0,      for i = 1...n
 where:
   nu = nuu - nul

What is unique about this interface is that the vector objects are hidden behind abstact interfaces. Clients can create vectors from the various vector spaces using the VectorSpace objects returned from this->space_x() (dim n), and this->space_c() (dim m). In this sense, an NLP object acks as an "Abstract Factory" to create the vectors needed by an optimization algorithm. This allows optimization software to be written in a way that is completly independent from the linear algebra components.

General Inequalities and Slacks

The underlying NLP may containe general inequality constraints which where converted to equalities using slack variables. The actual underlying NLP may take the form (in mathematical and ASCII notation):

\[ \begin{array}{lcl} \mbox{min} & & \hat{f}(\hat{x}) \\ \mbox{s.t.} & & \hat{c}(\hat{x}) = 0 \\ & & \hat{h}^L \leq \hat{h}(\hat{x}) \leq \hat{h}^U \\ & & \hat{x}^L \leq \hat{x} \leq \hat{x}^U \end{array} \]

where:

    min    f_breve(x_breve)
    s.t.   c_breve(x_breve) = 0
           hl_breve <= h_breve(x_breve) <= hu_breve
         xl_breve <= x_breve <= xu_breve
 where:
         x_breve    <: space_x
           f_breve(x_breve) <: space_x_breve -> R^1
         c_breve(x_breve) <: space_x_breve -> space_c_breve 
         h_breve(x_breve) <: space_x_breve -> space_h_breve
           space_x_breve <: R^n_breve
           space_c_breve <: R^n -> R^m_breve
           space_h_breve <: R^n -> R^mI_breve

ToDo: Finish!

Client Usage:

Before an NLP object can be used, the initialize() method must be called to make sure that all of the initializations needed for the NLP have been performed. This method also resets counters an other information. Before calling initialize() the client can specify whether the initial point for x must be in bounds by calling force_xinit_in_bounds(bool).

Smart reference counted pointers to the three vector spaces for x, c(x) and h(x) are returned by the methods space_x(), space_c() and space_h() respectively. The vector space objects returned by these methods are ment to be more than transient. In fact, it is expected that these vector space objects should remain valid for the entire run of an NLP algorithm. Only if the underlying NLP is changed in a fundamental way (i.e. n, or m changes) should the vector space objects returned from these function become invalid. In this case the client must call these methods again to get updated vector space objects.

The dimensionality of the NLP is returned by the methods n() and m() but they have default implementations based on space_x() and space_c() respectively.

The number of variables x with finite bounds is returned by the method num_bounded_x(). The methods xl() and xu() return references to vector objects representing these bounds. A lower bound is considered infinite if xl().get_ele(i) == -infinite_bound() and an upper bound is considered infinite if xu().get_ele(i) == +infinite_bound().

If ns() > 0, the methods hl_breve() and hu_breve() return references to the upper and lower bounds to the general inequality constraints h_breve(x_breve). While it is expected that hl_breve().get_ele(j) != -infinite_bound() || hu_breve().get_ele(j) != +infinite_bound for j = 1...ns() this is not required by this interface. On the other hand it seems silly to define general inequality constriants that are not bounded but there may be some reason to include these that makes things easier for the implementor of the NLP subclass.

The initial guess for the unknowns x (primal variables) is returned by the method xinit(). Vectors containing the initial guesses for the Lagrange multipliers can be obtained by calling the method get_init_lagrange_mult().

The bread and butter of an NLP interface is the calculation of the functions that define the objective and constraints and various points x using the methods calc_f() and calc_c(). The quantities that these functions update must be set prior by calling the methods set_f() and set_c() respectively. It may seem strange not to pass these quantities directly to the calculation functions but there is a good reason for this. The reason is that this interface supports the efficient update of mutiple quantities at a given point x. For example, in many NLPs some of the same terms are shared between the constriants functions and objective function. Therefore it is more efficient to compute f(x) and c(x) simultaneously rather than computing them separately. In order to allow for this possibility, the client can set the desired quantities (i.e. set_f() and set_c()) prior to calling calc_f() and calc_c(). Then, whatever quantities that have been set will be computed my any call to calc_?() method.

Once an optimization algorithm has the solution (or gives up with a suboptimal point), it should report this solution to the NLP object using the method report_final_solution().

Finally, the client can get the counts for the number of function evaluations since initialize() was called using the methods num_f_evals() and num_c_evals(). These counts do not include any function evaluations that may have been used internally for finite difference evaluations or anything of that nature. Also, if the client calls calc_info(x,false) (where info = f or c) several times then the default implementation will increment the count for each call even though the actual quantity may not actually be recalculated each time (i.e. if newx==false). This information is not known here in this base class but the subclasses can overide this behavior if desired.

Subclass developer's notes:

The calcuation methods calc_f() and calc_c() have default implementations in this base class that should meet the needs of all the subclasses. These method implementations call the protected pure virtual methods imp_calc_f() and imp_calc_c() to compute the actual quantities. Subclasses must override these methods (in addition to several methods from the public interface). Pointers to the quantities to be updated are passed to these methods to the subclasses in the form of an aggregate ZeroOrderInfo object that is returned from the protected method zero_order_info(). This ensures that the only interaction between an NLP base object and its subclass objects is through member functions and never through pulic or protected data members.

The following methods must be overridden by a subclass in order to create a concrete NLP object: force_xinit_in_bounds(bool), force_xinit_in_bounds(), is_initialized(), space_x(), space_c(), num_bounded_x(), xl(), xu(), xinit(),] scale_f(value_type), scale_f().

The following methods should be overridden by most subclasses but do not have to be: initialize(), get_init_lagrange_mult(), report_final_solution().

The following methods should never have to be overridden by most subclasses except in some very strange situations: set_f(), get_f(), f(), set_c(), get_c(), c(), calc_f(), calc_c(), num_f_evals(), num_c_evals().

Additional notes:

The bounds on the variables can play a very critical role in many optimization algorithms. It is desirable for the functions f(x) and c(x) to be defined and relatively well behaved (i.e. smooth, continous, differentiable etc.) in the region xl <= x <= xu. While this is not always possible, an NLP can often be reformulated to have this properly. For example, suppose there are constraints has the form:

log(x(1) - x(5)) - 4 == 0
x(1) - x(5) >= 1e-8

where x(1) and x(5) are unbounded. It is clear that if x(1) < x(5) that this constraint will be undefined and will return NaN on most computers (if you are lucky). The constraint x(1) < x(5) is very hard to enforce at every iteration in an NLP solver so that this will not happen A better approach would be to add an extra varaible (say x(51) for an NLP with n == 50) and add an extra constraint:

log(x(51)) == 0
x(51) - (x(1) - x(5)) == 0
x(51) >= 1e-8

In the above expanded formulation the simple bound x(51) >= 1e-8 is easy to inforce and these undefined regions can be avoided. While the property that f(x), c(x) and h(x) being bounded for all x <: { x | xl <= x <= x} is a desireable properly, this is not required by this interface. As a result the client should be prepaired to deal with return values of NaN or Inf for f, c and h.

Definition at line 303 of file NLPInterfacePack_NLP.hpp.

Member Typedef Documentation

Definition at line 310 of file NLPInterfacePack_NLP.hpp.

Definition at line 314 of file NLPInterfacePack_NLP.hpp.

Constructor & Destructor Documentation

NLPInterfacePack::NLP::NLP ( )

Initialize to no reference set to calculation quanities.

Definition at line 69 of file NLPInterfacePack_NLP.cpp.

NLPInterfacePack::NLP::~NLP ( )
virtual

Destructor that cleans all the memory it owns.

Definition at line 74 of file NLPInterfacePack_NLP.cpp.

Member Function Documentation

value_type NLPInterfacePack::NLP::infinite_bound ( )
static

Value for an infinite bound.

Definition at line 61 of file NLPInterfacePack_NLP.cpp.

virtual void NLPInterfacePack::NLP::force_xinit_in_bounds ( bool  force_xinit_in_bounds)
pure virtual

Set if the initial point must be within the bounds.

This method must be called before this->initialize() is called.

Postconditions:

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual bool NLPInterfacePack::NLP::force_xinit_in_bounds ( ) const
pure virtual

Returns if the initial point must be within the bounds.

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

void NLPInterfacePack::NLP::set_options ( const options_ptr_t options)
virtual

Set the options that this NLP may be interested in.

Note that it is allowed for the client to alter *options.get() after this method is called so this had better read the options inside of the this->initialize() method.

The default implementation is to just ignore these options.

Note that if the subclass overrides this method then it must also override the get_options() method.

Definition at line 77 of file NLPInterfacePack_NLP.cpp.

const NLP::options_ptr_t & NLPInterfacePack::NLP::get_options ( ) const
virtual

Get the OptionsFromStream object being used to extract the options from.

The default implementation returns return.get() == NULL.

Reimplemented in NLPInterfacePack::NLPSerialPreprocessExplJac.

Definition at line 81 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::initialize ( bool  test_setup = false)
virtual

Initialize the NLP before it is used.

Postconditions:

Note that subclasses must call this function to reset what needs to be reset in this base object.

Reimplemented in NLPInterfacePack::NLPDirect, NLPInterfacePack::NLPSerialPreprocess, NLPInterfacePack::NLPSerialPreprocessExplJac, NLPInterfacePack::NLPFirstOrder, NLPInterfacePack::NLPSecondOrder, NLPInterfacePack::NLPObjGrad, and NLPInterfacePack::NLPBarrier.

Definition at line 86 of file NLPInterfacePack_NLP.cpp.

virtual bool NLPInterfacePack::NLP::is_initialized ( ) const
pure virtual
size_type NLPInterfacePack::NLP::n ( ) const
virtual

Return the number of variables.

Preconditions:

Default implementation returns this-space_x()->dim().

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 94 of file NLPInterfacePack_NLP.cpp.

size_type NLPInterfacePack::NLP::m ( ) const
virtual

Return the number of general equality constraints.

Preconditions:

Default implementation returns ( this->space_c().get() != NULL ? this-space_c()->dim() : 0 ).

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 100 of file NLPInterfacePack_NLP.cpp.

virtual vec_space_ptr_t NLPInterfacePack::NLP::space_x ( ) const
pure virtual

Vector space object for unknown variables x (dimension n).

Preconditions:

Postconditions:

  • return.get() != NULL

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual vec_space_ptr_t NLPInterfacePack::NLP::space_c ( ) const
pure virtual

Vector space object for general equality constraints c(x) (dimension m).

Preconditions:

Postconditions:

  • [this->m() > 0] return.get() != NULL
  • [this->m() == 0] return.get() == NULL

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual size_type NLPInterfacePack::NLP::num_bounded_x ( ) const
pure virtual

Returns the number of variables in x(i) for which xl(i)> -infinite_bound() or xu(i) < +infinite_bound().

Preconditions:

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual const Vector& NLPInterfacePack::NLP::xl ( ) const
pure virtual

Returns the lower bounds on the variables x.

Preconditions:

Any bounds that are non-existant will return this->xl().get_ele(i) == -NLP::infinite_bound().

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual const Vector& NLPInterfacePack::NLP::xu ( ) const
pure virtual

Returns a reference to the vector of upper bounds on the variables x.

Preconditions:

Any bounds that are non-existant will return this->xu().get_ele(i) == +NLP::infinite_bound().

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual value_type NLPInterfacePack::NLP::max_var_bounds_viol ( ) const
pure virtual

Set the maximum absolute value for which the variable bounds may be violated by when computing function and gradient values.

In other words the client should never never call on the NLP to compute a function and gradient evaluation outside of:

  xl - max_var_bounds_viol <= x <= xu + max_var_bounds_viol

Implemented in NLPInterfacePack::NLPBarrier.

virtual const Vector& NLPInterfacePack::NLP::xinit ( ) const
pure virtual

Returns a reference to the vector of the initial guess for the solution x.

Preconditions:

Postconditions:

  • return.space().is_compatible(*this->space_x()) == true)

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

void NLPInterfacePack::NLP::get_init_lagrange_mult ( VectorMutable lambda,
VectorMutable nu 
) const
virtual

Get the initial value of the Lagrange multipliers lambda.

By default this function just sets them to zero.

Parameters
lambda[out] Pointer to lagrange multipliers for equalities. lambda == NULL is allowed in which case it will not be set. Must have been created by this->space_c()->create_member(). Must be NULL if m() == 0.
nu[out] Pointer to lagrange multipliers for bounds. nu == NULL is allowed in which case it will not be set. Must have been created by this->space_x()->create_member(). Must be NULL if num_bounded_x() == 0.

Preconditions:

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 108 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::set_f ( value_type *  f)
virtual

Set a pointer to an value to be updated when this->calc_f() is called.

Parameters
f[in] Pointer to objective function value. May be NULL.

Preconditions:

Postconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 133 of file NLPInterfacePack_NLP.cpp.

value_type * NLPInterfacePack::NLP::get_f ( )
virtual

Return pointer passed to this->set_f().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 138 of file NLPInterfacePack_NLP.cpp.

value_type & NLPInterfacePack::NLP::f ( )
virtual

Returns non-const *this->get_f().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 143 of file NLPInterfacePack_NLP.cpp.

const value_type & NLPInterfacePack::NLP::f ( ) const
virtual

Returns const *this->get_f().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 148 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::set_c ( VectorMutable c)
virtual

Set a pointer to a vector to be updated when this->calc_c() is called.

Parameters
c[in] Pointer to constraint residual vector. May be NULL.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • [c != NULL] c->space().is_compatible(*this->space_c()) == true (throw VectorSpace::IncompatibleVectorSpaces)

Postconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 155 of file NLPInterfacePack_NLP.cpp.

VectorMutable * NLPInterfacePack::NLP::get_c ( )
virtual

Return pointer passed to this->set_c().

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 164 of file NLPInterfacePack_NLP.cpp.

VectorMutable & NLPInterfacePack::NLP::c ( )
virtual

Returns non-const *this->get_c().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 172 of file NLPInterfacePack_NLP.cpp.

const Vector & NLPInterfacePack::NLP::c ( ) const
virtual

Returns const *this->get_c().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 180 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::unset_quantities ( )
virtual

Call to unset all storage quantities (both in this class and all subclasses).

Preconditions:

Postconditions:

This method must be called by all subclasses that override it.

Reimplemented in NLPInterfacePack::NLPFirstOrder, NLPInterfacePack::NLPObjGrad, and NLPInterfacePack::NLPSecondOrder.

Definition at line 188 of file NLPInterfacePack_NLP.cpp.

virtual void NLPInterfacePack::NLP::scale_f ( value_type  scale_f)
pure virtual

Set the scaling of the objective function.

Preconditions:

Postconditions:

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual value_type NLPInterfacePack::NLP::scale_f ( ) const
pure virtual

Return the scaling being used for the objective function.

Preconditions:

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

void NLPInterfacePack::NLP::calc_f ( const Vector x,
bool  newx = true 
) const
virtual

Update the value for the objective f at the point x and put it in the stored reference.

Parameters
x[in] Point at which to calculate the object function f.
newx[in] (default true) If false, the values in x are assumed to be the same as the last call to a this->calc_*(x,newx) member. If true, the values in x are assumed to not be the same as the last call to a this->calc_*(x,newx) member.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • x.space().is_compatible(*this->space_x()) == true (throw VectorSpace::IncompatibleVectorSpaces)
  • this->get_f() != NULL (throw NoRefSet)

Postconditions:

  • this->f() is updated to f(x)

The storage reference for c may also be updated at this point (if get_c() != NULL) but is not guarentied to be. But no other quanities from possible subclasses are allowed to be updated as a side effect.

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 198 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::calc_c ( const Vector x,
bool  newx = true 
) const
virtual

Update the constraint residual vector for c at the point x and put it in the stored reference.

Parameters
x[in] Point at which to calculate residual to the equality constraints c.
newx[in] (default true) If false, the values in x are assumed to be the same as the last call to a this->calc_*(x,newx) member. If true, the values in x are assumed to not be the same as the last call to a this->calc_*(x,newx) member.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • x.space().is_compatible(*this->space_x()) == true (throw VectorSpace::IncompatibleVectorSpaces)
  • this->get_c() != NULL (throw NoRefSet)

Postconditions:

  • this->c() is updated to c(x)

The storage reference for f may also be updated at this point (if get_f() != NULL) but is not guarentied to be. But no other quanities from possible subclasses are allowed to be updated as a side effect.

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 205 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::report_final_solution ( const Vector x,
const Vector lambda,
const Vector nu,
bool  is_optimal 
)
virtual

Used by the solver to report the final solution and multipliers.

Call this function to report the final solution of the unknows x and the Lagrange multipliers for the equality constriants lambda and the varaible bounds nu. If any of the Lagrange multipliers are not known then you can pass NULL in for them.

Preconditions:

The default behavior is to just ignore this.

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 215 of file NLPInterfacePack_NLP.cpp.

size_type NLPInterfacePack::NLP::num_f_evals ( ) const
virtual

Gives the number of object function f(x) evaluations called by the solver since initialize() was called.

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 225 of file NLPInterfacePack_NLP.cpp.

size_type NLPInterfacePack::NLP::num_c_evals ( ) const
virtual

Gives the number of constraint function c(x) evaluations called by the solver since initialize() was called. Throws exception if this->m() == 0.

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 230 of file NLPInterfacePack_NLP.cpp.

size_type NLPInterfacePack::NLP::ns ( ) const
virtual

Return the number of slack variables (i.e. number of general inequalities).

Default implementation returns (this->space_h_breve().get() ? this->space_h_breve()->dim() : 0).

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 240 of file NLPInterfacePack_NLP.cpp.

virtual vec_space_ptr_t NLPInterfacePack::NLP::space_c_breve ( ) const
virtual

Vector space object for the original equalities c_breve(x_breve)

Preconditions:

Postconditions:

  • [this->m() - this->ns() > 0] return.get() != NULL
  • [this->m() - this->ns() == 0] return.get() == NULL

The default implementation returns this->space_c().

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual vec_space_ptr_t NLPInterfacePack::NLP::space_h_breve ( ) const
virtual

Vector space object for the original inequalities h_breve(x_breve)

Preconditions:

Postconditions:

  • [this->ns() > 0] return.get() != NULL
  • [this->ns() == 0] return.get() == NULL

The default implementation returns return.get() == NULL.

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

const Vector & NLPInterfacePack::NLP::hl_breve ( ) const
virtual

Returns a reference to the vector of lower bounds on the general inequality constraints h_breve(x_breve).

Preconditions:

Any bounds that are non-existant will return this->hl_breve().get_ele(i) == -NLP::infinite_bound().

The default implementation throws an exception.

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 256 of file NLPInterfacePack_NLP.cpp.

const Vector & NLPInterfacePack::NLP::hu_breve ( ) const
virtual

Returns a reference to the vector of upper bounds on the general inequality constraints h_breve(x_breve).

Preconditions:

Any bounds that are non-existant will return this->hu_breve().get_ele(i) == +NLP::infinite_bound().

The default implementation throws an exception.

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 271 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::set_c_breve ( VectorMutable c_breve)
virtual

Set a pointer to a vector to be updated when this->calc_c_breve() is called.

Parameters
c_breve[in] Pointer to constraint residual vector. May be NULL.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • [c != NULL] c->space().is_compatible(*this->space_c_breve()) == true (throw VectorSpace::IncompatibleVectorSpaces)

Postconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 286 of file NLPInterfacePack_NLP.cpp.

VectorMutable * NLPInterfacePack::NLP::get_c_breve ( )
virtual

Return pointer passed to this->set_c_breve().

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 295 of file NLPInterfacePack_NLP.cpp.

VectorMutable & NLPInterfacePack::NLP::c_breve ( )
virtual

Returns non-const *this->get_c_breve().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 303 of file NLPInterfacePack_NLP.cpp.

const Vector & NLPInterfacePack::NLP::c_breve ( ) const
virtual

Returns const *this->get_c_breve().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 311 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::set_h_breve ( VectorMutable h_breve)
virtual

Set a pointer to a vector to be updated when this->calc_h_breve() is called.

Parameters
h_breve[in] Pointer to constraint residual vector. May be NULL.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • [c != NULL] c->space().is_compatible(*this->space_h_breve()) == true (throw VectorSpace::IncompatibleVectorSpaces)

Postconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 319 of file NLPInterfacePack_NLP.cpp.

VectorMutable * NLPInterfacePack::NLP::get_h_breve ( )
virtual

Return pointer passed to this->set_h_breve().

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 328 of file NLPInterfacePack_NLP.cpp.

VectorMutable & NLPInterfacePack::NLP::h_breve ( )
virtual

Returns non-const *this->get_h_breve().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 336 of file NLPInterfacePack_NLP.cpp.

const Vector & NLPInterfacePack::NLP::h_breve ( ) const
virtual

Returns const *this->get_h_breve().

Preconditions:

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 344 of file NLPInterfacePack_NLP.cpp.

virtual const Permutation& NLPInterfacePack::NLP::P_var ( ) const
virtual

Return the permutation object for the variables.

Preconditions:

Postconditions:

  • return.space()is_compatible(*this->space_x()) == true

The default returns return.is_identity() == true

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual const Permutation& NLPInterfacePack::NLP::P_equ ( ) const
virtual

Return the permutation object for the constraints.

Preconditions:

Postconditions:

  • return.space()is_compatible(*this->space_c()) == true

The default returns return.is_identity() == true

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

void NLPInterfacePack::NLP::calc_c_breve ( const Vector x,
bool  newx = true 
) const
virtual

Update the constraint residual vector for c_breve at the point x and put it in the stored reference.

Parameters
x[in] Point at which to calculate residual to the equality constraints c_breve.
newx[in] (default true) If true, the values in x are the same as the last call to a this->calc_*(x,newx) member. If false, the values in x are not the same as the last call to a this->calc_*(x,newx) member.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • x.space().is_compatible(*this->space_x()) == true (throw VectorSpace::IncompatibleVectorSpaces)
  • this->get_c_breve() != NULL (throw NoRefSet)

Postconditions:

  • this->c_breve() is updated to c_breve(x_breve)

The storage reference for f and/or h_breve may also be updated at this point (if get_f() != NULL and/or get_h_breve() != NULL) but is not guarentied to be. But no other quanities from possible subclasses are allowed to be updated as a side effect.

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 366 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::calc_h_breve ( const Vector x,
bool  newx = true 
) const
virtual

Update the constraint residual vector for h_breve at the point x and put it in the stored reference.

Parameters
x[in] Point at which to calculate residual to the equality constraints h_breve.
newx[in] (default true) If true, the values in x are the same as the last call to a this->calc_*(x,newx) member. If false, the values in x are not the same as the last call to a this->calc_*(x,newx) member.

Preconditions:

  • this->is_initialized() == true (throw NotInitialized)
  • x.space().is_compatible(*this->space_x()) == true (throw VectorSpace::IncompatibleVectorSpaces)
  • this->get_h_breve() != NULL (throw NoRefSet)

Postconditions:

  • this->h_breve() is updated to h_breve(x_breve)

The storage reference for f and/or c_breve may also be updated at this point (if get_f() != NULL and/or get_c_breve() != NULL) but is not guarentied to be. But no other quanities from possible subclasses are allowed to be updated as a side effect.

Reimplemented in NLPInterfacePack::NLPBarrier.

Definition at line 376 of file NLPInterfacePack_NLP.cpp.

const NLP::ZeroOrderInfo NLPInterfacePack::NLP::zero_order_info ( ) const
inline

Return pointer to set quantities.

Definition at line 1128 of file NLPInterfacePack_NLP.hpp.

const NLP::ZeroOrderInfo NLPInterfacePack::NLP::zero_order_info_breve ( ) const
inline

Return pointer to set hat quantities.

Definition at line 1134 of file NLPInterfacePack_NLP.hpp.

virtual void NLPInterfacePack::NLP::imp_calc_f ( const Vector x,
bool  newx,
const ZeroOrderInfo zero_order_info 
) const
protectedpure virtual

Overridden to compute f(x) (and perhaps other quantities if set).

Preconditions:

  • x.space().is_compatible(*this->space_x()) (throw IncompatibleType)
  • zero_order_info.f != NULL (throw std::invalid_argument)

Postconditions:

  • *zero_order_info.f is updated to f(x).
Parameters
x[in] Unknown vector (size n).
newx[in] True if is a new point.
zero_order_info[out] Pointers to f, c and h. On output, *zero_order_info.f is updated to f(x) If this->multi_calc() == true then any of the other quantities pointed to in zero_order_info may be set on output, but are not guaranteed to be.

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

virtual void NLPInterfacePack::NLP::imp_calc_c ( const Vector x,
bool  newx,
const ZeroOrderInfo zero_order_info 
) const
protectedpure virtual

Overridden to compute c(x) and perhaps f(x) and/or h(x) (if multiple calculaiton = true).

Preconditions:

  • x.space().is_compatible(*this->space_x()) (throw IncompatibleType)
  • zero_order_info.c != NULL (throw std::invalid_argument)

Postconditions:

  • *zero_order_info.c is updated to c(x).
Parameters
x[in] Unknown vector (size n).
newx[in] True if is a new point.
zero_order_info[out] Pointers to f, c and h. On output, *zero_order_info.c is updated to c(x) If this->multi_calc() == true then any of the other quantities pointed to in zero_order_info may be set on output, but are not guaranteed to be.

Implemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

void NLPInterfacePack::NLP::imp_calc_c_breve ( const Vector x,
bool  newx,
const ZeroOrderInfo zero_order_info_breve 
) const
protectedvirtual

Overridden to compute c_breve(x_breve) and perhaps f(x) and/or h_breve(x_breve)

Preconditions:

  • x.space().is_compatible(*this->space_x()) (throw IncompatibleType)
  • zero_order_info.c != NULL (throw std::invalid_argument)

Postconditions:

  • *zero_order_info.c is updated to c_breve(x_breve).
Parameters
x[in] Unknown vector (size n).
newx[in] True if is a new point.
zero_order_info_breve[out] Pointers to f, c_breve and h_breve. On output, *zero_order_info.c is updated to c_breve(x_breve)

The default implementation calls this->imp_calc_c().

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 388 of file NLPInterfacePack_NLP.cpp.

void NLPInterfacePack::NLP::imp_calc_h_breve ( const Vector x,
bool  newx,
const ZeroOrderInfo zero_order_info_breve 
) const
protectedvirtual

Overridden to compute h_breve(x_breve) and perhaps f(x) and/or c_breve(x_breve).

Preconditions:

  • x.space().is_compatible(*this->space_x()) (throw IncompatibleType)
  • zero_order_info.h != NULL (throw std::invalid_argument)

Postconditions:

  • *zero_order_info.h is updated to h_breve(x_breve).
Parameters
x[in] Unknown vector (size n).
newx[in] True if is a new point.
zero_order_info_breve[out] Pointers to f, c_breve and h_breve. On output, *zero_order_info.h is updated to h_breve(x_breve)

The default implementation throws an exception.

Reimplemented in NLPInterfacePack::NLPSerialPreprocess, and NLPInterfacePack::NLPBarrier.

Definition at line 397 of file NLPInterfacePack_NLP.cpp.

template<class T >
void NLPInterfacePack::NLP::assert_ref_set ( T *  p,
std::string  info 
) const
inlineprotected

Assert referece has been set for a quanity.

Definition at line 1097 of file NLPInterfacePack_NLP.hpp.


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