Logo ROOT   6.12/07
Reference Guide
TQpSolverBase.h
Go to the documentation of this file.
1 // @(#)root/quadp:$Id$
2 // Author: Eddy Offermann May 2004
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /*************************************************************************
13  * Parts of this file are copied from the OOQP distribution and *
14  * are subject to the following license: *
15  * *
16  * COPYRIGHT 2001 UNIVERSITY OF CHICAGO *
17  * *
18  * The copyright holder hereby grants you royalty-free rights to use, *
19  * reproduce, prepare derivative works, and to redistribute this software*
20  * to others, provided that any changes are clearly documented. This *
21  * software was authored by: *
22  * *
23  * E. MICHAEL GERTZ gertz@mcs.anl.gov *
24  * Mathematics and Computer Science Division *
25  * Argonne National Laboratory *
26  * 9700 S. Cass Avenue *
27  * Argonne, IL 60439-4844 *
28  * *
29  * STEPHEN J. WRIGHT swright@cs.wisc.edu *
30  * Computer Sciences Department *
31  * University of Wisconsin *
32  * 1210 West Dayton Street *
33  * Madison, WI 53706 FAX: (608)262-9777 *
34  * *
35  * Any questions or comments may be directed to one of the authors. *
36  * *
37  * ARGONNE NATIONAL LABORATORY (ANL), WITH FACILITIES IN THE STATES OF *
38  * ILLINOIS AND IDAHO, IS OWNED BY THE UNITED STATES GOVERNMENT, AND *
39  * OPERATED BY THE UNIVERSITY OF CHICAGO UNDER PROVISION OF A CONTRACT *
40  * WITH THE DEPARTMENT OF ENERGY. *
41  *************************************************************************/
42 
43 #ifndef ROOT_TQpSolverBase
44 #define ROOT_TQpSolverBase
45 
46 #include "TError.h"
47 
48 #include "TObject.h"
49 #include "TQpVar.h"
50 #include "TQpDataBase.h"
51 #include "TQpResidual.h"
52 #include "TQpProbBase.h"
53 #include "TQpLinSolverBase.h"
54 
55 #include "TMatrixD.h"
56 
57 ///////////////////////////////////////////////////////////////////////////
58 // //
59 // Abstract base class for QP solver using interior-point //
60 // //
61 ///////////////////////////////////////////////////////////////////////////
62 
64 {
70 };
71 
72 class TQpSolverBase : public TObject
73 {
74 
75 protected:
77 
78  Double_t fDnorm; // norm of problem data
79 
80  Double_t fMutol; // termination parameters
82 
83  Double_t fGamma_f; // parameters associated with the step length heuristic
85  Double_t fPhi; // merit function, defined as the sum of the complementarity gap
86  // the residual norms, divided by (1+norm of problem data)
87  Int_t fMaxit; // maximum number of iterations allowed
88 
89  Double_t *fMu_history; //[fMaxit] history of values of mu obtained on all iterations to date
90  Double_t *fRnorm_history; //[fMaxit] history of values of residual norm obtained on all iterations to date
91  Double_t *fPhi_history; //[fMaxit] history of values of phi obtained on all iterations to date
92 
93  Double_t *fPhi_min_history; //[fMaxit] the i-th entry of this array contains the minimum value of phi
94  // encountered by the algorithm on or before iteration i
95 
96 public:
97  Int_t fIter; // iteration counter
98 
99  TQpSolverBase();
100  TQpSolverBase(const TQpSolverBase &another);
101 
102  virtual ~TQpSolverBase();
103 
104  virtual void Start (TQpProbBase *formulation,
105  TQpVar *iterate,TQpDataBase *prob,
106  TQpResidual *resid,TQpVar *step);
107  // starting point heuristic
108  virtual void DefStart (TQpProbBase *formulation,
109  TQpVar *iterate,TQpDataBase *prob,
110  TQpResidual *resid,TQpVar *step);
111  // default starting point heuristic
112  virtual void SteveStart (TQpProbBase *formulation,
113  TQpVar *iterate,TQpDataBase *prob,
114  TQpResidual *resid,TQpVar *step);
115  // alternative starting point heuristic
116  virtual void DumbStart (TQpProbBase *formulation,
117  TQpVar *iterate,TQpDataBase *prob,
118  TQpResidual *resid,TQpVar *step);
119  // alternative starting point heuristic: sets the
120  // "complementary" variables to a large positive value
121  // (based on the norm of the problem data) and the
122  // remaining variables to zero
123 
124  virtual Int_t Solve (TQpDataBase *prob,TQpVar *iterate,
125  TQpResidual *resids) = 0;
126  // implements the interior-point method for solving the QP
127 
128  virtual Double_t FinalStepLength
129  (TQpVar *iterate,TQpVar *step);
130  // Mehrotra's heuristic to calculate the final step
131 
132  virtual void DoMonitor (TQpDataBase *data,TQpVar *vars,
133  TQpResidual *resids,Double_t alpha,
135  Int_t stop_code,Int_t level);
136  // perform monitor operation at each interior-point iteration
137  virtual void DefMonitor (TQpDataBase *data,TQpVar *vars,
138  TQpResidual *resids,Double_t alpha,
140  Int_t stop_code,Int_t level) = 0;
141  // default monitor: prints out one line of information
142  // on each interior-point iteration
143  virtual Int_t DoStatus (TQpDataBase *data,TQpVar *vars,
144  TQpResidual *resids,Int_t i,Double_t mu,
145  Int_t level);
146  // this method called to test for convergence status at
147  // at the end of each interior-point iteration
148  virtual Int_t DefStatus (TQpDataBase *data,TQpVar *vars,
149  TQpResidual *resids,Int_t i,Double_t mu,
150  Int_t level);
151  // default method for checking status. May be replaced
152  // by a user-defined method
153 
155  void SetMuTol (Double_t m) { fMutol = m; }
156  Double_t GetMuTol () { return fMutol; }
157 
158  void SetArTol (Double_t ar) { fArtol = ar; }
159  Double_t GetArTol () { return fArtol; }
160  Double_t DataNorm () { return fDnorm; }
161 
162  TQpSolverBase &operator= (const TQpSolverBase &source);
163 
164  ClassDef(TQpSolverBase,1) // Qp Solver class
165 };
166 #endif
Double_t * fRnorm_history
Definition: TQpSolverBase.h:90
int iterate(rng_state_t *X)
Definition: mixmax.icc:34
TQpSolverBase()
Default constructor.
virtual void DoMonitor(TQpDataBase *data, TQpVar *vars, TQpResidual *resids, Double_t alpha, Double_t sigma, Int_t i, Double_t mu, Int_t stop_code, Int_t level)
Monitor progress / convergence aat each interior-point iteration.
auto * m
Definition: textangle.C:8
Double_t GetMuTol()
virtual void SteveStart(TQpProbBase *formulation, TQpVar *iterate, TQpDataBase *prob, TQpResidual *resid, TQpVar *step)
Starting point algoritm according to Stephen Wright.
virtual void DefStart(TQpProbBase *formulation, TQpVar *iterate, TQpDataBase *prob, TQpResidual *resid, TQpVar *step)
Default starting point.
int Int_t
Definition: RtypesCore.h:41
virtual Double_t FinalStepLength(TQpVar *iterate, TQpVar *step)
Implements a version of Mehrotra starting point heuristic, modified to ensure identical steps in the ...
virtual Int_t DoStatus(TQpDataBase *data, TQpVar *vars, TQpResidual *resids, Int_t i, Double_t mu, Int_t level)
Tests for termination.
TQpSolverBase & operator=(const TQpSolverBase &source)
Assignment operator.
ETerminationCode
Definition: TQpSolverBase.h:63
#define ClassDef(name, id)
Definition: Rtypes.h:320
const Double_t sigma
virtual ~TQpSolverBase()
Deconstructor.
Double_t fGamma_a
Definition: TQpSolverBase.h:84
virtual void Start(TQpProbBase *formulation, TQpVar *iterate, TQpDataBase *prob, TQpResidual *resid, TQpVar *step)
Implements a default starting-point heuristic.
virtual void DumbStart(TQpProbBase *formulation, TQpVar *iterate, TQpDataBase *prob, TQpResidual *resid, TQpVar *step)
Alternative starting point heuristic: sets the "complementary" variables to a large positive value (b...
void SetMuTol(Double_t m)
Double_t fGamma_f
Definition: TQpSolverBase.h:83
Double_t fDnorm
Definition: TQpSolverBase.h:78
Double_t DataNorm()
virtual Int_t DefStatus(TQpDataBase *data, TQpVar *vars, TQpResidual *resids, Int_t i, Double_t mu, Int_t level)
Default status method.
Double_t * fMu_history
Definition: TQpSolverBase.h:89
TQpLinSolverBase * fSys
Definition: TQpSolverBase.h:76
void SetArTol(Double_t ar)
virtual void DefMonitor(TQpDataBase *data, TQpVar *vars, TQpResidual *resids, Double_t alpha, Double_t sigma, Int_t i, Double_t mu, Int_t stop_code, Int_t level)=0
double Double_t
Definition: RtypesCore.h:55
virtual Int_t Solve(TQpDataBase *prob, TQpVar *iterate, TQpResidual *resids)=0
Double_t fArtol
Definition: TQpSolverBase.h:81
Mother of all ROOT objects.
Definition: TObject.h:37
TQpLinSolverBase * GetLinearSystem()
Definition: TQpVar.h:59
Double_t fPhi
Definition: TQpSolverBase.h:85
Double_t fMutol
Definition: TQpSolverBase.h:80
Double_t * fPhi_history
Definition: TQpSolverBase.h:91
Double_t * fPhi_min_history
Definition: TQpSolverBase.h:93
Double_t GetArTol()