MOOCHO (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ConstrainedOptPack_QPSolverRelaxedQPSchur.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 QP_SOLVER_RELAXED_QP_SCHUR_H
43 #define QP_SOLVER_RELAXED_QP_SCHUR_H
44 
54 
55 namespace ConstrainedOptPack {
56 
64 public:
65 
71  class InitKKTSystem {
72  public:
74  typedef std::vector<size_type> i_x_free_t;
76  typedef std::vector<size_type> i_x_fixed_t;
78  typedef std::vector<EBounds> bnd_fixed_t;
80  typedef std::vector<size_type> j_f_decomp_t;
85  virtual ~InitKKTSystem() {}
159  virtual void initialize_kkt_system(
160  const Vector &g
161  ,const MatrixOp &G
162  ,value_type etaL
163  ,const Vector *dL
164  ,const Vector *dU
165  ,const MatrixOp *F
166  ,BLAS_Cpp::Transp trans_F
167  ,const Vector *f
168  ,const Vector *d
169  ,const Vector *nu
170  ,size_type *n_R
171  ,i_x_free_t *i_x_free
172  ,i_x_fixed_t *i_x_fixed
173  ,bnd_fixed_t *bnd_fixed
174  ,j_f_decomp_t *j_f_decomp
175  ,DVector *b_X
176  ,Ko_ptr_t *Ko
177  ,DVector *fo
178  ) const = 0;
179 
180  }; // end class InitKKTSystem
181 
192  public:
193  // ToDo: Create method reinitailze_kkt_system(...)
194  }; // end class ReinitKKTSystem
195 
199 
203 
206  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, max_qp_iter_frac );
207 
210  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, max_real_runtime );
211 
216  ,inequality_pick_policy
217  );
218 
221  USE_INPUT_ARG = -1 // Use the value input to solve_qp(...)
222  ,NO_OUTPUT = 0 //
223  ,OUTPUT_BASIC_INFO = 1 // values sent to QPSchur::solve_qp(...)
228  };
229 
233 
237 
241 
245 
249 
254 
259  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, huge_primal_step );
260 
266 
270 
274 
278 
282  STANDARD_MEMBER_COMPOSITION_MEMBERS( size_type, iter_refine_min_iter );
283 
287  STANDARD_MEMBER_COMPOSITION_MEMBERS( size_type, iter_refine_max_iter );
288 
292  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, iter_refine_opt_tol );
293 
297  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, iter_refine_feas_tol );
298 
302  STANDARD_MEMBER_COMPOSITION_MEMBERS( bool, iter_refine_at_solution );
303 
308  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, pivot_warning_tol );
309 
314  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, pivot_singular_tol );
315 
320  STANDARD_MEMBER_COMPOSITION_MEMBERS( value_type, pivot_wrong_inertia_tol );
321 
325  STANDARD_MEMBER_COMPOSITION_MEMBERS( bool, add_equalities_initially );
326 
329  const init_kkt_sys_ptr_t& init_kkt_sys = Teuchos::null
330  ,const constraints_ptr_t& constraints = Teuchos::rcp(new QPSchurPack::ConstraintsRelaxedStd)
331  ,value_type max_qp_iter_frac = 10.0
332  ,value_type max_real_runtime = 1e+20
334  inequality_pick_policy
336  ,ELocalOutputLevel print_level = USE_INPUT_ARG // Deduce from input arguments
337  ,value_type bounds_tol = -1.0 // use default
338  ,value_type inequality_tol = -1.0 // use default
339  ,value_type equality_tol = -1.0 // use default
340  ,value_type loose_feas_tol = -1.0 // use default
341  ,value_type dual_infeas_tol = -1.0 // use default
342  ,value_type huge_primal_step = -1.0 // use defalut
343  ,value_type huge_dual_step = -1.0 // use default
344  ,value_type bigM = 1e+10
345  ,value_type warning_tol = 1e-10
346  ,value_type error_tol = 1e-5
347  ,size_type iter_refine_min_iter = 1
348  ,size_type iter_refine_max_iter = 3
349  ,value_type iter_refine_opt_tol = 1e-12
350  ,value_type iter_refine_feas_tol = 1e-12
351  ,bool iter_refine_at_solution = true
352  ,value_type pivot_warning_tol = 1e-8
353  ,value_type pivot_singular_tol = 1e-11
354  ,value_type pivot_wrong_inertia_tol = 1e-11
355  ,bool add_equalities_initially= true
356  );
357 
360 
363 
365  QPSolverStats get_qp_stats() const;
367  void release_memory();
368 
370 
371 protected:
372 
375 
378  std::ostream* out, EOutputLevel olevel, ERunTests test_what
379  ,const Vector& g, const MatrixSymOp& G
380  ,value_type etaL
381  ,const Vector* dL, const Vector* dU
382  ,const MatrixOp* E, BLAS_Cpp::Transp trans_E, const Vector* b
383  ,const Vector* eL, const Vector* eU
384  ,const MatrixOp* F, BLAS_Cpp::Transp trans_F, const Vector* f
385  ,value_type* obj_d
386  ,value_type* eta, VectorMutable* d
387  ,VectorMutable* nu
388  ,VectorMutable* mu, VectorMutable* Ed
389  ,VectorMutable* lambda, VectorMutable* Fd
390  );
391 
393 
394 private:
395 
396  // ////////////////////////////
397  // Private data members
398 
403  VectorMutableDense bigM_vec_;
409 
410 }; // end class QPSolverRelaxedQPSchur
411 
412 } // end namespace ConstrainedOptPack
413 
414 #endif // QP_SOLVER_RELAXED_QP_SCHUR_H
AbstractLinAlgPack::size_type size_type
ERunTests
Enumeration for if to run internal tests or not.
QPSolverStats::ESolutionType imp_solve_qp(std::ostream *out, EOutputLevel olevel, ERunTests test_what, const Vector &g, const MatrixSymOp &G, value_type etaL, const Vector *dL, const Vector *dU, const MatrixOp *E, BLAS_Cpp::Transp trans_E, const Vector *b, const Vector *eL, const Vector *eU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const Vector *f, value_type *obj_d, value_type *eta, VectorMutable *d, VectorMutable *nu, VectorMutable *mu, VectorMutable *Ed, VectorMutable *lambda, VectorMutable *Fd)
Solves Quadratic Programs (QPs) of several different forms while allowing a relaxation of the constra...
virtual void initialize_kkt_system(const Vector &g, const MatrixOp &G, value_type etaL, const Vector *dL, const Vector *dU, const MatrixOp *F, BLAS_Cpp::Transp trans_F, const Vector *f, const Vector *d, const Vector *nu, size_type *n_R, i_x_free_t *i_x_free, i_x_fixed_t *i_x_fixed, bnd_fixed_t *bnd_fixed, j_f_decomp_t *j_f_decomp, DVector *b_X, Ko_ptr_t *Ko, DVector *fo) const =0
Initializes the KKT system.
Solves Quadratic Programming (QP) problems using QPSchur.
This class maintains the factorization of symmetric indefinite matrix using a Bunch & Kaufman factori...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Interface for the object that can reform an initial KKT system dynamically {abstract}.
std::ostream * out
Matrix class for non-singular Hessian matrix augmented with a terms for "Big M" relaxation variables...
Class for storing statistics about a run of a (active set?) QP solver.
Constraints subclass that is used to represent generic varaible bounds, and general inequality and eq...
EOutputLevel
Enumeration for the amount of output to create from solve_qp().
STANDARD_MEMBER_COMPOSITION_MEMBERS(value_type, max_qp_iter_frac)
Set the maximum number of QP iterations as max_itr = max_qp_iter_frac * n.
AbstractLinAlgPack::value_type value_type
General (and flexible) implementation class for a QPSchur QP problem.
STANDARD_COMPOSITION_MEMBERS(InitKKTSystem, init_kkt_sys)
Strategy object that sets up the initial KKT system.
Transp
TRANS.
Solves a Quadratic Program with a dual QP method using a schur complement factorization.
ESolutionType
Enumeration for the type of point returned from solve_qp(...).
Interface for the object that forms the initial KKT system {abstract}.
QPSolverRelaxedQPSchur(const init_kkt_sys_ptr_t &init_kkt_sys=Teuchos::null, const constraints_ptr_t &constraints=Teuchos::rcp(new QPSchurPack::ConstraintsRelaxedStd), value_type max_qp_iter_frac=10.0, value_type max_real_runtime=1e+20, QPSchurPack::ConstraintsRelaxedStd::EInequalityPickPolicy inequality_pick_policy=QPSchurPack::ConstraintsRelaxedStd::ADD_BOUNDS_THEN_MOST_VIOLATED_INEQUALITY, ELocalOutputLevel print_level=USE_INPUT_ARG, value_type bounds_tol=-1.0, value_type inequality_tol=-1.0, value_type equality_tol=-1.0, value_type loose_feas_tol=-1.0, value_type dual_infeas_tol=-1.0, value_type huge_primal_step=-1.0, value_type huge_dual_step=-1.0, value_type bigM=1e+10, value_type warning_tol=1e-10, value_type error_tol=1e-5, size_type iter_refine_min_iter=1, size_type iter_refine_max_iter=3, value_type iter_refine_opt_tol=1e-12, value_type iter_refine_feas_tol=1e-12, bool iter_refine_at_solution=true, value_type pivot_warning_tol=1e-8, value_type pivot_singular_tol=1e-11, value_type pivot_wrong_inertia_tol=1e-11, bool add_equalities_initially=true)