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