]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TPC/AliPerformanceObject.cxx
Fixed cmake files (Constantin)
[u/mrichter/AliRoot.git] / PWGPP / TPC / AliPerformanceObject.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERF, All rights reserved. *
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 
21 // Changes by M.Knichel 15/10/2010
22 //------------------------------------------------------------------------------
23
24 #include <iostream>
25
26 #include "TCanvas.h"
27 #include "TH1.h"
28 #include "TH2.h"
29 #include "TH3.h"
30 #include "TAxis.h"
31 #include "TPostScript.h"
32 #include "TList.h"
33 #include "TMath.h"
34
35 #include "AliLog.h" 
36 #include "AliESDVertex.h" 
37 #include "AliPerformanceObject.h" 
38
39 using namespace std;
40
41 ClassImp(AliPerformanceObject)
42
43 //_____________________________________________________________________________
44 AliPerformanceObject::AliPerformanceObject():
45   TNamed("AliPerformanceObject","AliPerformanceObject"),
46   fMergeTHnSparseObj(kFALSE),
47   fAnalysisMode(-1),
48   fRunNumber(-1),
49   fHptGenerator(kFALSE),
50   fTriggerClass(0),
51   fUseTrackVertex(kFALSE),
52   fHighMultiplicity(kFALSE),
53   fUseKinkDaughters(kTRUE),
54   fUseCentralityBin(0),
55   fUseTOFBunchCrossing(kTRUE)
56 {
57   // constructor
58 }
59
60 //_____________________________________________________________________________
61 AliPerformanceObject::AliPerformanceObject(const char* name, const char* title, Int_t run, Bool_t highMult):
62   TNamed(name,title),
63   fMergeTHnSparseObj(kFALSE),
64   fAnalysisMode(-1),
65   fRunNumber(run),
66   fHptGenerator(kFALSE),
67   fTriggerClass(0),
68   fUseTrackVertex(kFALSE),
69   fHighMultiplicity(highMult),
70   fUseKinkDaughters(kTRUE),
71   fUseCentralityBin(0),
72   fUseTOFBunchCrossing(kTRUE)
73 {
74   // constructor
75 }
76
77 //_____________________________________________________________________________
78 AliPerformanceObject::~AliPerformanceObject(){
79   // destructor 
80 }
81
82 //_____________________________________________________________________________
83 void AliPerformanceObject::PrintHisto(Bool_t logz, const Char_t * outFileName) {
84   // draw all histograms from the folder 
85   // and store them in the output *.ps file
86  
87   // use this folder
88   TFolder *folder = this->GetAnalysisFolder();
89   if (!folder) {
90      AliDebug(AliLog::kError, "folder not available");
91      return;
92   } 
93
94   TCanvas *can = new TCanvas("can");
95   can->Divide(2,2);
96
97   char fname[256];
98   const char* suffix=".ps"; 
99
100   if(outFileName) snprintf(fname,256,"%s",outFileName);
101   else snprintf(fname,256,"%s%s",folder->GetName(),suffix);
102   
103   TPostScript *ps = new TPostScript(fname,112);
104   Printf("Histograms are stored in %s", fname); 
105   TIter iter(folder->GetListOfFolders());
106
107   TH1 *obj = 0;
108   Int_t count = 0;
109   Int_t pad_count = 0;
110   while ((obj = (TH1*)iter()) !=0) {
111
112     TString name(obj->ClassName());
113  
114
115     // 4 figures per page
116     if((count%4) == 0) {
117       pad_count = 0;
118       ps->NewPage();
119     }
120
121     pad_count++; 
122     can->cd(pad_count);
123     
124     if(obj->TestBit(TH1::kLogX)) 
125        gPad->SetLogx(1);
126     else    
127        gPad->SetLogx(0);
128
129     if (obj->GetYaxis() && obj->GetZaxis()) {
130       if(logz) gPad->SetLogz();
131       if ( name.CompareTo("TH3D") )
132         obj->Draw("colz");
133     }
134     else { 
135       obj->SetMarkerStyle(24);
136       obj->SetMarkerSize(1.0);
137       if ( name.CompareTo("TH3D") )
138         obj->Draw();
139     }
140
141     if ((pad_count%4) == 0)  { 
142       can->Update();
143     }
144
145   count++;
146   }
147   ps->Close();
148 }
149  
150
151 //_____________________________________________________________________________
152 Double_t * AliPerformanceObject::CreateLogAxis(Int_t nbins, Double_t xmin, Double_t xmax) {
153   // retun pointer to the array with log axis
154   // it is user responsibility to delete the array
155  
156   Double_t logxmin = TMath::Log10(xmin);
157   Double_t logxmax = TMath::Log10(xmax);
158   Double_t binwidth = (logxmax-logxmin)/nbins;
159   
160   Double_t *xbins =  new Double_t[nbins+1];
161
162   xbins[0] = xmin;
163   for (Int_t i=1;i<=nbins;i++) {
164     xbins[i] = xmin + TMath::Power(10,logxmin+i*binwidth);
165   }
166
167 return xbins;
168 }
169
170 //_____________________________________________________________________________
171 void AliPerformanceObject::InitHighMult() {
172
173   fHighMultiplicity = kTRUE;
174   Init();
175
176 }
177
178
179 //_____________________________________________________________________________
180 void AliPerformanceObject::AddProjection(TObjArray* aFolderObj, TString nameSparse, THnSparse* hSparse, Int_t xDim, TString* selString) 
181 {
182   TH1 *h1=0;  
183   TString name = "h_tpc_" + nameSparse + '_';
184   if (selString) { name += *selString + '_'; }
185   name.ToLower();
186   name += xDim;
187   TString title = hSparse->GetAxis(xDim)->GetTitle();  
188   if (selString) { title += " (" + *selString + ")"; }
189   h1 = hSparse->Projection(xDim);
190   h1->SetName(name.Data());
191   h1->GetXaxis()->SetTitle(hSparse->GetAxis(xDim)->GetTitle());
192   h1->SetTitle(title.Data());  
193   aFolderObj->Add(h1);
194 }
195
196
197 //_____________________________________________________________________________
198 void AliPerformanceObject::AddProjection(TObjArray* aFolderObj, TString nameSparse, THnSparse *hSparse, Int_t yDim, Int_t xDim, TString* selString)
199 {
200   TH2 *h2=0;  
201   TString name = "h_tpc_" + nameSparse + '_';
202   if (selString) { name += *selString + '_'; }
203   name.ToLower();
204   name += yDim;
205   name += '_';
206   name += xDim;
207   TString title = hSparse->GetAxis(yDim)->GetTitle();
208   title += " vs ";
209   title += hSparse->GetAxis(xDim)->GetTitle();
210   if (selString) { title += " (" + *selString + ")"; }  
211   h2 = hSparse->Projection(yDim,xDim);
212   h2->SetName(name.Data());
213   h2->GetXaxis()->SetTitle(hSparse->GetAxis(xDim)->GetTitle());
214   h2->GetYaxis()->SetTitle(hSparse->GetAxis(yDim)->GetTitle());
215   h2->SetTitle(title.Data());  
216   aFolderObj->Add(h2);
217 }
218
219
220 //_____________________________________________________________________________
221 void AliPerformanceObject::AddProjection(TObjArray* aFolderObj, TString nameSparse, THnSparse *hSparse, Int_t xDim, Int_t yDim, Int_t zDim, TString* selString)
222 {
223   TH3 *h3=0;  
224   TString name = "h_tpc_" + nameSparse + '_';
225   if (selString) { name += *selString + '_'; }
226   name.ToLower();
227   name += xDim;
228   name += '_';
229   name += yDim;
230   name += '_';
231   name += zDim;
232   TString title = hSparse->GetAxis(xDim)->GetTitle();
233   title += " vs ";
234   title += hSparse->GetAxis(yDim)->GetTitle();
235   title += " vs ";
236   title += hSparse->GetAxis(zDim)->GetTitle();
237   if (selString) { title += " (" + *selString + ")"; }
238   h3 = hSparse->Projection(xDim,yDim,zDim);
239   h3->SetName(name.Data());
240   h3->GetXaxis()->SetTitle(hSparse->GetAxis(xDim)->GetTitle());
241   h3->GetYaxis()->SetTitle(hSparse->GetAxis(yDim)->GetTitle());
242   h3->GetZaxis()->SetTitle(hSparse->GetAxis(zDim)->GetTitle());
243   h3->SetTitle(title.Data());  
244   aFolderObj->Add(h3);
245 }