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 | List of all members
NLPInterfacePack::NLPVarReductPerm Class Referenceabstract

NLP interface class that adds variable and constriant permutations for variable reduction basis selections. More...

#include <NLPInterfacePack_NLPVarReductPerm.hpp>

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

Classes

class  InvalidBasis
 Thrown if an invalid basis selection is made. More...
 

Public types

typedef Teuchos::RCP< const
Teuchos::AbstractFactory
< Permutation > > 
perm_fcty_ptr_t
 

Abstract factories for Permutation objects

virtual const perm_fcty_ptr_t factory_P_var () const =0
 
virtual const perm_fcty_ptr_t factory_P_equ () const =0
 

Return ranges for the partitioning of variables and constraints

virtual Range1D var_dep () const =0
 
virtual Range1D var_indep () const =0
 
virtual Range1D equ_decomp () const =0
 
virtual Range1D equ_undecomp () const =0
 

Basis manipulation functions

virtual bool nlp_selects_basis () const =0
 Returns true if the NLP can suggest one or more basis selections. More...
 
virtual bool get_next_basis (Permutation *P_var, Range1D *var_dep, Permutation *P_equ, Range1D *equ_decomp)=0
 Returns the next basis the NLP has and sets the NLP to the returned basis. More...
 
virtual void set_basis (const Permutation &P_var, const Range1D &var_dep, const Permutation *P_equ, const Range1D *equ_decomp)=0
 Sets the basis the that the NLP will use to permute the problem. More...
 
virtual void get_basis (Permutation *P_var, Range1D *var_dep, Permutation *P_equ, Range1D *equ_decomp) const =0
 Returns the basis selection currently being used by the NLP. More...
 

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 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 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...
 
virtual size_type n () const
 Return the number of variables. More...
 
virtual size_type m () const
 Return the number of general equality constraints. More...
 
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...
 
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...
 
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...
 
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 unset_quantities ()
 Call to unset all storage quantities (both in this class and all subclasses). More...
 
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...
 
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...
 
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 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...
 
- Static Public Member Functions inherited from NLPInterfacePack::NLP
static value_type infinite_bound ()
 Value for an infinite bound. 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...
 
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 that adds variable and constriant permutations for variable reduction basis selections.

This class adds basis selection and manipulation. This functionality is needed by many optimization algorithms that categorize variables and constraints into specific sets according to a basis selection. To understand what sets these are, consider the following equality constraints (from the NLP interface).

           c(x)  = 0,    c(x) <: R^n -> R^m

This interface allows x and c(x) to be partitioned into different sets. The variables x are partitioned into a dependent set x(var_dep) and an independent set x(var_dep) by the permutation P_var. The equality constraints c(x) are partitioned into decomposed c(equ_decomp) and undecomposed c(equ_undecomp) sets by the permutation P_equ. These permutations permute from an original order to a new ordering. For example:

 Original Ordering      Permutation to new ordering                Partitioning
 -----------------    -------------------------------   -------------------------------------
      x_orig          P_var.permute(trans,x_orig,x)     -> x(var_dep),      x(var_indep)
      c_orig          P_equ.permute(trans,c_orig,c)     -> c(equ_decomp),   c(equ_undecomp)

Because of this partitioning, it is expected that the following vector sub-spaces will be non-null: space_x()->sub_space(var_indep), space_x()->sub_space(var_dep), space_c()->sub_space(equ_decomp), space_c()->sub_space(equ_undecomp). Other subspaces may be non-null also but these are the only ones that are required to be.

After initialization, the NLP subclass will be initialized to the first basis. This basis may be the original ordering if P_var and P_equ all return xxx_perm.is_identity(). If the concrete NLP is selecting the basis (nlp_selects_basis() == true) this basis will be that first basis. The first time that this->get_next_basis() is called it will return this initial basis (which may not be the original ordering).

The client can always see what this first basis is by calling this->get_basis(). If a basis goes singular the client can request other basis selections from the NLP by calling this->get_next_basis() (which will return true if more basis selections are available). The client can also select a basis itself and then set that basis by calling this->set_basis() to force the use of that basis selection. In this way a valid basis is automatically selected after initialization so that clients using another interface (NLP, NLPFirstOrder, or NLPSecondOrder) will be able to use the NLP object without even knowing about a basis selection.

Below are some obviouls assertions about the basis selection:

Definition at line 105 of file NLPInterfacePack_NLPVarReductPerm.hpp.

Member Typedef Documentation

Definition at line 114 of file NLPInterfacePack_NLPVarReductPerm.hpp.

Member Function Documentation

virtual const perm_fcty_ptr_t NLPInterfacePack::NLPVarReductPerm::factory_P_var ( ) const
pure virtual
virtual const perm_fcty_ptr_t NLPInterfacePack::NLPVarReductPerm::factory_P_equ ( ) const
pure virtual
virtual Range1D NLPInterfacePack::NLPVarReductPerm::var_dep ( ) const
pure virtual
virtual Range1D NLPInterfacePack::NLPVarReductPerm::var_indep ( ) const
pure virtual
virtual Range1D NLPInterfacePack::NLPVarReductPerm::equ_decomp ( ) const
pure virtual
virtual Range1D NLPInterfacePack::NLPVarReductPerm::equ_undecomp ( ) const
pure virtual
virtual bool NLPInterfacePack::NLPVarReductPerm::nlp_selects_basis ( ) const
pure virtual

Returns true if the NLP can suggest one or more basis selections.

Implemented in NLPInterfacePack::NLPSerialPreprocess.

virtual bool NLPInterfacePack::NLPVarReductPerm::get_next_basis ( Permutation *  P_var,
Range1D *  var_dep,
Permutation *  P_equ,
Range1D *  equ_decomp 
)
pure virtual

Returns the next basis the NLP has and sets the NLP to the returned basis.

Parameters
P_var[out] Variable permutations defined as P_var'*x_old -> x_new = [ x(var_dep); x(var_indep) ]
var_dep[out] Range of dependent variables in x_new(var_dep)
P_equ[out] Equality constraint permutations defined as P_equ'*c_old -> c_new = [ c(equ_decomp); c(equ_undecomp) ]
equ_decomp[out] Range of decomposed equalities in c_new(equ_decomp)

Postconditions: The NLP is set to the basis returned in the arguments and this->get_basis() will return the same basis.

This member returns true if the NLP has another basis to select, and is false if not. If false is returned the client has the option of selecting another basis on its own and passing it to the NLP by calling this->set_basis().

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

virtual void NLPInterfacePack::NLPVarReductPerm::set_basis ( const Permutation &  P_var,
const Range1D &  var_dep,
const Permutation *  P_equ,
const Range1D *  equ_decomp 
)
pure virtual

Sets the basis the that the NLP will use to permute the problem.

Parameters
P_var[in] Variable permutations defined as P_var'*x_old -> x_new = [ x(var_dep); x(var_indep) ]
var_dep[in] Range of dependent variables in x_new(var_dep)
P_equ[in] Equality constraint permutations defined as P_equ'*c_old -> c_new = [ c(equ_decomp); c(equ_undecomp) ]
equ_decomp[in] Range of decomposed equalities in c_new(equ_decomp)

Preconditions: The input basis meets the basis assertions stated above or an InvalidBasis exceptin is thrown.

Postconditions: The NLP is set to the basis given in the arguments and this->get_basis() will return this same basis.

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

virtual void NLPInterfacePack::NLPVarReductPerm::get_basis ( Permutation *  P_var,
Range1D *  var_dep,
Permutation *  P_equ,
Range1D *  equ_decomp 
) const
pure virtual

Returns the basis selection currently being used by the NLP.

Parameters
P_var[out] Variable permutations defined as P_var'*x_old -> x_new = [ x(var_dep); x(var_indep) ]
var_dep[out] Range of dependent variables in x_new(var_dep)
P_equ[out] Equality constraint permutations defined as P_equ'*c_old -> c_new = [ c(equ_decomp); c(equ_undecomp) ]
equ_decomp[out] Range of decomposed equalities in c_new(equ_decomp)

Implemented in NLPInterfacePack::NLPSerialPreprocess.


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