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