]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTClusterQA.cxx
New classes after fixing some bugs
[u/mrichter/AliRoot.git] / MFT / AliMFTClusterQA.cxx
1 #include "TObject.h"
2 #include "AliRunLoader.h"
3 #include "AliRun.h"
4 #include "AliLoader.h"
5 #include "AliMFT.h"
6 #include "TClonesArray.h"
7 #include "AliMFTCluster.h"
8 #include "AliMFTSegmentation.h"
9 #include "TFile.h"
10 #include "TH1D.h"
11 #include "TH2D.h"
12 #include "AliLog.h"
13 #include "TString.h"
14
15 #include "AliMFTClusterQA.h"
16
17 //====================================================================================================================================================
18 //
19 // Class for the analysis of the MFT clusters (a.k.a. rec points). Few QA histograms are created
20 //
21 // Contact author: antonio.uras@cern.ch
22 //
23 //====================================================================================================================================================
24
25 ClassImp(AliMFTClusterQA)
26
27 //====================================================================================================================================================
28
29 AliMFTClusterQA::AliMFTClusterQA():
30   TObject(),
31   fMFTLoader(0),
32   fRunLoader(0),
33   fMFT(0),
34   fNPlanes(0),
35   fNEvents(0),
36   fEv(0),
37   fFileOut(0),
38   fReadDir(0),
39   fOutDir(0)
40 {
41   
42   // default constructor
43
44   for (Int_t iPlane=0; iPlane<fMaxNPlanesMFT; iPlane++) {
45     fHistNClustersPerEvent[iPlane] = 0;
46     fHistNPixelsPerCluster[iPlane] = 0;
47     fHistClusterSizeX[iPlane] = 0; 
48     fHistClusterSizeY[iPlane] = 0;
49     fClusterScatterPlotXY[iPlane] = 0;
50   }
51
52 }
53
54 //====================================================================================================================================================
55
56 void AliMFTClusterQA::Init(Char_t *readDir, Char_t *outDir, Int_t nEventsToAnalyze) {
57
58   fReadDir = readDir;
59   fOutDir  = outDir;
60
61   fRunLoader = AliRunLoader::Open(Form("%s/galice.root", fReadDir.Data()));
62   gAlice = fRunLoader->GetAliRun();
63   if (!gAlice) fRunLoader->LoadgAlice();
64   fNEvents = fRunLoader->GetNumberOfEvents();
65   if (nEventsToAnalyze>0) fNEvents = TMath::Min(fNEvents, nEventsToAnalyze);
66
67   fMFT = (AliMFT*) gAlice->GetDetector("MFT"); 
68   fNPlanes = fMFT->GetSegmentation()->GetNPlanes();
69
70   BookHistos();
71
72   fMFTLoader = fRunLoader->GetDetectorLoader("MFT");
73   fMFTLoader -> LoadRecPoints("READ");
74
75 }
76
77 //====================================================================================================================================================
78
79 Bool_t AliMFTClusterQA::LoadNextEvent() {
80
81   if (fEv>=fNEvents) return kFALSE;
82   AliDebug(1, Form("event %5d",fEv));
83   
84   fRunLoader->GetEvent(fEv);
85   fEv++;
86
87   if (!fMFTLoader->TreeR()->GetEvent()) return kTRUE;
88
89   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
90     Int_t nClusters = fMFT->GetRecPointsList(iPlane)->GetEntries();
91     fHistNClustersPerEvent[iPlane] -> Fill(nClusters);
92     fClusterScatterPlotXY[iPlane]  -> Fill(0., 0.);    // "scaler" bin
93     AliDebug(1,Form("nClusters = %5d", nClusters));
94     for (Int_t iCluster=0; iCluster<nClusters; iCluster++) {
95       AliMFTCluster *cluster = (AliMFTCluster*) (fMFT->GetRecPointsList(iPlane))->At(iCluster);
96       fHistNPixelsPerCluster[iPlane] -> Fill(cluster->GetSize());
97       fHistClusterSizeX[iPlane]      -> Fill(cluster->GetErrX()*1.e4);   // converted in microns
98       fHistClusterSizeY[iPlane]      -> Fill(cluster->GetErrY()*1.e4);   // converted in microns
99       fClusterScatterPlotXY[iPlane]  -> Fill(cluster->GetX(), cluster->GetY());
100     }
101   }
102
103   return kTRUE;
104
105 }  
106
107 //====================================================================================================================================================
108
109 void AliMFTClusterQA::BookHistos() {
110
111   fFileOut = new TFile(Form("%s/MFT.RecPoints.QA.root",fOutDir.Data()), "recreate");
112
113   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
114
115     fHistNClustersPerEvent[iPlane] = new TH1D(Form("fHistNClustersPerEvent_Pl%02d",iPlane), 
116                                               Form("Number of clusters per event in Plane%02d",iPlane),
117                                               25000, -0.5, 24999.5);
118
119     fHistNPixelsPerCluster[iPlane] = new TH1D(Form("fHistNPixelsPerCluster_Pl%02d",iPlane), 
120                                               Form("Number of pixels per cluster in Plane%02d",iPlane),
121                                               15, -0.5, 14.5);
122
123     fHistClusterSizeX[iPlane]      = new TH1D(Form("fHistClusterSizeX_Pl%02d",iPlane), 
124                                               Form("#Deltax for clusters in Plane%02d",iPlane), 
125                                               100, 0., 100.);
126     
127     fHistClusterSizeY[iPlane]      = new TH1D(Form("fHistClusterSizeY_Pl%02d",iPlane), 
128                                               Form("#Deltay for clusters in Plane%02d",iPlane), 
129                                               100, 0., 100.);
130     
131     fHistNClustersPerEvent[iPlane] -> SetXTitle("N_{clusters} per Event");
132     fHistNPixelsPerCluster[iPlane] -> SetXTitle("N_{pixels} per Cluster");
133     fHistClusterSizeX[iPlane]      -> SetXTitle("#Deltax  [#mum]");
134     fHistClusterSizeY[iPlane]      -> SetXTitle("#Deltay  [#mum]");
135
136     fHistNClustersPerEvent[iPlane] -> Sumw2();
137     fHistNPixelsPerCluster[iPlane] -> Sumw2();
138     fHistClusterSizeX[iPlane]      -> Sumw2();
139     fHistClusterSizeY[iPlane]      -> Sumw2();
140
141     //------------------------------------------------------------
142
143     Int_t rMax = Int_t(10.*(fMFT->GetSegmentation()->GetPlane(iPlane)->GetRMaxSupport()));
144     fClusterScatterPlotXY[iPlane] = new TH2D(Form("fClusterScatterPlotXY_Pl%02d",iPlane),
145                                              Form("Cluster scatter plot (Plane%02d)",iPlane),
146                                              2*rMax+1, (-rMax-0.5)/10., (rMax+0.5)/10., 2*rMax+1, (-rMax-0.5)/10., (rMax+0.5)/10.);
147     
148     fClusterScatterPlotXY[iPlane] -> Sumw2();
149     
150   }
151   
152 }
153
154 //====================================================================================================================================================
155
156 void AliMFTClusterQA::Terminate() {
157
158   AliInfo("Writing QA histos...");
159
160   fFileOut->cd();
161
162   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
163     fHistNClustersPerEvent[iPlane] -> Write();
164     fHistNPixelsPerCluster[iPlane] -> Write();
165     fHistClusterSizeX[iPlane]      -> Write();
166     fHistClusterSizeY[iPlane]      -> Write();
167     fClusterScatterPlotXY[iPlane]  -> Write();
168   }
169
170   fFileOut -> Close();
171
172 }
173
174 //====================================================================================================================================================