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