MOOCHO (Single Doxygen Collection)
Version of the Day
|
Simple example NLP subclass to illustrate how to implement the NLPDirect
interface for a specialized NLP
.
More...
#include <NLPInterfacePack_ExampleNLPDirect.hpp>
Public Member Functions | |
ExampleNLPDirect (const VectorSpace::space_ptr_t &vec_space, value_type xo, bool has_bounds, bool dep_bounded) | |
Constructor. More... | |
Public Member Functions inherited from NLPInterfacePack::NLPDirect | |
void | set_factories (const mat_sym_fcty_ptr_t &factory_transDtD, const mat_sym_nonsing_fcty_ptr_t &factory_S) |
Initialize the factory objects for the special matrices for D'*D and S = I + D'*D . More... | |
virtual size_type | r () const |
Returns the number of decomposed equality constraints (r <= m ). More... | |
virtual Range1D | con_decomp () const |
Return the range of decomposed equality constraints. More... | |
virtual Range1D | con_undecomp () const |
Return the range of undecomposed equality constraints. More... | |
virtual const mat_fcty_ptr_t | factory_GcU () const |
Return a matrix factory object for creating GcU . More... | |
virtual const mat_fcty_ptr_t | factory_Uz () const |
Return a matrix factory object for Uz = F + E * D . More... | |
virtual const mat_fcty_ptr_t | factory_GcUD () const |
Return a matrix factory object for a mutable matrix compatible with GcU(var_dep) . More... | |
virtual const mat_sym_fcty_ptr_t | factory_transDtD () const |
Returns a matrix factory for the result of J = D'*D More... | |
virtual const mat_sym_nonsing_fcty_ptr_t | factory_S () const |
Returns a matrix factory for the result of S = I + D'*D 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 VectorMutable * | get_Gf () |
Return pointer passed to this->set_Gf() . More... | |
virtual VectorMutable & | Gf () |
Returns non-const *this->get_Gf() . More... | |
virtual const Vector & | Gf () 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_t & | get_options () const |
Get the OptionsFromStream object being used to extract the options from. 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 VectorMutable * | get_c () |
Return pointer passed to this->set_c() . More... | |
virtual VectorMutable & | c () |
Returns non-const *this->get_c() . More... | |
virtual const Vector & | c () 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 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 Vector & | hl_breve () const |
Returns a reference to the vector of lower bounds on the general inequality constraints h_breve(x_breve) . More... | |
virtual const Vector & | hu_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 VectorMutable * | get_c_breve () |
Return pointer passed to this->set_c_breve() . More... | |
virtual VectorMutable & | c_breve () |
Returns non-const *this->get_c_breve() . More... | |
virtual const Vector & | c_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 VectorMutable * | get_h_breve () |
Return pointer passed to this->set_h_breve() . More... | |
virtual VectorMutable & | h_breve () |
Returns non-const *this->get_h_breve() . More... | |
virtual const Vector & | h_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... | |
Public Member Functions inherited from Teuchos::VerboseObject< NLP > | |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT RCP< const ParameterList > | getValidVerboseObjectSublist () |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | setupVerboseObjectSublist (ParameterList *paramList) |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel) |
void | readVerboseObjectSublist (ParameterList *paramList, VerboseObject< NLP > *verboseObject) |
VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual void | initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
virtual const VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const |
virtual const VerboseObject & | setOverridingVerbLevel (const EVerbosityLevel verbLevel) const |
virtual EVerbosityLevel | getVerbLevel () const |
Public Member Functions inherited from NLPInterfacePack::ExampleNLPObjGrad | |
ExampleNLPObjGrad (const VectorSpace::space_ptr_t &vec_space, value_type xo, bool has_bounds, bool dep_bounded) | |
Constructor. More... | |
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 |
void | force_xinit_in_bounds (bool force_xinit_in_bounds) |
bool | force_xinit_in_bounds () const |
const Vector & | xinit () const |
const Vector & | xl () const |
const Vector & | xu () const |
value_type | max_var_bounds_viol () 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 optimal) |
Private Member Functions | |
void | assert_is_initialized () const |
Private Attributes | |
mat_fcty_ptr_t | factory_D_ |
bool | initialized_ |
Overridden public members from NLP | |
void | initialize (bool test_setup) |
bool | is_initialized () const |
Overridden public members from NLPDirect | |
Range1D | var_dep () const |
Range1D | var_indep () const |
const mat_fcty_ptr_t | factory_D () const |
void | calc_point (const Vector &x, value_type *f, VectorMutable *c, bool recalc_c, VectorMutable *Gf, VectorMutable *py, VectorMutable *rGf, MatrixOp *GcU, MatrixOp *D, MatrixOp *Uz) const |
void | calc_semi_newton_step (const Vector &x, VectorMutable *c, bool recalc_c, VectorMutable *py) const |
Additional Inherited Members | |
Public Types inherited from NLPInterfacePack::NLPDirect | |
typedef Teuchos::RCP< const Teuchos::AbstractFactory < MatrixOp > > | mat_fcty_ptr_t |
typedef Teuchos::RCP< const Teuchos::AbstractFactory < MatrixSymOp > > | mat_sym_fcty_ptr_t |
typedef Teuchos::RCP< const Teuchos::AbstractFactory < MatrixSymOpNonsing > > | mat_sym_nonsing_fcty_ptr_t |
Public Types inherited from NLPInterfacePack::NLP | |
typedef AbstractLinAlgPack::Vector | Vector |
typedef AbstractLinAlgPack::VectorMutable | VectorMutable |
typedef Teuchos::RCP< const VectorSpace > | vec_space_ptr_t |
typedef Teuchos::RCP< const OptionsFromStreamPack::OptionsFromStream > | options_ptr_t |
Static Public Member Functions inherited from NLPInterfacePack::NLP | |
static value_type | infinite_bound () |
Value for an infinite bound. More... | |
Static Public Member Functions inherited from Teuchos::VerboseObject< NLP > | |
static void | setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel) |
static EVerbosityLevel | getDefaultVerbLevel () |
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... | |
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... | |
Protected Member Functions inherited from NLPInterfacePack::ExampleNLPObjGrad | |
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_h (const Vector &x, bool newx, const ZeroOrderInfo &zero_order_info) const |
This implementation does nothing (should never be called though). More... | |
void | imp_calc_Gf (const Vector &x, bool newx, const ObjGradInfo &obj_grad_info) const |
Simple example NLP subclass to illustrate how to implement the NLPDirect
interface for a specialized NLP
.
For the NLP formulation see ExampleNLPObjGrad
.
In this NLP we will select the first n/2
variables as the dependent or basis variables. Note that the variable bounds are on the dependent or independent variables variables if has_bounds = true. Otherwise there are no variable bounds. Also the starting point can be varied.
The implementation of this NLP subclass is actually independent from the vector space for the dependent and independent variables (same vector space). This implementation is defined entirely based on an arbitrary VectorSpace
object that is passed to the constructor ExampleNLPDirect()
. This NLP subclass uses a VectorSpaceBlocked
object to represent the space for [ x_dep; x_indep ]
The quantities computed by this subclass include:
py = -inv(C)*c D = -inv(C)*N where: Gc' = [ C , N ] [ x(m+1) - 1 ] [ x(m+2) - 1 ] C = [ . ] [ . ] [ x(m+m) - 1 ] [ x(1) - 10 ] [ x(2) - 10 ] N = [ . ] [ . ] [ x(m) - 10 ]
Here Gc
is never computed explicitly.
To make this possible this subclass relies on some specialized RTOp operators which are implemented in C (for portability). These operator classes are declared in the header file NLPInterfacePack_ExampleNLPDirect.hpp
and are documented here.
Definition at line 103 of file NLPInterfacePack_ExampleNLPDirect.hpp.
NLPInterfacePack::ExampleNLPDirect::ExampleNLPDirect | ( | const VectorSpace::space_ptr_t & | vec_space, |
value_type | xo, | ||
bool | has_bounds, | ||
bool | dep_bounded | ||
) |
Constructor.
vec_space | [in] Smart pointer to a vector space object that will be used to define the spaces of dependent and independent variables. |
xo | [in] The initial starting guess for x. |
has_bounds | [in] If true , then the NLP will have bounds. If false then it will not have bounds. |
dep_bouned | [in] If true , then the bounds will be on the dependent variables. If false , then the bounds will be on the independent variable. This argument is ignored if has_bounds == false . |
Definition at line 75 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Reimplemented from NLPInterfacePack::NLPDirect.
Definition at line 97 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Implements NLPInterfacePack::NLP.
Definition at line 112 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Reimplemented from NLPInterfacePack::NLPDirect.
Definition at line 119 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Reimplemented from NLPInterfacePack::NLPDirect.
Definition at line 124 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Implements NLPInterfacePack::NLPDirect.
Definition at line 130 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Implements NLPInterfacePack::NLPDirect.
Definition at line 135 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
virtual |
Implements NLPInterfacePack::NLPDirect.
Definition at line 291 of file NLPInterfacePack_ExampleNLPDirect.cpp.
|
inlineprivate |
Definition at line 192 of file NLPInterfacePack_ExampleNLPDirect.hpp.
|
private |
Definition at line 176 of file NLPInterfacePack_ExampleNLPDirect.hpp.
|
private |
Definition at line 178 of file NLPInterfacePack_ExampleNLPDirect.hpp.