|
MOOCHO (Single Doxygen Collection)
Version of the Day
|
Orthogonal variable reduction subclass. More...
#include <ConstrainedOptPack_DecompositionSystemOrthogonal.hpp>

Private Types | |
| typedef Teuchos::RCP < MatrixSymOpNonsing > | S_ptr_t |
Private Member Functions | |
| DecompositionSystemOrthogonal (const DecompositionSystemOrthogonal &) | |
| DecompositionSystemOrthogonal & | operator= (const DecompositionSystemOrthogonal &) |
Private Attributes | |
| S_ptr_t | S_ptr_ |
Constructors / initializers | |
| DecompositionSystemOrthogonal (const VectorSpace::space_ptr_t &space_x=Teuchos::null, const VectorSpace::space_ptr_t &space_c=Teuchos::null, const basis_sys_ptr_t &basis_sys=Teuchos::null, const basis_sys_tester_ptr_t &basis_sys_tester=Teuchos::null, EExplicitImplicit D_imp=MAT_IMP_EXPLICIT, EExplicitImplicit Uz_imp=MAT_IMP_EXPLICIT) | |
Overridden from DecompositionSystem | |
| const mat_fcty_ptr_t | factory_Y () const |
| const mat_nonsing_fcty_ptr_t | factory_R () const |
| const mat_fcty_ptr_t | factory_Uy () const |
Overridden from DecompositionSystemVarReductImp | |
| void | update_D_imp_used (EExplicitImplicit *D_imp_used) const |
| mat_nonsing_fcty_ptr_t::element_type::obj_ptr_t | uninitialize_matrices (std::ostream *out, EOutputLevel olevel, MatrixOp *Y, MatrixOpNonsing *R, MatrixOp *Uy) const |
| void | initialize_matrices (std::ostream *out, EOutputLevel olevel, const mat_nonsing_fcty_ptr_t::element_type::obj_ptr_t &C, const mat_fcty_ptr_t::element_type::obj_ptr_t &D, MatrixOp *Y, MatrixOpNonsing *R, MatrixOp *Uy, EMatRelations mat_rel) const |
| void | print_update_matrices (std::ostream &out, const std::string &leading_str) const |
Additional Inherited Members | |
Public Types inherited from ConstrainedOptPack::DecompositionSystemVarReductImp | |
| typedef DecompositionSystem | inherited |
| typedef Teuchos::RCP< const BasisSystem > | basis_sys_ptr_t |
Public Types inherited from ConstrainedOptPack::DecompositionSystemVarReduct | |
| enum | EExplicitImplicit { MAT_IMP_EXPLICIT, MAT_IMP_IMPLICIT, MAT_IMP_AUTO } |
Public Types inherited from ConstrainedOptPack::DecompositionSystem | |
| enum | EOutputLevel { PRINT_NONE = 0, PRINT_BASIC_INFO = 1, PRINT_MORE_INFO = 2, PRINT_VECTORS = 3, PRINT_EVERY_THING = 4 } |
Enumeration for the amount of output to create from update_decomp(). More... | |
| enum | ERunTests { RUN_TESTS, NO_TESTS } |
| Enumeration for if to run internal tests or not. More... | |
| enum | EMatRelations { MATRICES_INDEP_IMPS, MATRICES_ALLOW_DEP_IMPS } |
| 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::DecompositionSystemVarReductImp | |
| void | initialize (const VectorSpace::space_ptr_t &space_x, const VectorSpace::space_ptr_t &space_c, const basis_sys_ptr_t &basis_sys) |
| Initialize. More... | |
| STANDARD_COMPOSITION_MEMBERS (BasisSystemTester, basis_sys_tester) | |
| Set the BasisSystem tester object. More... | |
| DecompositionSystemVarReductImp (const VectorSpace::space_ptr_t &space_x, const VectorSpace::space_ptr_t &space_c, const basis_sys_ptr_t &basis_sys, const basis_sys_tester_ptr_t &basis_sys_tester, EExplicitImplicit D_imp, EExplicitImplicit Uz_imp) | |
| Construct a variable reduction decomposition. More... | |
| const VectorSpace::space_ptr_t & | space_x () const |
| const VectorSpace::space_ptr_t & | space_c () const |
| const basis_sys_ptr_t & | basis_sys () const |
| void | get_basis_matrices (std::ostream *out, EOutputLevel olevel, ERunTests test_what, MatrixOp *Z, MatrixOp *Y, MatrixOpNonsing *R, MatrixOp *Uz, MatrixOp *Uy, Teuchos::RCP< MatrixOpNonsing > *C_ptr, Teuchos::RCP< MatrixOp > *D_ptr) |
| Called by client to uninitialize decomposition matrices in prepairation for selecting a different basis. More... | |
| void | set_basis_matrices (std::ostream *out, EOutputLevel olevel, ERunTests test_what, const Teuchos::RCP< MatrixOpNonsing > &C_ptr, const Teuchos::RCP< MatrixOp > &D_ptr, MatrixOp *Uz, const basis_sys_ptr_t &basis_sys=Teuchos::null) |
| Set updated basis matrices along with a possibly updated basis system object. More... | |
| EExplicitImplicit | D_imp_used () const |
| Get the type of D matrix to be used or is being used (returns MAT_IMP_EXPLICIT or MAT_IMP_IMPLICIT only). More... | |
| size_type | n () const |
Returns this->space_x()->dim(). More... | |
| size_type | m () const |
Returns this->space_c()->dim(). More... | |
| size_type | r () const |
Returns this->basis_sys()->equ_decomp().size(). More... | |
| const VectorSpace::space_ptr_t | space_range () const |
Returns this->space_x()->sub_space(var_dep) More... | |
| const VectorSpace::space_ptr_t | space_null () const |
Returns this->space_x()->sub_space(var_indep) More... | |
| const mat_fcty_ptr_t | factory_Z () const |
| const mat_fcty_ptr_t | factory_Uz () const |
| void | update_decomp (std::ostream *out, EOutputLevel olevel, ERunTests test_what, const MatrixOp &Gc, MatrixOp *Z, MatrixOp *Y, MatrixOpNonsing *R, MatrixOp *Uz, MatrixOp *Uy, EMatRelations mat_rel) const |
| Preconditions: More... | |
| void | print_update_decomp (std::ostream &out, const std::string &leading_str) const |
| Range1D | var_indep () const |
| Range1D | var_dep () const |
Public Member Functions inherited from ConstrainedOptPack::DecompositionSystemVarReduct | |
| 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... | |
| DecompositionSystemVarReduct (EExplicitImplicit D_imp=MAT_IMP_AUTO, EExplicitImplicit Uz_imp=MAT_IMP_AUTO) | |
Public Member Functions inherited from ConstrainedOptPack::DecompositionSystem | |
| virtual | ~DecompositionSystem () |
| 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... | |
Orthogonal variable reduction subclass.
This is the interface for the coordinate variable reduction decomposition where:
Y = [ I ] = [ I ] (class MatrixIdentConcatStd)
[ N'*inv(C') ] [ -D' ]
R = Gc(:,con_decomp)'*Y = [ C N ] * [ I ] = (C + N*N'*inv(C'))
[ N'*inv(C') ]
= C*(I + inv(C)*N*N'*inv(C'))
= C*(I + D*D')
Uy = Gc(:,con_undecomp)'*Y = [ E F ] * [ I ]
[ -D' ]
= E - F*D' See the matrix subclass MatrixDecompRangeOrthog for details on how linear systems with R are solved for.
For now the copy constructor and the assignment operator are not defined.
Definition at line 74 of file ConstrainedOptPack_DecompositionSystemOrthogonal.hpp.
|
private |
Definition at line 141 of file ConstrainedOptPack_DecompositionSystemOrthogonal.hpp.
| ConstrainedOptPack::DecompositionSystemOrthogonal::DecompositionSystemOrthogonal | ( | const VectorSpace::space_ptr_t & | space_x = Teuchos::null, |
| const VectorSpace::space_ptr_t & | space_c = Teuchos::null, |
||
| const basis_sys_ptr_t & | basis_sys = Teuchos::null, |
||
| const basis_sys_tester_ptr_t & | basis_sys_tester = Teuchos::null, |
||
| EExplicitImplicit | D_imp = MAT_IMP_EXPLICIT, |
||
| EExplicitImplicit | Uz_imp = MAT_IMP_EXPLICIT |
||
| ) |
Definition at line 60 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
private |
|
virtual |
Implements ConstrainedOptPack::DecompositionSystem.
Definition at line 76 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
virtual |
Implements ConstrainedOptPack::DecompositionSystem.
Definition at line 85 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
virtual |
Implements ConstrainedOptPack::DecompositionSystem.
Definition at line 94 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
protectedvirtual |
Reimplemented from ConstrainedOptPack::DecompositionSystemVarReductImp.
Definition at line 101 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
protectedvirtual |
Implements ConstrainedOptPack::DecompositionSystemVarReductImp.
Definition at line 107 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
protectedvirtual |
Implements ConstrainedOptPack::DecompositionSystemVarReductImp.
Definition at line 161 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
protectedvirtual |
Implements ConstrainedOptPack::DecompositionSystemVarReductImp.
Definition at line 259 of file ConstrainedOptPack_DecompositionSystemOrthogonal.cpp.
|
private |
|
mutableprivate |
Definition at line 146 of file ConstrainedOptPack_DecompositionSystemOrthogonal.hpp.
1.8.6