00001 #include "Riostream.h"
00002 #include "TCanvas.h"
00003 #include "TFile.h"
00004 #include "TMath.h"
00005 #include "TTree.h"
00006 #include "TArrayF.h"
00007 #include "TH1.h"
00008 #include "TF1.h"
00009 #include "TLegend.h"
00010 #include "TSystem.h"
00011
00012 #include "TMatrixD.h"
00013 #include "TMatrixDSym.h"
00014 #include "TVectorD.h"
00015 #include "TQpProbDens.h"
00016 #include "TGondzioSolver.h"
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115 const Int_t nrStocks = 10;
00116 static const Char_t *stocks[] =
00117 {"GE","SUNW","QCOM","BRCM","TYC","IBM","AMAT","C","PFE","HD"};
00118
00119 class TStockDaily {
00120 public:
00121 Int_t fDate;
00122 Int_t fOpen;
00123 Int_t fHigh;
00124 Int_t fLow;
00125 Int_t fClose;
00126 Int_t fVol;
00127 Int_t fCloseAdj;
00128
00129 TStockDaily() {
00130 fDate = fVol = fOpen = fHigh = fLow = fClose = fCloseAdj = 0;
00131 }
00132 virtual ~TStockDaily() {}
00133
00134 ClassDef(TStockDaily,1)
00135 };
00136
00137
00138 Double_t RiskProfile(Double_t *x, Double_t *par) {
00139 Double_t riskFactor = par[0];
00140 return 1-TMath::Exp(-riskFactor*x[0]);
00141 }
00142
00143
00144 TArrayF &StockReturn(TFile *f,const TString &name,Int_t sDay,Int_t eDay)
00145 {
00146 TTree *tDaily = (TTree*)f->Get(name);
00147 TStockDaily *data = 0;
00148 tDaily->SetBranchAddress("daily",&data);
00149 TBranch *b_closeAdj = tDaily->GetBranch("fCloseAdj");
00150 TBranch *b_date = tDaily->GetBranch("fDate");
00151
00152
00153 const Int_t nrEntries = (Int_t)tDaily->GetEntries();
00154 TArrayF closeAdj(nrEntries);
00155 for (Int_t i = 0; i < nrEntries; i++) {
00156 b_date->GetEntry(i);
00157 b_closeAdj->GetEntry(i);
00158 if (data->fDate >= sDay && data->fDate <= eDay)
00159 #ifdef __CINT__
00160 closeAdj.AddAt(data->fCloseAdj/100. , i );
00161 #else
00162 closeAdj[i] = data->fCloseAdj/100.;
00163 #endif
00164 }
00165
00166 TArrayF *r = new TArrayF(nrEntries-1);
00167 for (Int_t i = 1; i < nrEntries; i++)
00168
00169 #ifdef __CINT__
00170 r->AddAt(closeAdj[i]/closeAdj[i-1],1);
00171 #else
00172 (*r)[i-1] = closeAdj[i]/closeAdj[i-1];
00173 #endif
00174
00175 return *r;
00176 }
00177
00178 #ifndef __MAKECINT__
00179
00180 TVectorD OptimalInvest(Double_t riskFactor,TVectorD r,TMatrixDSym Covar)
00181 {
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198 const Int_t nrVar = nrStocks;
00199 const Int_t nrEqual = 1;
00200 const Int_t nrInEqual = 0;
00201
00202
00203 TVectorD c = -1.*r;
00204 TMatrixDSym Q = riskFactor*Covar;
00205
00206
00207 TMatrixD A(nrEqual,nrVar); A = 1;
00208 TVectorD b(nrEqual); b = 1;
00209
00210
00211
00212
00213
00214
00215
00216
00217 TMatrixD C (nrInEqual,nrVar);
00218 TVectorD clo (nrInEqual);
00219 TVectorD cup (nrInEqual);
00220 TVectorD iclo(nrInEqual);
00221 TVectorD icup(nrInEqual);
00222
00223
00224
00225
00226
00227
00228 TVectorD xlo (nrVar); xlo = 0;
00229 TVectorD xup (nrVar); xup = 0;
00230 TVectorD ixlo(nrVar); ixlo = 1;
00231 TVectorD ixup(nrVar); ixup = 0;
00232
00233
00234
00235
00236
00237
00238 TQpProbDens *qp = new TQpProbDens(nrVar,nrEqual,nrInEqual);
00239
00240
00241
00242 TQpDataDens *prob = (TQpDataDens *)qp->MakeData(c,Q,xlo,ixlo,xup,ixup,A,b,C,clo,iclo,cup,icup);
00243
00244
00245
00246 TQpVar *vars = qp->MakeVariables(prob);
00247 TQpResidual *resid = qp->MakeResiduals(prob);
00248
00249
00250
00251
00252
00253 TGondzioSolver *s = new TGondzioSolver(qp,prob);
00254 const Int_t status = s->Solve(prob,vars,resid);
00255
00256 const TVectorD weight = vars->fX;
00257
00258 delete qp; delete prob; delete vars; delete resid; delete s;
00259 if (status != 0) {
00260 cout << "Could not solve this problem." <<endl;
00261 return TVectorD(nrStocks);
00262 }
00263
00264 return weight;
00265 }
00266 #endif
00267
00268
00269 void portfolio()
00270 {
00271 const Int_t sDay = 20000809;
00272 const Int_t eDay = 20040602;
00273
00274 const char *fname = "stock.root";
00275 TFile *f = 0;
00276 if (!gSystem->AccessPathName(fname)) {
00277 f = TFile::Open(fname);
00278 } else {
00279 printf("accessing %s file from http://root.cern.ch/files\n",fname);
00280 f = TFile::Open(Form("http://root.cern.ch/files/%s",fname));
00281 }
00282 if (!f) return;
00283
00284 TArrayF *data = new TArrayF[nrStocks];
00285 for (Int_t i = 0; i < nrStocks; i++) {
00286 const TString symbol = stocks[i];
00287 data[i] = StockReturn(f,symbol,sDay,eDay);
00288 }
00289
00290 const Int_t nrData = data[0].GetSize();
00291
00292 TVectorD r(nrStocks);
00293 for (Int_t i = 0; i < nrStocks; i++)
00294 r[i] = data[i].GetSum()/nrData;
00295
00296 TMatrixDSym Covar(nrStocks);
00297 for (Int_t i = 0; i < nrStocks; i++) {
00298 for (Int_t j = 0; j <= i; j++) {
00299 Double_t sum = 0.;
00300 for (Int_t k = 0; k < nrData; k++)
00301 sum += (data[i][k]-r[i])*(data[j][k]-r[j]);
00302 Covar(i,j) = Covar(j,i) = sum/nrData;
00303 }
00304 }
00305
00306 const TVectorD weight1 = OptimalInvest(2.0,r,Covar);
00307 const TVectorD weight2 = OptimalInvest(10.,r,Covar);
00308
00309 cout << "stock daily daily w1 w2" <<endl;
00310 cout << "symb return sdv " <<endl;
00311 for (Int_t i = 0; i < nrStocks; i++)
00312 printf("%s\t: %.3f %.3f %.3f %.3f\n",stocks[i],r[i],TMath::Sqrt(Covar[i][i]),weight1[i],weight2[i]);
00313
00314 TCanvas *c1 = new TCanvas("c1","Portfolio Optimizations",10,10,800,900);
00315 c1->Divide(1,2);
00316
00317
00318
00319 c1->cd(1);
00320 gPad->SetGridx();
00321 gPad->SetGridy();
00322
00323 TF1 *f1 = new TF1("f1",RiskProfile,0,2.5,1);
00324 f1->SetParameter(0,2.0);
00325 f1->SetLineColor(49);
00326 f1->Draw("AC");
00327 f1->GetHistogram()->SetXTitle("dollar");
00328 f1->GetHistogram()->SetYTitle("utility");
00329 f1->GetHistogram()->SetMinimum(0.0);
00330 f1->GetHistogram()->SetMaximum(1.0);
00331 TF1 *f2 = new TF1("f2",RiskProfile,0,2.5,1);
00332 f2->SetParameter(0,10.);
00333 f2->SetLineColor(50);
00334 f2->Draw("CSAME");
00335
00336 TLegend *legend1 = new TLegend(0.50,0.65,0.70,0.82);
00337 legend1->AddEntry(f1,"1-exp(-2.0*x)","l");
00338 legend1->AddEntry(f2,"1-exp(-10.*x)","l");
00339 legend1->Draw();
00340
00341
00342
00343 c1->cd(2);
00344 TH1F *h1 = new TH1F("h1","Portfolio Distribution",nrStocks,0,0);
00345 TH1F *h2 = new TH1F("h2","Portfolio Distribution",nrStocks,0,0);
00346 h1->SetStats(0);
00347 h1->SetFillColor(49);
00348 h2->SetFillColor(50);
00349 h1->SetBarWidth(0.45);
00350 h1->SetBarOffset(0.1);
00351 h2->SetBarWidth(0.4);
00352 h2->SetBarOffset(0.55);
00353 for (Int_t i = 0; i < nrStocks; i++) {
00354 h1->Fill(stocks[i],weight1[i]);
00355 h2->Fill(stocks[i],weight2[i]);
00356 }
00357
00358 h1->Draw("BAR2");
00359 h2->Draw("BAR2SAME");
00360
00361 TLegend *legend2 = new TLegend(0.50,0.65,0.70,0.82);
00362 legend2->AddEntry(h1,"high risk","f");
00363 legend2->AddEntry(h2,"low risk","f");
00364 legend2->Draw();
00365 }