42 #ifndef MATRIX_SYM_POS_DEF_LBFGS_H 
   43 #define MATRIX_SYM_POS_DEF_LBFGS_H 
   49 #include "AbstractLinAlgPack/src/MatrixSymWithOpFactorized.hpp" 
   53 namespace ConstrainedOptPack {
 
  119 class MatrixSymPosDefLBFGS
 
  120   : 
public MatrixSymWithOpFactorized
 
  121   , 
public MatrixSymSecant
 
  122   , 
public MatrixSymAddDelUpdateable
 
  135     ,
bool       auto_rescaling     = 
false 
  186      ,
bool       auto_rescaling     = 
false 
  222   std::ostream& 
output(std::ostream& 
out) 
const;
 
  273     ,
bool              force_factorization
 
  275     ,PivotTolerances   pivot_tols
 
  294     ,
bool              force_refactorization
 
  295     ,EEigenValType     add_eigen_val
 
  296     ,PivotTolerances   pivot_tols
 
  301     ,
bool              force_refactorization
 
  302     ,EEigenValType     drop_eigen_val
 
  303     ,PivotTolerances   pivot_tols
 
  442 #endif  // MATRIX_SYM_POS_DEF_LBFGS_H 
void V_invQtV(DVectorSlice *y, const DVectorSlice &x) const 
y = inv(Q) * x 
const DMatrixSliceTri Lb() const 
Strictly lower triangular part of L. 
void init_identity(const VectorSpace &space_diag, value_type alpha)
AbstractLinAlgPack::size_type size_type
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. 
bool original_is_updated_
void V_InvMtV(VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const 
void update_Q() const 
Update Q. 
void initialize(value_type alpha, size_type max_size)
This is fine as long as alpha > 0.0. 
MatrixSymPosDefLBFGS(size_type m=10, bool maintain_original=true, bool maintain_inverse=true, bool auto_rescaling=false)
Calls this->initial_setup() 
void delete_update(size_type jd, bool force_refactorization, EEigenValType drop_eigen_val, PivotTolerances pivot_tols)
Should always succeed unless user gives a wrong value for drop_eigen_val. 
void init_diagonal(const Vector &diag)
Actually this calls init_identity( diag.space(), diag.norm_inf() ). 
size_type num_secant_updates_
const multi_vec_ptr_t S() const 
STANDARD_MEMBER_COMPOSITION_MEMBERS(bool, auto_rescaling)
Set whether automatic rescaling is used or not. 
const DMatrixSliceTri R() const 
const multi_vec_ptr_t Y() const 
std::ostream & output(std::ostream &out) const 
Inertia inertia() const 
Returns (0,0,rows()) 
void assert_initialized() const 
DenseLinAlgPack::VectorSliceTmpl< value_type > DVectorSlice
void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta) const 
AbstractLinAlgPack::value_type value_type
DenseLinAlgPack::DMatrixSliceTri DMatrixSliceTri
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)
DenseLinAlgPack::DMatrixSlice DMatrixSlice
DenseLinAlgPack::DMatrixSliceSym DMatrixSliceSym
value_type gamma_k() const 
void augment_update(const DVectorSlice *t, value_type alpha, bool force_refactorization, EEigenValType add_eigen_val, PivotTolerances pivot_tols)
Augment the matrix to add a row and column. 
void set_uninitialized()
Will set rows() == 0. 
size_type max_size() const 
void Vp_DtV(DVectorSlice *y, const DVectorSlice &x) const 
y += D * x