1 #include "AliHBTReaderESD.h"
2 //____________________________________________________________________
3 //////////////////////////////////////////////////////////////////////
5 // class AliHBTReaderESD //
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>
28 #include "AliHBTRun.h"
29 #include "AliHBTEvent.h"
30 #include "AliHBTParticle.h"
31 #include "AliHBTParticleCut.h"
32 #include "AliHBTTrackPoints.h"
33 #include "AliHBTClusterMap.h"
35 ClassImp(AliHBTReaderESD)
37 AliHBTReaderESD::AliHBTReaderESD(const Char_t* esdfilename, const Char_t* galfilename):
38 fESDFileName(esdfilename),
39 fGAlFileName(galfilename),
43 fReadParticles(kFALSE),
44 fCheckParticlePID(kFALSE),
50 fTPCChi2PerClustMin(0.0),
51 fTPCChi2PerClustMax(10e5),
77 if ( ((Int_t)kNSpecies) != ((Int_t)AliESDtrack::kSPECIES))
78 Fatal("AliHBTReaderESD","ESD defintions probobly changed. Ask Youra.");
80 /********************************************************************/
82 AliHBTReaderESD::AliHBTReaderESD(TObjArray* dirs,const Char_t* esdfilename, const Char_t* galfilename):
84 fESDFileName(esdfilename),
85 fGAlFileName(galfilename),
89 fReadParticles(kFALSE),
90 fCheckParticlePID(kFALSE),
96 fTPCChi2PerClustMin(0.0),
97 fTPCChi2PerClustMax(10e5),
122 if ( ((Int_t)kNSpecies) != ((Int_t)AliESDtrack::kSPECIES))
123 Fatal("AliHBTReaderESD","ESD defintions probobly changed. Ask Youra.");
125 /********************************************************************/
127 AliHBTReaderESD::~AliHBTReaderESD()
134 /**********************************************************/
135 Int_t AliHBTReaderESD::ReadNext()
137 //reads next event from fFile
138 //fRunLoader is for reading Kine
140 if (AliHBTParticle::GetDebug())
141 Info("ReadNext","Entered");
143 if (fParticlesEvent == 0x0) fParticlesEvent = new AliHBTEvent();
144 if (fTracksEvent == 0x0) fTracksEvent = new AliHBTEvent();
146 fParticlesEvent->Reset();
147 fTracksEvent->Reset();
149 do //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
153 fFile = OpenFile(fCurrentDir);//rl is opened here
156 Error("ReadNext","Cannot get fFile for dir no. %d",fCurrentDir);
161 fKeyIterator = new TIter(fFile->GetListOfKeys());
163 // fFile->GetListOfKeys()->Print();
165 TKey* key = (TKey*)fKeyIterator->Next();
168 if (AliHBTParticle::GetDebug() > 2 )
170 Info("ReadNext","No more keys.");
175 delete fFile;//we have to assume there is no more ESD objects in the fFile
184 // TObject* esdobj = key->ReadObj();
185 // if (esdobj == 0x0)
187 // if (AliHBTParticle::GetDebug() > 2 )
189 // Info("ReadNext","Key read NULL. Key Name is %s",key->GetName());
195 // AliESD* esd = dynamic_cast<AliESD*>(esdobj);
197 TString esdname = "ESD";
198 esdname+=fCurrentEvent;
199 AliESD* esd = dynamic_cast<AliESD*>(fFile->Get(esdname));
202 // if (AliHBTParticle::GetDebug() > 2 )
204 // Info("ReadNext","This key is not an AliESD object %s",key->GetName());
206 if (AliHBTParticle::GetDebug() > 2 )
208 Info("ReadNext","Can not find AliESD object named %s",esdname.Data());
213 delete fFile;//we have to assume there is no more ESD objects in the fFile
225 return 0;//success -> read one event
226 }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array
228 return 1; //no more directories to read
230 /**********************************************************/
232 Int_t AliHBTReaderESD::ReadESD(AliESD* esd)
234 //****** Tentative particle type "concentrations"
235 static const Double_t concentr[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
245 Error("ReadESD","ESD is NULL");
249 TDatabasePDG* pdgdb = TDatabasePDG::Instance();
252 Error("ReadESD","Can not get PDG Database Instance.");
256 Float_t mf = esd->GetMagneticField();
258 if ( (mf == 0.0) && (fNTrackPoints > 0) )
260 Error("ReadESD","Magnetic Field is 0 and Track Points Demended. Skipping to next event.");
264 AliStack* stack = 0x0;
265 if (fReadParticles && fRunLoader)
267 fRunLoader->GetEvent(fCurrentEvent);
268 stack = fRunLoader->Stack();
271 const AliESDVertex* vertex = esd->GetVertex();
274 Info("ReadESD","ESD returned NULL pointer to vertex - assuming (0.0,0.0,0.0)");
281 vertex->GetXYZ(vertexpos);
284 if (AliHBTParticle::GetDebug() > 0)
286 Info("ReadESD","Primary Vertex is (%f,%f,%f)",vertexpos[0],vertexpos[1],vertexpos[2]);
289 Info("ReadESD","Reading Event %d",fCurrentEvent);
291 Int_t ntr = esd->GetNumberOfTracks();
292 Info("ReadESD","Found %d tracks.",ntr);
293 for (Int_t i = 0;i<ntr; i++)
295 AliESDtrack *esdtrack = esd->GetTrack(i);
298 Error("Next","Can not get track %d", i);
302 //if (esdtrack->HasVertexParameters() == kFALSE)
303 if ((esdtrack->GetStatus() & AliESDtrack::kITSrefit) == kFALSE)
305 if (AliHBTParticle::GetDebug() > 2)
306 Info("ReadNext","Particle skipped: Data at vertex not available.");
310 if ((esdtrack->GetStatus() & AliESDtrack::kESDpid) == kFALSE)
312 if (AliHBTParticle::GetDebug() > 2)
313 Info("ReadNext","Particle skipped: PID BIT is not set.");
320 esdtrack->GetConstrainedExternalParameters(extx,extp);
323 if (AliHBTParticle::GetDebug() > 2)
324 Info("ReadNext","Track has 0 contrianed curvature -> Probobly parameters never updated. Skipping.");
327 esdtrack->GetESDpid(pidtable);
328 esdtrack->GetConstrainedPxPyPz(mom);
329 esdtrack->GetConstrainedXYZ(pos);
330 pos[0] -= vertexpos[0];//we are interested only in relative position to Primary vertex at this point
331 pos[1] -= vertexpos[1];
332 pos[2] -= vertexpos[2];
334 Int_t charge = (extp[4] > 0)?1:-1;//if curvature=charg/Pt is positive charge is positive
336 //Particle from kinematics
337 AliHBTParticle* particle = 0;
338 Bool_t keeppart = kFALSE;
339 if ( fReadParticles && stack )
341 if (esdtrack->GetLabel() < 0) continue;//this is fake - we are not able to match any track
342 TParticle *p = stack->Particle(esdtrack->GetLabel());
345 Error("ReadNext","Can not find track with such label.");
348 if(Pass(p->GetPdgCode()))
350 if ( AliHBTParticle::GetDebug() > 5 )
351 Info("ReadNext","Simulated Particle PID (%d) did not pass the cut.",p->GetPdgCode());
352 continue; //check if we are intersted with particles of this type
354 // if(p->GetPdgCode()<0) charge = -1;
355 particle = new AliHBTParticle(*p,i);
359 if(CheckTrack(esdtrack)) continue;
361 //Here we apply Bayes' formula
363 for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) rc+=concentr[s]*pidtable[s];
366 if (AliHBTParticle::GetDebug() > 2)
367 Info("ReadNext","Particle rejected since total bayessian PID probab. is zero.");
371 for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) w[s]=concentr[s]*pidtable[s]/rc;
373 if (AliHBTParticle::GetDebug() > 4)
375 Info("ReadNext","###########################################################################");
376 Info("ReadNext","Momentum: %f %f %f",mom[0],mom[1],mom[2]);
377 Info("ReadNext","Position: %f %f %f",pos[0],pos[1],pos[2]);
378 TString msg("Pid list got from track:");
379 for (Int_t s = 0;s<kNSpecies;s++)
384 msg+=charge*GetSpeciesPdgCode((ESpecies)s);
391 Info("ReadNext","%s",msg.Data());
392 }//if (AliHBTParticle::GetDebug()>4)
394 AliHBTTrackPoints* tpts = 0x0;
395 if (fNTrackPoints > 0)
397 tpts = new AliHBTTrackPoints(fNTrackPoints,esdtrack,mf,fdR);
400 AliHBTClusterMap* cmap = 0x0;
403 cmap = new AliHBTClusterMap(esdtrack);
406 for (Int_t s = 0; s<kNSpecies; s++)
408 Int_t pdgcode = charge*GetSpeciesPdgCode((ESpecies)s);
412 if ( AliHBTParticle::GetDebug() > 5 )
413 Info("ReadNext","Probability of being PID %d is zero. Continuing.",pdgcode);
419 if ( AliHBTParticle::GetDebug() > 5 )
420 Info("ReadNext","PID (%d) did not pass the cut.",pdgcode);
421 continue; //check if we are intersted with particles of this type
424 Double_t mass = pdgdb->GetParticle(pdgcode)->Mass();
425 Double_t tEtot = TMath::Sqrt( mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2] + mass*mass);//total energy of the track
427 AliHBTParticle* track = new AliHBTParticle(pdgcode, w[s],i,
428 mom[0], mom[1], mom[2], tEtot,
429 pos[0], pos[1], pos[2], 0.);
430 //copy probabilitis of other species (if not zero)
431 for (Int_t k = 0; k<kNSpecies; k++)
433 if (k == s) continue;
434 if (w[k] == 0.0) continue;
435 track->SetPIDprobability(charge*GetSpeciesPdgCode( (ESpecies)k ),w[k]);
438 if(Pass(track))//check if meets all criteria of any of our cuts
439 //if it does not delete it and take next good track
441 if ( AliHBTParticle::GetDebug() > 4 )
442 Info("ReadNext","Track did not pass the cut");
447 //Single Particle cuts on cluster map and track points rather do not have sense
450 track->SetTrackPoints(tpts);
455 track->SetClusterMap(cmap);
458 fTracksEvent->AddParticle(track);
459 if (particle) fParticlesEvent->AddParticle(particle);
462 if (AliHBTParticle::GetDebug() > 4 )
464 Info("ReadNext","\n\nAdding Particle with incarnation %d",pdgcode);
466 if (particle) particle->Print();
467 Info("ReadNext","\n----------------------------------------------\n");
469 }//for (Int_t s = 0; s<kNSpecies; s++)
471 if (keeppart == kFALSE)
473 delete particle;//particle was not stored in event
478 }//for (Int_t i = 0;i<ntr; i++) -- loop over tracks
480 Info("ReadNext","Read %d tracks and %d particles from event %d (event %d in dir %d).",
481 fTracksEvent->GetNumberOfParticles(), fParticlesEvent->GetNumberOfParticles(),
482 fNEventsRead,fCurrentEvent,fCurrentDir);
483 fTrackCounter->Fill(fTracksEvent->GetNumberOfParticles());
487 /**********************************************************/
489 void AliHBTReaderESD::Rewind()
500 if (fTrackCounter) fTrackCounter->Reset();
502 /**********************************************************/
504 TFile* AliHBTReaderESD::OpenFile(Int_t n)
506 //opens fFile with kine tree
508 const TString& dirname = GetDirName(n);
511 Error("OpenFiles","Can not get directory name");
514 TString filename = dirname +"/"+ fESDFileName;
515 TFile *ret = TFile::Open(filename.Data());
519 Error("OpenFiles","Can't open fFile %s",filename.Data());
524 Error("OpenFiles","Can't open fFile %s",filename.Data());
530 fRunLoader = AliRunLoader::Open(dirname +"/"+ fGAlFileName);
531 if (fRunLoader == 0x0)
533 Error("OpenFiles","Can't get RunLoader for directory %s",dirname.Data());
538 fRunLoader->LoadHeader();
539 if (fRunLoader->LoadKinematics())
541 Error("Next","Error occured while loading kinematics.");
550 /**********************************************************/
552 Int_t AliHBTReaderESD::GetSpeciesPdgCode(ESpecies spec)//skowron
554 //returns pdg code from the PID index
555 //ask jura about charge
574 ::Warning("GetSpeciesPdgCode","Specie with number %d is not defined.",(Int_t)spec);
579 /********************************************************************/
580 Bool_t AliHBTReaderESD::CheckTrack(AliESDtrack* t) const
582 //Performs check of the track
584 if ( (t->GetConstrainedChi2() < fChi2Min) || (t->GetConstrainedChi2() > fChi2Max) ) return kTRUE;
586 if ( (t->GetTPCclusters(0x0) < fNTPCClustMin) || (t->GetTPCclusters(0x0) > fNTPCClustMax) ) return kTRUE;
588 if (t->GetTPCclusters(0x0) > 0)
590 Float_t chisqpercl = t->GetTPCchi2()/((Double_t)t->GetTPCclusters(0x0));
591 if ( (chisqpercl < fTPCChi2PerClustMin) || (chisqpercl > fTPCChi2PerClustMax) ) return kTRUE;
595 t->GetConstrainedExternalCovariance(cc);
597 if ( (cc[0] < fC00Min) || (cc[0] > fC00Max) ) return kTRUE;
598 if ( (cc[2] < fC11Min) || (cc[2] > fC11Max) ) return kTRUE;
599 if ( (cc[5] < fC22Min) || (cc[5] > fC22Max) ) return kTRUE;
600 if ( (cc[9] < fC33Min) || (cc[9] > fC33Max) ) return kTRUE;
601 if ( (cc[14] < fC44Min) || (cc[14] > fC44Max) ) return kTRUE;
604 t->GetInnerExternalCovariance(cc);
606 if ( (cc[0] < fTPCC00Min) || (cc[0] > fTPCC00Max) ) return kTRUE;
607 if ( (cc[2] < fTPCC11Min) || (cc[2] > fTPCC11Max) ) return kTRUE;
608 if ( (cc[5] < fTPCC22Min) || (cc[5] > fTPCC22Max) ) return kTRUE;
609 if ( (cc[9] < fTPCC33Min) || (cc[9] > fTPCC33Max) ) return kTRUE;
610 if ( (cc[14] < fTPCC44Min) || (cc[14] > fTPCC44Max) ) return kTRUE;
615 /********************************************************************/
617 void AliHBTReaderESD::SetChi2(Float_t min, Float_t max)
619 //sets range of Chi2 per Cluster
623 /********************************************************************/
625 void AliHBTReaderESD::SetTPCNClustersRange(Int_t min,Int_t max)
627 //sets range of Number Of Clusters that tracks have to have
631 /********************************************************************/
633 void AliHBTReaderESD::SetTPCChi2PerCluserRange(Float_t min, Float_t max)
635 //sets range of Chi2 per Cluster
636 fTPCChi2PerClustMin = min;
637 fTPCChi2PerClustMax = max;
639 /********************************************************************/
641 void AliHBTReaderESD::SetC00Range(Float_t min, Float_t max)
643 //Sets range of C00 parameter of covariance matrix of the track
644 //it defines uncertainty of the momentum
648 /********************************************************************/
650 void AliHBTReaderESD::SetC11Range(Float_t min, Float_t max)
652 //Sets range of C11 parameter of covariance matrix of the track
653 //it defines uncertainty of the momentum
657 /********************************************************************/
659 void AliHBTReaderESD::SetC22Range(Float_t min, Float_t max)
661 //Sets range of C22 parameter of covariance matrix of the track
662 //it defines uncertainty of the momentum
666 /********************************************************************/
668 void AliHBTReaderESD::SetC33Range(Float_t min, Float_t max)
670 //Sets range of C33 parameter of covariance matrix of the track
671 //it defines uncertainty of the momentum
675 /********************************************************************/
677 void AliHBTReaderESD::SetC44Range(Float_t min, Float_t max)
679 //Sets range of C44 parameter of covariance matrix of the track
680 //it defines uncertainty of the momentum