]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG1/TPC/AliPerformanceObject.cxx
reduce memory consumption and add high multiplicity flag
[u/mrichter/AliRoot.git] / PWG1 / TPC / AliPerformanceObject.cxx
CommitLineData
7cc34f08 1/**************************************************************************
36ace53b 2* Copyright(c) 1998-1999, ALICE Experiment at CERF, All rights reserved. *
7cc34f08 3* *
4* Author: The ALICE Off-line Project. *
5* Contributors are mentioned in the code where appropriate. *
6* *
7* Permission to use, copy, modify and distribute this software and its *
8* documentation strictly for non-commercial purposes is hereby granted *
9* without fee, provided that the above copyright notice appears in all *
10* copies and that both the copyright notice and this permission notice *
11* appear in the supporting documentation. The authors make no claims *
12* about the suitability of this software for any purpose. It is *
13* provided "as is" without express or implied warranty. *
14**************************************************************************/
15
16//------------------------------------------------------------------------------
17// Implementation of abstract AliPerformanceObject class. It keeps information from
18// comparison of reconstructed and MC particle tracks.
19//
20// Author: J.Otwinowski 14/04/2008
814d192f 21// Changes by M.Knichel 15/10/2010
7cc34f08 22//------------------------------------------------------------------------------
23
24#include <iostream>
25
26#include "TCanvas.h"
27#include "TH1.h"
28#include "TH2.h"
28bb9d1f 29#include "TH3.h"
7cc34f08 30#include "TAxis.h"
31#include "TPostScript.h"
32#include "TList.h"
33#include "TMath.h"
34
35#include "AliLog.h"
758320f7 36#include "AliESDVertex.h"
7cc34f08 37#include "AliPerformanceObject.h"
38
39using namespace std;
40
41ClassImp(AliPerformanceObject)
42
43//_____________________________________________________________________________
44AliPerformanceObject::AliPerformanceObject():
45 TNamed("AliPerformanceObject","AliPerformanceObject"),
46 fAnalysisMode(-1),
814d192f 47 fRunNumber(-1),
e6a60a90 48 fHptGenerator(kFALSE),
758320f7 49 fTriggerClass(0),
28bb9d1f 50 fUseTrackVertex(kFALSE),
c1b69b58 51 fHighMultiplicity(kFALSE),
36ace53b 52 fUseKinkDaughters(kTRUE),
53 fUseCentralityBin(0)
7cc34f08 54{
55 // constructor
56}
57
58//_____________________________________________________________________________
c11cd0fa 59AliPerformanceObject::AliPerformanceObject(const char* name, const char* title, Int_t run, Bool_t highMult):
7cc34f08 60 TNamed(name,title),
61 fAnalysisMode(-1),
814d192f 62 fRunNumber(run),
e6a60a90 63 fHptGenerator(kFALSE),
758320f7 64 fTriggerClass(0),
28bb9d1f 65 fUseTrackVertex(kFALSE),
c11cd0fa 66 fHighMultiplicity(highMult),
36ace53b 67 fUseKinkDaughters(kTRUE),
68 fUseCentralityBin(0)
7cc34f08 69{
70 // constructor
71}
72
73//_____________________________________________________________________________
74AliPerformanceObject::~AliPerformanceObject(){
75 // destructor
76}
77
78//_____________________________________________________________________________
79void AliPerformanceObject::PrintHisto(Bool_t logz, Char_t * outFileName) {
80 // draw all histograms from the folder
81 // and store them in the output *.ps file
82
83 // use this folder
84 TFolder *folder = this->GetAnalysisFolder();
85 if (!folder) {
86 AliDebug(AliLog::kError, "folder not available");
87 return;
88 }
89
90 TCanvas *can = new TCanvas("can");
91 can->Divide(2,2);
92
93 char fname[256];
94 const char* suffix=".ps";
95
96 if(outFileName) sprintf(fname,"%s",outFileName);
97 else sprintf(fname,"%s%s",folder->GetName(),suffix);
98
99 TPostScript *ps = new TPostScript(fname,112);
100 Printf("Histograms are stored in %s", fname);
101 TIter iter(folder->GetListOfFolders());
102
103 TH1 *obj = 0;
104 Int_t count = 0;
105 Int_t pad_count = 0;
106 while ((obj = (TH1*)iter()) !=0) {
107
c1b69b58 108 TString name(obj->ClassName());
109
110
7cc34f08 111 // 4 figures per page
112 if((count%4) == 0) {
113 pad_count = 0;
114 ps->NewPage();
115 }
116
117 pad_count++;
118 can->cd(pad_count);
119
120 if(obj->TestBit(TH1::kLogX))
121 gPad->SetLogx(1);
122 else
123 gPad->SetLogx(0);
124
125 if (obj->GetYaxis() && obj->GetZaxis()) {
126 if(logz) gPad->SetLogz();
c1b69b58 127 if ( name.CompareTo("TH3D") )
128 obj->Draw("colz");
7cc34f08 129 }
130 else {
131 obj->SetMarkerStyle(24);
132 obj->SetMarkerSize(1.0);
c1b69b58 133 if ( name.CompareTo("TH3D") )
134 obj->Draw();
7cc34f08 135 }
136
137 if ((pad_count%4) == 0) {
138 can->Update();
139 }
140
7cc34f08 141 count++;
142 }
143 ps->Close();
144}
145
146
147//_____________________________________________________________________________
148Double_t * AliPerformanceObject::CreateLogAxis(Int_t nbins, Double_t xmin, Double_t xmax) {
149 // retun pointer to the array with log axis
150 // it is user responsibility to delete the array
151
152 Double_t logxmin = TMath::Log10(xmin);
153 Double_t logxmax = TMath::Log10(xmax);
154 Double_t binwidth = (logxmax-logxmin)/nbins;
155
156 Double_t *xbins = new Double_t[nbins+1];
157
158 xbins[0] = xmin;
159 for (Int_t i=1;i<=nbins;i++) {
160 xbins[i] = xmin + TMath::Power(10,logxmin+i*binwidth);
161 }
162
163return xbins;
164}
28bb9d1f 165
166//_____________________________________________________________________________
167void AliPerformanceObject::InitHighMult() {
168
169 fHighMultiplicity = kTRUE;
170 Init();
171
172}
173
174
175//_____________________________________________________________________________
176void AliPerformanceObject::AddProjection(TObjArray* aFolderObj, TString nameSparse, THnSparse* hSparse, Int_t xDim, TString* selString)
177{
178 TH1 *h1=0;
179 TString name = "h_tpc_" + nameSparse + '_';
180 if (selString) { name += *selString + '_'; }
181 name.ToLower();
182 name += xDim;
183 TString title = hSparse->GetAxis(xDim)->GetTitle();
184 if (selString) { title += " (" + *selString + ")"; }
185 h1 = hSparse->Projection(xDim);
186 h1->SetName(name.Data());
187 h1->GetXaxis()->SetTitle(hSparse->GetAxis(xDim)->GetTitle());
188 h1->SetTitle(title.Data());
189 aFolderObj->Add(h1);
190}
191
192
193//_____________________________________________________________________________
194void AliPerformanceObject::AddProjection(TObjArray* aFolderObj, TString nameSparse, THnSparse *hSparse, Int_t yDim, Int_t xDim, TString* selString)
195{
196 TH2 *h2=0;
197 TString name = "h_tpc_" + nameSparse + '_';
198 if (selString) { name += *selString + '_'; }
199 name.ToLower();
200 name += yDim;
201 name += '_';
202 name += xDim;
203 TString title = hSparse->GetAxis(yDim)->GetTitle();
204 title += " vs ";
205 title += hSparse->GetAxis(xDim)->GetTitle();
206 if (selString) { title += " (" + *selString + ")"; }
207 h2 = hSparse->Projection(yDim,xDim);
208 h2->SetName(name.Data());
209 h2->GetXaxis()->SetTitle(hSparse->GetAxis(xDim)->GetTitle());
210 h2->GetYaxis()->SetTitle(hSparse->GetAxis(yDim)->GetTitle());
211 h2->SetTitle(title.Data());
212 aFolderObj->Add(h2);
213}
214
215
216//_____________________________________________________________________________
217void AliPerformanceObject::AddProjection(TObjArray* aFolderObj, TString nameSparse, THnSparse *hSparse, Int_t xDim, Int_t yDim, Int_t zDim, TString* selString)
218{
219 TH3 *h3=0;
220 TString name = "h_tpc_" + nameSparse + '_';
221 if (selString) { name += *selString + '_'; }
222 name.ToLower();
223 name += xDim;
224 name += '_';
225 name += yDim;
226 name += '_';
227 name += zDim;
228 TString title = hSparse->GetAxis(xDim)->GetTitle();
229 title += " vs ";
230 title += hSparse->GetAxis(yDim)->GetTitle();
231 title += " vs ";
232 title += hSparse->GetAxis(zDim)->GetTitle();
233 if (selString) { title += " (" + *selString + ")"; }
234 h3 = hSparse->Projection(xDim,yDim,zDim);
235 h3->SetName(name.Data());
236 h3->GetXaxis()->SetTitle(hSparse->GetAxis(xDim)->GetTitle());
237 h3->GetYaxis()->SetTitle(hSparse->GetAxis(yDim)->GetTitle());
238 h3->GetZaxis()->SetTitle(hSparse->GetAxis(zDim)->GetTitle());
239 h3->SetTitle(title.Data());
240 aFolderObj->Add(h3);
241}