1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //____________________________________________________________________
19 //////////////////////////////////////////////////////////////////////
21 // class AliReaderESD //
23 // reader for ALICE Event Summary Data (ESD). //
25 // Piotr.Skowronski@cern.ch //
27 //////////////////////////////////////////////////////////////////////
30 #include <TGliteXmlEventlist.h>
33 #include <TParticle.h>
37 #include "AliAODParticle.h"
38 #include "AliClusterMap.h"
40 #include "AliESDtrack.h"
42 #include "AliReaderESD.h"
43 #include "AliRunLoader.h"
46 ClassImp(AliReaderESD)
48 AliReaderESD::AliReaderESD(const Char_t* esdfilename, const Char_t* galfilename):
49 fESDFileName(esdfilename),
50 fGAlFileName(galfilename),
55 fCheckParticlePID(kFALSE),
56 fReadMostProbableOnly(kFALSE),
60 fITSTrackPoints(kFALSE),
61 fITSTrackPointsType(AliTrackPoints::kITS),
63 fReadCentralBarrel(kTRUE),
68 fTPCChi2PerClustMin(0.0),
69 fTPCChi2PerClustMax(10e5),
95 if ( ((Int_t)kNSpecies) != ((Int_t)AliPID::kSPECIES))
96 Fatal("AliReaderESD","ESD defintions probobly changed. Ask Youra.");
98 /********************************************************************/
100 AliReaderESD::AliReaderESD(TObjArray* dirs,const Char_t* esdfilename, const Char_t* galfilename):
102 fESDFileName(esdfilename),
103 fGAlFileName(galfilename),
108 fCheckParticlePID(kFALSE),
109 fReadMostProbableOnly(kFALSE),
113 fITSTrackPoints(kFALSE),
114 fITSTrackPointsType(AliTrackPoints::kITS),
116 fReadCentralBarrel(kTRUE),
121 fTPCChi2PerClustMin(0.0),
122 fTPCChi2PerClustMax(10e5),
147 if ( ((Int_t)kNSpecies) != ((Int_t)AliPID::kSPECIES))
148 Fatal("AliReaderESD","ESD defintions probobly changed. Ask Youra.");
150 /********************************************************************/
152 AliReaderESD::~AliReaderESD()
160 /**********************************************************/
161 Int_t AliReaderESD::ReadNext()
163 //reads next event from fFile
164 //fRunLoader is for reading Kine
166 AliDebug(1,"Entered");
168 if (fEventSim == 0x0) fEventSim = new AliAOD();
169 if (fEventRec == 0x0) fEventRec = new AliAOD();
174 do //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
178 fFile = OpenFile(fCurrentDir);//rl is opened here
181 Error("ReadNext","Cannot get fFile for dir no. %d",fCurrentDir);
187 TString esdname = "ESD";
188 esdname+=fCurrentEvent;
189 AliESD* esd = dynamic_cast<AliESD*>(fFile->Get(esdname));
192 AliDebug(3,Form("Can not find AliESD object named %s",esdname.Data()));
194 delete fFile;//we have to assume there is no more ESD objects in the fFile
205 return 0;//success -> read one event
206 }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array
208 return 1; //no more directories to read
211 /**********************************************************/
212 Int_t AliReaderESD::ReadESD(AliESD* esd)
217 Error("ReadESD","ESD is NULL");
221 // seperate each method
222 if (fReadCentralBarrel) ReadESDCentral(esd);
224 if (fReadMuon) ReadESDMuon(esd);
226 if (fReadPHOS) ReadESDPHOS(esd);
231 /**********************************************************/
232 Int_t AliReaderESD::ReadESDCentral(AliESD* esd)
234 //****** Tentative particle type "concentrations"
235 static const Double_t kConcentr[5]={0.05, 0., 0.85, 0.10, 0.05};
237 Double_t pidtable[kNSpecies];//array used for reading pid probabilities from ESD track
238 Double_t w[kNSpecies];
239 Double_t mom[3];//momentum
240 Double_t pos[3];//position
241 Double_t vertexpos[3];//vertex position
244 TDatabasePDG* pdgdb = TDatabasePDG::Instance();
247 Error("ReadESD","Can not get PDG Database Instance.");
251 Float_t mf = esd->GetMagneticField()/10.; //AliESD::GetMagnField returns mf in kG
253 if ( (mf == 0.0) && ((fNTrackPoints > 0) || fITSTrackPoints) )
255 Error("ReadESD","Magnetic Field is 0 and Track Points Demended. Skipping to next event.");
261 Info("ReadESD","Magnetic Field is %f",mf);
262 //AliKalmanTrack::SetMagneticField(mf);
265 AliStack* stack = 0x0;
266 if (fReadSim && fRunLoader)
268 fRunLoader->GetEvent(fCurrentEvent);
269 stack = fRunLoader->Stack();
272 const AliESDVertex* vertex = esd->GetVertex();
275 Info("ReadESD","ESD returned NULL pointer to vertex - assuming (0.0,0.0,0.0)");
282 vertex->GetXYZ(vertexpos);
285 AliDebug(1,Form("Primary Vertex is (%f,%f,%f)",vertexpos[0],vertexpos[1],vertexpos[2]));
287 Info("ReadESD","Reading Event %d",fCurrentEvent);
289 Int_t ntr = esd->GetNumberOfTracks();
290 Info("ReadESD","Found %d tracks.",ntr);
291 for (Int_t i = 0;i<ntr; i++)
293 AliESDtrack *esdtrack = esd->GetTrack(i);
296 Error("Next","Can not get track %d", i);
300 //if (esdtrack->HasVertexParameters() == kFALSE)
301 if ((esdtrack->GetStatus() & AliESDtrack::kITSrefit) == kFALSE)
303 AliDebug(3,Form("Particle skipped: Data at vertex not available."));
309 if ((esdtrack->GetStatus() & AliESDtrack::kTPCin) == kFALSE)
311 AliDebug(3,"Particle skipped: Was not reconstructed in TPC.");
316 if ((esdtrack->GetStatus() & AliESDtrack::kESDpid) == kFALSE)
318 AliDebug(3,"Particle skipped: PID BIT is not set.");
325 esdtrack->GetConstrainedExternalParameters(alpha,extx,extp);
328 AliDebug(3,"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 AliDebug(6,Form("Simulated Particle PID (%d) did not pass the cut.",p->GetPdgCode()));
358 continue; //check if we are intersted with particles of this type
361 // if(p->GetPdgCode()<0) charge = -1;
362 particle = new AliAODParticle(*p,i);
366 if(CheckTrack(esdtrack)) continue;
368 //Here we apply Bayes' formula
370 for (Int_t s=0; s<AliPID::kSPECIES; s++) rc+=kConcentr[s]*pidtable[s];
373 AliDebug(3,"Particle rejected since total bayessian PID probab. is zero.");
377 for (Int_t s=0; s<AliPID::kSPECIES; s++) w[s]=kConcentr[s]*pidtable[s]/rc;
379 if (AliDebugLevel() > 4)
381 Info("ReadNext","###########################################################################");
382 Info("ReadNext","Momentum: %f %f %f",mom[0],mom[1],mom[2]);
383 Info("ReadNext","Position: %f %f %f",pos[0],pos[1],pos[2]);
384 TString msg("Pid list got from track:");
385 for (Int_t s = 0;s<kNSpecies;s++)
390 msg+=charge*GetSpeciesPdgCode((ESpecies)s);
397 Info("ReadNext","%s",msg.Data());
398 }//if (AliDebugLevel()>4)
400 AliTrackPoints* tpts = 0x0;
401 if (fNTrackPoints > 0)
403 tpts = new AliTrackPoints(fNTrackPoints,esdtrack,mf*10.0,fdR);
404 // tpts->Move(-vertexpos[0],-vertexpos[1],-vertexpos[2]);
407 AliTrackPoints* itstpts = 0x0;
410 itstpts = new AliTrackPoints(fITSTrackPointsType,esdtrack,mf*10.0);
411 // itstpts->Move(-vertexpos[0],-vertexpos[1],-vertexpos[2]);
414 AliClusterMap* cmap = 0x0;
417 cmap = new AliClusterMap(esdtrack);
420 //If this flag fReadMostProbableOnly is false the
421 //loop over species (see "LOOPPIDS") is over all possible PIDs
422 //in other case the most probablle PID is searched
423 //and the loop is limited to that PID only
425 Int_t firstspecie = 0;
426 Int_t lastspecie = kNSpecies;
428 if (fReadMostProbableOnly)
430 //find the most probable PID
432 Float_t maxprob = w[0];
433 for (Int_t s=1; s<AliPID::kSPECIES; s++)
442 lastspecie = spec + 1;
445 for (Int_t s = firstspecie; s<lastspecie; s++)//LOOPPIDS
447 Int_t pdgcode = charge*GetSpeciesPdgCode((ESpecies)s);
451 AliDebug(6,Form("Probability of being PID %d is zero. Continuing.",pdgcode));
455 if(Rejected(pdgcode))
457 AliDebug(6,Form("PID (%d) did not pass the cut.",pdgcode));
458 continue; //check if we are intersted with particles of this type
461 Double_t mass = pdgdb->GetParticle(pdgcode)->Mass();
462 Double_t tEtot = TMath::Sqrt( mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2] + mass*mass);//total energy of the track
464 AliAODParticle* track = new AliAODParticle(pdgcode, w[s],i,
465 mom[0], mom[1], mom[2], tEtot,
466 pos[0], pos[1], pos[2], 0.);
467 //copy probabilitis of other species (if not zero)
468 for (Int_t k = 0; k<kNSpecies; k++)
470 if (k == s) continue;
471 if (w[k] == 0.0) continue;
472 track->SetPIDprobability(charge*GetSpeciesPdgCode( (ESpecies)k ),w[k]);
475 if(Rejected(track))//check if meets all criteria of any of our cuts
476 //if it does not delete it and take next good track
478 AliDebug(5,"Track did not pass the cut");
483 //Single Particle cuts on cluster map and track points rather do not have sense
486 track->SetTPCTrackPoints(tpts);
491 track->SetITSTrackPoints(itstpts);
496 track->SetClusterMap(cmap);
499 fEventRec->AddParticle(track);
500 if (particle) fEventSim->AddParticle(particle);
503 if (AliDebugLevel() > 4 )
505 Info("ReadNext","\n\nAdding Particle with incarnation %d",pdgcode);
507 if (particle) particle->Print();
508 Info("ReadNext","\n----------------------------------------------\n");
510 }//for (Int_t s = 0; s<kNSpecies; s++)
512 if (keeppart == kFALSE)
514 delete particle;//particle was not stored in event
521 if ( fReadSim && stack )
523 if (particle->P() < 0.00001)
525 Info("ReadNext","###################################");
526 Info("ReadNext","###################################");
527 Info("ReadNext","Track Label %d",esdtrack->GetLabel());
528 TParticle *p = stack->Particle(esdtrack->GetLabel());
537 }//for (Int_t i = 0;i<ntr; i++) -- loop over tracks
539 Info("ReadNext","Read %d tracks and %d particles from event %d (event %d in dir %d).",
540 fEventRec->GetNumberOfParticles(), fEventSim->GetNumberOfParticles(),
541 fNEventsRead,fCurrentEvent,fCurrentDir);
542 fTrackCounter->Fill(fEventRec->GetNumberOfParticles());
544 /******************************************************/
545 /****** Setting glevet properties *************/
546 /******************************************************/
548 if (fEventRec->GetNumberOfParticles() > 0)
550 fEventRec->SetPrimaryVertex(vertexpos[0],vertexpos[1],vertexpos[2]);
555 /**********************************************************/
556 Int_t AliReaderESD::ReadESDMuon(AliESD* esd)
558 // Reads the muon tracks from the ESD
559 Double_t vertexpos[3];//vertex position, assuming no secondary decay
561 const AliESDVertex* vertex = esd->GetVertex();
564 Info("ReadESD","ESD returned NULL pointer to vertex - assuming (0.0,0.0,0.0)");
569 vertex->GetXYZ(vertexpos);
572 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
574 AliDebug(1,Form("Reading Event %d \nFound %d tracks.",fCurrentEvent,nTracks));
577 Float_t chi2Cut = 100.;
578 Float_t ptCutMin = 1.;
579 Float_t ptCutMax = 10000.;
580 Float_t muonMass = 0.105658389;
582 Double_t thetaX, thetaY, pYZ;
583 Double_t pxRec1, pyRec1, pzRec1, e1;
591 for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
593 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
595 thetaX = muonTrack->GetThetaX();
596 thetaY = muonTrack->GetThetaY();
598 pYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
599 pzRec1 = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaX));
600 pxRec1 = pzRec1 * TMath::Tan(thetaX);
601 pyRec1 = pzRec1 * TMath::Tan(thetaY);
602 charge = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
603 e1 = TMath::Sqrt(muonMass * muonMass + pxRec1 * pxRec1 + pyRec1 * pyRec1 + pzRec1 * pzRec1);
604 fV1.SetPxPyPzE(pxRec1, pyRec1, pzRec1, e1);
606 ntrackhits = muonTrack->GetNHit();
607 fitfmin = muonTrack->GetChi2();
609 // transverse momentum
610 Float_t pt1 = fV1.Pt();
613 Float_t ch1 = fitfmin / (2.0 * ntrackhits - 5);
615 if ((ch1 < chi2Cut) && (pt1 > ptCutMin) && (pt1 < ptCutMax)) {
616 AliAODParticle* track = new AliAODParticle(pdgcode*charge,1,iTrack,
617 pxRec1, pyRec1,pzRec1, e1,
618 vertexpos[0], vertexpos[1], vertexpos[2], 0.);
619 fEventRec->AddParticle(track);
623 fTrackCounter->Fill(fEventRec->GetNumberOfParticles());
627 /**********************************************************/
629 void AliReaderESD::Rewind()
632 // delete fKeyIterator;
635 // fKeyIterator = 0x0;
640 if (fEventList) fEventList->Reset();
641 if (fTrackCounter) fTrackCounter->Reset();
643 /**********************************************************/
645 TFile* AliReaderESD::OpenFile(Int_t n)
647 //opens fFile with tree
656 while (fCurrentDir < n)
664 const TString& dirname = GetDirName(n);
667 Error("OpenFiles","Can not get directory name");
674 filename = fEventList->GetURL(fESDFileName);
678 filename = dirname +"/"+ fESDFileName;
680 Info("OpenFile","%s ==> %s",fESDFileName.Data(),filename.Data());
681 TFile *ret = TFile::Open(filename.Data());
685 Error("OpenFiles","Can't open fFile %s",filename.Data());
690 Error("OpenFiles","Can't open fFile %s",filename.Data());
699 gafilename = fEventList->GetURL(fGAlFileName);
703 gafilename = dirname +"/"+ fGAlFileName;
705 Info("OpenFile","%s ==> %s",fGAlFileName.Data(),gafilename.Data());
707 fRunLoader = AliRunLoader::Open(gafilename);
709 if (fRunLoader == 0x0)
711 Error("OpenFiles","Can't get RunLoader for directory %s",dirname.Data());
716 fRunLoader->LoadHeader();
720 TString kinefilename = fEventList->GetURL("Kinematics.root");
721 fRunLoader->SetKineFileName(kinefilename);
722 Info("OpenFile","%s ==> %s","Kinematics.root",kinefilename.Data());
725 if (fRunLoader->LoadKinematics())
727 Error("Next","Error occured while loading kinematics.");
736 /**********************************************************/
738 Int_t AliReaderESD::GetSpeciesPdgCode(ESpecies spec)//skowron
740 //returns pdg code from the PID index
741 //ask jura about charge
760 ::Warning("GetSpeciesPdgCode","Specie with number %d is not defined.",(Int_t)spec);
765 /********************************************************************/
766 Bool_t AliReaderESD::CheckTrack(AliESDtrack* t) const
768 //Performs check of the track
770 if ( (t->GetConstrainedChi2() < fChi2Min) || (t->GetConstrainedChi2() > fChi2Max) ) return kTRUE;
772 if ( (t->GetTPCclusters(0x0) < fNTPCClustMin) || (t->GetTPCclusters(0x0) > fNTPCClustMax) ) return kTRUE;
774 if (t->GetTPCclusters(0x0) > 0)
776 Float_t chisqpercl = t->GetTPCchi2()/((Double_t)t->GetTPCclusters(0x0));
777 if ( (chisqpercl < fTPCChi2PerClustMin) || (chisqpercl > fTPCChi2PerClustMax) ) return kTRUE;
781 t->GetConstrainedExternalCovariance(cc);
783 if ( (cc[0] < fC00Min) || (cc[0] > fC00Max) ) return kTRUE;
784 if ( (cc[2] < fC11Min) || (cc[2] > fC11Max) ) return kTRUE;
785 if ( (cc[5] < fC22Min) || (cc[5] > fC22Max) ) return kTRUE;
786 if ( (cc[9] < fC33Min) || (cc[9] > fC33Max) ) return kTRUE;
787 if ( (cc[14] < fC44Min) || (cc[14] > fC44Max) ) return kTRUE;
790 t->GetInnerExternalCovariance(cc);
792 if ( (cc[0] < fTPCC00Min) || (cc[0] > fTPCC00Max) ) return kTRUE;
793 if ( (cc[2] < fTPCC11Min) || (cc[2] > fTPCC11Max) ) return kTRUE;
794 if ( (cc[5] < fTPCC22Min) || (cc[5] > fTPCC22Max) ) return kTRUE;
795 if ( (cc[9] < fTPCC33Min) || (cc[9] > fTPCC33Max) ) return kTRUE;
796 if ( (cc[14] < fTPCC44Min) || (cc[14] > fTPCC44Max) ) return kTRUE;
801 /********************************************************************/
803 void AliReaderESD::SetChi2Range(Float_t min, Float_t max)
805 //sets range of Chi2 per Cluster
809 /********************************************************************/
811 void AliReaderESD::SetTPCNClustersRange(Int_t min,Int_t max)
813 //sets range of Number Of Clusters that tracks have to have
817 /********************************************************************/
819 void AliReaderESD::SetTPCChi2PerCluserRange(Float_t min, Float_t max)
821 //sets range of Chi2 per Cluster
822 fTPCChi2PerClustMin = min;
823 fTPCChi2PerClustMax = max;
825 /********************************************************************/
827 void AliReaderESD::SetC00Range(Float_t min, Float_t max)
829 //Sets range of C00 parameter of covariance matrix of the track
830 //it defines uncertainty of the momentum
834 /********************************************************************/
836 void AliReaderESD::SetC11Range(Float_t min, Float_t max)
838 //Sets range of C11 parameter of covariance matrix of the track
839 //it defines uncertainty of the momentum
843 /********************************************************************/
845 void AliReaderESD::SetC22Range(Float_t min, Float_t max)
847 //Sets range of C22 parameter of covariance matrix of the track
848 //it defines uncertainty of the momentum
852 /********************************************************************/
854 void AliReaderESD::SetC33Range(Float_t min, Float_t max)
856 //Sets range of C33 parameter of covariance matrix of the track
857 //it defines uncertainty of the momentum
861 /********************************************************************/
863 void AliReaderESD::SetC44Range(Float_t min, Float_t max)
865 //Sets range of C44 parameter of covariance matrix of the track
866 //it defines uncertainty of the momentum