Logo ROOT   6.12/07
Reference Guide
TQpVar.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_TQpVar
44 #define ROOT_TQpVar
45 
46 #include "TError.h"
47 
48 #include "TMatrixD.h"
49 #include "TVectorD.h"
50 
51 ///////////////////////////////////////////////////////////////////////////
52 // //
53 // Class containing the variables for the general QP formulation //
54 // In terms of in our abstract problem formulation, these variables are //
55 // the vectors x, y, z and s. //
56 // //
57 ///////////////////////////////////////////////////////////////////////////
58 
59 class TQpVar : public TObject
60 {
61 
62 protected:
70 
71  // these variables will be "Used" not copied
76 
77  static Double_t StepBound (TVectorD &v,TVectorD &dir,Double_t maxStep);
78  static Double_t FindBlocking (TVectorD &w,TVectorD &wstep,TVectorD &u,TVectorD &ustep,
79  Double_t maxStep,Double_t &w_elt,Double_t &wstep_elt,Double_t &u_elt,
80  Double_t &ustep_elt,int& first_or_second);
81  static Double_t FindBlockingSub(Int_t n,Double_t *w,Int_t incw,Double_t *wstep,Int_t incwstep,
82  Double_t *u,Int_t incu,Double_t *ustep,Int_t incustep,
83  Double_t maxStep,Double_t &w_elt,Double_t &wstep_elt,
84  Double_t &u_elt,Double_t &ustep_elt,Int_t &first_or_second);
85 
86 public:
87 
88  Int_t fNComplementaryVariables; // number of complementary primal-dual variables.
89 
90  // these variables will be "Used" not copied
95 
98 
101 
104 
107 
108  TQpVar();
109  // constructor in which the data and variable pointers are set to point to the given arguments
110  TQpVar(TVectorD &x_in,TVectorD &s_in,TVectorD &y_in,TVectorD &z_in,
111  TVectorD &v_in,TVectorD &gamma_in,TVectorD &w_in,TVectorD &phi_in,
112  TVectorD &t_in,TVectorD &lambda_in,TVectorD &u_in,TVectorD &pi_in,
113  TVectorD &ixlow_in,TVectorD &ixupp_in,TVectorD &iclow_in,TVectorD &icupp_in);
114 
115  // constructor that creates variables objects of specified dimensions.
116  TQpVar(Int_t nx,Int_t my,Int_t mz,
117  TVectorD &ixlow,TVectorD &ixupp,TVectorD &iclow,TVectorD &icupp);
118  TQpVar(const TQpVar &another);
119 
120  virtual ~TQpVar() {}
121 
122  // Indicates what type is the blocking variable in the step length determination. If kt_block,
123  // then the blocking variable is one of the slack variables t for a general lower bound,
124  // and so on. Special value kno_block is for the case in which a step length of 1 can be
125  // taken without hitting the bound.
126 
129 
130  virtual Double_t GetMu (); // compute complementarity gap, obtained by taking the
131  // inner product of the complementary vectors and dividing
132  // by the total number of components
133  // computes mu = (t'lambda +u'pi + v'gamma + w'phi)/
134  // (mclow+mcupp+nxlow+nxupp)
135  virtual Double_t MuStep (TQpVar *step,Double_t alpha);
136  // compute the complementarity gap resulting from a step
137  // of length "alpha" along direction "step"
138  virtual void Saxpy (TQpVar *b,Double_t alpha);
139  // given variables b, compute a <- a + alpha b,
140  // where a are the variables in this class
141 
142  virtual void Negate (); // negate the value of all the variables in this structure
143 
144  virtual Double_t StepBound (TQpVar *b); // calculate the largest alpha in (0,1] such that the
145  // nonnegative variables stay nonnegative in the given
146  // search direction. In the general QP problem formulation
147  // this is the largest value of alpha such that
148  // (t,u,v,w,lambda,pi,phi,gamma) + alpha * (b->t,b->u,
149  // b->v,b->w,b->lambda,b->pi,b->phi,b->gamma) >= 0.
150 
151  virtual Double_t FindBlocking(TQpVar *step,Double_t &primalValue,Double_t &primalStep,Double_t &dualValue,
152  Double_t &dualStep,Int_t &firstOrSecond);
153  // Performs the same function as StepBound, and supplies
154  // additional information about which component of the
155  // nonnegative variables is responsible for restricting
156  // alpha. In terms of the abstract formulation, the
157  // components have the following meanings.
158  //
159  // primalValue: the value of the blocking component of the
160  // primal variables (u,t,v,w).
161  // primalStep: the corresponding value of the blocking
162  // component of the primal step variables (b->u,b->t,
163  // b->v,b->w)
164  // dualValue: the value of the blocking component of the
165  // dual variables (lambda,pi,phi,gamma).
166  // dualStep: the corresponding value of the blocking
167  // component of the dual step variables (b->lambda,b->pi,
168  // b->phi,b->gamma)
169  // firstOrSecond: 1 if the primal step is blocking, 2
170  // if the dual step is block, 0 if no step is blocking.
171 
172  virtual void InteriorPoint(Double_t alpha,Double_t beta);
173  // sets components of (u,t,v,w) to alpha and of
174  // (lambda,pi,phi,gamma) to beta
175  virtual void ShiftBoundVariables
176  (Double_t alpha,Double_t beta);
177  // add alpha to components of (u,t,v,w) and beta to
178  // components of (lambda,pi,phi,gamma)
179  virtual Bool_t IsInteriorPoint(); // check whether this is an interior point. Useful as a
180  // sanity check.
181  virtual Double_t Violation (); // The amount by which the current variables violate the
182  // non-negativity constraints.
183  virtual void Print (Option_t *option="") const;
184  virtual Double_t Norm1 (); // compute the 1-norm of the variables
185  virtual Double_t NormInf (); // compute the inf-norm of the variables
186  virtual Bool_t ValidNonZeroPattern();
187 
188  TQpVar &operator= (const TQpVar &source);
189 
190  ClassDef(TQpVar,1) // Qp Variables class
191 };
192 #endif
static Double_t FindBlockingSub(Int_t n, Double_t *w, Int_t incw, Double_t *wstep, Int_t incwstep, Double_t *u, Int_t incu, Double_t *ustep, Int_t incustep, Double_t maxStep, Double_t &w_elt, Double_t &wstep_elt, Double_t &u_elt, Double_t &ustep_elt, Int_t &first_or_second)
See FindBlocking function.
Definition: TQpVar.cxx:464
TVectorD fCloIndex
Definition: TQpVar.h:75
const char Option_t
Definition: RtypesCore.h:62
TVectorD fCupIndex
Definition: TQpVar.h:74
virtual ~TQpVar()
Definition: TQpVar.h:120
virtual Double_t MuStep(TQpVar *step, Double_t alpha)
Compute the complementarity gap resulting from a step of length "alpha" along direction "step"...
Definition: TQpVar.cxx:210
Int_t fMclo
Definition: TQpVar.h:69
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TVectorD fX
Definition: TQpVar.h:91
Int_t fNxup
Definition: TQpVar.h:66
double beta(double x, double y)
Calculates the beta function.
#define ClassDef(name, id)
Definition: Rtypes.h:320
virtual Double_t Violation()
The amount by which the current variables violate the non-negativity constraints. ...
Definition: TQpVar.cxx:573
virtual Double_t GetMu()
compute complementarity gap, obtained by taking the inner product of the complementary vectors and di...
Definition: TQpVar.cxx:191
static Double_t StepBound(TVectorD &v, TVectorD &dir, Double_t maxStep)
Find the maximum stepsize of v in direction dir before violating the nonnegativity constraints...
Definition: TQpVar.cxx:347
TVectorD fXupIndex
Definition: TQpVar.h:73
TVectorD fW
Definition: TQpVar.h:99
Int_t fMy
Definition: TQpVar.h:64
virtual void ShiftBoundVariables(Double_t alpha, Double_t beta)
Add alpha to components of (u,t,v,w) and beta to components of (lambda,pi,phi,gamma) ...
Definition: TQpVar.cxx:614
virtual Bool_t ValidNonZeroPattern()
Check that the variables conform to the non-zero indices.
Definition: TQpVar.cxx:739
EVarBlock
Definition: TQpVar.h:127
SVector< double, 2 > v
Definition: Dict.h:5
TVectorD fV
Definition: TQpVar.h:96
Int_t fMz
Definition: TQpVar.h:65
virtual void Negate()
Perform a "negate" operation on all data vectors : x = -x.
Definition: TQpVar.cxx:271
virtual Bool_t IsInteriorPoint()
Is the current position an interior point ?
Definition: TQpVar.cxx:374
Int_t fNxlo
Definition: TQpVar.h:67
Int_t fMcup
Definition: TQpVar.h:68
double Double_t
Definition: RtypesCore.h:55
static Double_t FindBlocking(TVectorD &w, TVectorD &wstep, TVectorD &u, TVectorD &ustep, Double_t maxStep, Double_t &w_elt, Double_t &wstep_elt, Double_t &u_elt, Double_t &ustep_elt, int &first_or_second)
See other FindBlocking function.
Definition: TQpVar.cxx:445
virtual void Saxpy(TQpVar *b, Double_t alpha)
Perform a "saxpy" operation on all data vectors : x += alpha*y.
Definition: TQpVar.cxx:231
TVectorD fY
Definition: TQpVar.h:93
TVectorD fLambda
Definition: TQpVar.h:103
TVectorD fU
Definition: TQpVar.h:105
TVectorD fPhi
Definition: TQpVar.h:97
TVectorD fGamma
Definition: TQpVar.h:100
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void InteriorPoint(Double_t alpha, Double_t beta)
Sets components of (u,t,v,w) to alpha and of (lambda,pi,phi,gamma) to beta.
Definition: TQpVar.cxx:533
virtual Double_t Norm1()
Return the sum of the vector-norm1&#39;s.
Definition: TQpVar.cxx:675
TQpVar()
Default constructor.
Definition: TQpVar.cxx:58
Definition: TQpVar.h:59
TVectorD fS
Definition: TQpVar.h:92
TVectorD fPi
Definition: TQpVar.h:106
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
TVectorD fT
Definition: TQpVar.h:102
virtual void Print(Option_t *option="") const
Print class members.
Definition: TQpVar.cxx:638
virtual Double_t NormInf()
Return the sum of the vector-normInf&#39;s.
Definition: TQpVar.cxx:699
TQpVar & operator=(const TQpVar &source)
Assignment operator.
Definition: TQpVar.cxx:771
Int_t fNComplementaryVariables
Definition: TQpVar.h:88
Int_t fNx
Definition: TQpVar.h:63
const Int_t n
Definition: legend1.C:16
TVectorD fXloIndex
Definition: TQpVar.h:72
TVectorD fZ
Definition: TQpVar.h:94