]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaAnalysis/AliTRDqaBasic.cxx
Coverity report problems
[u/mrichter/AliRoot.git] / TRD / qaAnalysis / AliTRDqaBasic.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, 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 /* $Id: AliTRDqaBasic.cxx 23387 2008-01-17 17:25:16Z cblume $ */
17
18 //
19 // This class is a part of a package of high level QA monitoring for TRD.
20 //
21 // S. Radomski
22 // radomski@physi.uni-heidelberg.de
23 // March 2008
24 //
25
26 #include "AliTRDqaBasic.h"
27 #include "AliTRDqaAT.h"
28
29 #include "TH1D.h"
30 #include "TH2D.h"
31 #include "TFile.h"
32 #include "TChain.h"
33
34 #include "AliESDEvent.h"
35 #include "AliESDtrack.h"
36 #include "AliESDtrackCuts.h"
37
38 //______________________________________________________________________________
39
40 AliTRDqaBasic::AliTRDqaBasic() 
41   : AliAnalysisTask("",""),  
42     fChain(0),
43     fESD(0),
44     fOutputContainer(0),
45     fTrackCuts(0),
46     fStatus(0),
47     fnTracks(0),
48     fPtIn(0),
49     fPtOut(0),
50     fPtVtx(0),
51     fPtVtxSec(0),
52     fPtPt(0)
53 {
54   //
55   // default constructor
56   //
57  
58 }
59 //______________________________________________________________________________
60
61 AliTRDqaBasic:: AliTRDqaBasic(const AliTRDqaBasic & /*trd*/)
62   : AliAnalysisTask("",""),  
63     fChain(0),
64     fESD(0),
65     fOutputContainer(0),
66     fTrackCuts(0),
67     fStatus(0),
68     fnTracks(0),
69     fPtIn(0),
70     fPtOut(0),
71     fPtVtx(0),
72     fPtVtxSec(0),
73     fPtPt(0)
74 {
75   //
76   // Copy constructor
77   //
78
79 }
80
81 //______________________________________________________________________________
82 AliTRDqaBasic::AliTRDqaBasic(const char *name) 
83   : AliAnalysisTask(name,""),  
84     fChain(0),
85     fESD(0),
86     fOutputContainer(0),
87     fTrackCuts(0),
88     fStatus(0),
89     fnTracks(0),
90     fPtIn(0),
91     fPtOut(0),
92     fPtVtx(0),
93     fPtVtxSec(0),
94     fPtPt(0)
95 {
96   // Constructor.
97   // Input slot #0 works with an Ntuple
98   DefineInput(0, TChain::Class());
99   // Output slot #0 writes into a TH1 container
100   DefineOutput(0,  TObjArray::Class()) ; 
101 }
102
103 //______________________________________________________________________________
104 void AliTRDqaBasic::ConnectInputData(const Option_t *)
105 {
106   // Initialisation of branch container and histograms 
107
108   //AliInfo(Form("*** Initialization of %s", GetName())) ; 
109
110   fChain = (TChain*)GetInputData(0);
111   fESD = new AliESDEvent();
112   fESD->ReadFromTree(fChain);
113 }
114
115 //________________________________________________________________________
116 void AliTRDqaBasic::CreateOutputObjects()
117 {
118   // build histograms
119  
120   fStatus   = new TH1D("status", ";status bit", 32, -0.5, 31.5);
121   fPtIn     = new TH1D("ptIn", ";p_{T}, GeV/c", 100, 0, 15);
122   fPtOut    = new TH1D("ptOut", ";p_{T}, GeV/c", 100, 0, 15);
123   fPtVtx    = new TH1D("ptVtx", ";p_{T}, GeV/c", 100, 0, 15);
124   fPtVtxSec = new TH1D("ptVtxSec", ";p_{T}, GeV/c", 100, 0, 15);
125   fnTracks  = new TH1D("nTracks", "", 200, -0.5, 199.5);
126   
127   fPtPt     = new TH2D("ptpt", ";p_{T} from inner plane, GeV/c;p_{T} at vertex, GeV/c", 
128                        100, 0, 10, 100, 0, 10);
129
130   Int_t c = 0;
131   fOutputContainer = new TObjArray(50);
132
133   fOutputContainer->AddAt(fStatus, c++);
134   fOutputContainer->AddAt(fPtIn, c++);
135   fOutputContainer->AddAt(fPtOut, c++);
136   fOutputContainer->AddAt(fPtVtx, c++);
137   fOutputContainer->AddAt(fPtVtxSec, c++);
138   fOutputContainer->AddAt(fnTracks, c++);
139   fOutputContainer->AddAt(fPtPt, c++);  
140
141   printf("n hist = %d\n", c);
142   
143   // initialize cuts
144   fTrackCuts =  new AliESDtrackCuts("AliESDtrackCuts");
145   fTrackCuts->DefineHistograms(0);
146
147   // default cuts for ITS+TPC
148   Double_t cov1 = 2;
149   Double_t cov2 = 2;
150   Double_t cov3 = 0.5;
151   Double_t cov4 = 0.5;
152   Double_t cov5 = 2;
153   Double_t nSigma = 3;
154
155   TString tag("Global tracking");
156
157   fTrackCuts->SetMaxCovDiagonalElements(cov1, cov2, cov3, cov4, cov5);
158
159   fTrackCuts->SetMaxNsigmaToVertex(nSigma);
160   fTrackCuts->SetRequireSigmaToVertex(kTRUE);
161
162   fTrackCuts->SetRequireTPCRefit(kTRUE);
163   fTrackCuts->SetAcceptKinkDaughters(kFALSE);
164
165   fTrackCuts->SetMinNClustersTPC(50);
166   fTrackCuts->SetMaxChi2PerClusterTPC(3.5);
167 }
168 //______________________________________________________________________________
169 void AliTRDqaBasic::Exec(Option_t *) 
170 {
171   // Process one event
172   Long64_t entry = fChain->GetReadEntry() ;
173   if (!(entry%100)) Info("Exec", "Entry = %lld", entry);
174
175   // Processing of one event 
176    
177   if (!fESD) {
178     //AliError("fESD is not connected to the input!") ; 
179     return ; 
180   }
181   
182   Int_t nTracks = fESD->GetNumberOfTracks();
183   fnTracks->Fill(nTracks); 
184
185   // track loop
186   for(Int_t i=0; i<nTracks; i++) {
187     
188     //
189     // track selection 
190     //
191     // param in and Out
192     // TRDrefit and TRDPid bit
193     //
194  
195     AliESDtrack *track = fESD->GetTrack(i);
196     const AliExternalTrackParam *paramOut = track->GetOuterParam();
197     const AliExternalTrackParam *paramIn = track->GetInnerParam();
198
199     // long track ..
200     if (!paramIn) continue;
201     if (!paramOut) continue;
202     
203     UInt_t status = track->GetStatus();
204     if (!(status & AliESDtrack::kTPCrefit)) continue;
205
206     //if (!fTrackCuts->AcceptTrack(track)) continue;
207     
208     AliTRDqaAT::FillStatus(fStatus, status);
209
210     fPtIn->Fill(paramIn->Pt());
211     fPtOut->Fill(paramOut->Pt());
212     fPtVtx->Fill(track->Pt());
213     if (track->GetConstrainedParam())
214       fPtVtxSec->Fill(track->GetConstrainedParam()->Pt());
215     
216     fPtPt->Fill(paramIn->Pt(), track->Pt());
217   }
218
219   PostData(0, fOutputContainer);
220 }
221
222 //______________________________________________________________________________
223 void AliTRDqaBasic::Terminate(Option_t *)
224 {
225   // save histograms
226   fOutputContainer = (TObjArray*)GetOutputData(0);
227   
228   TFile *file = new TFile("outBasic.root", "RECREATE");
229   fOutputContainer->Write();
230  
231   file->Flush();
232   file->Close();
233   delete file;
234
235   //for(Int_t i=0; i<fOutputContainer->GetEntries(); i++) {
236   //  TObject *obj = fOu
237   // }
238 }
239
240 //______________________________________________________________________________