]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTReaderESD.cxx
Default constructor added to ITS separation cut
[u/mrichter/AliRoot.git] / HBTAN / AliHBTReaderESD.cxx
1 #include "AliHBTReaderESD.h"
2 //____________________________________________________________________
3 //////////////////////////////////////////////////////////////////////
4 //                                                                  //
5 // class AliHBTReaderESD                                            //
6 //                                                                  //
7 // reader for ALICE Event Summary Data (ESD).                       //
8 //                                                                  //
9 // Piotr.Skowronski@cern.ch                                         //
10 //                                                                  //
11 //////////////////////////////////////////////////////////////////////
12
13 #include <Riostream.h>
14
15 #include <TPDGCode.h>
16 #include <TString.h>
17 #include <TObjString.h>
18 #include <TTree.h>
19 #include <TFile.h>
20 #include <TKey.h>
21 #include <TParticle.h>
22 #include <TH1.h>
23
24 #include <AliRun.h>
25 #include <AliRunLoader.h>
26 #include <AliStack.h>
27 #include <AliESDtrack.h>
28 #include <AliKalmanTrack.h>
29 #include <AliESD.h>
30
31 #include "AliHBTRun.h"
32 #include "AliHBTEvent.h"
33 #include "AliHBTParticle.h"
34 #include "AliHBTParticleCut.h"
35 #include "AliHBTTrackPoints.h"
36 #include "AliHBTClusterMap.h"
37
38 ClassImp(AliHBTReaderESD)
39
40 AliHBTReaderESD::AliHBTReaderESD(const Char_t* esdfilename, const Char_t* galfilename):
41  fESDFileName(esdfilename),
42  fGAlFileName(galfilename),
43  fFile(0x0),
44  fRunLoader(0x0),
45  fKeyIterator(0x0),
46  fReadParticles(kFALSE),
47  fCheckParticlePID(kFALSE),
48  fNTrackPoints(0),
49  fdR(0.0),
50  fClusterMap(kFALSE),
51  fITSTrackPoints(kFALSE),
52  fMustTPC(kFALSE),
53  fNTPCClustMin(0),
54  fNTPCClustMax(150),
55  fTPCChi2PerClustMin(0.0),
56  fTPCChi2PerClustMax(10e5),
57  fChi2Min(0.0),
58  fChi2Max(10e5),
59  fC00Min(0.0),
60  fC00Max(10e5),
61  fC11Min(0.0),
62  fC11Max(10e5),
63  fC22Min(0.0),
64  fC22Max(10e5),
65  fC33Min(0.0),
66  fC33Max(10e5),
67  fC44Min(0.0),
68  fC44Max(10e5),
69  fTPCC00Min(0.0),
70  fTPCC00Max(10e5),
71  fTPCC11Min(0.0),
72  fTPCC11Max(10e5),
73  fTPCC22Min(0.0),
74  fTPCC22Max(10e5),
75  fTPCC33Min(0.0),
76  fTPCC33Max(10e5),
77  fTPCC44Min(0.0),
78  fTPCC44Max(10e5)
79
80 {
81   //cosntructor
82   if ( ((Int_t)kNSpecies) != ((Int_t)AliESDtrack::kSPECIES))
83     Fatal("AliHBTReaderESD","ESD defintions probobly changed. Ask Youra.");
84 }
85 /********************************************************************/
86   
87 AliHBTReaderESD::AliHBTReaderESD(TObjArray* dirs,const Char_t* esdfilename, const Char_t* galfilename):
88  AliHBTReader(dirs), 
89  fESDFileName(esdfilename),
90  fGAlFileName(galfilename),
91  fFile(0x0),
92  fRunLoader(0x0),
93  fKeyIterator(0x0),
94  fReadParticles(kFALSE),
95  fCheckParticlePID(kFALSE),
96  fNTrackPoints(0),
97  fdR(0.0),
98  fClusterMap(kFALSE),
99  fITSTrackPoints(kFALSE),
100  fMustTPC(kFALSE),
101  fNTPCClustMin(0),
102  fNTPCClustMax(150),
103  fTPCChi2PerClustMin(0.0),
104  fTPCChi2PerClustMax(10e5),
105  fChi2Min(0.0),
106  fChi2Max(10e5),
107  fC00Min(0.0),
108  fC00Max(10e5),
109  fC11Min(0.0),
110  fC11Max(10e5),
111  fC22Min(0.0),
112  fC22Max(10e5),
113  fC33Min(0.0),
114  fC33Max(10e5),
115  fC44Min(0.0),
116  fC44Max(10e5),
117  fTPCC00Min(0.0),
118  fTPCC00Max(10e5),
119  fTPCC11Min(0.0),
120  fTPCC11Max(10e5),
121  fTPCC22Min(0.0),
122  fTPCC22Max(10e5),
123  fTPCC33Min(0.0),
124  fTPCC33Max(10e5),
125  fTPCC44Min(0.0),
126  fTPCC44Max(10e5)
127 {
128   //cosntructor
129   if ( ((Int_t)kNSpecies) != ((Int_t)AliESDtrack::kSPECIES))
130     Fatal("AliHBTReaderESD","ESD defintions probobly changed. Ask Youra.");
131 }
132 /********************************************************************/
133
134 AliHBTReaderESD::~AliHBTReaderESD()
135 {
136  //desctructor
137   delete fRunLoader;
138   delete fKeyIterator;
139   delete fFile;
140 }
141 /**********************************************************/
142 Int_t AliHBTReaderESD::ReadNext()
143 {
144 //reads next event from fFile
145   //fRunLoader is for reading Kine
146   
147   if (AliHBTParticle::GetDebug())
148     Info("ReadNext","Entered");
149     
150   if (fParticlesEvent == 0x0)  fParticlesEvent = new AliHBTEvent();
151   if (fTracksEvent == 0x0)  fTracksEvent = new AliHBTEvent();
152   
153   fParticlesEvent->Reset();
154   fTracksEvent->Reset();
155         
156   do  //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
157    {
158      if (fFile == 0x0)
159       {
160        fFile = OpenFile(fCurrentDir);//rl is opened here
161        if (fFile == 0x0)
162          {
163            Error("ReadNext","Cannot get fFile for dir no. %d",fCurrentDir);
164            fCurrentDir++;
165            continue;
166          }
167        fCurrentEvent = 0;
168        fKeyIterator = new TIter(fFile->GetListOfKeys());  
169 //       fFile->Dump();
170 //       fFile->GetListOfKeys()->Print();
171       } 
172      TKey* key = (TKey*)fKeyIterator->Next();
173      if (key == 0x0)
174       {
175         if (AliHBTParticle::GetDebug() > 2 )
176           {
177             Info("ReadNext","No more keys.");
178           }
179         fCurrentDir++;
180         delete fKeyIterator;
181         fKeyIterator = 0x0;
182         delete fFile;//we have to assume there is no more ESD objects in the fFile
183         fFile = 0x0;
184         delete fRunLoader;
185         fRunLoader = 0x0;
186         continue;
187       }
188      //try to read
189      
190      
191 //     TObject* esdobj = key->ReadObj();
192 //     if (esdobj == 0x0)
193 //      {
194 //        if (AliHBTParticle::GetDebug() > 2 )
195 //          {
196 //            Info("ReadNext","Key read NULL. Key Name is %s",key->GetName());
197 //            key->Dump();
198 //          }
199 //        continue;
200 //      }
201 //     esdobj->Dump();
202 //     AliESD* esd = dynamic_cast<AliESD*>(esdobj);
203      
204      TString esdname = "ESD";
205      esdname+=fCurrentEvent;
206      AliESD* esd = dynamic_cast<AliESD*>(fFile->Get(esdname));
207      if (esd == 0x0)
208       {
209 //        if (AliHBTParticle::GetDebug() > 2 )
210 //          {
211 //            Info("ReadNext","This key is not an AliESD object %s",key->GetName());
212 //          }
213         if (AliHBTParticle::GetDebug() > 2 )
214           {
215             Info("ReadNext","Can not find AliESD object named %s",esdname.Data());
216           }
217         fCurrentDir++;
218         delete fKeyIterator;
219         fKeyIterator = 0x0;
220         delete fFile;//we have to assume there is no more ESD objects in the fFile
221         fFile = 0x0;
222         delete fRunLoader;
223         fRunLoader = 0x0;
224         continue;
225       }
226      
227      ReadESD(esd);
228       
229      fCurrentEvent++;
230      fNEventsRead++;
231      delete esd;
232      return 0;//success -> read one event
233    }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array  
234    
235   return 1; //no more directories to read
236 }
237 /**********************************************************/
238
239 Int_t AliHBTReaderESD::ReadESD(AliESD* esd)
240 {
241   //****** Tentative particle type "concentrations"
242   static const Double_t concentr[5]={0.05, 0., 0.85, 0.10, 0.05};
243   
244   Double_t pidtable[kNSpecies];//array used for reading pid probabilities from ESD track
245   Double_t w[kNSpecies];
246   Double_t mom[3];//momentum
247   Double_t pos[3];//position
248   Double_t vertexpos[3];//vertex position
249   //Reads one ESD
250   if (esd == 0x0)
251    {
252      Error("ReadESD","ESD is NULL");
253      return 1;
254    }
255
256   TDatabasePDG* pdgdb = TDatabasePDG::Instance();
257   if (pdgdb == 0x0)
258    {
259      Error("ReadESD","Can not get PDG Database Instance.");
260      return 1;
261    }
262    
263   Float_t mf = esd->GetMagneticField(); 
264
265   if ( (mf == 0.0) && ((fNTrackPoints > 0) || fITSTrackPoints) )
266    {
267       Error("ReadESD","Magnetic Field is 0 and Track Points Demended. Skipping to next event.");
268       return 1;
269    }
270   if (fITSTrackPoints)
271    {
272      Info("ReadESD","Magnetic Field is %f",mf/10.);
273      AliKalmanTrack::SetMagneticField(mf/10.);
274    }
275    
276   AliStack* stack = 0x0;
277   if (fReadParticles && fRunLoader)
278    {
279      fRunLoader->GetEvent(fCurrentEvent);
280      stack = fRunLoader->Stack();
281    }
282
283   const AliESDVertex* vertex = esd->GetVertex();
284   if (vertex == 0x0)
285    {
286      Info("ReadESD","ESD returned NULL pointer to vertex - assuming (0.0,0.0,0.0)");
287      vertexpos[0] = 0.0;
288      vertexpos[1] = 0.0;
289      vertexpos[2] = 0.0;
290    }
291   else
292    {
293      vertex->GetXYZ(vertexpos);
294    }
295    
296   if (AliHBTParticle::GetDebug() > 0)
297    {
298      Info("ReadESD","Primary Vertex is (%f,%f,%f)",vertexpos[0],vertexpos[1],vertexpos[2]);
299    }
300    
301   Info("ReadESD","Reading Event %d",fCurrentEvent);
302
303   Int_t ntr = esd->GetNumberOfTracks();
304   Info("ReadESD","Found %d tracks.",ntr);
305   for (Int_t i = 0;i<ntr; i++)
306    {
307      AliESDtrack *esdtrack = esd->GetTrack(i);
308      if (esdtrack == 0x0)
309       {
310         Error("Next","Can not get track %d", i);
311         continue;
312       }
313
314      //if (esdtrack->HasVertexParameters() == kFALSE) 
315      if ((esdtrack->GetStatus() & AliESDtrack::kITSrefit) == kFALSE)
316       {
317         if (AliHBTParticle::GetDebug() > 2) 
318           Info("ReadNext","Particle skipped: Data at vertex not available.");
319         continue;
320       }
321
322      if (fMustTPC)
323       {
324        if ((esdtrack->GetStatus() & AliESDtrack::kTPCin) == kFALSE)
325         {
326           if (AliHBTParticle::GetDebug() > 2) 
327             Info("ReadNext","Particle skipped: Was not reconstructed in TPC.");
328           continue;
329         }
330       }     
331      if ((esdtrack->GetStatus() & AliESDtrack::kESDpid) == kFALSE) 
332       {
333         if (AliHBTParticle::GetDebug() > 2) 
334           Info("ReadNext","Particle skipped: PID BIT is not set.");
335         continue;
336       }
337
338
339      Double_t extx;
340      Double_t extp[5];
341      esdtrack->GetConstrainedExternalParameters(extx,extp);
342      if (extp[4] == 0.0)
343       {
344         if (AliHBTParticle::GetDebug() > 2) 
345           Info("ReadNext","Track has 0 contrianed curvature -> Probobly parameters never updated. Skipping.");
346         continue;
347       } 
348      esdtrack->GetESDpid(pidtable);
349      esdtrack->GetConstrainedPxPyPz(mom);
350      esdtrack->GetConstrainedXYZ(pos);
351      pos[0] -= vertexpos[0];//we are interested only in relative position to Primary vertex at this point
352      pos[1] -= vertexpos[1];
353      pos[2] -= vertexpos[2];
354
355      Int_t charge = (extp[4] > 0)?1:-1;//if curvature=charg/Pt is positive charge is positive
356
357      //Particle from kinematics
358      AliHBTParticle* particle = 0;
359      Bool_t keeppart = kFALSE;
360      if ( fReadParticles && stack  )
361       {
362         if (esdtrack->GetLabel() < 0) continue;//this is fake -  we are not able to match any track 
363         TParticle *p = stack->Particle(esdtrack->GetLabel());
364         if (p==0x0) 
365          {
366            Error("ReadNext","Can not find track with such label.");
367            continue;
368          }
369         if (fCheckParticlePID)
370          {
371            if(Pass(p->GetPdgCode())) 
372             {
373               if ( AliHBTParticle::GetDebug() > 5 )
374                 Info("ReadNext","Simulated Particle PID (%d) did not pass the cut.",p->GetPdgCode());
375               continue; //check if we are intersted with particles of this type 
376             }
377          }
378 //           if(p->GetPdgCode()<0) charge = -1;
379         particle = new AliHBTParticle(*p,i);
380
381       }
382       
383      if(CheckTrack(esdtrack)) continue;
384       
385      //Here we apply Bayes' formula
386      Double_t rc=0.;
387      for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) rc+=concentr[s]*pidtable[s];
388      if (rc==0.0) 
389       {
390         if (AliHBTParticle::GetDebug() > 2) 
391           Info("ReadNext","Particle rejected since total bayessian PID probab. is zero.");
392         continue;
393       }
394
395      for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) w[s]=concentr[s]*pidtable[s]/rc;
396
397      if (AliHBTParticle::GetDebug() > 4)
398       { 
399         Info("ReadNext","###########################################################################");
400         Info("ReadNext","Momentum: %f %f %f",mom[0],mom[1],mom[2]);
401         Info("ReadNext","Position: %f %f %f",pos[0],pos[1],pos[2]);
402         TString msg("Pid list got from track:");
403         for (Int_t s = 0;s<kNSpecies;s++)
404          {
405            msg+="\n    ";
406            msg+=s;
407            msg+="(";
408            msg+=charge*GetSpeciesPdgCode((ESpecies)s);
409            msg+="): ";
410            msg+=w[s];
411            msg+=" (";
412            msg+=pidtable[s];
413            msg+=")";
414          }
415         Info("ReadNext","%s",msg.Data());
416       }//if (AliHBTParticle::GetDebug()>4)
417
418       AliHBTTrackPoints* tpts = 0x0;
419       if (fNTrackPoints > 0) 
420        {
421          tpts = new AliHBTTrackPoints(fNTrackPoints,esdtrack,mf,fdR);
422          tpts->Move(-vertexpos[0],-vertexpos[1],-vertexpos[2]);
423        }
424       
425       AliHBTTrackPoints* itstpts = 0x0;
426       if (fITSTrackPoints) 
427        {
428          itstpts = new AliHBTTrackPoints(AliHBTTrackPoints::kITS,esdtrack);
429          itstpts->Move(-vertexpos[0],-vertexpos[1],-vertexpos[2]);
430        }
431       
432
433       AliHBTClusterMap* cmap = 0x0; 
434       if (  fClusterMap ) 
435        {
436          cmap = new AliHBTClusterMap(esdtrack);
437        }
438
439      for (Int_t s = 0; s<kNSpecies; s++)
440       {
441         Int_t pdgcode = charge*GetSpeciesPdgCode((ESpecies)s);
442         Float_t pp = w[s];
443         if (pp == 0.0) 
444          {
445            if ( AliHBTParticle::GetDebug() > 5 )
446              Info("ReadNext","Probability of being PID %d is zero. Continuing.",pdgcode);
447            continue;
448          }
449
450         if(Pass(pdgcode)) 
451          {
452            if ( AliHBTParticle::GetDebug() > 5 )
453              Info("ReadNext","PID (%d) did not pass the cut.",pdgcode);
454            continue; //check if we are intersted with particles of this type 
455          }
456
457         Double_t mass = pdgdb->GetParticle(pdgcode)->Mass();
458         Double_t tEtot = TMath::Sqrt( mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2] + mass*mass);//total energy of the track
459
460         AliHBTParticle* track = new AliHBTParticle(pdgcode, w[s],i, 
461                                                    mom[0], mom[1], mom[2], tEtot,
462                                                    pos[0], pos[1], pos[2], 0.);
463         //copy probabilitis of other species (if not zero)
464         for (Int_t k = 0; k<kNSpecies; k++)
465          {
466            if (k == s) continue;
467            if (w[k] == 0.0) continue;
468            track->SetPIDprobability(charge*GetSpeciesPdgCode( (ESpecies)k ),w[k]);
469          }
470
471         if(Pass(track))//check if meets all criteria of any of our cuts
472                        //if it does not delete it and take next good track
473          { 
474            if ( AliHBTParticle::GetDebug() > 4 )
475              Info("ReadNext","Track did not pass the cut");
476            delete track;
477            continue;
478          }
479
480         //Single Particle cuts on cluster map and track points rather do not have sense
481         if (tpts)
482          {
483            track->SetTrackPoints(tpts); 
484          }
485
486         if (itstpts)
487          {
488            track->SetITSTrackPoints(itstpts); 
489          }
490
491         if (cmap) 
492          { 
493            track->SetClusterMap(cmap);
494          }
495
496         fTracksEvent->AddParticle(track);
497         if (particle) fParticlesEvent->AddParticle(particle);
498         keeppart = kTRUE;
499
500         if (AliHBTParticle::GetDebug() > 4 )
501          {
502            Info("ReadNext","\n\nAdding Particle with incarnation %d",pdgcode);
503            track->Print();
504            if (particle) particle->Print();
505            Info("ReadNext","\n----------------------------------------------\n");
506          }
507       }//for (Int_t s = 0; s<kNSpecies; s++)
508
509      if (keeppart == kFALSE) 
510       {
511         delete particle;//particle was not stored in event
512         delete tpts;
513         delete itstpts;
514         delete cmap;
515       }
516      else
517       {
518         if (particle->P() < 0.00001)
519          {
520            Info("ReadNext","###################################");
521            Info("ReadNext","###################################");
522            Info("ReadNext","Track Label %d",esdtrack->GetLabel());
523            TParticle *p = stack->Particle(esdtrack->GetLabel());
524            Info("ReadNext","");
525            p->Print();
526            Info("ReadNext","");
527            particle->Print();
528            TString command("touch BadInput");
529            ofstream sfile("BadEvent",ios::out);
530            sfile<<fCurrentDir<<endl;
531          }
532       } 
533
534    }//for (Int_t i = 0;i<ntr; i++)  -- loop over tracks
535
536   Info("ReadNext","Read %d tracks and %d particles from event %d (event %d in dir %d).",
537          fTracksEvent->GetNumberOfParticles(), fParticlesEvent->GetNumberOfParticles(),
538          fNEventsRead,fCurrentEvent,fCurrentDir);
539   fTrackCounter->Fill(fTracksEvent->GetNumberOfParticles());
540   return 0;
541 }
542
543 /**********************************************************/
544
545 void AliHBTReaderESD::Rewind()
546 {
547   //rewinds reading 
548   delete fKeyIterator;
549   delete fFile;
550   fFile = 0x0;
551   fKeyIterator = 0x0;
552   delete fRunLoader;
553   fRunLoader = 0x0;
554   fCurrentDir = 0;
555   fNEventsRead = 0;
556   if (fTrackCounter) fTrackCounter->Reset();
557 }
558 /**********************************************************/
559
560 TFile* AliHBTReaderESD::OpenFile(Int_t n)
561 {
562 //opens fFile with kine tree
563
564  const TString& dirname = GetDirName(n);
565  if (dirname == "")
566   {
567    Error("OpenFiles","Can not get directory name");
568    return 0x0;
569   }
570  TString filename = dirname +"/"+ fESDFileName;
571  TFile *ret = TFile::Open(filename.Data()); 
572
573  if ( ret == 0x0)
574   {
575     Error("OpenFiles","Can't open fFile %s",filename.Data());
576     return 0x0;
577   }
578  if (!ret->IsOpen())
579   {
580     Error("OpenFiles","Can't open fFile  %s",filename.Data());
581     return 0x0;
582   }
583  
584  if (fReadParticles )
585   {
586    fRunLoader = AliRunLoader::Open(dirname +"/"+ fGAlFileName);
587    if (fRunLoader == 0x0)
588     {
589       Error("OpenFiles","Can't get RunLoader for directory %s",dirname.Data());
590       delete ret;
591       return 0x0;
592     }
593     
594    fRunLoader->LoadHeader();
595    if (fRunLoader->LoadKinematics())
596     {
597       Error("Next","Error occured while loading kinematics.");
598       delete fRunLoader;
599       delete ret;
600       return 0x0;
601     }
602   }
603    
604  return ret;
605 }
606 /**********************************************************/
607
608 Int_t AliHBTReaderESD::GetSpeciesPdgCode(ESpecies spec)//skowron
609 {
610   //returns pdg code from the PID index
611   //ask jura about charge
612   switch (spec)
613    {
614      case kESDElectron:
615        return kPositron;
616        break;
617      case kESDMuon:
618        return kMuonPlus;
619        break;
620      case kESDPion:
621        return kPiPlus;
622        break;
623      case kESDKaon:
624        return kKPlus;
625        break;
626      case kESDProton:
627        return kProton;
628        break;
629      default:
630        ::Warning("GetSpeciesPdgCode","Specie with number %d is not defined.",(Int_t)spec);
631        break;
632    }
633   return 0;
634 }
635 /********************************************************************/
636 Bool_t AliHBTReaderESD::CheckTrack(AliESDtrack* t) const
637 {
638   //Performs check of the track
639   
640   if ( (t->GetConstrainedChi2() < fChi2Min) || (t->GetConstrainedChi2() > fChi2Max) ) return kTRUE;
641   
642   if ( (t->GetTPCclusters(0x0) < fNTPCClustMin) || (t->GetTPCclusters(0x0) > fNTPCClustMax) ) return kTRUE;
643
644   if (t->GetTPCclusters(0x0) > 0)
645    {
646      Float_t chisqpercl = t->GetTPCchi2()/((Double_t)t->GetTPCclusters(0x0));
647      if ( (chisqpercl < fTPCChi2PerClustMin) || (chisqpercl > fTPCChi2PerClustMax) ) return kTRUE;
648    }
649
650   Double_t cc[15];
651   t->GetConstrainedExternalCovariance(cc);
652
653   if ( (cc[0]  < fC00Min) || (cc[0]  > fC00Max) ) return kTRUE;
654   if ( (cc[2]  < fC11Min) || (cc[2]  > fC11Max) ) return kTRUE;
655   if ( (cc[5]  < fC22Min) || (cc[5]  > fC22Max) ) return kTRUE;
656   if ( (cc[9]  < fC33Min) || (cc[9]  > fC33Max) ) return kTRUE;
657   if ( (cc[14] < fC44Min) || (cc[14] > fC44Max) ) return kTRUE;
658
659
660   t->GetInnerExternalCovariance(cc);
661
662   if ( (cc[0]  < fTPCC00Min) || (cc[0]  > fTPCC00Max) ) return kTRUE;
663   if ( (cc[2]  < fTPCC11Min) || (cc[2]  > fTPCC11Max) ) return kTRUE;
664   if ( (cc[5]  < fTPCC22Min) || (cc[5]  > fTPCC22Max) ) return kTRUE;
665   if ( (cc[9]  < fTPCC33Min) || (cc[9]  > fTPCC33Max) ) return kTRUE;
666   if ( (cc[14] < fTPCC44Min) || (cc[14] > fTPCC44Max) ) return kTRUE;
667
668   return kFALSE;
669
670 }
671 /********************************************************************/
672
673 void AliHBTReaderESD::SetChi2Range(Float_t min, Float_t max)
674 {
675   //sets range of Chi2 per Cluster
676   fChi2Min = min;
677   fChi2Max = max;
678 }
679 /********************************************************************/
680
681 void AliHBTReaderESD::SetTPCNClustersRange(Int_t min,Int_t max)
682 {
683  //sets range of Number Of Clusters that tracks have to have
684  fNTPCClustMin = min;
685  fNTPCClustMax = max;
686 }
687 /********************************************************************/
688
689 void AliHBTReaderESD::SetTPCChi2PerCluserRange(Float_t min, Float_t max)
690 {
691   //sets range of Chi2 per Cluster
692   fTPCChi2PerClustMin = min;
693   fTPCChi2PerClustMax = max;
694 }
695 /********************************************************************/
696
697 void AliHBTReaderESD::SetC00Range(Float_t min, Float_t max)
698 {
699  //Sets range of C00 parameter of covariance matrix of the track
700  //it defines uncertainty of the momentum
701   fC00Min = min;
702   fC00Max = max;
703 }
704 /********************************************************************/
705
706 void AliHBTReaderESD::SetC11Range(Float_t min, Float_t max)
707 {
708  //Sets range of C11 parameter of covariance matrix of the track
709  //it defines uncertainty of the momentum
710   fC11Min = min;
711   fC11Max = max;
712 }
713 /********************************************************************/
714
715 void AliHBTReaderESD::SetC22Range(Float_t min, Float_t max)
716 {
717  //Sets range of C22 parameter of covariance matrix of the track
718  //it defines uncertainty of the momentum
719   fC22Min = min;
720   fC22Max = max;
721 }
722 /********************************************************************/
723
724 void AliHBTReaderESD::SetC33Range(Float_t min, Float_t max)
725 {
726  //Sets range of C33 parameter of covariance matrix of the track
727  //it defines uncertainty of the momentum
728   fC33Min = min;
729   fC33Max = max;
730 }
731 /********************************************************************/
732
733 void AliHBTReaderESD::SetC44Range(Float_t min, Float_t max)
734 {
735  //Sets range of C44 parameter of covariance matrix of the track
736  //it defines uncertainty of the momentum
737   fC44Min = min;
738   fC44Max = max;
739 }