AliPhysics  720d1f3 (720d1f3)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliCaloTrackAODReader.cxx
Go to the documentation of this file.
1 
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * *
5  * Author: The ALICE Off-line Project. *
6  * Contributors are mentioned in the code where appropriate. *
7  * *
8  * Permission to use, copy, modify and distribute this software and its *
9  * documentation strictly for non-commercial purposes is hereby granted *
10  * without fee, provided that the above copyright notice appears in all *
11  * copies and that both the copyright notice and this permission notice *
12  * appear in the supporting documentation. The authors make no claims *
13  * about the suitability of this software for any purpose. It is *
14  * provided "as is" without express or implied warranty. *
15  **************************************************************************/
16 
17 //---- ANALYSIS system ----
18 #include "AliCaloTrackAODReader.h"
19 #include "AliAODInputHandler.h"
20 #include "AliMultiEventInputHandler.h"
21 #include "AliAnalysisManager.h"
22 #include "AliMixedEvent.h"
23 #include "AliAODEvent.h"
24 #include "AliLog.h"
25 
29 
30 //______________________________________________
32 //______________________________________________
34  AliCaloTrackReader(), fOrgInputEvent(0x0),
35  fSelectHybridTracks(0), fSelectPrimaryTracks(0),
36  fTrackFilterMask(0), fTrackFilterMaskComplementary(0),
37  fSelectFractionTPCSharedClusters(0), fCutTPCSharedClustersFraction(0)
38 {
39  fDataType = kAOD;
40 
41  fReadStack = kTRUE;
42  fReadAODMCParticles = kFALSE;
43 
44  fTrackFilterMask = 128;
45  fTrackFilterMaskComplementary = 0; // in case of hybrid tracks, without using the standard method
46 
49 }
50 
51 //_________________________________________________________
54 //_________________________________________________________
56 {
57  AliAODEvent * aodevent = dynamic_cast<AliAODEvent*>(fInputEvent);
58  if(!aodevent) return kFALSE;
59 
60  if (aodevent->GetPrimaryVertex() != NULL)
61  {
62  if(aodevent->GetPrimaryVertex()->GetNContributors() > 0)
63  {
64  return kTRUE;
65  }
66  }
67 
68  if(aodevent->GetPrimaryVertexSPD() != NULL)
69  {
70  if(aodevent->GetPrimaryVertexSPD()->GetNContributors() > 0)
71  {
72  return kTRUE;
73  }
74  else
75  {
76  AliWarning(Form("Number of contributors from bad vertex type:: %s",
77  aodevent->GetPrimaryVertex()->GetName()));
78  return kFALSE;
79  }
80  }
81 
82  return kFALSE;
83 }
84 
85 
86 //___________________________________________________
91 //___________________________________________________
93 {
95 
96  if(fFillCTS)
97  {
98  for(Int_t i = 0; i < 4; i++)
99  {
100  TString names[] = {"FilterBit_Hybrid", "SPDHit", "SharedCluster", "Primary"};
101 
102  fhCTSAODTrackCutsPt[i] = new TH1F(Form("hCTSReaderAODClusterCuts_%d_%s",i,names[i].Data()),
103  Form("AOD CTS Cut %d, %s",i,names[i].Data()),
105  fhCTSAODTrackCutsPt[i]->SetYTitle("# tracks");
106  fhCTSAODTrackCutsPt[i]->SetXTitle("#it{p}_{T} (GeV)");
108  }
109  }
110 
111  return fOutputContainer ;
112 }
113 
114 //________________________________________________________
116 //________________________________________________________
118 {
119  // Recover the string from the mother class
120  TString parList = (AliCaloTrackReader::GetListOfParameters())->GetString();
121 
122  const Int_t buffersize = 255;
123  char onePar[buffersize] ;
124 
125  snprintf(onePar,buffersize,"AOD Track: Hybrid %d, Filter bit %d, Complementary bit %d, Primary %d; ",
127  parList+=onePar ;
128 
130  {
131  snprintf(onePar,buffersize,"Fraction of TPC shared clusters ON: %2.2f ", fCutTPCSharedClustersFraction) ;
132  parList+=onePar ;
133  }
134 
135  return new TObjString(parList) ;
136 
137 }
138 
139 //____________________________________________________________
141 //____________________________________________________________
143 {
144  TClonesArray * particles = NULL ;
145 
146  AliAODEvent * aod = dynamic_cast<AliAODEvent*> (fInputEvent) ;
147  if(aod) particles = (TClonesArray*) aod->FindListObject("mcparticles");
148 
149  return particles ;
150 }
151 
152 //___________________________________________________________
154 //___________________________________________________________
156 {
157  AliAODMCHeader *mch = NULL;
158 
159  AliAODEvent * aod = dynamic_cast<AliAODEvent*> (fInputEvent);
160  if(aod) mch = dynamic_cast<AliAODMCHeader*>(aod->FindListObject("mcHeader"));
161 
162  return mch;
163 }
164 
165 //_____________________________
176 //_____________________________
178 {
179  Int_t id = track->GetID();
180 
181  if( id < 0 ) id = TMath::Abs(id) - 1;
182 
183  return id;
184 }
185 
186 
187 //_____________________________________________________________________________
189 //_____________________________________________________________________________
191 {
192  AliAODTrack *aodtrack = dynamic_cast <AliAODTrack*>(track);
193 
194  if(!aodtrack) return kFALSE;
195 
196  AliDebug(2,Form("AOD track type: %d (primary %d), hybrid? %d",
197  aodtrack->GetType(),AliAODTrack::kPrimary,
198  aodtrack->IsHybridGlobalConstrainedGlobal()));
199 
200  // Hybrid?
202  {
203  if (!aodtrack->IsHybridGlobalConstrainedGlobal()) return kFALSE ;
204  }
205  else // Filter Bit?
206  {
207  Bool_t accept = aodtrack->TestFilterBit(fTrackFilterMask);
208 
209  if(!fSelectHybridTracks && !accept) return kFALSE ;
210 
211  if(fSelectHybridTracks) // Second filter bit for hybrids?
212  {
213  Bool_t acceptcomplement = aodtrack->TestFilterBit(fTrackFilterMaskComplementary);
214  if (!accept && !acceptcomplement) return kFALSE ;
215  }
216  }
217 
218  fhCTSAODTrackCutsPt[0]->Fill(aodtrack->Pt());
219 
220  //
222  { // Not much sense to use with TPC only or Hybrid tracks
223  if(!aodtrack->HasPointOnITSLayer(0) && !aodtrack->HasPointOnITSLayer(1)) return kFALSE ;
224  }
225 
226  fhCTSAODTrackCutsPt[1]->Fill(aodtrack->Pt());
227 
228  //
230  {
231  Double_t frac = 0;
232  Float_t ncls = Float_t(aodtrack->GetTPCncls ());
233  Float_t nclsS = Float_t(aodtrack->GetTPCnclsS());
234  if ( ncls> 0 ) frac = nclsS / ncls ;
235 
237  {
238  AliDebug(2,Form("\t Reject track, shared cluster fraction %f > %f",frac, fCutTPCSharedClustersFraction));
239  return kFALSE ;
240  }
241  }
242 
243  fhCTSAODTrackCutsPt[2]->Fill(aodtrack->Pt());
244 
245  //
246  if ( fSelectPrimaryTracks )
247  {
248  if ( aodtrack->GetType()!= AliAODTrack::kPrimary )
249  {
250  AliDebug(2,"\t Remove not primary track");
251  return kFALSE ;
252  }
253  }
254 
255  AliDebug(2,"\t accepted track!");
256 
257  fhCTSAODTrackCutsPt[3]->Fill(aodtrack->Pt());
258 
259  track->GetPxPyPz(pTrack) ;
260 
261  return kTRUE;
262 }
263 
264 //_________________________________________________________________
267 //_________________________________________________________________
269  AliAODEvent* aod,
270  AliMCEvent* mc)
271 {
272  //printf("AODInputHandler %p, MergeEvents %d \n",aodIH, aodIH->GetMergeEvents());
273 
274  Bool_t tesd = kFALSE ;
275  Bool_t taod = kTRUE ;
276  if ( strcmp(input->GetName(), "AliMixedEvent") == 0 )
277  {
278  AliMultiEventInputHandler* multiEH = dynamic_cast<AliMultiEventInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
279  if(multiEH){
280  if (multiEH->GetFormat() == 0 )
281  {
282  tesd = kTRUE ;
283  } else if (multiEH->GetFormat() == 1)
284  {
285  taod = kTRUE ;
286  }
287  }
288  else
289  {
290  AliFatal("MultiEventHandler is NULL");
291  return;
292  }
293  }
294  if (strcmp(input->GetName(),"AliESDEvent") == 0)
295  {
296  tesd = kTRUE ;
297  } else if (strcmp(input->GetName(),"AliAODEvent") == 0)
298  {
299  taod = kTRUE ;
300  }
301 
302  if(tesd)
303  {
304  SetInputEvent(aod);
305  SetOutputEvent(aod);
306  fOrgInputEvent = input;
307  }
308  else if(taod)
309  {
310  AliAODInputHandler* aodIH = dynamic_cast<AliAODInputHandler*>((AliAnalysisManager::GetAnalysisManager())->GetInputEventHandler());
311 
312  if (aodIH && aodIH->GetMergeEvents())
313  {
314  //Merged events, use output AOD.
315  SetInputEvent(aod);
316  SetOutputEvent(aod);
317  fOrgInputEvent = input;
318  }
319  else
320  {
321  SetInputEvent(input);
322  SetOutputEvent(aod);
323  }
324  }
325  else
326  {
327  AliFatal(Form("STOP : Wrong data format: %s",input->GetName()));
328  }
329 
330  SetMC(mc);
331 }
332 
Bool_t fReadAODMCParticles
Access kine information from filtered AOD MC particles.
double Double_t
Definition: External.C:58
virtual void SetMC(AliMCEvent *const mc)
Bool_t fSelectHybridTracks
Select CTS tracks of type hybrid.
AliVEvent * fOrgInputEvent
! Original input event, not from filtering
AliVEvent * fInputEvent
! pointer to esd or aod input.
ULong_t fTrackFilterMaskComplementary
Complementary Track selection bit, for AODs in case hybrid option selected.
virtual void SetInputEvent(AliVEvent *input)
Bool_t SelectTrack(AliVTrack *track, Double_t *pTrack)
Select AOD track using the AOD filter bits or predefined selection methods.
TList * fOutputContainer
! Output container with cut control histograms.
Int_t fEnergyHistogramNbins
Binning of the control histograms, min and max window.
Bool_t fFillCTS
Use data from CTS.
AliAODMCHeader * GetAODMCHeader() const
ULong_t fTrackFilterMask
Track selection bit, for AODs (any difference with track status?)
Class for event, clusters and tracks filtering and preparation for the AOD analysis.
Float_t fEnergyHistogramLimit[2]
Binning of the control histograms, number of bins.
int Int_t
Definition: External.C:63
Bool_t fSelectSPDHitTracks
Ensure that track hits SPD layers.
float Float_t
Definition: External.C:68
virtual void SetOutputEvent(AliAODEvent *aod)
Bool_t fSelectPrimaryTracks
Select CTS tracks of type primary.
virtual TObjString * GetListOfParameters()
Save parameters used for analysis in a string.
Base class for event, clusters and tracks filtering and preparation for the analysis.
TClonesArray * GetAODMCParticles() const
Float_t fCutTPCSharedClustersFraction
Fraction of TPC shared clusters to be accepted.
Bool_t Data(TH1F *h, Double_t *rangefit, Bool_t writefit, Double_t &sgn, Double_t &errsgn, Double_t &bkg, Double_t &errbkg, Double_t &sgnf, Double_t &errsgnf, Double_t &sigmafit, Int_t &status)
TObjString * GetListOfParameters()
Save parameters used for analysis in a string.
void SetInputOutputMCEvent(AliVEvent *esd, AliAODEvent *aod, AliMCEvent *mc)
Int_t GetTrackID(AliVTrack *track)
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
Bool_t fReadStack
Access kine information from stack.
Bool_t fSelectFractionTPCSharedClusters
Accept only TPC tracks with over a given fraction of shared clusters.
AliCaloTrackAODReader()
Default constructor. Initialize parameters.
bool Bool_t
Definition: External.C:53
Bool_t CheckForPrimaryVertex() const
virtual TList * GetCreateControlHistograms()
Int_t fDataType
Select MC: Kinematics, Data: ESD/AOD, MCData: Both.
TH1F * fhCTSAODTrackCutsPt[4]
! control histogram on the different CTS tracks selection cuts, pT