MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ConstrainedOptPack_MatrixSymPosDefLBFGSSerial.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Moocho: Multi-functional Object-Oriented arCHitecture for Optimization
5 // Copyright (2003) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef MATRIX_SYM_POS_DEF_LBFGS_H
43 #define MATRIX_SYM_POS_DEF_LBFGS_H
44 
45 #include <vector>
46 
49 #include "AbstractLinAlgPack/src/MatrixSymWithOpFactorized.hpp"
52 
53 namespace ConstrainedOptPack {
54 
119 class MatrixSymPosDefLBFGS
120  : public MatrixSymWithOpFactorized
121  , public MatrixSymSecant
122  , public MatrixSymAddDelUpdateable
123 {
124 public:
125 
126  // //////////////////////////////////////////////
127  // Constructors and initializers
128 
131  size_type max_size = 0
132  ,size_type m = 10
133  ,bool maintain_original = true
134  ,bool maintain_inverse = true
135  ,bool auto_rescaling = false
136  );
137 
143  STANDARD_MEMBER_COMPOSITION_MEMBERS( bool, auto_rescaling );
144 
181  void initial_setup(
182  size_type max_size = 0
183  ,size_type m = 10
184  ,bool maintain_original = true
185  ,bool maintain_inverse = true
186  ,bool auto_rescaling = false
187  );
188 
189  // //////////////////////////////////
190  // Representation access
191 
193  size_type m() const;
195  size_type m_bar() const;
197  size_type k_bar() const;
199  value_type gamma_k() const;
201  const DMatrixSlice S() const;
203  const DMatrixSlice Y() const;
205  bool maintain_original() const;
207  bool maintain_inverse() const;
210 
211  // /////////////////////////////////////////////////////
212  // Overridden from Matrix
213 
215  size_type rows() const;
216 
217  // /////////////////////////////////////////////////////////
220 
222  std::ostream& output(std::ostream& out) const;
224  MatrixOp& operator=(const MatrixOp& m);
226  void Vp_StMtV(DVectorSlice* vs_lhs, value_type alpha, BLAS_Cpp::Transp trans_rhs1
227  , const DVectorSlice& vs_rhs2, value_type beta) const;
228 
230 
231  // ////////////////////////////////////////////////////////////
234 
236  void V_InvMtV( DVectorSlice* v_lhs, BLAS_Cpp::Transp trans_rhs1
237  , const DVectorSlice& vs_rhs2) const;
238 
240 
241  // ///////////////////////////////////////////////////////////
244 
246  void init_identity( size_type n, value_type alpha );
253  void init_diagonal( const DVectorSlice& diag );
256 
257  // end Overridden from MatrixSymSecant
259 
260  // ////////////////////////////////////////////////////////
263 
265  void initialize(
266  value_type alpha
268  );
270  void initialize(
271  const DMatrixSliceSym &A
273  ,bool force_factorization
274  ,Inertia inertia
275  ,PivotTolerances pivot_tols
276  );
278  size_type max_size() const;
280  Inertia inertia() const;
282  void set_uninitialized();
291  void augment_update(
292  const DVectorSlice *t
293  ,value_type alpha
294  ,bool force_refactorization
295  ,EEigenValType add_eigen_val
296  ,PivotTolerances pivot_tols
297  );
299  void delete_update(
300  size_type jd
301  ,bool force_refactorization
302  ,EEigenValType drop_eigen_val
303  ,PivotTolerances pivot_tols
304  );
305 
307 
308 private:
309 
310  // //////////////////////////////////
311  // Private types
312 
313  // //////////////////////////////////
314  // Private data members
315 
316  bool maintain_original_; // If true, qualities needed for Bk will be maintained
317  bool original_is_updated_;// If true, qualities needed for Bk are already updated
318  bool maintain_inverse_; // If true, quantities needed for Hk will be maintained
319  bool inverse_is_updated_; // If true, quantities needed for Hk are already updated
320 
321  size_type n_max_, // The maximum size the matrix is allowed to become.
322  n_, // Size of the matrix. If 0 then is uninitialized
323  m_, // Maximum number of update vectors that can be stored.
324  m_bar_, // Current number of update vectors being stored.
325  // 0 <= m_bar <= m
326  k_bar_, // Position of the most recently stored update vector in S & Y
327  // 1 <= k_bar <= m_bar
328  num_secant_updates_; // Records the number of secant updates performed
329  value_type gamma_k_;// Scaling factor for Bo = (1/gamma_k) * I.
330 
331  DMatrix S_, // (n_max x m) Matrix of stored update vectors = [ s1, ..., sm ]
332  // S(:,k_bar) is the most recently stored s update vector
333  Y_, // (n_max x m) Matrix of stored update vectors = [ y1, ..., ym ]
334  // Y(:,k_bar) is the most recently stored y update vector
335  STY_, // (m x m) The matrix S'Y
336  STSYTY_;// ((m+1) x (m+1)) The strictly upper triangular part stores the
337  // upper triangular part Y'Y and the strictly lower triangular
338  // part stores the lower triangular part of S'S. The diagonal
339  // can be used for workspace.
340 
341  mutable bool Q_updated_; // True if Q has been updated for the most current update.
342  mutable DMatrix QJ_; // Used to store factorization of the schur complement of Q.
343 
344  mutable DVector work_; // workspace for performing operations.
345 
346  // //////////////////////////////////
347  // Private member functions
348 
349  // Access to important matrices.
350 
352  const DMatrixSliceTri R() const;
354  const DMatrixSliceTri Lb() const;
356  DMatrixSlice STY();
358  const DMatrixSlice STY() const;
362  const DMatrixSliceSym STS() const;
366  const DMatrixSliceSym YTY() const;
368  void V_invQtV( DVectorSlice* y, const DVectorSlice& x ) const;
370  void Vp_DtV( DVectorSlice* y, const DVectorSlice& x ) const;
371 
372  // Updates
373 
375  void update_Q() const;
376 
378  void assert_initialized() const;
379 
380 }; // end class MatrixSymPosDefLBFGS
381 
382 // //////////////////////////////////////////////
383 // Inline member functions
384 
385 inline
387 {
388  return m_;
389 }
390 
391 inline
393 {
394  return m_bar_;
395 }
396 
397 inline
399 {
400  return k_bar_;
401 }
402 
403 inline
405 {
406  return gamma_k_;
407 }
408 
409 inline
411 {
412  return S_(1,n_,1,m_bar_);
413 }
414 
415 inline
417 {
418  return Y_(1,n_,1,m_bar_);
419 }
420 
421 inline
423 {
424  return maintain_original_;
425 }
426 
427 inline
429 {
430  return maintain_inverse_;
431 }
432 
433 inline
435 {
436  return num_secant_updates_;
437 }
438 
439 
440 } // end namespace ConstrainedOptPack
441 
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
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 V_InvMtV(VectorMutable *v_lhs, BLAS_Cpp::Transp trans_rhs1, const Vector &v_rhs2) const
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() ).
STANDARD_MEMBER_COMPOSITION_MEMBERS(bool, auto_rescaling)
Set whether automatic rescaling is used or not.
std::ostream * out
Inertia inertia() const
Returns (0,0,rows())
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
size_type num_secant_updates() const
Returns the total number of successful secant updates performed since this->init_identity() was calle...
DenseLinAlgPack::DMatrixSlice DMatrixSlice
DenseLinAlgPack::DMatrixSliceSym DMatrixSliceSym
Transp
TRANS.
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.
void Vp_DtV(DVectorSlice *y, const DVectorSlice &x) const
y += D * x