1 #include "AliReaderESD.h"
2 //____________________________________________________________________
3 //////////////////////////////////////////////////////////////////////
5 // class AliReaderESD //
7 // reader for ALICE Event Summary Data (ESD). //
9 // Piotr.Skowronski@cern.ch //
11 //////////////////////////////////////////////////////////////////////
15 #include <TObjString.h>
19 #include <TParticle.h>
23 #include <AliRunLoader.h>
25 #include <AliESDtrack.h>
26 #include <AliKalmanTrack.h>
29 #include "AliAnalysis.h"
30 #include "AliAODRun.h"
32 #include "AliAODParticle.h"
33 #include "AliAODParticleCut.h"
34 #include "AliTrackPoints.h"
35 #include "AliClusterMap.h"
37 ClassImp(AliReaderESD)
39 AliReaderESD::AliReaderESD(const Char_t* esdfilename, const Char_t* galfilename):
40 fESDFileName(esdfilename),
41 fGAlFileName(galfilename),
46 fCheckParticlePID(kFALSE),
47 fReadMostProbableOnly(kFALSE),
51 fITSTrackPoints(kFALSE),
53 fReadCentralBarrel(kTRUE),
58 fTPCChi2PerClustMin(0.0),
59 fTPCChi2PerClustMax(10e5),
85 if ( ((Int_t)kNSpecies) != ((Int_t)AliESDtrack::kSPECIES))
86 Fatal("AliReaderESD","ESD defintions probobly changed. Ask Youra.");
88 /********************************************************************/
90 AliReaderESD::AliReaderESD(TObjArray* dirs,const Char_t* esdfilename, const Char_t* galfilename):
92 fESDFileName(esdfilename),
93 fGAlFileName(galfilename),
98 fCheckParticlePID(kFALSE),
99 fReadMostProbableOnly(kFALSE),
103 fITSTrackPoints(kFALSE),
105 fReadCentralBarrel(kTRUE),
110 fTPCChi2PerClustMin(0.0),
111 fTPCChi2PerClustMax(10e5),
136 if ( ((Int_t)kNSpecies) != ((Int_t)AliESDtrack::kSPECIES))
137 Fatal("AliReaderESD","ESD defintions probobly changed. Ask Youra.");
139 /********************************************************************/
141 AliReaderESD::~AliReaderESD()
149 /**********************************************************/
150 Int_t AliReaderESD::ReadNext()
152 //reads next event from fFile
153 //fRunLoader is for reading Kine
155 if (AliVAODParticle::GetDebug())
156 Info("ReadNext","Entered");
158 if (fEventSim == 0x0) fEventSim = new AliAOD();
159 if (fEventRec == 0x0) fEventRec = new AliAOD();
164 do //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
168 fFile = OpenFile(fCurrentDir);//rl is opened here
171 Error("ReadNext","Cannot get fFile for dir no. %d",fCurrentDir);
177 TString esdname = "ESD";
178 esdname+=fCurrentEvent;
179 AliESD* esd = dynamic_cast<AliESD*>(fFile->Get(esdname));
182 if (AliVAODParticle::GetDebug() > 2 )
184 Info("ReadNext","Can not find AliESD object named %s",esdname.Data());
187 delete fFile;//we have to assume there is no more ESD objects in the fFile
198 return 0;//success -> read one event
199 }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array
201 return 1; //no more directories to read
204 /**********************************************************/
205 Int_t AliReaderESD::ReadESD(AliESD* esd)
210 Error("ReadESD","ESD is NULL");
214 // seperate each method
215 if (fReadCentralBarrel) ReadESDCentral(esd);
217 if (fReadMuon) ReadESDMuon(esd);
219 if (fReadPHOS) ReadESDPHOS(esd);
224 /**********************************************************/
225 Int_t AliReaderESD::ReadESDCentral(AliESD* esd)
227 //****** Tentative particle type "concentrations"
228 static const Double_t concentr[5]={0.05, 0., 0.85, 0.10, 0.05};
230 Double_t pidtable[kNSpecies];//array used for reading pid probabilities from ESD track
231 Double_t w[kNSpecies];
232 Double_t mom[3];//momentum
233 Double_t pos[3];//position
234 Double_t vertexpos[3];//vertex position
237 TDatabasePDG* pdgdb = TDatabasePDG::Instance();
240 Error("ReadESD","Can not get PDG Database Instance.");
244 Float_t mf = esd->GetMagneticField()/10.; //AliESD::GetMagnField returns mf in kG
246 if ( (mf == 0.0) && ((fNTrackPoints > 0) || fITSTrackPoints) )
248 Error("ReadESD","Magnetic Field is 0 and Track Points Demended. Skipping to next event.");
254 Info("ReadESD","Magnetic Field is %f",mf);
255 AliKalmanTrack::SetMagneticField(mf);
258 AliStack* stack = 0x0;
259 if (fReadSim && fRunLoader)
261 fRunLoader->GetEvent(fCurrentEvent);
262 stack = fRunLoader->Stack();
265 const AliESDVertex* vertex = esd->GetVertex();
268 Info("ReadESD","ESD returned NULL pointer to vertex - assuming (0.0,0.0,0.0)");
275 vertex->GetXYZ(vertexpos);
278 if (AliVAODParticle::GetDebug() > 0)
280 Info("ReadESD","Primary Vertex is (%f,%f,%f)",vertexpos[0],vertexpos[1],vertexpos[2]);
283 Info("ReadESD","Reading Event %d",fCurrentEvent);
285 Int_t ntr = esd->GetNumberOfTracks();
286 Info("ReadESD","Found %d tracks.",ntr);
287 for (Int_t i = 0;i<ntr; i++)
289 AliESDtrack *esdtrack = esd->GetTrack(i);
292 Error("Next","Can not get track %d", i);
296 //if (esdtrack->HasVertexParameters() == kFALSE)
297 if ((esdtrack->GetStatus() & AliESDtrack::kITSrefit) == kFALSE)
299 if (AliVAODParticle::GetDebug() > 2)
300 Info("ReadNext","Particle skipped: Data at vertex not available.");
306 if ((esdtrack->GetStatus() & AliESDtrack::kTPCin) == kFALSE)
308 if (AliVAODParticle::GetDebug() > 2)
309 Info("ReadNext","Particle skipped: Was not reconstructed in TPC.");
314 if ((esdtrack->GetStatus() & AliESDtrack::kESDpid) == kFALSE)
316 if (AliVAODParticle::GetDebug() > 2)
317 Info("ReadNext","Particle skipped: PID BIT is not set.");
324 esdtrack->GetConstrainedExternalParameters(extx,extp);
327 if (AliVAODParticle::GetDebug() > 2)
328 Info("ReadNext","Track has 0 contrianed curvature -> Probobly parameters never updated. Skipping.");
331 esdtrack->GetESDpid(pidtable);
332 esdtrack->GetConstrainedPxPyPz(mom);
333 esdtrack->GetConstrainedXYZ(pos);
334 pos[0] -= vertexpos[0];//we are interested only in relative position to Primary vertex at this point
335 pos[1] -= vertexpos[1];
336 pos[2] -= vertexpos[2];
338 Int_t charge = (extp[4] > 0)?1:-1;//if curvature=charg/Pt is positive charge is positive
340 //Particle from kinematics
341 AliAODParticle* particle = 0;
342 Bool_t keeppart = kFALSE;
343 if ( fReadSim && stack )
345 if (esdtrack->GetLabel() < 0) continue;//this is fake - we are not able to match any track
346 TParticle *p = stack->Particle(esdtrack->GetLabel());
349 Error("ReadNext","Can not find track with such label.");
353 if (fCheckParticlePID)
355 if(Rejected(p->GetPdgCode()))
357 if ( AliVAODParticle::GetDebug() > 5 )
358 Info("ReadNext","Simulated Particle PID (%d) did not pass the cut.",p->GetPdgCode());
359 continue; //check if we are intersted with particles of this type
362 // if(p->GetPdgCode()<0) charge = -1;
363 particle = new AliAODParticle(*p,i);
367 if(CheckTrack(esdtrack)) continue;
369 //Here we apply Bayes' formula
371 for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) rc+=concentr[s]*pidtable[s];
374 if (AliVAODParticle::GetDebug() > 2)
375 Info("ReadNext","Particle rejected since total bayessian PID probab. is zero.");
379 for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) w[s]=concentr[s]*pidtable[s]/rc;
381 if (AliVAODParticle::GetDebug() > 4)
383 Info("ReadNext","###########################################################################");
384 Info("ReadNext","Momentum: %f %f %f",mom[0],mom[1],mom[2]);
385 Info("ReadNext","Position: %f %f %f",pos[0],pos[1],pos[2]);
386 TString msg("Pid list got from track:");
387 for (Int_t s = 0;s<kNSpecies;s++)
392 msg+=charge*GetSpeciesPdgCode((ESpecies)s);
399 Info("ReadNext","%s",msg.Data());
400 }//if (AliVAODParticle::GetDebug()>4)
402 AliTrackPoints* tpts = 0x0;
403 if (fNTrackPoints > 0)
405 tpts = new AliTrackPoints(fNTrackPoints,esdtrack,mf*10.0,fdR);
406 // tpts->Move(-vertexpos[0],-vertexpos[1],-vertexpos[2]);
409 AliTrackPoints* itstpts = 0x0;
412 itstpts = new AliTrackPoints(AliTrackPoints::kITS,esdtrack,mf*10.0);
413 // itstpts->Move(-vertexpos[0],-vertexpos[1],-vertexpos[2]);
416 AliClusterMap* cmap = 0x0;
419 cmap = new AliClusterMap(esdtrack);
422 //If this flag fReadMostProbableOnly is false the
423 //loop over species (see "LOOPPIDS") is over all possible PIDs
424 //in other case the most probablle PID is searched
425 //and the loop is limited to that PID only
427 Int_t firstspecie = 0;
428 Int_t lastspecie = kNSpecies;
430 if (fReadMostProbableOnly)
432 //find the most probable PID
434 Float_t maxprob = w[0];
435 for (Int_t s=1; s<AliESDtrack::kSPECIES; s++)
444 lastspecie = spec + 1;
447 for (Int_t s = firstspecie; s<lastspecie; s++)//LOOPPIDS
449 Int_t pdgcode = charge*GetSpeciesPdgCode((ESpecies)s);
453 if ( AliVAODParticle::GetDebug() > 5 )
454 Info("ReadNext","Probability of being PID %d is zero. Continuing.",pdgcode);
458 if(Rejected(pdgcode))
460 if ( AliVAODParticle::GetDebug() > 5 )
461 Info("ReadNext","PID (%d) did not pass the cut.",pdgcode);
462 continue; //check if we are intersted with particles of this type
465 Double_t mass = pdgdb->GetParticle(pdgcode)->Mass();
466 Double_t tEtot = TMath::Sqrt( mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2] + mass*mass);//total energy of the track
468 AliAODParticle* track = new AliAODParticle(pdgcode, w[s],i,
469 mom[0], mom[1], mom[2], tEtot,
470 pos[0], pos[1], pos[2], 0.);
471 //copy probabilitis of other species (if not zero)
472 for (Int_t k = 0; k<kNSpecies; k++)
474 if (k == s) continue;
475 if (w[k] == 0.0) continue;
476 track->SetPIDprobability(charge*GetSpeciesPdgCode( (ESpecies)k ),w[k]);
479 if(Rejected(track))//check if meets all criteria of any of our cuts
480 //if it does not delete it and take next good track
482 if ( AliVAODParticle::GetDebug() > 4 )
483 Info("ReadNext","Track did not pass the cut");
488 //Single Particle cuts on cluster map and track points rather do not have sense
491 track->SetTPCTrackPoints(tpts);
496 track->SetITSTrackPoints(itstpts);
501 track->SetClusterMap(cmap);
504 fEventRec->AddParticle(track);
505 if (particle) fEventSim->AddParticle(particle);
508 if (AliVAODParticle::GetDebug() > 4 )
510 Info("ReadNext","\n\nAdding Particle with incarnation %d",pdgcode);
512 if (particle) particle->Print();
513 Info("ReadNext","\n----------------------------------------------\n");
515 }//for (Int_t s = 0; s<kNSpecies; s++)
517 if (keeppart == kFALSE)
519 delete particle;//particle was not stored in event
526 if ( fReadSim && stack )
528 if (particle->P() < 0.00001)
530 Info("ReadNext","###################################");
531 Info("ReadNext","###################################");
532 Info("ReadNext","Track Label %d",esdtrack->GetLabel());
533 TParticle *p = stack->Particle(esdtrack->GetLabel());
542 }//for (Int_t i = 0;i<ntr; i++) -- loop over tracks
544 Info("ReadNext","Read %d tracks and %d particles from event %d (event %d in dir %d).",
545 fEventRec->GetNumberOfParticles(), fEventSim->GetNumberOfParticles(),
546 fNEventsRead,fCurrentEvent,fCurrentDir);
547 fTrackCounter->Fill(fEventRec->GetNumberOfParticles());
549 /******************************************************/
550 /****** Setting glevet properties *************/
551 /******************************************************/
553 if (fEventRec->GetNumberOfParticles() > 0)
555 fEventRec->SetPrimaryVertex(vertexpos[0],vertexpos[1],vertexpos[2]);
560 /**********************************************************/
561 Int_t AliReaderESD::ReadESDMuon(AliESD* esd)
564 Double_t vertexpos[3];//vertex position, assuming no secondary decay
566 const AliESDVertex* vertex = esd->GetVertex();
569 Info("ReadESD","ESD returned NULL pointer to vertex - assuming (0.0,0.0,0.0)");
574 vertex->GetXYZ(vertexpos);
577 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
579 if (AliVAODParticle::GetDebug() > 0) {
580 Info("ReadESD","Reading Event %d",fCurrentEvent);
581 Info("ReadESD","Found %d tracks.",nTracks);
584 Float_t Chi2Cut = 100.;
585 Float_t PtCutMin = 1.;
586 Float_t PtCutMax = 10000.;
587 Float_t muonMass = 0.105658389;
589 Double_t thetaX, thetaY, pYZ;
590 Double_t pxRec1, pyRec1, pzRec1, E1;
598 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
600 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
602 thetaX = muonTrack->GetThetaX();
603 thetaY = muonTrack->GetThetaY();
605 pYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
606 pzRec1 = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
607 pxRec1 = pzRec1 * TMath::Tan(thetaX);
608 pyRec1 = pzRec1 * TMath::Tan(thetaY);
609 charge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
610 E1 = TMath::Sqrt(muonMass * muonMass + pxRec1 * pxRec1 + pyRec1 * pyRec1 + pzRec1 * pzRec1);
611 fV1.SetPxPyPzE(pxRec1, pyRec1, pzRec1, E1);
613 ntrackhits = muonTrack->GetNHit();
614 fitfmin = muonTrack->GetChi2();
616 // transverse momentum
617 Float_t pt1 = fV1.Pt();
620 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
622 if ((ch1 < Chi2Cut) && (pt1 > PtCutMin) && (pt1 < PtCutMax)) {
623 AliAODParticle* track = new AliAODParticle(pdgcode*charge,1,iTrack,
624 pxRec1, pyRec1,pzRec1, E1,
625 vertexpos[0], vertexpos[1], vertexpos[2], 0.);
626 fEventRec->AddParticle(track);
630 fTrackCounter->Fill(fEventRec->GetNumberOfParticles());
634 /**********************************************************/
636 void AliReaderESD::Rewind()
639 // delete fKeyIterator;
642 // fKeyIterator = 0x0;
647 if (fTrackCounter) fTrackCounter->Reset();
649 /**********************************************************/
651 TFile* AliReaderESD::OpenFile(Int_t n)
653 //opens fFile with tree
655 const TString& dirname = GetDirName(n);
658 Error("OpenFiles","Can not get directory name");
661 TString filename = dirname +"/"+ fESDFileName;
662 TFile *ret = TFile::Open(filename.Data());
666 Error("OpenFiles","Can't open fFile %s",filename.Data());
671 Error("OpenFiles","Can't open fFile %s",filename.Data());
677 fRunLoader = AliRunLoader::Open(dirname +"/"+ fGAlFileName);
678 if (fRunLoader == 0x0)
680 Error("OpenFiles","Can't get RunLoader for directory %s",dirname.Data());
685 fRunLoader->LoadHeader();
686 if (fRunLoader->LoadKinematics())
688 Error("Next","Error occured while loading kinematics.");
697 /**********************************************************/
699 Int_t AliReaderESD::GetSpeciesPdgCode(ESpecies spec)//skowron
701 //returns pdg code from the PID index
702 //ask jura about charge
721 ::Warning("GetSpeciesPdgCode","Specie with number %d is not defined.",(Int_t)spec);
726 /********************************************************************/
727 Bool_t AliReaderESD::CheckTrack(AliESDtrack* t) const
729 //Performs check of the track
731 if ( (t->GetConstrainedChi2() < fChi2Min) || (t->GetConstrainedChi2() > fChi2Max) ) return kTRUE;
733 if ( (t->GetTPCclusters(0x0) < fNTPCClustMin) || (t->GetTPCclusters(0x0) > fNTPCClustMax) ) return kTRUE;
735 if (t->GetTPCclusters(0x0) > 0)
737 Float_t chisqpercl = t->GetTPCchi2()/((Double_t)t->GetTPCclusters(0x0));
738 if ( (chisqpercl < fTPCChi2PerClustMin) || (chisqpercl > fTPCChi2PerClustMax) ) return kTRUE;
742 t->GetConstrainedExternalCovariance(cc);
744 if ( (cc[0] < fC00Min) || (cc[0] > fC00Max) ) return kTRUE;
745 if ( (cc[2] < fC11Min) || (cc[2] > fC11Max) ) return kTRUE;
746 if ( (cc[5] < fC22Min) || (cc[5] > fC22Max) ) return kTRUE;
747 if ( (cc[9] < fC33Min) || (cc[9] > fC33Max) ) return kTRUE;
748 if ( (cc[14] < fC44Min) || (cc[14] > fC44Max) ) return kTRUE;
751 t->GetInnerExternalCovariance(cc);
753 if ( (cc[0] < fTPCC00Min) || (cc[0] > fTPCC00Max) ) return kTRUE;
754 if ( (cc[2] < fTPCC11Min) || (cc[2] > fTPCC11Max) ) return kTRUE;
755 if ( (cc[5] < fTPCC22Min) || (cc[5] > fTPCC22Max) ) return kTRUE;
756 if ( (cc[9] < fTPCC33Min) || (cc[9] > fTPCC33Max) ) return kTRUE;
757 if ( (cc[14] < fTPCC44Min) || (cc[14] > fTPCC44Max) ) return kTRUE;
762 /********************************************************************/
764 void AliReaderESD::SetChi2Range(Float_t min, Float_t max)
766 //sets range of Chi2 per Cluster
770 /********************************************************************/
772 void AliReaderESD::SetTPCNClustersRange(Int_t min,Int_t max)
774 //sets range of Number Of Clusters that tracks have to have
778 /********************************************************************/
780 void AliReaderESD::SetTPCChi2PerCluserRange(Float_t min, Float_t max)
782 //sets range of Chi2 per Cluster
783 fTPCChi2PerClustMin = min;
784 fTPCChi2PerClustMax = max;
786 /********************************************************************/
788 void AliReaderESD::SetC00Range(Float_t min, Float_t max)
790 //Sets range of C00 parameter of covariance matrix of the track
791 //it defines uncertainty of the momentum
795 /********************************************************************/
797 void AliReaderESD::SetC11Range(Float_t min, Float_t max)
799 //Sets range of C11 parameter of covariance matrix of the track
800 //it defines uncertainty of the momentum
804 /********************************************************************/
806 void AliReaderESD::SetC22Range(Float_t min, Float_t max)
808 //Sets range of C22 parameter of covariance matrix of the track
809 //it defines uncertainty of the momentum
813 /********************************************************************/
815 void AliReaderESD::SetC33Range(Float_t min, Float_t max)
817 //Sets range of C33 parameter of covariance matrix of the track
818 //it defines uncertainty of the momentum
822 /********************************************************************/
824 void AliReaderESD::SetC44Range(Float_t min, Float_t max)
826 //Sets range of C44 parameter of covariance matrix of the track
827 //it defines uncertainty of the momentum