42 #ifndef MATRIX_SYM_POS_DEF_LBFGS_H
43 #define MATRIX_SYM_POS_DEF_LBFGS_H
47 #include "ConstrainedOptPack_Types.hpp"
48 #include "AbstractLinAlgPack_MatrixSymSecant.hpp"
49 #include "AbstractLinAlgPack_MatrixSymOpNonsing.hpp"
50 #include "AbstractLinAlgPack_MultiVectorMutable.hpp"
51 #include "AbstractLinAlgPack_VectorSpace.hpp"
52 #include "DenseLinAlgPack_DMatrixAsTriSym.hpp"
55 namespace ConstrainedOptPack {
112 ,
public MatrixSymSecant
131 ,
bool auto_rescaling =
false
136 ,auto_rescaling_(auto_rescaling)
141 p->
initial_setup(m_,maintain_original_,maintain_inverse_,auto_rescaling_);
145 bool maintain_original_;
146 bool maintain_inverse_;
147 bool auto_rescaling_;
160 ,
bool auto_rescaling =
false
207 ,
bool auto_rescaling =
false
240 std::ostream&
output(std::ostream& out)
const;
242 MatrixOp&
operator=(
const MatrixOp& mwo);
246 ,
const Vector& v_rhs2, value_type beta )
const;
256 ,
const Vector& v_rhs2 )
const;
264 void init_identity(
const VectorSpace& space_diag, value_type alpha );
286 typedef VectorSpace::multi_vec_mut_ptr_t multi_vec_mut_ptr_t;
291 bool maintain_original_;
292 bool original_is_updated_;
293 bool maintain_inverse_;
294 bool inverse_is_updated_;
296 VectorSpace::space_ptr_t
317 mutable bool Q_updated_;
326 const DMatrixSliceTri R()
const;
328 const DMatrixSliceTri Lb()
const;
332 const DMatrixSlice STY()
const;
334 DMatrixSliceSym STS();
336 const DMatrixSliceSym STS()
const;
338 DMatrixSliceSym YTY();
340 const DMatrixSliceSym YTY()
const;
342 void V_invQtV( DVectorSlice* y,
const DVectorSlice& x )
const;
344 void Vp_DtV( DVectorSlice* y,
const DVectorSlice& x )
const;
349 void update_Q()
const;
352 void assert_initialized()
const;
381 return S_->mv_sub_view(1,n_,1,m_bar_);
388 return Y_->mv_sub_view(1,n_,1,m_bar_);
394 return maintain_original_;
400 return maintain_inverse_;
406 return num_secant_updates_;
411 #endif // MATRIX_SYM_POS_DEF_LBFGS_H
void init_identity(const VectorSpace &space_diag, value_type alpha)
bool maintain_original() const
void secant_update(VectorMutable *s, VectorMutable *y, VectorMutable *Bs)
void initial_setup(size_type m=10, bool maintain_original=true, bool maintain_inverse=true, bool auto_rescaling=false)
Initial setup for the matrix.
PostMod(size_type m=10, bool maintain_original=true, bool maintain_inverse=true, bool auto_rescaling=false)
const VectorSpace & space_cols() const
void V_InvMtV(VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const
MatrixSymPosDefLBFGS(size_type m=10, bool maintain_original=true, bool maintain_inverse=true, bool auto_rescaling=false)
Calls this->initial_setup()
void init_diagonal(const Vector &diag)
Actually this calls init_identity( diag.space(), diag.norm_inf() ).
PostMod class to use with MemMngPack::AbstractFactorStd.
const multi_vec_ptr_t S() const
STANDARD_MEMBER_COMPOSITION_MEMBERS(bool, auto_rescaling)
Set whether automatic rescaling is used or not.
Teuchos::RCP< const MultiVector > multi_vec_ptr_t
const multi_vec_ptr_t Y() const
std::ostream & output(std::ostream &out) const
void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const
Implementation of limited Memory BFGS matrix for arbitrary vector spaces.
void initialize(MatrixSymPosDefLBFGS *p) const
bool maintain_inverse() const
size_type num_secant_updates() const
Returns the total number of successful secant updates performed since this->init_identity() was calle...
MatrixOp & operator=(const MatrixOp &mwo)
value_type gamma_k() const