ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
portfolio.C
Go to the documentation of this file.
1 #include "Riostream.h"
2 #include "TCanvas.h"
3 #include "TFile.h"
4 #include "TMath.h"
5 #include "TTree.h"
6 #include "TArrayF.h"
7 #include "TH1.h"
8 #include "TF1.h"
9 #include "TLegend.h"
10 #include "TSystem.h"
11 
12 #include "TMatrixD.h"
13 #include "TMatrixDSym.h"
14 #include "TVectorD.h"
15 #include "TQpProbDens.h"
16 #include "TGondzioSolver.h"
17 
18 //+ This macro shows in detail the use of the quadratic programming package quadp .
19 // Running this macro :
20 // .x portfolio.C+
21 // or gSystem->Load("libQuadp"); .L portFolio.C+; portfolio()
22 //
23 // Let's first review what we exactly mean by "quadratic programming" :
24 //
25 // We want to minimize the following objective function :
26 //
27 // c^T x + ( 1/2 ) x^T Q x wrt. the vector x
28 //
29 // c is a vector and Q a symmetric positive definite matrix
30 //
31 // You might wonder what is so special about this objective which is quadratic in
32 // the unknowns, that can not be done by Minuit/Fumili . Well, we have in addition
33 // the following boundary conditions on x:
34 //
35 // A x = b
36 // clo <= C x <= cup
37 // xlo <= x <= xup , where A and C are arbitray matrices and the rest are vectors
38 //
39 // Not all these constraints have to be defined . Our example will only use xlo, A and b
40 // Still, this could be handled by a general non-linear minimizer like Minuit by introducing
41 // so-called "slack" variables . However, quadp is tailored to objective functions not more
42 // complex than being quadratic . This allows usage of solving techniques which are even
43 // stable for problems involving for instance 500 variables, 100 inequality conditions
44 // and 50 equality conditions .
45 //
46 // Enough said about quadratic programming, let's return to our example .
47 // Suppose, after a long day of doing physics, you have a look at your investments and
48 // realize that an early retirement is not possible, given the returns of your stocks .
49 // So what now ? ROOT to the rescue ......
50 //
51 // In 1990 Harry Markowitz was awarded the Nobel prize for economics : " his work provided new tools
52 // for weighing the risks and rewards of different investments and for valuing corporate stocks and bonds" .
53 // In plain English, he developed the tools to balance greed and fear, we want the maximum return
54 // with the minimum amount of risk. Our stock portfolio should be at the "Efficient Frontier",
55 // see http://www.riskglossary.com/articles/efficient_frontier.htm .
56 // To quantify better the risk we are willing to take, we define a utility function U(x) . It describes
57 // as a function of our total assets x, our "satisfaction" . A common choice is 1-exp(-k*x) (the reason for
58 // the exponent will be clear later) . The parameter k is the risk-aversion factor . For small values of k
59 // the satisfaction is small for small values of x; by increasing x the satisfaction can still be increased
60 // significantly . For large values of k, U(x) increases rapidly to 1, there is no increase in satisfaction
61 // for additional dollars earned .
62 //
63 // In summary : small k ==> risk-loving investor
64 // large k ==> risk-averse investor
65 //
66 // Suppose we have for nrStocks the historical daily returns r = closing_price(n) - closing_price(n-1) .
67 // Define a vector x of length of nrStocks, which contains the fraction of our money invested in
68 // each stock . We can calculate the average daily return z of our portfolio and its variance using
69 // the portfolio covariance Covar :
70 //
71 // z = r^T x and var = x^T Covar x
72 //
73 // Assuming that the daily returns have a Normal distribution, N(x), so will z with mean r^T x
74 // and variance x^T Covar x
75 //
76 // The expected value of the utility function is : E(u(x)) = Int (1-exp(-k*x) N(x) dx
77 // = 1-exp(-k (r^T x - 0.5 k x^T Covar x) )
78 //
79 // Its value is maximized by maximizing r^T x -0.5 k x^T Covar x
80 // under the condition sum (x_i) = 1, meaning we want all our money invested and
81 // x_i >= 0 , we can not "short" a stock
82 //
83 // For 10 stocks we got the historical daily data for Sep-2000 to Jun-2004:
84 //
85 // GE : General Electric Co
86 // SUNW : Sun Microsystems Inc
87 // QCOM : Qualcomm Inc
88 // BRCM : Broadcom Corp
89 // TYC : Tyco International Ltd
90 // IBM : International Business Machines Corp
91 // AMAT : Applied Materials Inc
92 // C : Citigroup Inc
93 // PFE : Pfizer Inc
94 // HD : Home Depot Inc
95 //
96 // We calculate the optimal portfolio for 2.0 and 10.0 .
97 //
98 // Food for thought :
99 // - We assumed that the stock returns have a Normal distribution . Check this assumption by
100 // histogramming the stock returns !
101 // - We used for the expected return in the objective function, the flat average over a time
102 // period . Investment firms will put significant resources in improving the return predicton .
103 // - If you want to trade significant number of shares, several other considerations have
104 // to be taken into account :
105 // + If you are going to buy, you will drive the price up (so-called "slippage") .
106 // This can be taken into account by adding terms to the objective
107 // (Google for "slippage optimization")
108 // + FTC regulations might have to be added to the inequality constraints
109 // - Investment firms do not want to be exposed to the "market" as defined by a broad
110 // index like the S&P and "hedge" this exposure away . A perfect hedge this can be added
111 // as an equality constrain, otherwise add an inequality constrain .
112 //
113 //Author: Eddy Offermann
114 
115 const Int_t nrStocks = 10;
116 static const Char_t *stocks[] =
117  {"GE","SUNW","QCOM","BRCM","TYC","IBM","AMAT","C","PFE","HD"};
118 
119 class TStockDaily {
120 public:
121  Int_t fDate;
122  Int_t fOpen; // 100*open_price
123  Int_t fHigh; // 100*high_price
124  Int_t fLow; // 100*low_price
125  Int_t fClose; // 100*close_price
126  Int_t fVol;
127  Int_t fCloseAdj; // 100*close_price adjusted for splits and dividend
128 
129  TStockDaily() {
130  fDate = fVol = fOpen = fHigh = fLow = fClose = fCloseAdj = 0;
131  }
132  virtual ~TStockDaily() {}
133 
134  ClassDef(TStockDaily,1)
135 };
136 
137 //---------------------------------------------------------------------------
139  Double_t riskFactor = par[0];
140  return 1-TMath::Exp(-riskFactor*x[0]);
141 }
142 
143 //---------------------------------------------------------------------------
145 {
146  TTree *tDaily = (TTree*)f->Get(name);
147  TStockDaily *data = 0;
148  tDaily->SetBranchAddress("daily",&data);
149  TBranch *b_closeAdj = tDaily->GetBranch("fCloseAdj");
150  TBranch *b_date = tDaily->GetBranch("fDate");
151 
152  //read only the "adjusted close" branch for all entries
153  const Int_t nrEntries = (Int_t)tDaily->GetEntries();
154  TArrayF closeAdj(nrEntries);
155  for (Int_t i = 0; i < nrEntries; i++) {
156  b_date->GetEntry(i);
157  b_closeAdj->GetEntry(i);
158  if (data->fDate >= sDay && data->fDate <= eDay)
159 #ifdef __CINT__
160  closeAdj.AddAt(data->fCloseAdj/100. , i );
161 #else
162  closeAdj[i] = data->fCloseAdj/100.;
163 #endif
164  }
165 
166  TArrayF *r = new TArrayF(nrEntries-1);
167  for (Int_t i = 1; i < nrEntries; i++)
168 // (*r)[i-1] = closeAdj[i]-closeAdj[i-1];
169 #ifdef __CINT__
170  r->AddAt(closeAdj[i]/closeAdj[i-1],1);
171 #else
172  (*r)[i-1] = closeAdj[i]/closeAdj[i-1];
173 #endif
174 
175  return *r;
176 }
177 
178 #ifndef __MAKECINT__
179 //---------------------------------------------------------------------------
181 {
182 // what the quadratic programming package will do:
183 //
184 // minimize c^T x + ( 1/2 ) x^T Q x
185 // subject to A x = b
186 // clo <= C x <= cup
187 // xlo <= x <= xup
188 // what we want :
189 //
190 // maximize c^T x - k ( 1/2 ) x^T Q x
191 // subject to sum_x x_i = 1
192 // 0 <= x_i
193 
194  // We have nrStocks weights to determine,
195  // 1 equality- and 0 inequality- equations (the simple square boundary
196  // condition (xlo <= x <= xup) does not count)
197 
198  const Int_t nrVar = nrStocks;
199  const Int_t nrEqual = 1;
200  const Int_t nrInEqual = 0;
201 
202  // flip the sign of the objective function because we want to maximize
203  TVectorD c = -1.*r;
204  TMatrixDSym Q = riskFactor*Covar;
205 
206  // equality equation
207  TMatrixD A(nrEqual,nrVar); A = 1;
208  TVectorD b(nrEqual); b = 1;
209 
210  // inequality equation
211  //
212  // - although not applicable in the current situatio since nrInEqual = 0, one
213  // has to specify not only clo and cup but also an index vector iclo and icup,
214  // whose values are either 0 or 1 . If iclo[j] = 1, the lower boundary condition
215  // is active on x[j], etc. ...
216 
217  TMatrixD C (nrInEqual,nrVar);
218  TVectorD clo (nrInEqual);
219  TVectorD cup (nrInEqual);
220  TVectorD iclo(nrInEqual);
221  TVectorD icup(nrInEqual);
222 
223  // simple square boundary condition : 0 <= x_i, so only xlo is relevant .
224  // Like for clo and cup above, we have to define an index vector ixlo and ixup .
225  // Since each variable has the lower boundary, we can set the whole vector
226  // ixlo = 1
227 
228  TVectorD xlo (nrVar); xlo = 0;
229  TVectorD xup (nrVar); xup = 0;
230  TVectorD ixlo(nrVar); ixlo = 1;
231  TVectorD ixup(nrVar); ixup = 0;
232 
233  // setup the quadratic programming problem . Since a small number of variables are
234  // involved and "Q" has everywhere entries, we chose the dense version "TQpProbDens" .
235  // In case of a sparse formulation, simply replace all "Dens" by "Sparse" below and
236  // use TMatrixDSparse instead of TMatrixDSym and TMatrixD
237 
238  TQpProbDens *qp = new TQpProbDens(nrVar,nrEqual,nrInEqual);
239 
240  // stuff all the matrices/vectors defined above in the proper places
241 
242  TQpDataDens *prob = (TQpDataDens *)qp->MakeData(c,Q,xlo,ixlo,xup,ixup,A,b,C,clo,iclo,cup,icup);
243 
244  // setup the nrStock variables, vars->fX will contain the final solution
245 
246  TQpVar *vars = qp->MakeVariables(prob);
247  TQpResidual *resid = qp->MakeResiduals(prob);
248 
249  // Now we have to choose the method of solving, either TGondzioSolver or TMehrotraSolver
250  // The Gondzio method is more sophisticated and therefore numerically more involved
251  // If one want the Mehrotra method, simply replace "Gondzio" by "Mehrotra" .
252 
253  TGondzioSolver *s = new TGondzioSolver(qp,prob);
254  const Int_t status = s->Solve(prob,vars,resid);
255 
256  const TVectorD weight = vars->fX;
257 
258  delete qp; delete prob; delete vars; delete resid; delete s;
259  if (status != 0) {
260  cout << "Could not solve this problem." <<endl;
261  return TVectorD(nrStocks);
262  }
263 
264  return weight;
265 }
266 #endif
267 
268 //---------------------------------------------------------------------------
269 void portfolio()
270 {
271  const Int_t sDay = 20000809;
272  const Int_t eDay = 20040602;
273 
274  const char *fname = "stock.root";
275  TFile *f = 0;
276  if (!gSystem->AccessPathName(fname)) {
277  f = TFile::Open(fname);
278  } else {
279  printf("accessing %s file from http://root.cern.ch/files\n",fname);
280  f = TFile::Open(Form("http://root.cern.ch/files/%s",fname));
281  }
282  if (!f) return;
283 
284  TArrayF *data = new TArrayF[nrStocks];
285  for (Int_t i = 0; i < nrStocks; i++) {
286  const TString symbol = stocks[i];
287  data[i] = StockReturn(f,symbol,sDay,eDay);
288  }
289 
290  const Int_t nrData = data[0].GetSize();
291 
292  TVectorD r(nrStocks);
293  for (Int_t i = 0; i < nrStocks; i++)
294  r[i] = data[i].GetSum()/nrData;
295 
296  TMatrixDSym Covar(nrStocks);
297  for (Int_t i = 0; i < nrStocks; i++) {
298  for (Int_t j = 0; j <= i; j++) {
299  Double_t sum = 0.;
300  for (Int_t k = 0; k < nrData; k++)
301  sum += (data[i][k]-r[i])*(data[j][k]-r[j]);
302  Covar(i,j) = Covar(j,i) = sum/nrData;
303  }
304  }
305 
306  const TVectorD weight1 = OptimalInvest(2.0,r,Covar);
307  const TVectorD weight2 = OptimalInvest(10.,r,Covar);
308 
309  cout << "stock daily daily w1 w2" <<endl;
310  cout << "symb return sdv " <<endl;
311  for (Int_t i = 0; i < nrStocks; i++)
312  printf("%s\t: %.3f %.3f %.3f %.3f\n",stocks[i],r[i],TMath::Sqrt(Covar[i][i]),weight1[i],weight2[i]);
313 
314  TCanvas *c1 = new TCanvas("c1","Portfolio Optimizations",10,10,800,900);
315  c1->Divide(1,2);
316 
317  // utility function / risk profile
318 
319  c1->cd(1);
320  gPad->SetGridx();
321  gPad->SetGridy();
322 
323  TF1 *f1 = new TF1("f1",RiskProfile,0,2.5,1);
324  f1->SetParameter(0,2.0);
325  f1->SetLineColor(49);
326  f1->Draw("AC");
327  f1->GetHistogram()->SetXTitle("dollar");
328  f1->GetHistogram()->SetYTitle("utility");
329  f1->GetHistogram()->SetMinimum(0.0);
330  f1->GetHistogram()->SetMaximum(1.0);
331  TF1 *f2 = new TF1("f2",RiskProfile,0,2.5,1);
332  f2->SetParameter(0,10.);
333  f2->SetLineColor(50);
334  f2->Draw("CSAME");
335 
336  TLegend *legend1 = new TLegend(0.50,0.65,0.70,0.82);
337  legend1->AddEntry(f1,"1-exp(-2.0*x)","l");
338  legend1->AddEntry(f2,"1-exp(-10.*x)","l");
339  legend1->Draw();
340 
341  // vertical bar chart of portfolio distribution
342 
343  c1->cd(2);
344  TH1F *h1 = new TH1F("h1","Portfolio Distribution",nrStocks,0,0);
345  TH1F *h2 = new TH1F("h2","Portfolio Distribution",nrStocks,0,0);
346  h1->SetStats(0);
347  h1->SetFillColor(49);
348  h2->SetFillColor(50);
349  h1->SetBarWidth(0.45);
350  h1->SetBarOffset(0.1);
351  h2->SetBarWidth(0.4);
352  h2->SetBarOffset(0.55);
353  for (Int_t i = 0; i < nrStocks; i++) {
354  h1->Fill(stocks[i],weight1[i]);
355  h2->Fill(stocks[i],weight2[i]);
356  }
357 
358  h1->Draw("BAR2 HIST");
359  h2->Draw("BAR2SAME HIST");
360 
361  TLegend *legend2 = new TLegend(0.50,0.65,0.70,0.82);
362  legend2->AddEntry(h1,"high risk","f");
363  legend2->AddEntry(h2,"low risk","f");
364  legend2->Draw();
365 }
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1213
double par[1]
Definition: unuranDistr.cxx:38
virtual void SetBarOffset(Float_t offset=0.25)
Definition: TH1.h:356
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3159
virtual void SetMaximum(Double_t maximum=-1111)
Definition: TH1.h:394
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
void AddAt(Float_t c, Int_t i)
Add float c at position i. Check for out of bounds.
Definition: TArrayF.cxx:93
return c
TCanvas * c1
Definition: legend1.C:2
tuple f2
Definition: surfaces.py:24
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format...
Definition: TFile.h:45
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
virtual TObject * Get(const char *namecycle)
Return pointer to object identified by namecycle.
virtual void SetMinimum(Double_t minimum=-1111)
Definition: TH1.h:395
Basic string class.
Definition: TString.h:137
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
void portfolio()
Definition: portfolio.C:269
Array of floats (32 bits per element).
Definition: TArrayF.h:29
virtual TQpResidual * MakeResiduals(const TQpDataBase *data)
Setup the residuals.
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
virtual void SetYTitle(const char *title)
Definition: TH1.h:409
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1052
TVectorD fX
Definition: TQpVar.h:97
static double A[]
TFile * f
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to vusualize the function.
Definition: TF1.cxx:1274
virtual void SetBarWidth(Float_t width=0.5)
Definition: TH1.h:357
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=1, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3851
Double_t x[n]
Definition: legend1.C:17
#define ClassDef(name, id)
Definition: Rtypes.h:254
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=0)
Change branch address, dealing with clone trees properly.
Definition: TTree.cxx:7510
TH2D * h2
Definition: fit2dHist.C:45
TH1F * h1
Definition: legend1.C:5
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:24
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:4803
static const Char_t * stocks[]
Definition: portfolio.C:116
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual Int_t Solve(TQpDataBase *prob, TQpVar *iterate, TQpResidual *resid)
Solve the quadratic programming problem as formulated through prob, store the final solution in itera...
static double C[]
ROOT::R::TRInterface & r
Definition: Object.C:4
R__EXTERN TSystem * gSystem
Definition: TSystem.h:545
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
char * Form(const char *fmt,...)
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all leaves of entry and return total number of bytes read.
Definition: TBranch.cxx:1199
Double_t RiskProfile(Double_t *x, Double_t *par)
Definition: portfolio.C:138
Int_t GetSize() const
Definition: TArray.h:49
TVectorD OptimalInvest(Double_t riskFactor, TVectorD r, TMatrixDSym Covar)
Definition: portfolio.C:180
The Canvas class.
Definition: TCanvas.h:48
Double_t Exp(Double_t x)
Definition: TMath.h:495
double Double_t
Definition: RtypesCore.h:55
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
const Int_t nrStocks
Definition: portfolio.C:115
TArrayF & StockReturn(TFile *f, const TString &name, Int_t sDay, Int_t eDay)
Definition: portfolio.C:144
#define name(a, b)
Definition: linkTestLib0.cpp:5
char Char_t
Definition: RtypesCore.h:29
virtual void SetXTitle(const char *title)
Definition: TH1.h:408
Definition: TQpVar.h:65
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1073
1-Dim function class
Definition: TF1.h:149
TF1 * f1
Definition: legend1.C:11
#define gPad
Definition: TVirtualPad.h:288
virtual Long64_t GetEntries() const
Definition: TTree.h:386
A TTree object has a header with a name and a title.
Definition: TTree.h:98
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:424
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual TQpDataBase * MakeData(Double_t *c, Double_t *Q, Double_t *xlo, Bool_t *ixlo, Double_t *xup, Bool_t *ixup, Double_t *A, Double_t *bA, Double_t *C, Double_t *clo, Bool_t *iclo, Double_t *cup, Bool_t *icup)
Setup the data.
Definition: TQpProbDens.cxx:80
A TTree is a list of TBranches.
Definition: TBranch.h:58
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320
static double Q[]
virtual TQpVar * MakeVariables(const TQpDataBase *data)
Setup the variables.