Logo ROOT  
Reference Guide
TQpLinSolverBase.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_TQpLinSolverBase
44#define ROOT_TQpLinSolverBase
45
46#include "TError.h"
47
48#include "TQpVar.h"
49#include "TQpDataBase.h"
50#include "TQpResidual.h"
51#include "TQpProbBase.h"
52
53#include "TMatrixD.h"
54
55///////////////////////////////////////////////////////////////////////////
56// //
57// Implements the main solver for linear systems that arise in //
58// primal-dual interior-point methods for QP . This class contains //
59// definitions of methods and data common to the sparse and dense //
60// special cases of the general formulation. The derived classes contain //
61// the aspects that are specific to the sparse and dense forms. //
62// //
63///////////////////////////////////////////////////////////////////////////
64
65class TQpProbBase;
67{
68
69protected:
70
71 TVectorD fNomegaInv; // stores a critical diagonal matrix as a vector
72 TVectorD fRhs; // right-hand side of the system
73
74 Int_t fNx; // dimensions of the vectors in the general QP formulation
77
78 TVectorD fDd; // temporary storage vectors
80
81 TVectorD fXupIndex; // index matrices for the upper and lower bounds on x and Cx
85
86 Int_t fNxup; // dimensions of the upper and lower bound vectors
90
92
93public:
96 TQpLinSolverBase(const TQpLinSolverBase &another);
97
98 virtual ~TQpLinSolverBase() {}
99
100 virtual void Factor (TQpDataBase *prob,TQpVar *vars);
101 // sets up the matrix for the main linear system in
102 // "augmented system" form. The actual factorization is
103 // performed by a routine specific to either the sparse
104 // or dense case
105 virtual void Solve (TQpDataBase *prob,TQpVar *vars,TQpResidual *resids,TQpVar *step);
106 // solves the system for a given set of residuals.
107 // Assembles the right-hand side appropriate to the
108 // matrix factored in factor, solves the system using
109 // the factorization produced there, partitions the
110 // solution vector into step components, then recovers
111 // the step components eliminated during the block
112 // elimination that produced the augmented system form
113 virtual void JoinRHS (TVectorD &rhs, TVectorD &rhs1,TVectorD &rhs2,TVectorD &rhs3);
114 // assembles a single vector object from three given vectors
115 // rhs (output) final joined vector
116 // rhs1 (input) first part of rhs
117 // rhs2 (input) middle part of rhs
118 // rhs3 (input) last part of rhs
119 virtual void SeparateVars (TVectorD &vars1,TVectorD &vars2,TVectorD &vars3,TVectorD &vars);
120 // extracts three component vectors from a given aggregated
121 // vector.
122 // vars (input) aggregated vector
123 // vars1 (output) first part of vars
124 // vars2 (output) middle part of vars
125 // vars3 (output) last part of vars
126
127 virtual void SolveXYZS (TVectorD &stepx,TVectorD &stepy,TVectorD &stepz,TVectorD &steps,
128 TVectorD &ztemp,TQpDataBase *data);
129 // assemble right-hand side of augmented system and call
130 // SolveCompressed to solve it
131
132 virtual void SolveCompressed (TVectorD &rhs) = 0;
133 // perform the actual solve using the factors produced in
134 // factor.
135 // rhs on input contains the aggregated right-hand side of
136 // the augmented system; on output contains the solution in
137 // aggregated form
138
139 virtual void PutXDiagonal (TVectorD &xdiag) = 0;
140 // places the diagonal resulting from the bounds on x into
141 // the augmented system matrix
142 virtual void PutZDiagonal (TVectorD& zdiag) = 0;
143 // places the diagonal resulting from the bounds on Cx into
144 // the augmented system matrix
145 virtual void ComputeDiagonals(TVectorD &dd,TVectorD &omega,TVectorD &t, TVectorD &lambda,
147 TVectorD &w, TVectorD &phi);
148 // computes the diagonal matrices in the augmented system
149 // from the current set of variables
150
152
153 ClassDef(TQpLinSolverBase,1) // Qp linear solver base class
154};
155#endif
#define ClassDef(name, id)
Definition: Rtypes.h:322
Mother of all ROOT objects.
Definition: TObject.h:37
Data for the general QP formulation.
Definition: TQpDataBase.h:61
Implementation of main solver for linear systems.
virtual ~TQpLinSolverBase()
virtual void JoinRHS(TVectorD &rhs, TVectorD &rhs1, TVectorD &rhs2, TVectorD &rhs3)
Assembles a single vector object from three given vectors .
virtual void Solve(TQpDataBase *prob, TQpVar *vars, TQpResidual *resids, TQpVar *step)
Solves the system for a given set of residuals.
TQpLinSolverBase & operator=(const TQpLinSolverBase &source)
Assignment opeartor.
virtual void SolveCompressed(TVectorD &rhs)=0
virtual void Factor(TQpDataBase *prob, TQpVar *vars)
Sets up the matrix for the main linear system in "augmented system" form.
virtual void ComputeDiagonals(TVectorD &dd, TVectorD &omega, TVectorD &t, TVectorD &lambda, TVectorD &u, TVectorD &pi, TVectorD &v, TVectorD &gamma, TVectorD &w, TVectorD &phi)
Computes the diagonal matrices in the augmented system from the current set of variables.
virtual void PutXDiagonal(TVectorD &xdiag)=0
TQpLinSolverBase()
Default constructor.
TQpProbBase * fFactory
virtual void SeparateVars(TVectorD &vars1, TVectorD &vars2, TVectorD &vars3, TVectorD &vars)
Extracts three component vectors from a given aggregated vector.
virtual void SolveXYZS(TVectorD &stepx, TVectorD &stepy, TVectorD &stepz, TVectorD &steps, TVectorD &ztemp, TQpDataBase *data)
Assemble right-hand side of augmented system and call SolveCompressed to solve it.
virtual void PutZDiagonal(TVectorD &zdiag)=0
default general problem formulation:
Definition: TQpProbBase.h:89
The Residuals class calculates and stores the quantities that appear on the right-hand side of the li...
Definition: TQpResidual.h:62
Class containing the variables for the general QP formulation.
Definition: TQpVar.h:60
double gamma(double x)
static constexpr double pi