| 
    ConstrainedOptPack: C++ Tools for Constrained (and Unconstrained) Optimization
    Version of the Day
    
   | 
 
Specialization of DecompositionSystem for variable reduction decompositions.  
 More...
#include <ConstrainedOptPack_DecompositionSystemVarReduct.hpp>

Public types | |
| enum | EExplicitImplicit | 
Matrix representations | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (EExplicitImplicit, D_imp) | |
Set whether to use explicit or implicit D = -inv(C)*N matrix.  More... | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (EExplicitImplicit, Uz_imp) | |
Set whether to use explicit or implicit Uz = F + E * D matrix.  More... | |
Constructors / initializers | |
| DecompositionSystemVarReduct (EExplicitImplicit D_imp=MAT_IMP_AUTO, EExplicitImplicit Uz_imp=MAT_IMP_AUTO) | |
Variable partitions. | |
| virtual Range1D | var_indep () const =0 | 
| virtual Range1D | var_dep () const =0 | 
Additional Inherited Members | |
  Public Types inherited from ConstrainedOptPack::DecompositionSystem | |
| 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 | 
  Public Member Functions inherited from ConstrainedOptPack::DecompositionSystem | |
| virtual | ~DecompositionSystem () | 
| 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... | |
| 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... | |
| 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... | |
| 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... | |
Specialization of DecompositionSystem for variable reduction decompositions. 
This interface abstracts a variable reduction decomposition where:
 Gc' = [ C  N ] 
       [ E  F ]
 Z   = [ D ]
       [ I ]
 Uz  = F + E * D
     where:
          C = Gc(var_dep,con_decomp)'     [nonsingular]
          N = Gc(var_indep,con_decomp)'
          E = Gc(var_dep,con_undecomp)'
          F = Gc(var_indep,con_undecomp)'
          D = -inv(C) * N
This interface simply allows clients to determine how D and Uz are implemented (implicitly or explicity). 
Definition at line 76 of file ConstrainedOptPack_DecompositionSystemVarReduct.hpp.
Definition at line 83 of file ConstrainedOptPack_DecompositionSystemVarReduct.hpp.
      
  | 
  inline | 
Definition at line 107 of file ConstrainedOptPack_DecompositionSystemVarReduct.hpp.
| ConstrainedOptPack::DecompositionSystemVarReduct::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | EExplicitImplicit | , | 
| D_imp | |||
| ) | 
Set whether to use explicit or implicit D = -inv(C)*N matrix. 
| ConstrainedOptPack::DecompositionSystemVarReduct::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | EExplicitImplicit | , | 
| Uz_imp | |||
| ) | 
Set whether to use explicit or implicit Uz = F + E * D matrix. 
      
  | 
  pure virtual | 
      
  | 
  pure virtual | 
 1.8.6