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 Member Functions | Static Public Member Functions | List of all members
NLPInterfacePack::NLPSerialPreprocess Class Referenceabstract

NLP node implementation subclass for preprocessing and basis manipulation. More...

#include <NLPInterfacePack_NLPSerialPreprocess.hpp>

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

Classes

class  InconsistantBounds
 Thrown if xl(i) > xu(i) More...
 
struct  ObjGradInfoSerial
 Struct for serial gradient (objective), objective and constriants (pointers) More...
 
struct  ZeroOrderInfoSerial
 Struct for objective and constriants (pointer) as serial vectors. More...
 

Public Member Functions

 NLPSerialPreprocess ()
 Default Constructor. More...
 
- Public Member Functions inherited from NLPInterfacePack::NLPObjGrad
 NLPObjGrad ()
 Initialize to no reference set to calculation quanities. More...
 
virtual bool supports_Gf () const
 Determine if the objective gradient is supported or not. More...
 
virtual bool supports_Gf_prod () const
 Determine if the objective gradient product is supported or not. More...
 
virtual void set_Gf (VectorMutable *Gf)
 Set a pointer to a vector to be updated when this->calc_Gf() is called. More...
 
virtual VectorMutableget_Gf ()
 Return pointer passed to this->set_Gf(). More...
 
virtual VectorMutableGf ()
 Returns non-const *this->get_Gf(). More...
 
virtual const VectorGf () const
 Returns const *this->get_Gf(). More...
 
void unset_quantities ()
 Call to unset all storage quantities (both in this class and all subclasses). More...
 
virtual void calc_Gf (const Vector &x, bool newx=true) const
 Update the vector for Gf at the point x and put it in the stored reference. More...
 
virtual value_type calc_Gf_prod (const Vector &x, const Vector &d, bool newx=true) const
 Calculate the inner product Gf(x)'*d at the point x and put it in the stored reference. More...
 
virtual size_type num_Gf_evals () const
 Objective gradient evaluations count. More...
 
- Public Member Functions inherited from NLPInterfacePack::NLP
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...
 
 NLP ()
 Initialize to no reference set to calculation quanities. More...
 
virtual ~NLP ()
 Destructor that cleans all the memory it owns. 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 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...
 
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...
 
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...
 
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...
 
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...
 
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 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...
 
- Public Member Functions inherited from NLPInterfacePack::NLPVarReductPerm

Static Public Member Functions

static value_type fixed_var_mult ()
 Gives the value of a Lagrange multipler for a fixed variable bound .that has been preprocessed out of the problem. More...
 
- Static Public Member Functions inherited from NLPInterfacePack::NLP
static value_type infinite_bound ()
 Value for an infinite bound. More...
 

Overridden public members from NLP

void force_xinit_in_bounds (bool force_xinit_in_bounds)
 
bool force_xinit_in_bounds () const
 
void initialize (bool test_setup)
 
bool is_initialized () const
 
size_type n () const
 
size_type m () const
 
vec_space_ptr_t space_x () const
 
vec_space_ptr_t space_c () const
 
size_type num_bounded_x () const
 
const Vectorxl () const
 
const Vectorxu () const
 
const Vectorxinit () const
 
void get_init_lagrange_mult (VectorMutable *lambda, VectorMutable *nu) const
 
void scale_f (value_type scale_f)
 
value_type scale_f () const
 
void report_final_solution (const Vector &x, const Vector *lambda, const Vector *nu, bool is_optimal)
 Overridden to permute the variables back into an order that is natural to the subclass. More...
 
virtual size_type ns () const
 
vec_space_ptr_t space_c_breve () const
 
vec_space_ptr_t space_h_breve () const
 
const Vectorhl_breve () const
 
const Vectorhu_breve () const
 
const Permutation & P_var () const
 
const Permutation & P_equ () const
 

Overridden public members from NLPVarReductPerm

const perm_fcty_ptr_t factory_P_var () const
 
const perm_fcty_ptr_t factory_P_equ () const
 
Range1D var_dep () const
 
Range1D var_indep () const
 
Range1D equ_decomp () const
 
Range1D equ_undecomp () const
 
bool nlp_selects_basis () const
 
bool get_next_basis (Permutation *P_var, Range1D *var_dep, Permutation *P_equ, Range1D *equ_decomp)
 
void set_basis (const Permutation &P_var, const Range1D &var_dep, const Permutation *P_equ, const Range1D *equ_decomp)
 
void get_basis (Permutation *P_var, Range1D *var_dep, Permutation *P_equ, Range1D *equ_decomp) const
 

Overridden protected members from NLP

void imp_calc_f (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info) const
 
void imp_calc_c (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info) const
 
void imp_calc_c_breve (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info_breve) const
 
void imp_calc_h_breve (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info_breve) const
 

Overridden protected members from NLPObjGrad

void imp_calc_Gf (const Vector &x, bool newx, const ObjGradInfo &obj_grad_info) const
 

Pure virtual methods to be defined by subclasses

virtual bool imp_nlp_has_changed () const
 Return if the definition of the NLP has changed since the last call to initialize() More...
 
virtual size_type imp_n_orig () const =0
 Return the number of variables in the original problem (including those fixed by bounds) More...
 
virtual size_type imp_m_orig () const =0
 Return the number of general equality constraints in the original problem. More...
 
virtual size_type imp_mI_orig () const =0
 Return the number of general inequality constraints in the original problem. More...
 
virtual const DVectorSlice imp_xinit_orig () const =0
 Return the original initial point (size imp_n_orig()). More...
 
virtual bool imp_has_var_bounds () const =0
 Return if the NLP has bounds. More...
 
virtual const DVectorSlice imp_xl_orig () const =0
 Return the original lower variable bounds (size imp_n_orig()). More...
 
virtual const DVectorSlice imp_xu_orig () const =0
 Return the original upper variable bounds (size imp_n_orig()). More...
 
virtual const DVectorSlice imp_hl_orig () const =0
 Return the original lower general inequality bounds (size imp_mI_orig()). More...
 
virtual const DVectorSlice imp_hu_orig () const =0
 Return the original upper general inequality bounds (size imp_mI_orig()). More...
 
virtual void imp_calc_f_orig (const DVectorSlice &x_full, bool newx, const ZeroOrderInfoSerial &zero_order_info) const =0
 Calculate the objective function for the original NLP. More...
 
virtual void imp_calc_c_orig (const DVectorSlice &x_full, bool newx, const ZeroOrderInfoSerial &zero_order_info) const =0
 Calculate the vector for all of the general equality constaints in the original NLP. More...
 
virtual void imp_calc_h_orig (const DVectorSlice &x_full, bool newx, const ZeroOrderInfoSerial &zero_order_info) const =0
 Calculate the vector for all of the general inequality constaints in the original NLP. More...
 
virtual void imp_calc_Gf_orig (const DVectorSlice &x_full, bool newx, const ObjGradInfoSerial &obj_grad_info) const =0
 Calculate the vector for the gradient of the objective in the original NLP. More...
 
virtual bool imp_get_next_basis (IVector *var_perm_full, IVector *equ_perm_full, size_type *rank_full, size_type *rank)
 Return the next basis selection (default returns false). More...
 
virtual void imp_report_orig_final_solution (const DVectorSlice &x_full, const DVectorSlice *lambda_orig, const DVectorSlice *lambdaI_orig, const DVectorSlice *nu_orig, bool optimal)
 To be overridden by subclasses to report the final solution in the original ordering natural to the subclass. More...
 

Other protected implementation functions for subclasses to call

void set_not_initialized ()
 Used by subclasses to set the state of the NLP to not initialized. More...
 
void assert_initialized () const
 Assert if we have been initizlized (throws UnInitialized) More...
 
void set_x_full (const DVectorSlice &x, bool newx, DVectorSlice *x_full) const
 Set the full x vector if newx == true More...
 
DVectorSlice x_full () const
 Give reference to current x_full. More...
 
const ZeroOrderInfoSerial zero_order_orig_info () const
 
const ObjGradInfoSerial obj_grad_orig_info () const
 
const IVector & var_remove_fixed_to_full () const
 Permutation vector for partitioning free and fixed variables. More...
 
const IVector & var_full_to_remove_fixed () const
 Inverse permutation vector of var_remove_fixed_to_full(). More...
 
const IVector & var_perm () const
 Permutes from the compated variable vector (removing fixed variables) to the current basis selection. More...
 
const IVector & equ_perm () const
 Permutes from the original constriant ordering to the current basis selection. More...
 
const IVector & inv_equ_perm () const
 Inverse of equ_perm() More...
 
void var_from_full (DVectorSlice::const_iterator vec_full, DVectorSlice::iterator vec) const
 
void var_to_full (DVectorSlice::const_iterator vec, DVectorSlice::iterator vec_full) const
 
void equ_from_full (const DVectorSlice &c_orig, const DVectorSlice &h_orig, const DVectorSlice &s_orig, DVectorSlice *c_full) const
 

Additional Inherited Members

- Public Types inherited from NLPInterfacePack::NLP
typedef Teuchos::RCP< const
VectorSpace > 
vec_space_ptr_t
 
typedef Teuchos::RCP< const
OptionsFromStreamPack::OptionsFromStream
options_ptr_t
 
- Public Types inherited from NLPInterfacePack::NLPVarReductPerm
typedef Teuchos::RCP< const
Teuchos::AbstractFactory
< Permutation > > 
perm_fcty_ptr_t
 
- Protected Member Functions inherited from NLPInterfacePack::NLPObjGrad
const ObjGradInfo obj_grad_info () const
 Return objective gradient and zero order information. More...
 
- Protected Member Functions inherited from NLPInterfacePack::NLP
template<class T >
void assert_ref_set (T *p, std::string info) const
 Assert referece has been set for a quanity. More...
 

Detailed Description

NLP node implementation subclass for preprocessing and basis manipulation.

This is an implementation node class that takes an original NLP and transforms it by:

Original NLP formulation

The original NLP (as specified by the subclass) takes the form:

    min    f_orig(x_orig)
    s.t.   c_orig(x_orig) = 0
           hl_orig <= h(x_orig) <= hu_orig
           xl_orig <= x_orig <= xu_orig
 where:
         x_orig         <: REAL^n_orig
         f_orig(x_orig) <: REAL^n_orig -> REAL
         c_orig(x_orig) <: REAL^n_orig -> REAL^m_orig
         h_orig(x_orig) <: REAL^n_orig -> REAL^mI_orig

Conversion of general inequalities to equalities using slack variables

The original NLP formulation above is transformed by adding slack variables s_orig <: REAL^mI_orig, defining a new x_full = [ x_orig; s_orig ] and forming the new NLP:

    min    f_full(x_full)
    s.t.   c_full(x_full) = 0
           xl_full <= x_full <= xu_full

 where:

         x_full           = [ x_orig ]  n_orig
                            [ s_orig ]  mI_orig

         f_full(x_full)   = f_orig(x_orig)

         c_full(x_full)   = [ c_orig(x_orig)          ]  m_orig
                            [ h_orig(x_orig) - s_orig ]  mI_orig

         xl_full          = [ xl_orig ]  n_orig
                            [ hl_orig ]  mI_orig

         xu_full          = [ xu_orig ]  n_orig
                            [ hu_orig ]  mI_orig

Note that in this case, the Jacobian of the new equality constraints and the gradient of the new objective become:

   Gc_full = [  Gc_orig    Gh_orig   ]  n_orig
             [     0          -I     ]  mI_orig

                 m_orig     mI_orig

   Gf_full = [ Gf_orig ]  n_orig
             [    0    ]  mI_orig

It is up to the subclass to implement imp_calc_Gc() and imp_calc_Gh() in a way that is consistent with the above transformation while also considering basis permutations (see NLPSerialPreprocessExplJac). As for the gradient Gc_full, the subclass can actually include terms for the slack variables in the objective function but the most common behavior will be to just ignore slack variables in the subclass.

Preprocessing and basis manipulation

The initial basis selection is the original order (x_full = [ x_orig; s_orig ]) with the variables fixed by bounds being removed, and assumes there are no dependent equations (r == m).

The implementations of the Jacobian matrices Gc and Gh are not determined here and must be defined by an NLP subclass (see NLPSerialPreprocessExplJac for example).

This class stores the variable permutations and processing information in two parts. In the first state, the fixed variables are removed as:

   var_remove_fixed_to_full =  [ not fixed by bounds |  fixed by bounds  ]
                               [1 ..                n|n+1 ..        n_full]

The mapping i_full = var_remove_fixed_to_full()(i_free_fixed) gives the index of the original variable (i_full) for the sets of variables not fixed and fixed by bounds.

The inverse mapping i_free_fixed = var_full_to_remove_fixed()(i_full) can be used to determine if a variable is fixed by bounds or not..

On top of this partitioning of free and fixed variables, there is a second stage which is a permutation of the free variables into dependent and independent sets that is needed by the client.

   var_perm = [ dependent variables | independent variables ]
              [1..               n-r|n-r+1...              n]

The mapping i_free_fixed = var_perm()(i_perm) is used to determine the index of a free variable in var_remove_fixed_to_full() given its index (i_perm) for the current basis selection.

For example, if x is the vector of variables for the current basis selection and x_full is the vector of variables in the original order including the fixed variables then the following is true:

x(i) == x_full(var_remove_fixed_to_full()(var_perm()(i))), for i = 1...n

The permutation equ_perm() gives the partitioning of the equality constraints into decomposed and undecomposed equalities. Decomposed inequality constraints are not supported currently.

Subclass developers notes

Handling of multiple updates by subclasses: Here we discuss the protocol for the handling of multiple updates to quantities during the calculation of other quantities. In order to simplify the implementation of subclasses as much as possible, storage for all iteration quantities will be passed to the subclass in the methods imp_calc_f_orig(), imp_calc_c_orig(), imp_calc_h_orig() and imp_calc_Gf_orig() regardless of what quantities where set by the user in the NLP interface. The subclass can always find out what was set by the client by calling get_f(), get_c(), get_Gf() etc. Therefore, in general, clients should just only compute what is required in each call to imp_calc_xxx_orig() and only update other quantities if it is absolutely free to do so (e.g. computing a function value when a gradient is computed using AD) or is required to do so (e.g.an external interface that forces both f_orig(x_orig), c_orig(x_orig) and h_orig(x_orig) be computed at the same time). It is up to the subclass to remember when a quantity has already been computed so that it will not be computed again unnecessarily. It is always safe for the subclass to ignore these issues and just do what is easiest. More careful implementations can be handled by the subclass by keeping track of get_xxx() and newx and remembering when quantities are computed.

The following methods from the NLP interface must be overridden by the NLP subclass: max_var_bounds_viol(), set_multi_calc(), multi_calc().

The following methods from the NLPVarReductPerm interface must be overridden by the NLP subclass: nlp_selects_basis().

In addition, the methods from this interface that must be overridden are: imp_n_orig(), imp_m_orig(), imp_mI_orig(), imp_xinit_orig(), imp_has_var_bounds(), imp_xl_orig(), imp_xu_orig(), imp_hl_orig(), imp_hu_orig(), imp_calc_f_orig(), imp_calc_c_orig(), imp_calc_h_orig() and imp_calc_Gf_orig().

The NLP method initialize() should also be overridden by all of the subclasses (and call initialize() on its direct subclass).

The following methods (with default implementations) may also be overridden by a subclass: imp_get_next_basis() and imp_report_orig_final_solution().

Definition at line 221 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

Constructor & Destructor Documentation

NLPInterfacePack::NLPSerialPreprocess::NLPSerialPreprocess ( )

Default Constructor.

This initalizes the basis to the first basis if the subclass specifies one and if not picks to first r variables as the dependent variables and the last n-r variables as the independent variables. Also the default behavior is to force the initial point in bounds.

Definition at line 85 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

Member Function Documentation

value_type NLPInterfacePack::NLPSerialPreprocess::fixed_var_mult ( )
static

Gives the value of a Lagrange multipler for a fixed variable bound .that has been preprocessed out of the problem.

Definition at line 78 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

void NLPInterfacePack::NLPSerialPreprocess::force_xinit_in_bounds ( bool  force_xinit_in_bounds)
virtual

Implements NLPInterfacePack::NLP.

Definition at line 96 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

bool NLPInterfacePack::NLPSerialPreprocess::force_xinit_in_bounds ( ) const
virtual

Implements NLPInterfacePack::NLP.

Definition at line 101 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

void NLPInterfacePack::NLPSerialPreprocess::initialize ( bool  test_setup)
virtual
bool NLPInterfacePack::NLPSerialPreprocess::is_initialized ( ) const
virtual
size_type NLPInterfacePack::NLPSerialPreprocess::n ( ) const
virtual

Reimplemented from NLPInterfacePack::NLP.

Definition at line 338 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 344 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Implements NLPInterfacePack::NLP.

Definition at line 350 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Implements NLPInterfacePack::NLP.

Definition at line 356 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

size_type NLPInterfacePack::NLPSerialPreprocess::num_bounded_x ( ) const
virtual

Implements NLPInterfacePack::NLP.

Definition at line 362 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

const Vector & NLPInterfacePack::NLPSerialPreprocess::xl ( ) const
virtual

Implements NLPInterfacePack::NLP.

Definition at line 367 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

const Vector & NLPInterfacePack::NLPSerialPreprocess::xu ( ) const
virtual

Implements NLPInterfacePack::NLP.

Definition at line 373 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

const Vector & NLPInterfacePack::NLPSerialPreprocess::xinit ( ) const
virtual

Implements NLPInterfacePack::NLP.

Definition at line 379 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 385 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

void NLPInterfacePack::NLPSerialPreprocess::scale_f ( value_type  scale_f)
virtual

Implements NLPInterfacePack::NLP.

Definition at line 398 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

value_type NLPInterfacePack::NLPSerialPreprocess::scale_f ( ) const
virtual

Implements NLPInterfacePack::NLP.

Definition at line 404 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Overridden to permute the variables back into an order that is natural to the subclass.

The default implementation of this function is to call the method imp_report_full_final_solution(x_full,lambda_full,nu_full). This function translates from x, lambda and nu into the original order with fixed variables added back to form x_full, lambda_full, lambdaI_full and nu_full.

Reimplemented from NLPInterfacePack::NLP.

Definition at line 410 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 456 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 463 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 470 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 477 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 483 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 489 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 495 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

const NLPVarReductPerm::perm_fcty_ptr_t NLPInterfacePack::NLPSerialPreprocess::factory_P_var ( ) const
virtual
const NLPVarReductPerm::perm_fcty_ptr_t NLPInterfacePack::NLPSerialPreprocess::factory_P_equ ( ) const
virtual
Range1D NLPInterfacePack::NLPSerialPreprocess::var_dep ( ) const
virtual
Range1D NLPInterfacePack::NLPSerialPreprocess::var_indep ( ) const
virtual
Range1D NLPInterfacePack::NLPSerialPreprocess::equ_decomp ( ) const
virtual
Range1D NLPInterfacePack::NLPSerialPreprocess::equ_undecomp ( ) const
virtual
bool NLPInterfacePack::NLPSerialPreprocess::nlp_selects_basis ( ) const
virtual
bool NLPInterfacePack::NLPSerialPreprocess::get_next_basis ( Permutation *  P_var,
Range1D *  var_dep,
Permutation *  P_equ,
Range1D *  equ_decomp 
)
virtual
void NLPInterfacePack::NLPSerialPreprocess::set_basis ( const Permutation &  P_var,
const Range1D &  var_dep,
const Permutation *  P_equ,
const Range1D *  equ_decomp 
)
virtual
void NLPInterfacePack::NLPSerialPreprocess::get_basis ( Permutation *  P_var,
Range1D *  var_dep,
Permutation *  P_equ,
Range1D *  equ_decomp 
) const
virtual
void NLPInterfacePack::NLPSerialPreprocess::imp_calc_f ( const Vector x,
bool  newx,
const ZeroOrderInfo zero_order_info 
) const
protectedvirtual

Implements NLPInterfacePack::NLP.

Definition at line 650 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

void NLPInterfacePack::NLPSerialPreprocess::imp_calc_c ( const Vector x,
bool  newx,
const ZeroOrderInfo zero_order_info 
) const
protectedvirtual

Implements NLPInterfacePack::NLP.

Definition at line 663 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 685 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

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

Reimplemented from NLPInterfacePack::NLP.

Definition at line 702 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

void NLPInterfacePack::NLPSerialPreprocess::imp_calc_Gf ( const Vector x,
bool  newx,
const ObjGradInfo obj_grad_info 
) const
protectedvirtual
virtual bool NLPInterfacePack::NLPSerialPreprocess::imp_nlp_has_changed ( ) const
inlineprotectedvirtual

Return if the definition of the NLP has changed since the last call to initialize()

The default return is true. This function is present in order to avoid preprocessing when initialize() is called but nothing has changed.

Definition at line 460 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

virtual size_type NLPInterfacePack::NLPSerialPreprocess::imp_n_orig ( ) const
protectedpure virtual

Return the number of variables in the original problem (including those fixed by bounds)

virtual size_type NLPInterfacePack::NLPSerialPreprocess::imp_m_orig ( ) const
protectedpure virtual

Return the number of general equality constraints in the original problem.

virtual size_type NLPInterfacePack::NLPSerialPreprocess::imp_mI_orig ( ) const
protectedpure virtual

Return the number of general inequality constraints in the original problem.

virtual const DVectorSlice NLPInterfacePack::NLPSerialPreprocess::imp_xinit_orig ( ) const
protectedpure virtual

Return the original initial point (size imp_n_orig()).

virtual bool NLPInterfacePack::NLPSerialPreprocess::imp_has_var_bounds ( ) const
protectedpure virtual

Return if the NLP has bounds.

virtual const DVectorSlice NLPInterfacePack::NLPSerialPreprocess::imp_xl_orig ( ) const
protectedpure virtual

Return the original lower variable bounds (size imp_n_orig()).

Only to be called if this->imp_has_var_bounds() == true. A lower bound is considered free if it is less than or equal to:

<tt>-NLP::infinite_bound()</tt>
virtual const DVectorSlice NLPInterfacePack::NLPSerialPreprocess::imp_xu_orig ( ) const
protectedpure virtual

Return the original upper variable bounds (size imp_n_orig()).

Only to be called if this->imp_has_var_bounds() == true. An upper bound is considered free if it is greater than or equal to:

<tt>+NLP::infinite_bound()</tt>
virtual const DVectorSlice NLPInterfacePack::NLPSerialPreprocess::imp_hl_orig ( ) const
protectedpure virtual

Return the original lower general inequality bounds (size imp_mI_orig()).

Only to be called if this->imp_mI_orig() == true. A lower bound is considered free if it is equal to:

-NLP::infinite_bound()

virtual const DVectorSlice NLPInterfacePack::NLPSerialPreprocess::imp_hu_orig ( ) const
protectedpure virtual

Return the original upper general inequality bounds (size imp_mI_orig()).

Only to be called if this->imp_mI_orig() == true. An upper bound is considered free if it is equal to:

+NLP::infinite_bound()

virtual void NLPInterfacePack::NLPSerialPreprocess::imp_calc_f_orig ( const DVectorSlice &  x_full,
bool  newx,
const ZeroOrderInfoSerial zero_order_info 
) const
protectedpure virtual

Calculate the objective function for the original NLP.

virtual void NLPInterfacePack::NLPSerialPreprocess::imp_calc_c_orig ( const DVectorSlice &  x_full,
bool  newx,
const ZeroOrderInfoSerial zero_order_info 
) const
protectedpure virtual

Calculate the vector for all of the general equality constaints in the original NLP.

virtual void NLPInterfacePack::NLPSerialPreprocess::imp_calc_h_orig ( const DVectorSlice &  x_full,
bool  newx,
const ZeroOrderInfoSerial zero_order_info 
) const
protectedpure virtual

Calculate the vector for all of the general inequality constaints in the original NLP.

virtual void NLPInterfacePack::NLPSerialPreprocess::imp_calc_Gf_orig ( const DVectorSlice &  x_full,
bool  newx,
const ObjGradInfoSerial obj_grad_info 
) const
protectedpure virtual

Calculate the vector for the gradient of the objective in the original NLP.

Note that the dimension of obj_grad_info.Gf->dim() is n_orig + mI_orig.

On input, if mI_orig > 0 then (*obj_grad_info.Gf)(n_orig+1,n_orig+mI_orig) is initialized to 0.0 (since slacks do not ordinarily do not appear in the objective function). However, the subclass can assign (smooth) contributions for the slacks if desired.

bool NLPInterfacePack::NLPSerialPreprocess::imp_get_next_basis ( IVector *  var_perm_full,
IVector *  equ_perm_full,
size_type rank_full,
size_type rank 
)
protectedvirtual

Return the next basis selection (default returns false).

Parameters
var_perm_fullout. Contains the variable permutations (including slack variables possibly).
equ_perm_full[out] (size = m_orig + mI_orig)). Contains the constriant permutations (including general inequalities possibly).
rank_full[out] Returns the rank of the basis before fixed variables are taken out of var_perm_full and equ_perm.
rank[out] Returns the rank of the basis after fixed variables are taken out of var_perm_full and equ_perm. If there are no fixed variables then rank should be equal to rank_full.

Postconditions:

  • [return == true] var_perm_full(i) < var_perm_full(i+1), for i = 1...rank_full-1
  • [return == true] var_perm_full(i) < var_perm_full(i+1), for i = rank_full...n_full-1
  • [return == true] equ_perm_full(i) < equ_perm_full(i+1), for i = 1...rank-1
  • [return == true] equ_perm_full(i) < equ_perm_full(i+1), for i = rank...m_full-1

This method will only be called if this->nlp_selects_basis() == true.

The basis returned by the subclass must be sorted var_perm_full = [ dep | indep ] and equ_perm_full = [ equ_decomp | equ_undecomp ]. The subclass should not remove the variables fixed by bounds from var_perm_full as they will be removed by this class as they are translated. In addition, the subclass can also include slack variables in the basis (if mI_orig > 0>/tt>). Therefore, a nonsingular basis before fixed variables are removed may not be nonsingular once the fixed variables are removed. During the translation of var_perm_perm, the variables fixed by bounds are removed by compacting var_perm_full and adjusting the remaining indices. For this to be correct with variables fixed by bounds, it is assumed that the subclass knows which variables are fixed by bounds and can construct var_perm_full so that after the translated the basis will be nonsingular. The first rank entries in var_perm_full[1:rank_full] left after the fixed variables have been removed give the indices of the dependent (basic) variables and the remaining variables in var_perm_full[rank_full+1,n_full] are the indices for the independent (nonbasic) variables. To simplify things, it would be wise for the NLP subclass not to put fixed variables in the basis since this will greatly simplify selecting a nonsingular basis.

The first time this method is called, the subclass should return the first suggested basis selection (even if it happens to be identical to the original ordering).

The default implementation returns false which implies that the NLP subclass has no idea what a good basis selection should be..

Definition at line 739 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

virtual void NLPInterfacePack::NLPSerialPreprocess::imp_report_orig_final_solution ( const DVectorSlice &  x_full,
const DVectorSlice *  lambda_orig,
const DVectorSlice *  lambdaI_orig,
const DVectorSlice *  nu_orig,
bool  optimal 
)
inlineprotectedvirtual

To be overridden by subclasses to report the final solution in the original ordering natural to the subclass.

Note that the lagrange multipliers for fixed variables that have been preprocessed out of the problem are not computed by the optimization algorithm and are therefore not available. These multipliers are designated with the special value fixed_var_mult() but the numerical value is not significant.

The default implementation of this function is to do nothing.

Definition at line 623 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

void NLPInterfacePack::NLPSerialPreprocess::set_not_initialized ( )
inlineprotected

Used by subclasses to set the state of the NLP to not initialized.

Definition at line 845 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

void NLPInterfacePack::NLPSerialPreprocess::assert_initialized ( ) const
protected

Assert if we have been initizlized (throws UnInitialized)

Definition at line 749 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

void NLPInterfacePack::NLPSerialPreprocess::set_x_full ( const DVectorSlice &  x,
bool  newx,
DVectorSlice *  x_full 
) const
protected

Set the full x vector if newx == true

Definition at line 756 of file NLPInterfacePack_NLPSerialPreprocess.cpp.

DVectorSlice NLPInterfacePack::NLPSerialPreprocess::x_full ( ) const
inlineprotected

Give reference to current x_full.

Definition at line 851 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const NLPSerialPreprocess::ZeroOrderInfoSerial NLPInterfacePack::NLPSerialPreprocess::zero_order_orig_info ( ) const
inlineprotected

Definition at line 858 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const NLPSerialPreprocess::ObjGradInfoSerial NLPInterfacePack::NLPSerialPreprocess::obj_grad_orig_info ( ) const
inlineprotected

Definition at line 865 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const IVector & NLPInterfacePack::NLPSerialPreprocess::var_remove_fixed_to_full ( ) const
inlineprotected

Permutation vector for partitioning free and fixed variables.

var_remove_fixed_to_full =  [ not fixed by bounds |  fixed by bounds  ]
                            [1 ..                n|n + 1 ..     n_full]

The mapping i_full = var_remove_fixed_to_full()(i_free_fixed) gives the index of the original variable (i_full) for the sets of variables not fixed and fixed (upper and lower bounds where equal).

Definition at line 871 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const IVector & NLPInterfacePack::NLPSerialPreprocess::var_full_to_remove_fixed ( ) const
inlineprotected

Inverse permutation vector of var_remove_fixed_to_full().

The inverse mapping i_free_fixed = var_full_to_remove_fixed()(i_full) can be used to determine if a variable is free for fixed.

Definition at line 877 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const IVector & NLPInterfacePack::NLPSerialPreprocess::var_perm ( ) const
inlineprotected

Permutes from the compated variable vector (removing fixed variables) to the current basis selection.

On top of this partitioning of free and fixed variables, there is a permutation of the free variables into dependent and independent variables that is needed by the optimization algorithm.

var_perm = [ dependent variables | independent variables ]
           [1..                 r|r+1..                 n]

The mapping i_free_fixed = var_perm()(i_perm) is used to determine the index of a free variable in var_remove_fixed_to_full() given its index (i_perm) being used by the client.

Definition at line 883 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const IVector & NLPInterfacePack::NLPSerialPreprocess::equ_perm ( ) const
inlineprotected

Permutes from the original constriant ordering to the current basis selection.

equ_perm = [ decomposed equalities | undecomposed equalities ]
           [1..                   r|n-r+1...                n]

The mapping j_full = equ_perm()(j_perm) is used to determine the index of the constriant in c_full given its index i_perm being used by the NLP client.

Definition at line 889 of file NLPInterfacePack_NLPSerialPreprocess.hpp.

const IVector & NLPInterfacePack::NLPSerialPreprocess::inv_equ_perm ( ) const
inlineprotected

Inverse of equ_perm()

The mapping j_perm = inv_equ_perm()(j_full) is used to determine the index j_perm of the constriant c being used by the client given the index in c_full.

Definition at line 895 of file NLPInterfacePack_NLPSerialPreprocess.hpp.


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