]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTReaderESD.cxx
c477189d7b8be1b1e1cc4891d4bfa49e76ca2455
[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(Pass(p->GetPdgCode())) 
349          {
350            if ( AliHBTParticle::GetDebug() > 5 )
351              Info("ReadNext","Simulated Particle PID (%d) did not pass the cut.",p->GetPdgCode());
352            continue; //check if we are intersted with particles of this type 
353          }
354 //           if(p->GetPdgCode()<0) charge = -1;
355         particle = new AliHBTParticle(*p,i);
356
357       }
358       
359      if(CheckTrack(esdtrack)) continue;
360       
361      //Here we apply Bayes' formula
362      Double_t rc=0.;
363      for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) rc+=concentr[s]*pidtable[s];
364      if (rc==0.0) 
365       {
366         if (AliHBTParticle::GetDebug() > 2) 
367           Info("ReadNext","Particle rejected since total bayessian PID probab. is zero.");
368         continue;
369       }
370
371      for (Int_t s=0; s<AliESDtrack::kSPECIES; s++) w[s]=concentr[s]*pidtable[s]/rc;
372
373      if (AliHBTParticle::GetDebug() > 4)
374       { 
375         Info("ReadNext","###########################################################################");
376         Info("ReadNext","Momentum: %f %f %f",mom[0],mom[1],mom[2]);
377         Info("ReadNext","Position: %f %f %f",pos[0],pos[1],pos[2]);
378         TString msg("Pid list got from track:");
379         for (Int_t s = 0;s<kNSpecies;s++)
380          {
381            msg+="\n    ";
382            msg+=s;
383            msg+="(";
384            msg+=charge*GetSpeciesPdgCode((ESpecies)s);
385            msg+="): ";
386            msg+=w[s];
387            msg+=" (";
388            msg+=pidtable[s];
389            msg+=")";
390          }
391         Info("ReadNext","%s",msg.Data());
392       }//if (AliHBTParticle::GetDebug()>4)
393
394       AliHBTTrackPoints* tpts = 0x0;
395       if (fNTrackPoints > 0) 
396        {
397          tpts = new AliHBTTrackPoints(fNTrackPoints,esdtrack,mf,fdR);
398        }
399
400       AliHBTClusterMap* cmap = 0x0; 
401       if (  fClusterMap ) 
402        {
403          cmap = new AliHBTClusterMap(esdtrack);
404        }
405
406      for (Int_t s = 0; s<kNSpecies; s++)
407       {
408         Int_t pdgcode = charge*GetSpeciesPdgCode((ESpecies)s);
409         Float_t pp = w[s];
410         if (pp == 0.0) 
411          {
412            if ( AliHBTParticle::GetDebug() > 5 )
413              Info("ReadNext","Probability of being PID %d is zero. Continuing.",pdgcode);
414            continue;
415          }
416
417         if(Pass(pdgcode)) 
418          {
419            if ( AliHBTParticle::GetDebug() > 5 )
420              Info("ReadNext","PID (%d) did not pass the cut.",pdgcode);
421            continue; //check if we are intersted with particles of this type 
422          }
423
424         Double_t mass = pdgdb->GetParticle(pdgcode)->Mass();
425         Double_t tEtot = TMath::Sqrt( mom[0]*mom[0] + mom[1]*mom[1] + mom[2]*mom[2] + mass*mass);//total energy of the track
426
427         AliHBTParticle* track = new AliHBTParticle(pdgcode, w[s],i, 
428                                                    mom[0], mom[1], mom[2], tEtot,
429                                                    pos[0], pos[1], pos[2], 0.);
430         //copy probabilitis of other species (if not zero)
431         for (Int_t k = 0; k<kNSpecies; k++)
432          {
433            if (k == s) continue;
434            if (w[k] == 0.0) continue;
435            track->SetPIDprobability(charge*GetSpeciesPdgCode( (ESpecies)k ),w[k]);
436          }
437
438         if(Pass(track))//check if meets all criteria of any of our cuts
439                        //if it does not delete it and take next good track
440          { 
441            if ( AliHBTParticle::GetDebug() > 4 )
442              Info("ReadNext","Track did not pass the cut");
443            delete track;
444            continue;
445          }
446
447         //Single Particle cuts on cluster map and track points rather do not have sense
448         if (tpts)
449          {
450            track->SetTrackPoints(tpts); 
451          }
452
453         if (cmap) 
454          { 
455            track->SetClusterMap(cmap);
456          }
457
458         fTracksEvent->AddParticle(track);
459         if (particle) fParticlesEvent->AddParticle(particle);
460         keeppart = kTRUE;
461
462         if (AliHBTParticle::GetDebug() > 4 )
463          {
464            Info("ReadNext","\n\nAdding Particle with incarnation %d",pdgcode);
465            track->Print();
466            if (particle) particle->Print();
467            Info("ReadNext","\n----------------------------------------------\n");
468          }
469       }//for (Int_t s = 0; s<kNSpecies; s++)
470
471      if (keeppart == kFALSE) 
472       {
473         delete particle;//particle was not stored in event
474         delete tpts;
475         delete cmap;
476       }
477
478    }//for (Int_t i = 0;i<ntr; i++)  -- loop over tracks
479
480   Info("ReadNext","Read %d tracks and %d particles from event %d (event %d in dir %d).",
481          fTracksEvent->GetNumberOfParticles(), fParticlesEvent->GetNumberOfParticles(),
482          fNEventsRead,fCurrentEvent,fCurrentDir);
483   fTrackCounter->Fill(fTracksEvent->GetNumberOfParticles());
484   return 0;
485 }
486
487 /**********************************************************/
488
489 void AliHBTReaderESD::Rewind()
490 {
491   //rewinds reading 
492   delete fKeyIterator;
493   delete fFile;
494   fFile = 0x0;
495   fKeyIterator = 0x0;
496   delete fRunLoader;
497   fRunLoader = 0x0;
498   fCurrentDir = 0;
499   fNEventsRead = 0;
500   if (fTrackCounter) fTrackCounter->Reset();
501 }
502 /**********************************************************/
503
504 TFile* AliHBTReaderESD::OpenFile(Int_t n)
505 {
506 //opens fFile with kine tree
507
508  const TString& dirname = GetDirName(n);
509  if (dirname == "")
510   {
511    Error("OpenFiles","Can not get directory name");
512    return 0x0;
513   }
514  TString filename = dirname +"/"+ fESDFileName;
515  TFile *ret = TFile::Open(filename.Data()); 
516
517  if ( ret == 0x0)
518   {
519     Error("OpenFiles","Can't open fFile %s",filename.Data());
520     return 0x0;
521   }
522  if (!ret->IsOpen())
523   {
524     Error("OpenFiles","Can't open fFile  %s",filename.Data());
525     return 0x0;
526   }
527  
528  if (fReadParticles )
529   {
530    fRunLoader = AliRunLoader::Open(dirname +"/"+ fGAlFileName);
531    if (fRunLoader == 0x0)
532     {
533       Error("OpenFiles","Can't get RunLoader for directory %s",dirname.Data());
534       delete ret;
535       return 0x0;
536     }
537     
538    fRunLoader->LoadHeader();
539    if (fRunLoader->LoadKinematics())
540     {
541       Error("Next","Error occured while loading kinematics.");
542       delete fRunLoader;
543       delete ret;
544       return 0x0;
545     }
546   }
547    
548  return ret;
549 }
550 /**********************************************************/
551
552 Int_t AliHBTReaderESD::GetSpeciesPdgCode(ESpecies spec)//skowron
553 {
554   //returns pdg code from the PID index
555   //ask jura about charge
556   switch (spec)
557    {
558      case kESDElectron:
559        return kPositron;
560        break;
561      case kESDMuon:
562        return kMuonPlus;
563        break;
564      case kESDPion:
565        return kPiPlus;
566        break;
567      case kESDKaon:
568        return kKPlus;
569        break;
570      case kESDProton:
571        return kProton;
572        break;
573      default:
574        ::Warning("GetSpeciesPdgCode","Specie with number %d is not defined.",(Int_t)spec);
575        break;
576    }
577   return 0;
578 }
579 /********************************************************************/
580 Bool_t AliHBTReaderESD::CheckTrack(AliESDtrack* t) const
581 {
582   //Performs check of the track
583   
584   if ( (t->GetConstrainedChi2() < fChi2Min) || (t->GetConstrainedChi2() > fChi2Max) ) return kTRUE;
585   
586   if ( (t->GetTPCclusters(0x0) < fNTPCClustMin) || (t->GetTPCclusters(0x0) > fNTPCClustMax) ) return kTRUE;
587
588   if (t->GetTPCclusters(0x0) > 0)
589    {
590      Float_t chisqpercl = t->GetTPCchi2()/((Double_t)t->GetTPCclusters(0x0));
591      if ( (chisqpercl < fTPCChi2PerClustMin) || (chisqpercl > fTPCChi2PerClustMax) ) return kTRUE;
592    }
593
594   Double_t cc[15];
595   t->GetConstrainedExternalCovariance(cc);
596
597   if ( (cc[0]  < fC00Min) || (cc[0]  > fC00Max) ) return kTRUE;
598   if ( (cc[2]  < fC11Min) || (cc[2]  > fC11Max) ) return kTRUE;
599   if ( (cc[5]  < fC22Min) || (cc[5]  > fC22Max) ) return kTRUE;
600   if ( (cc[9]  < fC33Min) || (cc[9]  > fC33Max) ) return kTRUE;
601   if ( (cc[14] < fC44Min) || (cc[14] > fC44Max) ) return kTRUE;
602
603
604   t->GetInnerExternalCovariance(cc);
605
606   if ( (cc[0]  < fTPCC00Min) || (cc[0]  > fTPCC00Max) ) return kTRUE;
607   if ( (cc[2]  < fTPCC11Min) || (cc[2]  > fTPCC11Max) ) return kTRUE;
608   if ( (cc[5]  < fTPCC22Min) || (cc[5]  > fTPCC22Max) ) return kTRUE;
609   if ( (cc[9]  < fTPCC33Min) || (cc[9]  > fTPCC33Max) ) return kTRUE;
610   if ( (cc[14] < fTPCC44Min) || (cc[14] > fTPCC44Max) ) return kTRUE;
611
612   return kFALSE;
613
614 }
615 /********************************************************************/
616
617 void AliHBTReaderESD::SetChi2Range(Float_t min, Float_t max)
618 {
619   //sets range of Chi2 per Cluster
620   fChi2Min = min;
621   fChi2Max = max;
622 }
623 /********************************************************************/
624
625 void AliHBTReaderESD::SetTPCNClustersRange(Int_t min,Int_t max)
626 {
627  //sets range of Number Of Clusters that tracks have to have
628  fNTPCClustMin = min;
629  fNTPCClustMax = max;
630 }
631 /********************************************************************/
632
633 void AliHBTReaderESD::SetTPCChi2PerCluserRange(Float_t min, Float_t max)
634 {
635   //sets range of Chi2 per Cluster
636   fTPCChi2PerClustMin = min;
637   fTPCChi2PerClustMax = max;
638 }
639 /********************************************************************/
640
641 void AliHBTReaderESD::SetC00Range(Float_t min, Float_t max)
642 {
643  //Sets range of C00 parameter of covariance matrix of the track
644  //it defines uncertainty of the momentum
645   fC00Min = min;
646   fC00Max = max;
647 }
648 /********************************************************************/
649
650 void AliHBTReaderESD::SetC11Range(Float_t min, Float_t max)
651 {
652  //Sets range of C11 parameter of covariance matrix of the track
653  //it defines uncertainty of the momentum
654   fC11Min = min;
655   fC11Max = max;
656 }
657 /********************************************************************/
658
659 void AliHBTReaderESD::SetC22Range(Float_t min, Float_t max)
660 {
661  //Sets range of C22 parameter of covariance matrix of the track
662  //it defines uncertainty of the momentum
663   fC22Min = min;
664   fC22Max = max;
665 }
666 /********************************************************************/
667
668 void AliHBTReaderESD::SetC33Range(Float_t min, Float_t max)
669 {
670  //Sets range of C33 parameter of covariance matrix of the track
671  //it defines uncertainty of the momentum
672   fC33Min = min;
673   fC33Max = max;
674 }
675 /********************************************************************/
676
677 void AliHBTReaderESD::SetC44Range(Float_t min, Float_t max)
678 {
679  //Sets range of C44 parameter of covariance matrix of the track
680  //it defines uncertainty of the momentum
681   fC44Min = min;
682   fC44Max = max;
683 }