76 #include "ConstrainedOptPack_MatrixSymPosDefLBFGS.hpp"
77 #include "AbstractLinAlgPack_BFGS_helpers.hpp"
78 #include "AbstractLinAlgPack_VectorDenseEncap.hpp"
79 #include "AbstractLinAlgPack_MatrixOpOut.hpp"
80 #include "AbstractLinAlgPack_VectorStdOps.hpp"
81 #include "AbstractLinAlgPack_LinAlgOpPack.hpp"
82 #include "DenseLinAlgPack_LinAlgOpPack.hpp"
83 #include "DenseLinAlgPack_DMatrixOut.hpp"
84 #include "DenseLinAlgLAPack.hpp"
85 #include "Teuchos_Workspace.hpp"
86 #include "Teuchos_Assert.hpp"
90 using DenseLinAlgPack::DVectorSlice;
100 void comp_Cb(
const DMatrixSlice& Lb,
const DVectorSlice& Db_diag
101 , DMatrixSlice* Cb );
105 namespace ConstrainedOptPack {
111 const DMatrixSliceTri MatrixSymPosDefLBFGS::R()
const
113 return DenseLinAlgPack::tri( STY_(1,m_bar_,1,m_bar_), BLAS_Cpp::upper, BLAS_Cpp::nonunit );
117 const DMatrixSliceTri MatrixSymPosDefLBFGS::Lb()
const
119 return DenseLinAlgPack::tri( STY_(2,m_bar_,1,m_bar_-1), BLAS_Cpp::lower, BLAS_Cpp::nonunit );
123 DMatrixSlice MatrixSymPosDefLBFGS::STY()
125 return STY_(1,m_bar_,1,m_bar_);
129 const DMatrixSlice MatrixSymPosDefLBFGS::STY()
const
131 return STY_(1,m_bar_,1,m_bar_);
135 DMatrixSliceSym MatrixSymPosDefLBFGS::STS()
137 return DenseLinAlgPack::nonconst_sym( STSYTY_(2,m_bar_+1,1,m_bar_),BLAS_Cpp::lower );
141 const DMatrixSliceSym MatrixSymPosDefLBFGS::STS()
const
143 return DenseLinAlgPack::sym( STSYTY_(2,m_bar_+1,1,m_bar_),BLAS_Cpp::lower );
147 DMatrixSliceSym MatrixSymPosDefLBFGS::YTY()
149 return DenseLinAlgPack::nonconst_sym( STSYTY_(1,m_bar_,2,m_bar_+1),BLAS_Cpp::upper );
153 const DMatrixSliceSym MatrixSymPosDefLBFGS::YTY()
const
155 return DenseLinAlgPack::sym( STSYTY_(1,m_bar_,2,m_bar_+1),BLAS_Cpp::upper );
163 ,
bool maintain_original
164 ,
bool maintain_inverse
168 initial_setup(m,maintain_original,maintain_inverse,auto_rescaling);
173 bool maintain_original,
174 bool maintain_inverse,
180 !maintain_original && !maintain_inverse, std::invalid_argument
181 ,
"MatrixSymPosDefLBFGS::initial_setup(...) : "
182 "Error, both maintain_original and maintain_inverse can not both be false!" );
184 m < 1, std::invalid_argument
185 ,
"MatrixSymPosDefLBFGS::set_num_updates_stored(m) : "
186 "Error, the number of storage locations must be > 0" );
187 vec_spc_ = Teuchos::null;
192 num_secant_updates_= 0;
193 auto_rescaling_ = auto_rescaling;
205 assert_initialized();
206 out <<
"\n*** Limited Memory BFGS matrix\n";
207 out <<
"\nConversion to dense =\n";
208 MatrixOp::output(out);
209 out <<
"\nStored quantities\n"
212 <<
"\nm_bar = " << m_bar_
213 <<
"\ngamma_k = " << gamma_k_ << std::endl;
215 out <<
"\nS =\n" << *
S()
217 <<
"\nS'Y =\n" << STY_(1,m_bar_,1,m_bar_)
218 <<
"\nlower(S'S) \\ zero diagonal \\ upper(Y'Y) =\n"
219 << STSYTY_(1,m_bar_+1,1,m_bar_+1)
220 <<
"\nQ updated? = " << Q_updated_ << std::endl;
222 out <<
"\nCholesky of schur complement of Q, QJ =\n" << QJ_(1,m_bar_,1,m_bar_);
231 if( p_m ==
this )
return *
this;
233 auto_rescaling_ = p_m->auto_rescaling_;
234 maintain_original_ = p_m->maintain_original_;
235 original_is_updated_ = p_m->original_is_updated_;
236 maintain_inverse_ = p_m->maintain_inverse_;
237 inverse_is_updated_ = p_m->inverse_is_updated_;
238 vec_spc_ = p_m->vec_spc_.get() ? p_m->vec_spc_->clone() : Teuchos::null;
241 m_bar_ = p_m->m_bar_;
242 num_secant_updates_ = p_m->num_secant_updates_;
243 gamma_k_ = p_m->gamma_k_;
247 STSYTY_ = p_m->STSYTY_;
248 Q_updated_ = p_m->Q_updated_;
253 true,std::invalid_argument
254 ,
"MatrixSymPosDefLBFGS::operator=(const MatrixOp& mwo) : Error, "
255 "The concrete type of mwo \'" <<
typeName(mwo) <<
"\' is not "
256 "as subclass of MatrixSymPosDefLBFGS as required" );
265 ,
const Vector& x, value_type beta
271 using LinAlgOpPack::V_StMtV;
272 using LinAlgOpPack::V_MtV;
273 typedef VectorDenseEncap vde;
274 typedef VectorDenseMutableEncap vdme;
278 assert_initialized();
302 invgk = 1.0 / gamma_k_;
322 Workspace<value_type> t1_ws(wss,2*mb);
323 DVectorSlice t1(&t1_ws[0],t1_ws.size());
324 Workspace<value_type> t2_ws(wss,2*mb);
325 DVectorSlice t2(&t2_ws[0],t2_ws.size());
327 VectorSpace::vec_mut_ptr_t
328 t = S->space_rows().create_member();
334 t1(1,mb) = vde(*t)();
336 t1(mb+1,2*mb) = vde(*t)();
342 (vdme(*t)() = t2(1,mb));
346 (vdme(*t)() = t2(mb+1,2*mb));
363 using LinAlgOpPack::V_MtV;
364 using LinAlgOpPack::V_StMtV;
365 using LinAlgOpPack::Vp_MtV;
367 typedef VectorDenseEncap vde;
368 typedef VectorDenseMutableEncap vdme;
372 assert_initialized();
405 V_StV( y, gamma_k_, x );
419 Workspace<value_type> t1_ws(wss,2*mb);
420 DVectorSlice t1(&t1_ws[0],t1_ws.size());
421 Workspace<value_type> t2_ws(wss,mb);
422 DVectorSlice t2(&t2_ws[0],t2_ws.size());
423 Workspace<value_type> t3_ws(wss,mb);
424 DVectorSlice t3(&t3_ws[0],t3_ws.size());
425 Workspace<value_type> t4_ws(wss,mb);
426 DVectorSlice t4(&t4_ws[0],t4_ws.size());
427 Workspace<value_type> t5_ws(wss,mb);
428 DVectorSlice t5(&t5_ws[0],t5_ws.size());
430 VectorSpace::vec_mut_ptr_t
431 t = S->space_rows().create_member();
433 const DMatrixSliceTri
436 const DMatrixSliceSym
442 t1(1,mb) = vde(*t)();
444 t1(mb+1,2*mb) = vde(*t)();
450 V_mV( &t3, t1(mb+1,2*mb) );
451 V_InvMtV( &t3, R, BLAS_Cpp::trans, t3 );
460 V_InvMtV( &t5, R, BLAS_Cpp::trans, t4 );
481 alpha <= 0.0, std::invalid_argument
482 ,
"MatrixSymPosDefLBFGS::init_identity(n,alpha) : Error, "
483 "alpha = " << alpha <<
" <= 0 is not allowed!" );
486 vec_spc_ = space_diag.clone();
490 S_ = vec_spc_->create_members(m_);
491 Y_ = vec_spc_->create_members(m_);
494 STY_.resize( m_, m_ );
495 STSYTY_.resize( m_+1, m_+1 );
496 STSYTY_.diag(0) = 0.0;
498 gamma_k_ = 1.0/alpha;
503 n_ = vec_spc_->dim();
504 original_is_updated_ =
true;
505 inverse_is_updated_ =
true;
506 num_secant_updates_ = 0;
515 VectorMutable* s, VectorMutable* y, VectorMutable* Bs
518 using AbstractLinAlgPack::BFGS_sTy_suff_p_d;
520 using LinAlgOpPack::V_MtV;
524 assert_initialized();
529 std::ostringstream omsg;
530 if( !BFGS_sTy_suff_p_d(*s,*y,&sTy,&omsg,
"MatrixSymPosDefLBFGS::secant_update(...)" ) ) {
531 throw UpdateSkippedException( omsg.str() );
541 {
for(
size_type k = 1; k <= m_-1; ++k ) {
542 S_->col(k) = S_->col(k+1);
543 Y_->col(k) = Y_->col(k+1);
544 STY_.col(k)(1,m_-1) = STY_.col(k+1)(2,m_);
545 STSYTY_.col(k)(k+1,m_) = STSYTY_.col(k+1)(k+2,m_+1);
546 STSYTY_.col(k+1)(1,k) = STSYTY_.col(k+2)(2,k+1);
555 *S_->col(m_bar_) = *s;
556 *Y_->col(m_bar_) = *y;
581 VectorSpace::vec_mut_ptr_t
582 t = S->space_rows().create_member();
586 STY_.col(m_bar_)(1,m_bar_) = VectorDenseEncap(*t)();
590 STY_.row(m_bar_)(1,m_bar_) = VectorDenseEncap(*t)();
610 Workspace<value_type> work_ws(wss,m_bar_);
611 DVectorSlice work(&work_ws[0],work_ws.size());
615 work = VectorDenseEncap(*t)();
618 STSYTY_.row(m_bar_+1)(1,m_bar_) = work;
620 STSYTY_.col(m_bar_)(m_bar_+1,m_bar_+1) = work(m_bar_,m_bar_);
644 work = VectorDenseEncap(*t)();
647 STSYTY_.col(m_bar_+1)(1,m_bar_) = work;
649 STSYTY_.row(m_bar_)(m_bar_+1,m_bar_+1) = work(m_bar_,m_bar_);
656 gamma_k_ = STY_(m_bar_,m_bar_) / STSYTY_(m_bar_,m_bar_+1);
662 num_secant_updates_++;
676 void MatrixSymPosDefLBFGS::Vp_DtV( DVectorSlice* y,
const DVectorSlice& x )
const
678 DenseLinAlgPack::Vp_MtV_assert_sizes(
681 DVectorSlice::const_iterator
682 d_itr = STY_.diag(0).begin(),
684 DVectorSlice::iterator
687 while( y_itr != y->end() )
688 *y_itr++ += (*d_itr++) * (*x_itr++);
716 void MatrixSymPosDefLBFGS::update_Q()
const
718 using DenseLinAlgPack::tri;
719 using DenseLinAlgPack::tri_ele;
749 if( QJ_.rows() < m_ )
750 QJ_.resize( m_, m_ );
779 comp_Cb( STY_(2,mb,1,mb-1), STY_.diag(0)(1,mb-1), &C(2,mb,2,mb) );
783 const DMatrixSliceSym &STS = this->STS();
784 Mp_StM( &C, (1/gamma_k_), tri( STS.gms(), STS.uplo(), BLAS_Cpp::nonunit )
791 DMatrixSliceTriEle C_upper = tri_ele( C, BLAS_Cpp::upper );
793 DenseLinAlgLAPack::potrf( &C_upper );
797 true, UpdateFailedException
798 ,
"Error, the factorization of Q which should be s.p.d. failed with"
799 " the error message: {" << fe.what() <<
"}";
806 void MatrixSymPosDefLBFGS::V_invQtV( DVectorSlice* x,
const DVectorSlice& y )
const
808 using DenseLinAlgPack::sym;
809 using DenseLinAlgPack::tri;
814 using LinAlgOpPack::V_MtV;
832 x2 = (*x)(mb+1,2*mb);
844 DVectorSlice::const_iterator
845 d_itr = STY_.diag(0).begin(),
847 DVectorSlice::iterator
849 while( x1_itr != x1.end() )
850 *x1_itr++ = *y2_itr++ / *d_itr++;
867 std::copy( x1a.rbegin(), x1a.rend(), x1b.rbegin() );
868 V_MtV( &x1b, Lb(), BLAS_Cpp::no_trans, x1b );
876 const DMatrixSliceTri J = tri( QJ_(1,mb,1,mb), BLAS_Cpp::upper, BLAS_Cpp::nonunit );
880 V_InvMtV( &x1, J, BLAS_Cpp::no_trans, x1 );
903 DVectorSlice::const_iterator
904 d_itr = STY_.diag(0).begin();
905 DVectorSlice::iterator
907 while( x2_itr != x2.end() )
908 *x2_itr++ /= *d_itr++;
912 void MatrixSymPosDefLBFGS::assert_initialized()
const
915 throw std::logic_error(
"MatrixSymPosDefLBFGS::assert_initialized() : "
916 "Error, matrix not initialized" );
924 const DMatrixSlice& Lb,
const DVectorSlice& Db_diag
946 typedef DenseLinAlgPack::value_type value_type;
950 const size_type p = Db_diag.dim();
952 for( size_type i = 1; i <= p; ++i ) {
953 for( size_type j = i; j <= p; ++j ) {
954 value_type &c = (*Cb)(i,j) = 0.0;
955 for( size_type k = 1; k <= i; ++k ) {
956 c += Lb(i,k) * Lb(j,k) / Db_diag(k);
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.
void Vt_S(VectorMutable *v_lhs, const value_type &alpha)
const VectorSpace & space_cols() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void Vp_StV(VectorMutable *v_lhs, const value_type &alpha, const Vector &v_rhs)
void V_InvMtV(VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const
void V_StV(VectorMutable *v_lhs, value_type alpha, const V &V_rhs)
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() ).
void V_mV(VectorMutable *v_lhs, const V &V_rhs)
const multi_vec_ptr_t S() const
const multi_vec_ptr_t Y() const
void V_InvMtV(VectorMutable *v_lhs, const MatrixNonsing &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2)
void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
void Vp_StMtV(VectorMutable *v_lhs, value_type alpha, const MatrixOp &M_rhs1, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2, value_type beta=1.0)
std::ostream & output(std::ostream &out) const
friend friend void Mp_StM(MatrixOp *mwo_lhs, value_type alpha, const MatrixOp &M_rhs, BLAS_Cpp::Transp trans_rhs)
value_type dot(const Vector &v_rhs1, const Vector &v_rhs2)
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.
bool maintain_inverse() const
MatrixOp & operator=(const MatrixOp &mwo)
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
TEUCHOSCORE_LIB_DLL_EXPORT Teuchos::RCP< WorkspaceStore > get_default_workspace_store()
std::string typeName(const T &t)
void Vp_V(VectorMutable *v_lhs, const V &V_rhs)