1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Authors: Svein Lindal, Daniel Lohner *
8 * based on: on older version (see aliroot up to v5-04-42-AN) *
10 * Authors: Kathrin Koch, Kenneth Aamodt, Ana Marin *
12 * Permission to use, copy, modify and distribute this software and its *
13 * documentation strictly for non-commercial purposes is hereby granted *
14 * without fee, provided that the above copyright notice appears in all *
15 * copies and that both the copyright notice and this permission notice *
16 * appear in the supporting documentation. The authors make no claims *
17 * about the suitability of this software for any purpose. It is *
18 * provided "as is" without express or implied warranty. *
19 **************************************************************************/
21 ////////////////////////////////////////////////
22 //---------------------------------------------
23 // Class reconstructing conversion photons from V0s
24 //---------------------------------------------
25 ////////////////////////////////////////////////
27 #include "AliV0ReaderV1.h"
28 #include "AliKFParticle.h"
31 #include "AliAODEvent.h"
32 #include "AliESDEvent.h"
34 #include "AliMCEvent.h"
36 #include "AliMCEventHandler.h"
37 #include "AliESDpid.h"
38 #include "AliESDtrackCuts.h"
40 #include "AliGenCocktailEventHeader.h"
42 #include "AliKFConversionPhoton.h"
43 #include "AliAODConversionPhoton.h"
44 #include "AliConversionPhotonBase.h"
46 #include "AliKFVertex.h"
47 #include "AliAODTrack.h"
48 #include "AliESDtrack.h"
49 #include "AliAnalysisManager.h"
50 #include "AliInputEventHandler.h"
51 #include "AliAODHandler.h"
52 #include "AliPIDResponse.h"
56 #include "TObjArray.h"
63 ClassImp(AliV0ReaderV1)
65 //________________________________________________________________________
66 AliV0ReaderV1::AliV0ReaderV1(const char *name) : AliAnalysisTaskSE(name),
67 fConversionCuts(NULL),
68 fConversionGammas(NULL),
69 fUseImprovedVertex(kTRUE),
70 fUseOwnXYZCalculation(kTRUE),
71 fUseConstructGamma(kFALSE),
72 kUseAODConversionPhoton(kTRUE),
74 fDeltaAODBranchName("GammaConv"),
75 fDeltaAODFilename("AliAODGammaConversion.root"),
77 fEventIsSelected(kFALSE),
80 // Default constructor
82 DefineInput(0, TChain::Class());
85 //________________________________________________________________________
86 AliV0ReaderV1::~AliV0ReaderV1()
88 // default deconstructor
90 if(fConversionGammas){
91 fConversionGammas->Delete();// Clear Objects
92 delete fConversionGammas;
93 fConversionGammas=0x0;
97 //________________________________________________________________________
98 AliV0ReaderV1::AliV0ReaderV1(AliV0ReaderV1 &original) : AliAnalysisTaskSE(original),
99 fConversionCuts(NULL),
100 fConversionGammas(NULL),
101 fUseImprovedVertex(original.fUseImprovedVertex),
102 fUseOwnXYZCalculation(original.fUseOwnXYZCalculation),
103 fUseConstructGamma(original.fUseConstructGamma),
104 kUseAODConversionPhoton(original.kUseAODConversionPhoton),
105 fCreateAOD(original.fCreateAOD),
106 fDeltaAODBranchName(original.fDeltaAODBranchName),
107 fDeltaAODFilename(original.fDeltaAODFilename),
108 fEventIsSelected(original.fEventIsSelected)
110 // Default constructor
112 DefineInput(0, TChain::Class());
115 //____________________________________________________________
116 AliV0ReaderV1 &AliV0ReaderV1::operator=(const AliV0ReaderV1 &ref){
118 // Assignment operator
119 // Only copies pointers, object is not the owner of the references
122 AliAnalysisTaskSE::operator=(ref);
123 fUseImprovedVertex=ref.fUseImprovedVertex;
124 fUseOwnXYZCalculation=ref.fUseOwnXYZCalculation;
125 fUseConstructGamma=ref.fUseConstructGamma;
126 kUseAODConversionPhoton=ref.kUseAODConversionPhoton;
127 fCreateAOD=ref.fCreateAOD;
128 fDeltaAODBranchName=ref.fDeltaAODBranchName;
129 fDeltaAODFilename=ref.fDeltaAODFilename;
130 fEventIsSelected=ref.fEventIsSelected;
135 //________________________________________________________________________
136 void AliV0ReaderV1::Init()
138 // Initialize function to be called once before analysis
139 if(fConversionCuts==NULL){
140 if(fConversionCuts==NULL)AliError("No Cut Selection initialized");
142 if(fCreateAOD){kUseAODConversionPhoton=kTRUE;}
144 if(fConversionGammas != NULL){
145 delete fConversionGammas;
146 fConversionGammas=NULL;
149 if(fConversionGammas == NULL){
150 if(kUseAODConversionPhoton){
151 fConversionGammas = new TClonesArray("AliAODConversionPhoton",100);}
153 fConversionGammas = new TClonesArray("AliKFConversionPhoton",100);}
155 fConversionGammas->Delete();//Reset the TClonesArray
158 //________________________________________________________________________
159 void AliV0ReaderV1::UserCreateOutputObjects()
165 fDeltaAODBranchName.Append("_");
166 fDeltaAODBranchName.Append(fConversionCuts->GetCutNumber());
167 fDeltaAODBranchName.Append("_gamma");
169 fConversionGammas->SetName(fDeltaAODBranchName.Data());
171 AddAODBranch("TClonesArray", &fConversionGammas, fDeltaAODFilename.Data());
172 AliAnalysisManager::GetAnalysisManager()->RegisterExtraFile(fDeltaAODFilename.Data());
176 //________________________________________________________________________
177 Bool_t AliV0ReaderV1::Notify()
179 if (fPeriodName.CompareTo("") == 0){
180 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
182 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
184 TTree* tree = (TTree*) inputHandler->GetTree();
185 TFile* file = (TFile*) tree->GetCurrentFile();
186 TString fileName(file->GetName());
187 TObjArray *arr = fileName.Tokenize("/");
188 for (Int_t i = 0; i < arr->GetEntriesFast();i++ ){
189 TObjString* testObjString = (TObjString*)arr->At(i);
190 if (testObjString->GetString().Contains("LHC")){
191 fPeriodName = testObjString->GetString();
192 i = arr->GetEntriesFast();
197 if(!fConversionCuts->GetDoEtaShift()) return kTRUE; // No Eta Shift requested, continue
198 if(fConversionCuts->GetEtaShift() == 0.0){ // Eta Shift requested but not set, get shift automatically
199 fConversionCuts->GetCorrectEtaShiftFromPeriod(fPeriodName);
200 fConversionCuts->DoEtaShift(kFALSE); // Eta Shift Set, make sure that it is called only once
204 printf(" Gamma Conversion Reader %s :: Eta Shift Manually Set to %f \n\n",
205 (fConversionCuts->GetCutNumber()).Data(),fConversionCuts->GetEtaShift());
206 fConversionCuts->DoEtaShift(kFALSE); // Eta Shift Set, make sure that it is called only once
212 //________________________________________________________________________
213 void AliV0ReaderV1::UserExec(Option_t *option){
215 // Check if correctly initialized
216 if(!fConversionGammas)Init();
219 fEventIsSelected=ProcessEvent(fInputEvent,fMCEvent);
226 //________________________________________________________________________
227 Bool_t AliV0ReaderV1::ProcessEvent(AliVEvent *inputEvent,AliMCEvent *mcEvent)
230 //Reset the TClonesArray
231 fConversionGammas->Delete();
233 fInputEvent=inputEvent;
237 AliError("No Input event");
241 if(!fConversionCuts){AliError("No ConversionCuts");return kFALSE;}
244 if(!fConversionCuts->EventIsSelected(fInputEvent,fMCEvent))return kFALSE;
246 // Set Magnetic Field
247 AliKFParticle::SetField(fInputEvent->GetMagneticField());
249 if(fInputEvent->IsA()==AliESDEvent::Class()){
252 if(fInputEvent->IsA()==AliAODEvent::Class()){
253 GetAODConversionGammas();
258 ///________________________________________________________________________
259 void AliV0ReaderV1::FillAODOutput()
261 // Fill AOD Output with reconstructed Photons
263 if(fInputEvent->IsA()==AliESDEvent::Class()){
264 ///Make sure delta aod is filled if standard aod is filled (for synchronization when reading aod with standard aod)
266 AliAODHandler * aodhandler = dynamic_cast<AliAODHandler*>(AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler());
267 if (aodhandler && aodhandler->GetFillAOD()) {
268 AliAnalysisManager::GetAnalysisManager()->GetOutputEventHandler()->SetFillExtension(kTRUE);
269 //PostData(0, fConversionGammas);
276 ///________________________________________________________________________
277 const AliExternalTrackParam *AliV0ReaderV1::GetExternalTrackParam(AliESDv0 *fCurrentV0,Int_t &tracklabel,Int_t charge){
279 // Get External Track Parameter with given charge
281 if(!(charge==1||charge==-1)){AliError("Charge not defined");return 0x0;}
283 // Check for sign flip
285 if(!fCurrentV0->GetParamN()||!fCurrentV0->GetParamP())return 0x0;
286 if(!fConversionCuts->GetTrack(fInputEvent,fCurrentV0->GetNindex())||!fConversionCuts->GetTrack(fInputEvent,fCurrentV0->GetPindex()))return 0x0;
287 if((fConversionCuts->GetTrack(fInputEvent,fCurrentV0->GetPindex()))->Charge()==charge){
288 tracklabel=fCurrentV0->GetPindex();
289 return fCurrentV0->GetParamP();}
290 if((fConversionCuts->GetTrack(fInputEvent,fCurrentV0->GetNindex()))->Charge()==charge){
291 tracklabel=fCurrentV0->GetNindex();
292 return fCurrentV0->GetParamN();}
297 ///________________________________________________________________________
298 Bool_t AliV0ReaderV1::ProcessESDV0s()
300 // Process ESD V0s for conversion photon reconstruction
301 AliESDEvent *fESDEvent=dynamic_cast<AliESDEvent*>(fInputEvent);
303 AliKFConversionPhoton *fCurrentMotherKFCandidate=NULL;
307 for(Int_t currentV0Index=0;currentV0Index<fESDEvent->GetNumberOfV0s();currentV0Index++){
308 AliESDv0 *fCurrentV0=(AliESDv0*)(fESDEvent->GetV0(currentV0Index));
310 printf("Requested V0 does not exist");
314 fCurrentMotherKFCandidate=ReconstructV0(fCurrentV0,currentV0Index);
316 if(fCurrentMotherKFCandidate){
317 // Add Gamma to the TClonesArray
319 if(kUseAODConversionPhoton){
320 new((*fConversionGammas)[fConversionGammas->GetEntriesFast()]) AliAODConversionPhoton(fCurrentMotherKFCandidate);
321 AliAODConversionPhoton * currentConversionPhoton = (AliAODConversionPhoton*)(fConversionGammas->At(fConversionGammas->GetEntriesFast()-1));
322 currentConversionPhoton->SetMass(fCurrentMotherKFCandidate->M());
323 currentConversionPhoton->SetMassToZero();
327 new((*fConversionGammas)[fConversionGammas->GetEntriesFast()]) AliKFConversionPhoton(*fCurrentMotherKFCandidate);
330 delete fCurrentMotherKFCandidate;
331 fCurrentMotherKFCandidate=NULL;
338 ///________________________________________________________________________
339 AliKFConversionPhoton *AliV0ReaderV1::ReconstructV0(AliESDv0 *fCurrentV0,Int_t currentV0Index)
341 // Reconstruct conversion photon from ESD v0
342 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kPhotonIn);
344 //checks if on the fly mode is set
345 if(!fConversionCuts->SelectV0Finder(fCurrentV0->GetOnFlyStatus())){
346 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kOnFly);
351 Int_t currentTrackLabels[2]={-1,-1};
353 // Get Daughter KF Particles
355 const AliExternalTrackParam *fCurrentExternalTrackParamPositive=GetExternalTrackParamP(fCurrentV0,currentTrackLabels[0]);
356 const AliExternalTrackParam *fCurrentExternalTrackParamNegative=GetExternalTrackParamN(fCurrentV0,currentTrackLabels[1]);
358 if(!fCurrentExternalTrackParamPositive||!fCurrentExternalTrackParamNegative)return 0x0;
360 // Apply some Cuts before Reconstruction
362 AliVTrack * posTrack = fConversionCuts->GetTrack(fInputEvent,currentTrackLabels[0]);
363 AliVTrack * negTrack = fConversionCuts->GetTrack(fInputEvent,currentTrackLabels[1]);
365 if(!negTrack || !posTrack) {
366 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kNoTracks);
371 if(!fConversionCuts->TracksAreSelected(negTrack, posTrack)){
372 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kTrackCuts);
377 if(!fConversionCuts->dEdxCuts(negTrack) || !fConversionCuts->dEdxCuts(posTrack)) {
378 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kdEdxCuts);
382 // Reconstruct Photon
384 AliKFConversionPhoton *fCurrentMotherKF=NULL;
386 AliKFParticle fCurrentNegativeKFParticle(*(fCurrentExternalTrackParamNegative),11);
387 AliKFParticle fCurrentPositiveKFParticle(*(fCurrentExternalTrackParamPositive),-11);
391 if(fUseConstructGamma){
392 fCurrentMotherKF = new AliKFConversionPhoton();
393 fCurrentMotherKF->ConstructGamma(fCurrentNegativeKFParticle,fCurrentPositiveKFParticle);
395 fCurrentMotherKF = new AliKFConversionPhoton(fCurrentNegativeKFParticle,fCurrentPositiveKFParticle);
396 fCurrentMotherKF->SetMassConstraint(0,0.0001);
401 fCurrentMotherKF->SetTrackLabels(currentTrackLabels[0],currentTrackLabels[1]);
405 fCurrentMotherKF->SetV0Index(currentV0Index);
411 AliStack *fMCStack= fMCEvent->Stack();
413 Int_t labelp=TMath::Abs(fConversionCuts->GetTrack(fInputEvent,fCurrentMotherKF->GetTrackLabelPositive())->GetLabel());
414 Int_t labeln=TMath::Abs(fConversionCuts->GetTrack(fInputEvent,fCurrentMotherKF->GetTrackLabelNegative())->GetLabel());
416 TParticle *fNegativeMCParticle = fMCStack->Particle(labeln);
417 TParticle *fPositiveMCParticle = fMCStack->Particle(labelp);
419 if(fPositiveMCParticle&&fNegativeMCParticle){
420 fCurrentMotherKF->SetMCLabelPositive(labelp);
421 fCurrentMotherKF->SetMCLabelNegative(labeln);
425 // Update Vertex (moved for same eta compared to old)
426 if(fUseImprovedVertex == kTRUE){
427 AliKFVertex primaryVertexImproved(*fInputEvent->GetPrimaryVertex());
428 primaryVertexImproved+=*fCurrentMotherKF;
429 fCurrentMotherKF->SetProductionVertex(primaryVertexImproved);
434 Double_t PsiPair=GetPsiPair(fCurrentV0,fCurrentExternalTrackParamPositive,fCurrentExternalTrackParamNegative);
435 fCurrentMotherKF->SetPsiPair(PsiPair);
438 // Recalculate ConversionPoint
439 Double_t dca[2]={0,0};
440 if(fUseOwnXYZCalculation){
441 Double_t convpos[3]={0,0,0};
442 if(!GetConversionPoint(fCurrentExternalTrackParamPositive,fCurrentExternalTrackParamNegative,convpos,dca)){
443 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kConvPointFail);
444 delete fCurrentMotherKF;
445 fCurrentMotherKF=NULL;
449 fCurrentMotherKF->SetConversionPoint(convpos);
452 if(fCurrentMotherKF->GetNDF() > 0.)
453 fCurrentMotherKF->SetChi2perNDF(fCurrentMotherKF->GetChi2()/fCurrentMotherKF->GetNDF()); //->Photon is created before all chi2 relevant changes are performed, set it "by hand"
456 // Set Dilepton Mass (moved down for same eta compared to old)
457 fCurrentMotherKF->SetMass(fCurrentMotherKF->M());
461 if(!fConversionCuts->PhotonCuts(fCurrentMotherKF,fInputEvent)){
462 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kPhotonCuts);
463 delete fCurrentMotherKF;
464 fCurrentMotherKF=NULL;
468 fConversionCuts->FillPhotonCutIndex(AliConversionCuts::kPhotonOut);
469 return fCurrentMotherKF;
472 ///________________________________________________________________________
473 Double_t AliV0ReaderV1::GetPsiPair(const AliESDv0* v0, const AliExternalTrackParam *positiveparam,const AliExternalTrackParam *negativeparam) const {
475 // Angle between daughter momentum plane and plane
478 AliExternalTrackParam nt(*negativeparam);
479 AliExternalTrackParam pt(*positiveparam);
481 Float_t magField = fInputEvent->GetMagneticField();
483 Double_t xyz[3] = {0.,0.,0.};
484 v0->GetXYZ(xyz[0],xyz[1],xyz[2]);
486 // Double_t pPlus[3] = {pt.Px(),pt.Py(),pt.Pz()};
487 // Double_t pMinus[3] = {nt.Px(),nt.Py(),nt.Pz()};
489 // Double_t u[3] = {pPlus[0]+pMinus[0],pPlus[1]+pMinus[1],pPlus[2]+pMinus[2]};
490 // Double_t normu = sqrt( (u[0]*u[0]) + (u[1]*u[1]) + (u[2]*u[2]) );
492 // u[0] = u[0] / normu;
493 // u[1] = u[1] / normu;
494 // u[2] = u[2] / normu;
496 // Double_t normpPlus = sqrt( (pPlus[0]*pPlus[0]) + (pPlus[1]*pPlus[1]) + (pPlus[2]*pPlus[2]) );
497 // Double_t normpMinus = sqrt( (pMinus[0]*pMinus[0]) + (pMinus[1]*pMinus[1]) + (pMinus[2]*pMinus[2]) );
499 // pPlus[0] = pPlus[0] / normpPlus;
500 // pPlus[1] = pPlus[1] / normpPlus;
501 // pPlus[2] = pPlus[2] / normpPlus;
503 // pMinus[0] = pMinus[0] / normpMinus;
504 // pMinus[1] = pMinus[1] / normpMinus;
505 // pMinus[2] = pMinus[2] / normpMinus;
507 // Double_t v[3] = {0,0,0}; // pPlus X pMinus
508 // v[0] = (pPlus[1]*pMinus[2]) - (pPlus[2]*pMinus[1]);
509 // v[1] = (pPlus[2]*pMinus[0]) - (pPlus[0]*pMinus[2]);
510 // v[2] = (pPlus[0]*pMinus[1]) - (pPlus[1]*pMinus[0]);
512 // Double_t w[3] = {0,0,0}; // u X v
513 // w[0] = (u[1]*v[2]) - (u[2]*v[1]);
514 // w[1] = (u[2]*v[0]) - (u[0]*v[2]);
515 // w[2] = (u[0]*v[1]) - (u[1]*v[0]);
517 // Double_t z[3] = {0,0,1};
518 // Double_t wc[3] = {0,0,0}; // u X v
519 // wc[0] = (u[1]*z[2]) - (u[2]*z[1]);
520 // wc[1] = (u[2]*z[0]) - (u[0]*z[2]);
521 // wc[2] = (u[0]*z[1]) - (u[1]*z[0]);
523 // Double_t PhiV = TMath::ACos((w[0]*wc[0]) + (w[1]*wc[1]) + (w[2]*wc[2]));
527 // TVector3 pPlus(pt.Px(),pt.Py(),pt.Pz());
528 // TVector3 pMinus(nt.Px(),nt.Py(),nt.Pz());
530 // TVector3 u = pMinus + pPlus;
531 // u = u*(1/u.Mag());
533 // TVector3 pHPlus = pPlus*(1/pPlus.Mag());
534 // TVector3 pHMinus = pMinus*(1/pMinus.Mag());
536 // TVector3 v = pHPlus.Cross(pHMinus);
537 // TVector3 w = u.Cross(v);
538 // TVector3 z(0,0,1);
539 // TVector3 wc = u.Cross(z);
541 // Double_t PhiV = w * wc;
543 Double_t mn[3] = {0,0,0};
544 Double_t mp[3] = {0,0,0};
546 v0->GetNPxPyPz(mn[0],mn[1],mn[2]);//reconstructed cartesian momentum components of negative daughter;
547 v0->GetPPxPyPz(mp[0],mp[1],mp[2]);//reconstructed cartesian momentum components of positive daughter;
549 Double_t deltat = 1.;
550 deltat = TMath::ATan(mp[2]/(TMath::Sqrt(mp[0]*mp[0] + mp[1]*mp[1])+1.e-13)) - TMath::ATan(mn[2]/(TMath::Sqrt(mn[0]*mn[0] + mn[1]*mn[1])+1.e-13));//difference of angles of the two daughter tracks with z-axis
551 Double_t radiussum = TMath::Sqrt(xyz[0]*xyz[0] + xyz[1]*xyz[1]) + 50;//radius to which tracks shall be propagated
553 Double_t momPosProp[3] = {0,0,0};
554 Double_t momNegProp[3] = {0,0,0};
556 Double_t psiPair = 4.;
557 if(nt.PropagateTo(radiussum,magField) == 0) return psiPair; //propagate tracks to the outside -> Better Purity and Efficiency
558 if(pt.PropagateTo(radiussum,magField) == 0) return psiPair; //propagate tracks to the outside -> Better Purity and Efficiency
560 pt.GetPxPyPz(momPosProp);//Get momentum vectors of tracks after propagation
561 nt.GetPxPyPz(momNegProp);
564 TMath::Sqrt(momNegProp[0]*momNegProp[0]+momNegProp[1]*momNegProp[1]+momNegProp[2]*momNegProp[2]);//absolute momentum value of negative daughter
567 TMath::Sqrt(momPosProp[0]*momPosProp[0]+momPosProp[1]*momPosProp[1]+momPosProp[2]*momPosProp[2]);//absolute momentum value of positive daughter
569 Double_t scalarproduct =
570 momPosProp[0]*momNegProp[0]+momPosProp[1]*momNegProp[1]+momPosProp[2]*momNegProp[2];//scalar product of propagated positive and negative daughters' momenta
572 Double_t chipair = TMath::ACos(scalarproduct/(pEle*pPos));//Angle between propagated daughter tracks
574 psiPair = TMath::Abs(TMath::ASin(deltat/chipair));
579 ///________________________________________________________________________
580 Bool_t AliV0ReaderV1::GetHelixCenter(const AliExternalTrackParam *track,Double_t center[2]){
582 // Get Center of the helix track parametrization
584 Int_t charge=track->Charge();
585 Double_t b=fInputEvent->GetMagneticField();
588 track->GetHelixParameters(helix,b);
590 Double_t xpos = helix[5];
591 Double_t ypos = helix[0];
592 Double_t radius = TMath::Abs(1./helix[4]);
593 Double_t phi = helix[2];
596 phi = phi + 2*TMath::Pi();
599 phi -= TMath::Pi()/2.;
600 Double_t xpoint = radius * TMath::Cos(phi);
601 Double_t ypoint = radius * TMath::Sin(phi);
625 center[0] = xpos + xpoint;
626 center[1] = ypos + ypoint;
630 ///________________________________________________________________________
631 Bool_t AliV0ReaderV1::GetConversionPoint(const AliExternalTrackParam *pparam,const AliExternalTrackParam *nparam,Double_t convpos[3],Double_t dca[2]){
633 // Recalculate Conversion Point
635 if(!pparam||!nparam)return kFALSE;
637 Double_t helixcenterpos[2];
638 GetHelixCenter(pparam,helixcenterpos);
640 Double_t helixcenterneg[2];
641 GetHelixCenter(nparam,helixcenterneg);
643 Double_t helixpos[6];
644 pparam->GetHelixParameters(helixpos,fInputEvent->GetMagneticField());
645 Double_t posradius = TMath::Abs(1./helixpos[4]);
647 Double_t helixneg[6];
648 nparam->GetHelixParameters(helixneg,fInputEvent->GetMagneticField());
649 Double_t negradius = TMath::Abs(1./helixneg[4]);
651 // Calculate xy-position
653 Double_t xpos = helixcenterpos[0];
654 Double_t ypos = helixcenterpos[1];
655 Double_t xneg = helixcenterneg[0];
656 Double_t yneg = helixcenterneg[1];
658 convpos[0] = (xpos*negradius + xneg*posradius)/(negradius+posradius);
659 convpos[1] = (ypos*negradius+ yneg*posradius)/(negradius+posradius);
662 // Calculate Track XY vertex position
664 Double_t deltaXPos = convpos[0] - xpos;
665 Double_t deltaYPos = convpos[1] - ypos;
667 Double_t deltaXNeg = convpos[0] - xneg;
668 Double_t deltaYNeg = convpos[1] - yneg;
670 Double_t alphaPos = TMath::Pi() + TMath::ATan2(-deltaYPos,-deltaXPos);
671 Double_t alphaNeg = TMath::Pi() + TMath::ATan2(-deltaYNeg,-deltaXNeg);
673 Double_t vertexXNeg = xneg + TMath::Abs(negradius)*TMath::Cos(alphaNeg);
674 Double_t vertexYNeg = yneg + TMath::Abs(negradius)*TMath::Sin(alphaNeg);
676 Double_t vertexXPos = xpos + TMath::Abs(posradius)*TMath::Cos(alphaPos);
677 Double_t vertexYPos = ypos + TMath::Abs(posradius)*TMath::Sin(alphaPos);
679 AliExternalTrackParam p(*pparam);
680 AliExternalTrackParam n(*nparam);
682 TVector2 vertexPos(vertexXPos,vertexYPos);
683 TVector2 vertexNeg(vertexXNeg,vertexYNeg);
685 // Convert to local coordinate system
686 TVector2 vertexPosRot=vertexPos.Rotate(-p.GetAlpha());
687 TVector2 vertexNegRot=vertexNeg.Rotate(-n.GetAlpha());
689 // Propagate Track Params to Vertex
691 if(!p.PropagateTo(vertexPosRot.X(),fInputEvent->GetMagneticField()))return kFALSE;
692 if(!n.PropagateTo(vertexNegRot.X(),fInputEvent->GetMagneticField()))return kFALSE;
694 // Check whether propagation was sucessful
696 if(TMath::Abs(vertexPos.Mod()-TMath::Sqrt(p.GetX()*p.GetX()+p.GetY()*p.GetY()))>0.01)return kFALSE;
697 if(TMath::Abs(vertexNeg.Mod()-TMath::Sqrt(n.GetX()*n.GetX()+n.GetY()*n.GetY()))>0.01)return kFALSE;
699 // Calculate z position
701 convpos[2] = (p.GetZ()*negradius+n.GetZ()*posradius)/(negradius+posradius);
704 TVector2 vdca=vertexPos-vertexNeg;
706 dca[1]=TMath::Abs(n.GetZ()-p.GetZ());
710 //________________________________________________________________________
711 Bool_t AliV0ReaderV1::GetAODConversionGammas(){
713 // Get reconstructed conversion photons from satellite AOD file
715 AliAODEvent *fAODEvent=dynamic_cast<AliAODEvent*>(fInputEvent);
719 if(fConversionGammas == NULL){
720 fConversionGammas = new TClonesArray("AliAODConversionPhoton",100);
722 fConversionGammas->Delete();//Reset the TClonesArray
724 //Get Gammas from satellite AOD gamma branch
726 AliAODConversionPhoton *gamma=0x0;
728 TClonesArray *fInputGammas=dynamic_cast<TClonesArray*>(fAODEvent->FindListObject(fDeltaAODBranchName.Data()));
730 FindDeltaAODBranchName();
731 fInputGammas=dynamic_cast<TClonesArray*>(fAODEvent->FindListObject(fDeltaAODBranchName.Data()));}
732 if(!fInputGammas){AliError("No Gamma Satellites found");return kFALSE;}
733 // Apply Selection Cuts to Gammas and create local working copy
735 for(Int_t i=0;i<fInputGammas->GetEntriesFast();i++){
736 gamma=dynamic_cast<AliAODConversionPhoton*>(fInputGammas->At(i));
738 if(fRelabelAODs)RelabelAODPhotonCandidates(gamma);
739 if(fConversionCuts->PhotonIsSelected(gamma,fInputEvent)){
740 new((*fConversionGammas)[fConversionGammas->GetEntriesFast()]) AliAODConversionPhoton(*gamma);}
746 if(fConversionGammas->GetEntries()){return kTRUE;}
751 //________________________________________________________________________
752 void AliV0ReaderV1::FindDeltaAODBranchName(){
754 // Find delta AOD branchname containing reconstructed photons
756 TList *list=fInputEvent->GetList();
757 for(Int_t ii=0;ii<list->GetEntries();ii++){
758 TString name((list->At(ii))->GetName());
759 if(name.BeginsWith(fDeltaAODBranchName)&&name.EndsWith("gamma")){
760 fDeltaAODBranchName=name;
761 AliInfo(Form("Set DeltaAOD BranchName to: %s",fDeltaAODBranchName.Data()));
766 //________________________________________________________________________
767 void AliV0ReaderV1::RelabelAODPhotonCandidates(AliAODConversionPhoton *PhotonCandidate){
769 // Relabeling For AOD Event
771 // MCLabel -> AODMCLabel
772 Bool_t AODLabelPos = kFALSE;
773 Bool_t AODLabelNeg = kFALSE;
775 for(Int_t i = 0; i<fInputEvent->GetNumberOfTracks();i++){
776 AliAODTrack *tempDaughter = static_cast<AliAODTrack*>(fInputEvent->GetTrack(i));
778 if( tempDaughter->GetID() == PhotonCandidate->GetTrackLabelPositive() ){
779 PhotonCandidate->SetMCLabelPositive(abs(tempDaughter->GetLabel()));
780 PhotonCandidate->SetLabelPositive(i);
785 if( tempDaughter->GetID() == PhotonCandidate->GetTrackLabelNegative()){
786 PhotonCandidate->SetMCLabelNegative(abs(tempDaughter->GetLabel()));
787 PhotonCandidate->SetLabelNegative(i);
791 if(AODLabelNeg && AODLabelPos){
795 if(!AODLabelPos || !AODLabelNeg){
796 AliError(Form("NO AOD Daughters Found Pos: %i %i Neg: %i %i",AODLabelPos,PhotonCandidate->GetTrackLabelPositive(),AODLabelNeg,PhotonCandidate->GetTrackLabelNegative()));
801 //________________________________________________________________________
802 void AliV0ReaderV1::Terminate(Option_t *)