]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/esdTrackCuts/AliTestESDtrackCutsSelector.cxx
Added selector to check esd track cuts.
[u/mrichter/AliRoot.git] / PWG0 / esdTrackCuts / AliTestESDtrackCutsSelector.cxx
1 /* $Id$ */
2
3 #include "AliTestESDtrackCutsSelector.h"
4
5 #include <TStyle.h>
6 #include <TSystem.h>
7 #include <TCanvas.h>
8 #include <TParticle.h>
9 #include <TParticlePDG.h>
10 #include <TVector3.h>
11 #include <TChain.h>
12 #include <TFile.h>
13 #include <TH1F.h>
14
15 #include <TSelector.h>
16 #include <TFile.h>
17
18 #include <AliLog.h>
19 #include <AliESD.h>
20 #include <AliRunLoader.h>
21 #include <AliStack.h>
22
23 #include "esdTrackCuts/AliESDtrackCuts.h"
24 #include "AliPWG0Helper.h"
25
26
27 ClassImp(AliTestESDtrackCutsSelector)
28
29 AliTestESDtrackCutsSelector::AliTestESDtrackCutsSelector() :
30   AliSelectorRL(),
31   fEsdTrackCutsAll(0),
32   fEsdTrackCutsPri(0),
33   fEsdTrackCutsSec(0)
34 {
35   //
36   // Constructor. Initialization of pointers
37   //
38 }
39
40 AliTestESDtrackCutsSelector::~AliTestESDtrackCutsSelector()
41 {
42   //
43   // Destructor
44   //
45
46   // histograms are in the output list and deleted when the output
47   // list is deleted by the TSelector dtor
48 }
49
50 void AliTestESDtrackCutsSelector::Begin(TTree* tree)
51 {
52   // Begin function
53
54   AliSelectorRL::Begin(tree);
55
56   ReadUserObjects(tree);
57 }
58
59 void AliTestESDtrackCutsSelector::ReadUserObjects(TTree* tree)
60 {
61   // read the user objects, called from slavebegin and begin
62
63   if (!fEsdTrackCutsAll && fInput)
64     fEsdTrackCutsAll = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("esdTrackCutsAll"));
65   if (!fEsdTrackCutsPri && fInput)
66     fEsdTrackCutsPri = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("esdTrackCutsPri"));
67   if (!fEsdTrackCutsSec && fInput)
68     fEsdTrackCutsSec = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("esdTrackCutsSec"));
69
70   if (!fEsdTrackCutsAll && tree)
71     fEsdTrackCutsAll = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("esdTrackCutsAll"));
72   if (!fEsdTrackCutsPri && tree)
73     fEsdTrackCutsPri = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("esdTrackCutsPri"));
74   if (!fEsdTrackCutsSec && tree)
75     fEsdTrackCutsSec = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("esdTrackCutsSec"));
76
77   if (!fEsdTrackCutsAll || !fEsdTrackCutsPri || !fEsdTrackCutsSec)
78      AliDebug(AliLog::kError, "ERROR: Could not read esdTrackCutsXXX from input list.");
79 }
80
81 void AliTestESDtrackCutsSelector::SlaveBegin(TTree* tree)
82 {
83   // The SlaveBegin() function is called after the Begin() function.
84   // When running with PROOF SlaveBegin() is called on each slave server.
85   // The tree argument is deprecated (on PROOF 0 is passed).
86
87   AliSelectorRL::SlaveBegin(tree);
88
89   ReadUserObjects(tree);
90 }
91
92 Bool_t AliTestESDtrackCutsSelector::Process(Long64_t entry)
93 {
94   // The Process() function is called for each entry in the tree (or possibly
95   // keyed object in the case of PROOF) to be processed. The entry argument
96   // specifies which entry in the currently loaded tree is to be processed.
97   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
98   // to read either all or the required parts of the data. When processing
99   // keyed objects with PROOF, the object is already loaded and is available
100   // via the fObject pointer.
101   //
102   // This function should contain the "body" of the analysis. It can contain
103   // simple or elaborate selection criteria, run algorithms on the data
104   // of the event and typically fill histograms.
105
106   // WARNING when a selector is used with a TChain, you must use
107   //  the pointer to the current TTree to call GetEntry(entry).
108   //  The entry is always the local entry number in the current tree.
109   //  Assuming that fTree is the pointer to the TChain being processed,
110   //  use fTree->GetTree()->GetEntry(entry).
111
112   if (AliSelectorRL::Process(entry) == kFALSE)
113     return kFALSE;
114
115   // Check prerequisites
116   if (!fESD) {
117     AliDebug(AliLog::kError, "ESD branch not available");
118     return kFALSE;
119   }  
120
121   if (!AliPWG0Helper::IsVertexReconstructed(fESD)) {
122     AliDebug(AliLog::kDebug+5, "Vertex is not reconstructed");
123     return kFALSE;
124   }
125
126   // check if the esd track cut objects are there
127   if (!fEsdTrackCutsAll || !fEsdTrackCutsPri || !fEsdTrackCutsSec) {
128     AliDebug(AliLog::kError, "fEsdTrackCutsXXX not available");
129     return kFALSE;
130   }
131
132   // get particle stack
133   AliStack* stack = GetStack();
134   if (!stack) {
135     AliDebug(AliLog::kError, "Stack not available");
136     return kFALSE;
137   }
138   Int_t nPrim  = stack->GetNprimary();
139   
140   // ########################################################
141   // loop over esd tracks
142   Int_t nTracks = fESD->GetNumberOfTracks();
143
144   // count the number of "good" tracks as parameter for vertex reconstruction efficiency
145   for (Int_t t=0; t<nTracks; t++) {
146     AliDebug(AliLog::kDebug+1, Form("ESD Loop: Processing track %d.", t));
147
148     AliESDtrack* esdTrack = fESD->GetTrack(t);
149
150     fEsdTrackCutsAll->AcceptTrack(esdTrack);
151
152     // using the properties of the mc particle
153     Int_t label = TMath::Abs(esdTrack->GetLabel());
154     if (label == 0) {
155       AliDebug(AliLog::kWarning, Form("WARNING: cannot find corresponding mc part for track %d.", t));
156       continue;
157     }
158     TParticle* particle = stack->Particle(label);
159     if (!particle) {
160       AliDebug(AliLog::kError, Form("UNEXPECTED: part with label %d not found in stack (track loop).", label));
161       continue;
162     }
163     if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
164       fEsdTrackCutsPri->AcceptTrack(esdTrack);
165     else
166       fEsdTrackCutsSec->AcceptTrack(esdTrack);
167   }
168   
169   return kTRUE;
170 }
171
172 void AliTestESDtrackCutsSelector::SlaveTerminate()
173 {
174   // The SlaveTerminate() function is called after all entries or objects
175   // have been processed. When running with PROOF SlaveTerminate() is called
176   // on each slave server.
177
178   AliSelectorRL::SlaveTerminate();
179   
180   // Add the histograms to the output on each slave server
181   if (!fOutput)
182   {
183     AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
184     return;
185   }
186
187   fOutput->Add(fEsdTrackCutsAll);
188   fOutput->Add(fEsdTrackCutsPri);
189   fOutput->Add(fEsdTrackCutsSec);
190 }
191
192 void AliTestESDtrackCutsSelector::Terminate()
193 {
194   // The Terminate() function is the last function to be called during
195   // a query. It always runs on the client, it can be used to present
196   // the results graphically or save the results to file.
197
198   AliSelectorRL::Terminate();
199
200   TFile* file = TFile::Open("trackCuts.root", "RECREATE");
201   fEsdTrackCutsAll->SaveHistograms("esdTrackCutsAll");
202   fEsdTrackCutsPri->SaveHistograms("esdTrackCutsPri");
203   fEsdTrackCutsSec->SaveHistograms("esdTrackCutsSec");
204
205   file->Close();
206 }