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