AliPhysics  a56b849 (a56b849)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliV0ReaderStrange.cxx
Go to the documentation of this file.
1 #include <vector>
2 #include <TGeoGlobalMagField.h>
3 
4 #include "AliV0ReaderStrange.h"
5 #include "AliKFParticle.h"
6 #include "AliV0ParticleStrange.h"
7 #include "AliAODv0.h"
8 #include "AliESDv0.h"
9 #include "AliV0.h"
10 #include "AliAODEvent.h"
11 #include "AliESDEvent.h"
12 #include "AliVEvent.h"
13 #include "AliPID.h"
14 #include "AliMCEvent.h"
15 #include "AliMCEventHandler.h"
16 #include "AliESDpid.h"
17 #include "AliESDtrackCuts.h"
18 #include "TRandom3.h"
19 #include "AliGenCocktailEventHeader.h"
20 #include "TList.h"
21 #include "AliKFConversionPhoton.h"
22 #include "AliAODConversionPhoton.h"
24 #include "TVector.h"
25 #include "AliKFVertex.h"
26 #include "AliAODTrack.h"
27 #include "AliESDtrack.h"
28 #include "AliAnalysisManager.h"
29 #include "AliInputEventHandler.h"
30 #include "AliAODHandler.h"
31 #include "AliPIDResponse.h"
32 #include "TChain.h"
33 #include "TFile.h"
34 #include "TString.h"
35 #include "TObjArray.h"
36 #include "AliVTrack.h"
37 #include "AliKFParticle.h"
38 
39 class iostream;
40 
41 using namespace std;
42 
43 ClassImp(AliV0ReaderStrange)
44 
45 //________________________________________________________________________
47  fEventCuts(NULL),
48  fV0Cuts(NULL),
49  fConversionGammas(NULL),
50  fEventIsSelected(kFALSE),
51  fVectorFoundGammas(0)
52 {
53  // Default constructor
54 
55  DefineInput(0, TChain::Class());
56 }
57 
58 //________________________________________________________________________
60 {
61  // default deconstructor
62 
63  if(fConversionGammas){
64  fConversionGammas->Delete();// Clear Objects
65  delete fConversionGammas;
66  fConversionGammas=0x0;
67  }
68 }
69 
70 
71 //________________________________________________________________________
73 {
74  // Initialize function to be called once before analysis
75  if(fV0Cuts==NULL){
76  if(fV0Cuts==NULL)AliError("No V0 Cut Selection initialized");
77  }
78  if(fEventCuts==NULL){
79  if(fEventCuts==NULL)AliError("No Event Cut Selection initialized");
80  }
81 
82  if(fConversionGammas != NULL){
83  delete fConversionGammas;
84  fConversionGammas=NULL;
85  }
86 
87  if(fConversionGammas == NULL){
88  fConversionGammas = new TClonesArray("AliKFParticle",1000);
89 // if(kUseAODConversionPhoton){
90 // fConversionGammas = new TClonesArray("AliAODConversionPhoton",100);}
91 // else{
92 // fConversionGammas = new TClonesArray("AliKFConversionPhoton",100);}
93  }
94  fConversionGammas->Delete();//Reset the TClonesArray
95 }
96 
97 //________________________________________________________________________
99 {
100  fVectorFoundGammas.clear();
101 }
102 
103 //________________________________________________________________________
105 
106  AliESDEvent * esdEvent = dynamic_cast<AliESDEvent*>(fInputEvent);
107  if(esdEvent) {
108  if (!TGeoGlobalMagField::Instance()->GetField()) esdEvent->InitMagneticField();
109  }
110 
111  // Check if correctly initialized
112  if(!fConversionGammas)Init();
113 
114  // User Exec
115  fEventIsSelected=ProcessEvent(fInputEvent,fMCEvent);
116 }
117 
118 //________________________________________________________________________
119 Bool_t AliV0ReaderStrange::ProcessEvent(AliVEvent *inputEvent,AliMCEvent *mcEvent)
120 {
121  //Reset the TClonesArray
122  fConversionGammas->Delete();
123 
124  fInputEvent=inputEvent;
125  fMCEvent=mcEvent;
126 
127  if(!fInputEvent){
128  AliError("No Input event");
129  return kFALSE;
130  }
131  if(!fEventCuts){AliError("No EventCuts");return kFALSE;}
132  if(!fV0Cuts){AliError("No V0 Cuts");return kFALSE;}
133 
134  // Event Cuts
135  if(!fEventCuts->EventIsSelected(fInputEvent,fMCEvent))return kFALSE;
136 
137  // Set Magnetic Field
138  AliKFParticle::SetField(fInputEvent->GetMagneticField());
139 
140  if(fInputEvent->IsA()==AliESDEvent::Class()){
141  ProcessESDV0s();
142  }
143  if(fInputEvent->IsA()==AliAODEvent::Class()){
144  ProcessAODV0s();
145  }
146 
147  return kTRUE;
148 }
149 
152 {
153  // Process ESD V0s for V0 reconstruction
154 
155  AliV0ParticleStrange *fCurrentMotherLambdaCandidate=NULL;
156 
157  AliESDEvent *fESDEvent=dynamic_cast<AliESDEvent*>(fInputEvent);
158 
159  for(Int_t iV0=0; iV0<fESDEvent->GetNumberOfV0s(); ++iV0){
160  AliESDv0 *fCurrentV0=(AliESDv0*)(fESDEvent->GetV0(iV0));
161  if(!fCurrentV0){
162  printf("Requested V0 does not exist");
163  continue;
164  }
165 
166  fCurrentMotherLambdaCandidate = ReconstructV0(fESDEvent, fCurrentV0, iV0);
167  if(fCurrentMotherLambdaCandidate){
168  new((*fConversionGammas)[fConversionGammas->GetEntriesFast()]) AliV0ParticleStrange(*fCurrentMotherLambdaCandidate);
169  delete fCurrentMotherLambdaCandidate;
170  fCurrentMotherLambdaCandidate=NULL;
171  }
172  }
173 
174  return kTRUE;
175 }
176 
178 AliV0ParticleStrange *AliV0ReaderStrange::ReconstructV0(AliESDEvent *fESDEvent, AliESDv0 *fCurrentV0,Int_t currentV0Index)
179 {
180  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kV0In);
181 
182  //checks if on the fly mode is set
183  if(!fV0Cuts->SelectV0Finder(fCurrentV0->GetOnFlyStatus())){
184  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kOnFly);
185  return 0x0;
186  }
187 
188 // Double_t lV0CosineOfPointingAngle = fCurrentV0->GetV0CosineOfPointingAngle();
189 // GetV0CosineOfPointingAngle(lBestPrimaryVtxPos[0],lBestPrimaryVtxPos[1],lBestPrimaryVtxPos[2]);
190 
191  if(fCurrentV0->GetV0CosineOfPointingAngle() < 0.98) return 0x0;
192 
193  AliVTrack * pos = (AliVTrack*)fESDEvent->GetTrack(fCurrentV0->GetPindex());
194  AliVTrack * neg = (AliVTrack*)fESDEvent->GetTrack(fCurrentV0->GetNindex());
195  const AliExternalTrackParam * paramPos = fCurrentV0->GetParamP();
196  const AliExternalTrackParam * paramNeg = fCurrentV0->GetParamN();
197 
198  if(pos->GetSign() <0){//change tracks
199  pos=neg;
200  neg=fESDEvent->GetTrack(fCurrentV0->GetPindex());
201  paramPos=paramNeg;
202  paramNeg=fCurrentV0->GetParamP();
203  }
204 
205  if(!pos || !neg ) {
206  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kNoTracks);
207  return 0x0;
208 // return;
209  }
210  //cuts ------------ //remove like sign pairs
211  if(pos->GetSign() == neg->GetSign()){
212  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kSameSign);
213  return 0x0;
214  }
215 
216 // if(pos->Pt() < .16) return 0x0;
217 // if(neg->Pt() < .16) return 0x0;
218 // if(TMath::Abs(pos->Eta()) > 0.8) return 0x0;
219 // if(TMath::Abs(neg->Eta()) > 0.8) return 0x0;
220 
221 // if( neg->GetKinkIndex(0) > 0 || pos->GetKinkIndex(0) > 0) return 0x0 ;
222 // if( !(pos->GetStatus() & AliESDtrack::kTPCrefit) || !(neg->GetStatus() & AliESDtrack::kTPCrefit) ) return 0x0;
223 
224  Bool_t isPiPlus = fV0Cuts->GetPIDpion(pos);
225  Bool_t isPiMinus = fV0Cuts->GetPIDpion(neg);
226  Bool_t isProton = fV0Cuts->GetPIDproton(pos);
227  Bool_t isAntiProton = fV0Cuts->GetPIDproton(neg);
228 
229 
230  if(!((isProton && isPiMinus) || (isAntiProton && isPiPlus))) return 0x0;
231 
232  // pi plus: 211, pi minus: -211, proton: 2212, anti-proton: -2212
233 
234  Int_t pdgCodePos = 2211; //proton
235  Int_t pdgCodeNeg = -211; // pi minus
236  if(isPiPlus){
237  pdgCodePos = -211;
238  }
239  if(isAntiProton){
240  pdgCodeNeg = -2212;
241  }
242 
243 // AliVVertex *vertex = fESDEvent->GetPrimaryVertex();
244 // fCurrentV0->DecayLength(vertex);
245 
246  const AliKFParticle fCurrentPositiveKFParticle(*(paramPos), pdgCodePos);
247  const AliKFParticle fCurrentNegativeKFParticle(*(paramNeg), pdgCodeNeg);
248  const Bool_t gamma = kFALSE;
249  AliKFParticle *fCurrentMother = new AliKFParticle(fCurrentPositiveKFParticle, fCurrentNegativeKFParticle, gamma);
250 
251  AliV0ParticleStrange *fCurrentMotherV0 = new AliV0ParticleStrange(fCurrentMother);
252 
253  if(fMCEvent){
254 
255  Int_t labelp=TMath::Abs(fV0Cuts->GetTrack(fESDEvent,fCurrentMotherV0->GetTrackLabelPositive())->GetLabel());
256  Int_t labeln=TMath::Abs(fV0Cuts->GetTrack(fESDEvent,fCurrentMotherV0->GetTrackLabelNegative())->GetLabel());
257 
258 // cout << "rec: " << currentTrackLabels[0] << "\t" << currentTrackLabels[1] << endl;
259 // cout << "recProp: " << fCurrentMotherKF->GetTrackLabelPositive() << "\t" << fCurrentMotherKF->GetTrackLabelNegative() << endl;
260 // cout << "MC: " << labeln << "\t" << labelp << endl;
261 
262  TParticle *fNegativeMCParticle = 0x0;
263  if(labeln>-1) fNegativeMCParticle = fMCEvent->Particle(labeln);
264  TParticle *fPositiveMCParticle = 0x0;
265  if(labelp>-1) fPositiveMCParticle = fMCEvent->Particle(labelp);
266 
267  if(fPositiveMCParticle&&fNegativeMCParticle){
268  fCurrentMotherV0->SetMCLabelPositive(labelp);
269  fCurrentMotherV0->SetMCLabelNegative(labeln);
270  }
271  }
272 
273  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kV0Out);
274  return fCurrentMotherV0;
275 }
276 
277 
280 {
281  // Process ESD V0s for V0 reconstruction
282 
283  AliV0ParticleStrange *fCurrentMotherLambdaCandidate=NULL;
284 
285  AliAODEvent *fAODEvent=dynamic_cast<AliAODEvent*>(fInputEvent);
286 
287  for(Int_t iV0=0; iV0<fAODEvent->GetNumberOfV0s(); ++iV0){
288  AliAODv0 *fCurrentV0=(AliAODv0*)(fAODEvent->GetV0(iV0));
289  if(!fCurrentV0){
290  printf("Requested V0 does not exist");
291  continue;
292  }
293 
294  fCurrentMotherLambdaCandidate = ReconstructV0(fAODEvent, fCurrentV0, iV0);
295  if(fCurrentMotherLambdaCandidate){
296  new((*fConversionGammas)[fConversionGammas->GetEntriesFast()]) AliV0ParticleStrange(*fCurrentMotherLambdaCandidate);
297  delete fCurrentMotherLambdaCandidate;
298  fCurrentMotherLambdaCandidate=NULL;
299  }
300  }
301 
302  return kTRUE;
303 }
304 
306 AliV0ParticleStrange *AliV0ReaderStrange::ReconstructV0(AliAODEvent *fAODEvent, AliAODv0 *fCurrentV0,Int_t currentV0Index)
307 {
308  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kV0In);
309 
310  //checks if on the fly mode is set
311  if(!fV0Cuts->SelectV0Finder(fCurrentV0->GetOnFlyStatus())){
312  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kOnFly);
313  return 0x0;
314  }
315 
316 // Double_t lV0CosineOfPointingAngle = fCurrentV0->GetV0CosineOfPointingAngle();
317 // GetV0CosineOfPointingAngle(lBestPrimaryVtxPos[0],lBestPrimaryVtxPos[1],lBestPrimaryVtxPos[2]);
318 
319  AliAODVertex *vertex = fAODEvent->GetPrimaryVertex();
320  if(fCurrentV0->CosPointingAngle(vertex) < 0.98) return 0x0;
321 
322  AliVTrack* pos = (AliVTrack*)fCurrentV0->GetDaughter(0);
323  AliVTrack* neg = (AliVTrack*)fCurrentV0->GetDaughter(1);
324 
325  AliExternalTrackParam paramPos;
326  paramPos.CopyFromVTrack(pos);
327  AliExternalTrackParam paramNeg;
328  paramPos.CopyFromVTrack(neg);
329 
330  if(pos->GetSign() <0){//change tracks
331  pos=neg;
332  neg=(AliVTrack*)fCurrentV0->GetDaughter(0);
333  paramPos=paramNeg;
334  paramNeg.CopyFromVTrack(pos);
335  }
336 
337  if(!pos || !neg ) {
338  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kNoTracks);
339  return 0x0;
340  }
341 
342  //cuts ------------ //remove like sign pairs
343  if(pos->GetSign() == neg->GetSign()){
344  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kSameSign);
345  return 0x0;
346  }
347 
348 // if(pos->Pt() < .16) return 0x0;
349 // if(neg->Pt() < .16) return 0x0;
350 // if(TMath::Abs(pos->Eta()) > 0.8) return 0x0;
351 // if(TMath::Abs(neg->Eta()) > 0.8) return 0x0;
352 
353 // if( neg->GetKinkIndex(0) > 0 || pos->GetKinkIndex(0) > 0) return 0x0 ;
354 // if( !(pos->GetStatus() & AliESDtrack::kTPCrefit) || !(neg->GetStatus() & AliESDtrack::kTPCrefit) ) return 0x0;
355 
356  Bool_t isPiPlus = fV0Cuts->GetPIDpion(pos);
357  Bool_t isPiMinus = fV0Cuts->GetPIDpion(neg);
358  Bool_t isProton = fV0Cuts->GetPIDproton(pos);
359  Bool_t isAntiProton = fV0Cuts->GetPIDproton(neg);
360 
361 
362  if(!((isProton && isPiMinus) || (isAntiProton && isPiPlus))) return 0x0;
363 
364  // pi plus: 211, pi minus: -211, proton: 2212, anti-proton: -2212
365 
366  Int_t pdgCodePos = 2211; //proton
367  Int_t pdgCodeNeg = -211; // pi minus
368  if(isPiPlus){
369  pdgCodePos = -211;
370  }
371  if(isAntiProton){
372  pdgCodeNeg = -2212;
373  }
374 
375  const AliKFParticle fCurrentPositiveKFParticle(paramPos, pdgCodePos);
376  const AliKFParticle fCurrentNegativeKFParticle(paramNeg, pdgCodeNeg);
377  const Bool_t gamma = kFALSE;
378  AliKFParticle *fCurrentMother = new AliKFParticle(fCurrentPositiveKFParticle, fCurrentNegativeKFParticle, gamma);
379 
380  AliV0ParticleStrange *fCurrentMotherV0 = new AliV0ParticleStrange(fCurrentMother);
381 
382  if(fMCEvent){
383 
384  Int_t labelp=TMath::Abs(fV0Cuts->GetTrack(fAODEvent,fCurrentMotherV0->GetTrackLabelPositive())->GetLabel());
385  Int_t labeln=TMath::Abs(fV0Cuts->GetTrack(fAODEvent,fCurrentMotherV0->GetTrackLabelNegative())->GetLabel());
386 
387 // cout << "rec: " << currentTrackLabels[0] << "\t" << currentTrackLabels[1] << endl;
388 // cout << "recProp: " << fCurrentMotherKF->GetTrackLabelPositive() << "\t" << fCurrentMotherKF->GetTrackLabelNegative() << endl;
389 // cout << "MC: " << labeln << "\t" << labelp << endl;
390 
391  TParticle *fNegativeMCParticle = 0x0;
392  if(labeln>-1) fNegativeMCParticle = fMCEvent->Particle(labeln);
393  TParticle *fPositiveMCParticle = 0x0;
394  if(labelp>-1) fPositiveMCParticle = fMCEvent->Particle(labelp);
395 
396  if(fPositiveMCParticle&&fNegativeMCParticle){
397  fCurrentMotherV0->SetMCLabelPositive(labelp);
398  fCurrentMotherV0->SetMCLabelNegative(labeln);
399  }
400  }
401 
402 
403  fV0Cuts->FillV0CutIndex(AliV0CutsStrange::kV0Out);
404  return fCurrentMotherV0;
405 }
void SetMCLabelNegative(Int_t label)
int Int_t
Definition: External.C:63
void SetMCLabelPositive(Int_t label)
virtual void UserExec(Option_t *option)
Bool_t ProcessEvent(AliVEvent *inputEvent, AliMCEvent *mcEvent=NULL)
AliV0ParticleStrange * ReconstructV0(AliESDEvent *fESDEvent, AliESDv0 *fCurrentV0, Int_t currentV0Index)
const char Option_t
Definition: External.C:48
bool Bool_t
Definition: External.C:53