ConstrainedOptPack: C++ Tools for Constrained (and Unconstrained) Optimization  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Classes | Public Member Functions | List of all members
ConstrainedOptPack::DecompositionSystem Class Referenceabstract

This class abstracts a decomposition choice for the quasi-range space Y and null space Z matrices for a linearly independent set of columns of Gc. More...

#include <ConstrainedOptPack_DecompositionSystem.hpp>

Inheritance diagram for ConstrainedOptPack::DecompositionSystem:
Inheritance graph
[legend]

Classes

class  InvalidMatrixType
 
class  SingularDecomposition
 
class  TestFailed
 

Public Member Functions

virtual ~DecompositionSystem ()
 

Public types

enum  EOutputLevel
 Enumeration for the amount of output to create from update_decomp(). More...
 
enum  ERunTests
 Enumeration for if to run internal tests or not. More...
 
enum  EMatRelations
 
typedef Teuchos::RCP< const
Teuchos::AbstractFactory
< MatrixOpNonsing > > 
mat_nonsing_fcty_ptr_t
 
typedef Teuchos::RCP< const
Teuchos::AbstractFactory
< MatrixOp > > 
mat_fcty_ptr_t
 

Dimensionality of the decomposition

virtual size_type n () const
 Return the number of rows in Gc. More...
 
virtual size_type m () const =0
 Return the number of columns in Gc. More...
 
virtual size_type r () const
 Returns the rank of Gc(:,equ_decomp()). More...
 
virtual Range1D equ_decomp () const
 Returns the range of the decomposed equalities. More...
 
virtual Range1D equ_undecomp () const
 Returns the range of the undecomposed equalities. More...
 

Range and null vector spaces

virtual const
VectorSpace::space_ptr_t 
space_range () const =0
 Return a VectorSpace object for the range space. More...
 
virtual const
VectorSpace::space_ptr_t 
space_null () const =0
 Return a VectorSpace object for the range space. More...
 

Matrix factories

virtual const mat_fcty_ptr_t factory_Z () const =0
 Return a matrix factory object for Z More...
 
virtual const mat_fcty_ptr_t factory_Y () const =0
 Return a matrix factory object for Y More...
 
virtual const
mat_nonsing_fcty_ptr_t 
factory_R () const =0
 Return a matrix factory object for R. More...
 
virtual const mat_fcty_ptr_t factory_Uz () const =0
 Return a matrix factory object for Uz More...
 
virtual const mat_fcty_ptr_t factory_Uy () const =0
 Return a matrix factory object for Uy More...
 

Update range/null decomposition

virtual void update_decomp (std::ostream *out, EOutputLevel olevel, ERunTests test_what, const MatrixOp &Gc, MatrixOp *Z, MatrixOp *Y, MatrixOpNonsing *R, MatrixOp *Uz, MatrixOp *Vy, EMatRelations mat_rel=MATRICES_INDEP_IMPS) const =0
 Creates the range/null decomposition for Gc(:,equ_decomp)'. More...
 
virtual void print_update_decomp (std::ostream &out, const std::string &leading_str) const =0
 Print the sub-algorithm by which the decomposition is formed. More...
 

Detailed Description

This class abstracts a decomposition choice for the quasi-range space Y and null space Z matrices for a linearly independent set of columns of Gc.

Gc = [ Gc(:,equ_decomp), Gc(:,equ_undecomp) ]

where Gc is n x m, Gc(:,equ_decomp) is n x r and Gc(:,equ_undecomp) is n x (m - r).

Note that the columns in Gc(:,equ_undecomp) may be linearly dependent with the columns in Gc(:,equ_undecomp) or they may just be undecomposed linearly independent equality constraints.

The decomposition formed by subclasses must have the properties:

  Z s.t. Gc(:,equ_decomp)' * Z = 0
  Y s.t. [Z  Y] is nonsingular
  R = Gc(:,equ_decomp)' * Y is nonsingular
  Uz = Gc(:,equ_undecomp)' * Z
  Uy = Gc(:,equ_undecomp)' * Y

The matrix factory objects returned by ??? are ment to have a lifetime that is independent of this.

The decomposition matrices Z, Y, R, Uz and Uy which are updated in this->update_decomp() must be completely independent from this and from each other and Gc that they based on. For example, Once update_decomp() is called, this, Gc can be destroyed and the behaviors of the decomposition matrices must not be altered. In this respect the DecompositionSystem interface is really nothing more than a "Strategy" interface (with some state data of course) for computing range/null decompositions. This gives the client great flexibility in how the decomposition matrices are used.

ToDo: Finish documentation!

Definition at line 89 of file ConstrainedOptPack_DecompositionSystem.hpp.

Member Typedef Documentation

Definition at line 97 of file ConstrainedOptPack_DecompositionSystem.hpp.

Definition at line 100 of file ConstrainedOptPack_DecompositionSystem.hpp.

Member Enumeration Documentation

Enumeration for the amount of output to create from update_decomp().

Definition at line 111 of file ConstrainedOptPack_DecompositionSystem.hpp.

Enumeration for if to run internal tests or not.

Definition at line 119 of file ConstrainedOptPack_DecompositionSystem.hpp.

Definition at line 121 of file ConstrainedOptPack_DecompositionSystem.hpp.

Constructor & Destructor Documentation

virtual ConstrainedOptPack::DecompositionSystem::~DecompositionSystem ( )
inlinevirtual

Definition at line 126 of file ConstrainedOptPack_DecompositionSystem.hpp.

Member Function Documentation

size_type ConstrainedOptPack::DecompositionSystem::n ( ) const
virtual

Return the number of rows in Gc.

Postconditions:

  • n > m

The default implementation returns this->space_range()->dim() + this->space_null()->dim().

Reimplemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.

Definition at line 46 of file ConstrainedOptPack_DecompositionSystem.cpp.

virtual size_type ConstrainedOptPack::DecompositionSystem::m ( ) const
pure virtual

Return the number of columns in Gc.

Postconditions:

  • m > 0

Implemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.

size_type ConstrainedOptPack::DecompositionSystem::r ( ) const
virtual

Returns the rank of Gc(:,equ_decomp()).

Postconditions:

  • r =< m

The default implementation returns this->space_range()->dim().

Reimplemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.

Definition at line 51 of file ConstrainedOptPack_DecompositionSystem.cpp.

Range1D ConstrainedOptPack::DecompositionSystem::equ_decomp ( ) const
virtual

Returns the range of the decomposed equalities.

The default implementation returns Range1D(1,this->r()).

Reimplemented in ConstrainedOptPack::DecompositionSystemVarReductPermStd.

Definition at line 56 of file ConstrainedOptPack_DecompositionSystem.cpp.

Range1D ConstrainedOptPack::DecompositionSystem::equ_undecomp ( ) const
virtual

Returns the range of the undecomposed equalities.

The default implementation returns Range1D(this->r()+1,this->m()) or Range1D::Invalid if this->r() == this->m()

Reimplemented in ConstrainedOptPack::DecompositionSystemVarReductPermStd.

Definition at line 63 of file ConstrainedOptPack_DecompositionSystem.cpp.

virtual const VectorSpace::space_ptr_t ConstrainedOptPack::DecompositionSystem::space_range ( ) const
pure virtual

Return a VectorSpace object for the range space.

Postconditions:

  • return.get() != NULL
  • return->dim() == this->r()

Implemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.

virtual const VectorSpace::space_ptr_t ConstrainedOptPack::DecompositionSystem::space_null ( ) const
pure virtual

Return a VectorSpace object for the range space.

Postconditions:

  • return.get() != NULL
  • return->dim() == this->n() - this->r()

Implemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.

virtual const mat_fcty_ptr_t ConstrainedOptPack::DecompositionSystem::factory_Z ( ) const
pure virtual
virtual const mat_fcty_ptr_t ConstrainedOptPack::DecompositionSystem::factory_Y ( ) const
pure virtual
virtual const mat_nonsing_fcty_ptr_t ConstrainedOptPack::DecompositionSystem::factory_R ( ) const
pure virtual
virtual const mat_fcty_ptr_t ConstrainedOptPack::DecompositionSystem::factory_Uz ( ) const
pure virtual
virtual const mat_fcty_ptr_t ConstrainedOptPack::DecompositionSystem::factory_Uy ( ) const
pure virtual
virtual void ConstrainedOptPack::DecompositionSystem::update_decomp ( std::ostream *  out,
EOutputLevel  olevel,
ERunTests  test_what,
const MatrixOp &  Gc,
MatrixOp *  Z,
MatrixOp *  Y,
MatrixOpNonsing *  R,
MatrixOp *  Uz,
MatrixOp *  Vy,
EMatRelations  mat_rel = MATRICES_INDEP_IMPS 
) const
pure virtual

Creates the range/null decomposition for Gc(:,equ_decomp)'.

The decomposition is based on the linearly independent columns Gc(:,equ_decomp) of Gc

Gc = [ Gc(:,equ_decomp), Gc(:,equ_undecomp) ]

Specifically this operation finds the matrices:

Z s.t. Gc(:,equ_deomp)' * Z = 0
Y s.t. [Z  Y] is nonsingular
R = Gc(:,equ_decomp)' * Y is nonsingular
Uz = Gc(:,equ_undecomp)' * Z
Uy = Gc(:,equ_undecomp)' * Y

If there is some problem creating the decomposition then exceptions with the base class std::exception may be thrown. The meaning of these exceptions are more associated with the subclasses that implement this operation.

The concrete types for Gc, Z, Y, Uz and Uy must be compatable with the concrete implementation of this or an InvalidMatrixType exeption will be thrown..

Parameters
out[out] If out!=NULL then output is printed to this stream depending on the value of olevel.
olevel[in] Determines the amount of output to print to *out. The exact type of output is determined by the implementing subclass but here is the sugguested behavior:
  • PRINT_NONE : Don't print anything (same as out==NULL).
  • PRINT_BASIC_INFO : Only print basic information about how the decomposition is formed. Amount of output = O(1).
  • PRINT_VECTORS : Prints out important vectors computed durring the computations (usually only durring testing). This level is only useful for debugging. Amount of output = O(n).
  • PRINT_EVERY_THING : Print out nearly every important quantity that is computed (except for the output matrices themselves, clients can do that) while the matrices are being formed or tests are being conducted. This level is only useful for debugging. Amount of output = O(m*n)
test_what[in] Determines if internal validation tests are performed. The post conditions for the output matrices are not checked internally. This is something that client can (and should) do independently (see DecompositionSystemTester). Values:
  • RUN_TESTS : As many validation/consistency tests are performed internally as possible. If a test fails then a TestFailed execption will be thrown. The subclasses determine what the tests are and what failing a test means.
  • NO_TEST : No tests are performed internally. This is to allow the fastest possible execution.
If a test fails, then a TestFailed exception will be thrown with a helpful error message.
Gc[in] The matrix for which the range/null decomposition is defined.
Z[out] On output represents the n x (n-r) null space matrix such that Gc(:,equ_decomp) * Z == 0. This matrix object must have been created by this->factory_Z()->create().
Y[out] On output represents the n x r range space matrix such that [ Y Z ] is nonsingular. This matrix object must have been created by this->factory_Y()->create().
R[out] On output represents the nonsingular r x r matrix Gc(:,equ_decomp) * Y. This matrix object must have been created by this->factory_R()->create().
Uz[in/out] If Uz != NULL (this->m() > this->r() only) then on output *Uz represents the (m-r) x (n-r) matrix Gc(:,equ_undecomp) * Z. If this->m() == this->r() then Uz == NULL must be true. If Uz!=NULL, then this matrix object must have been created by this->factory_Uz()->create().
Uy[in/out] If Uy != NULL (this->m() > this->r() only) then on output *Uy represents the (m-r) x r matrix Gc(:,equ_undecomp) * Y. If this->m() == this->r() then Uy == NULL must be true. If Uy!=NULL, then this matrix object must have been created by this->factory_Uy()->create().
mat_rel[in] Determines if the ouput matrix objects must be completely independent or not.
  • MATRICES_INDEP_IMPS: The matrix objects must have independent implementations (default).
  • MATRICES_ALLOW_DEP_IMPS: The matrix objects can have implementation dependencies.

Preconditions:

  • Gc.rows() == this->n() (throw std::invalid_argument)
  • Gc.cols() == this->m() (throw std::invalid_argument)
  • [this->m() == this->r()] Uz == NULL (throw std::invalid_argument)
  • [this->m() == this->r()] Uy == NULL (throw std::invalid_argument)
  • Z!=NULL || Y!=NULL || R!=NULL || Uz!=NULL || Uy!=NULL (throw std::invalid_argument)

Postconditions:

  • [Z != NULL] Z.space_cols().is_compatible(Gc.space_cols()) == true)
  • [Z != NULL] Z.space_rows().is_compatible(*space_null()) == true
  • [Z != NULL] Gc(:,equ_decomp())' * Z == 0
  • [Y != NULL] Y.space_cols().is_compatible(Gc.space_cols()) == true)
  • [Y != NULL] Y.cols() == this->r()
  • [Y != NULL] [ Y Z ] is nonsingular
  • [R != NULL] R->space_cols().is_compatible(*Gc.space_cols()->sub_space(equ_decomp())) == true
  • [R != NULL] R->space_rows().is_compatible(*space_range()) == true
  • [R != NULL] R == Gc(:,equ_decomp())'*Y
  • [Uz != NULL] Uz.space_cols().is_compatible(*Gc.space_rows()->sub_space(equ_undecomp())) == true
  • [Uz != NULL] Uz.space_rows().is_compatible(*space_null()) == true
  • [Uz != NULL] Uz == Gc(:,equ_undecomp())'*Z
  • [Uy != NULL] Uy.space_cols().is_compatible(*Gc.space_rows()->sub_space(equ_undecomp())) == true
  • [Uy != NULL] Uy.space_rows().is_compatible(*space_range()) == true
  • [Uy != NULL] Uy == Gc(:,equ_undecomp())'*Y
  • [mat_rel == MATRICES_INDEP_IMPS] The behaviors of all of the participating output matrix objects must not be altered by changes to other matrix objects.
  • [mat_rel == MATRICES_ALLOW_DEP_IMPS] The behaviors of all of the participating output matrix objects may change when other matrix objects or this is altered.

Note that this method requires that all of the output matrix objects Z, Y, R, Uz, Uy, Vz and Vy must be independent of this and of each other. For example, the behavior of R must be be altered if this is destroyed or if Z is modified of destroyed. This requirment constrains the implementations somewhat but makes things much easier for the client and gives the client much more power.

Implemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.

virtual void ConstrainedOptPack::DecompositionSystem::print_update_decomp ( std::ostream &  out,
const std::string &  leading_str 
) const
pure virtual

Print the sub-algorithm by which the decomposition is formed.

Implemented in ConstrainedOptPack::DecompositionSystemVarReductImp, and ConstrainedOptPack::DecompositionSystemVarReductPermStd.


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