]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliUEHistograms.cxx
Rewritten spectrum2 task, old running is optional, various minor updates to AddTask...
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliUEHistograms.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id: AliUEHistograms.cxx 20164 2007-08-14 15:31:50Z morsch $ */
17
18 //
19 //
20 // encapsulates several AliUEHist objects for a full UE analysis plus additional control histograms
21 //
22 //
23 // Author: Jan Fiete Grosse-Oetringhaus, Sara Vallero
24
25 #include "AliUEHistograms.h"
26
27 #include "AliCFContainer.h"
28 #include "AliVParticle.h"
29
30 #include "TList.h"
31 #include "TH2F.h"
32 #include "TH1F.h"
33 #include "TH3F.h"
34 #include "TMath.h"
35
36 ClassImp(AliUEHistograms)
37
38 const Int_t AliUEHistograms::fgkUEHists = 3;
39
40 AliUEHistograms::AliUEHistograms(const char* name, const char* histograms) : 
41   TNamed(name, name),
42   fNumberDensitypT(0),
43   fSumpT(0),
44   fNumberDensityPhi(0),
45   fCorrelationpT(0),
46   fCorrelationEta(0),
47   fCorrelationPhi(0),
48   fCorrelationR(0),
49   fCorrelationLeading2Phi(0),
50   fCorrelationMultiplicity(0),
51   fEventCount(0),
52   fEventCountDifferential(0),
53   fVertexContributors(0),
54   fCentralityDistribution(0)
55 {
56   // Constructor
57   //
58   // the string histograms defines which histograms are created:
59   //    1 = NumberDensitypT
60   //    2 = SumpT
61   //    3 = NumberDensityPhi
62   //    4 = NumberDensityPhiCentrality (other multiplicity for Pb)
63   
64   TString histogramsStr(histograms);
65   
66   if (histogramsStr.Contains("1"))
67     fNumberDensitypT = new AliUEHist("NumberDensitypT");
68   if (histogramsStr.Contains("2"))
69     fSumpT = new AliUEHist("SumpT");
70   
71   if (histogramsStr.Contains("3"))
72     fNumberDensityPhi = new AliUEHist("NumberDensityPhi");
73   else if (histogramsStr.Contains("4"))
74     fNumberDensityPhi = new AliUEHist("NumberDensityPhiCentrality");
75   
76   // do not add this hists to the directory
77   Bool_t oldStatus = TH1::AddDirectoryStatus();
78   TH1::AddDirectory(kFALSE);
79   
80   fCorrelationpT  = new TH2F("fCorrelationpT", ";p_{T,lead} (MC);p_{T,lead} (RECO)", 200, 0, 50, 200, 0, 50);
81   fCorrelationEta = new TH2F("fCorrelationEta", ";#eta_{T,lead} (MC);#eta_{T,lead} (RECO)", 200, -1, 1, 200, -1, 1);
82   fCorrelationPhi = new TH2F("fCorrelationPhi", ";#phi_{T,lead} (MC);#phi_{T,lead} (RECO)", 200, 0, TMath::TwoPi(), 200, 0, TMath::TwoPi());
83   fCorrelationR =   new TH2F("fCorrelationR", ";R;p_{T,lead} (MC)", 200, 0, 2, 200, 0, 50);
84   fCorrelationLeading2Phi = new TH2F("fCorrelationLeading2Phi", ";#Delta #phi;p_{T,lead} (MC)", 200, -TMath::Pi(), TMath::Pi(), 200, 0, 50);
85   fCorrelationMultiplicity = new TH2F("fCorrelationMultiplicity", ";MC tracks;Reco tracks", 100, -0.5, 99.5, 100, -0.5, 99.5);
86   
87   fEventCount = new TH2F("fEventCount", ";step;event type;count", AliUEHist::fgkCFSteps+2, -2.5, -0.5 + AliUEHist::fgkCFSteps, 3, -0.5, 2.5);
88   fEventCount->GetYaxis()->SetBinLabel(1, "ND");
89   fEventCount->GetYaxis()->SetBinLabel(2, "SD");
90   fEventCount->GetYaxis()->SetBinLabel(3, "DD");
91   
92   fEventCountDifferential = new TH3F("fEventCountDifferential", ";p_{T,lead};step;event type", 100, 0, 50, AliUEHist::fgkCFSteps, -0.5, -0.5 + AliUEHist::fgkCFSteps, 3, -0.5, 2.5);
93   fEventCountDifferential->GetZaxis()->SetBinLabel(1, "ND");
94   fEventCountDifferential->GetZaxis()->SetBinLabel(2, "SD");
95   fEventCountDifferential->GetZaxis()->SetBinLabel(3, "DD");
96   
97   fVertexContributors = new TH1F("fVertexContributors", ";contributors;count", 100, -0.5, 99.5);
98   
99   fCentralityDistribution = new TH1F("fCentralityDistribution", ";;count", fNumberDensityPhi->GetEventHist()->GetNBins(1), fNumberDensityPhi->GetEventHist()->GetAxis(1, 0)->GetBinLowEdge(1), fNumberDensityPhi->GetEventHist()->GetAxis(1, 0)->GetBinUpEdge(fNumberDensityPhi->GetEventHist()->GetNBins(1)));
100   
101   TH1::AddDirectory(oldStatus);
102 }
103
104 //_____________________________________________________________________________
105 AliUEHistograms::AliUEHistograms(const AliUEHistograms &c) :
106   TNamed(fName, fTitle),
107   fNumberDensitypT(0),
108   fSumpT(0),
109   fNumberDensityPhi(0),
110   fCorrelationpT(0),
111   fCorrelationEta(0),
112   fCorrelationPhi(0),
113   fCorrelationR(0),
114   fCorrelationLeading2Phi(0),
115   fCorrelationMultiplicity(0),
116   fEventCount(0),
117   fEventCountDifferential(0),
118   fVertexContributors(0),
119   fCentralityDistribution(0)
120 {
121   //
122   // AliUEHistograms copy constructor
123   //
124
125   ((AliUEHistograms &) c).Copy(*this);
126 }
127
128 //____________________________________________________________________
129 AliUEHistograms::~AliUEHistograms()
130 {
131   // Destructor
132   
133   if (fNumberDensitypT)
134   {
135     delete fNumberDensitypT;
136     fNumberDensitypT = 0;
137   }
138   
139   if (fSumpT)
140   {
141     delete fSumpT;
142     fSumpT = 0;
143   }
144   
145   if (fNumberDensityPhi)
146   {
147     delete fNumberDensityPhi;
148     fNumberDensityPhi = 0;
149   }
150   
151   if (fCorrelationpT)
152   {
153     delete fCorrelationpT;
154     fCorrelationpT = 0;
155   }
156   
157   if (fCorrelationEta)
158   {
159     delete fCorrelationEta;
160     fCorrelationEta = 0;
161   }
162   
163   if (fCorrelationPhi)
164   {
165     delete fCorrelationPhi;
166     fCorrelationPhi = 0;
167   }
168   
169   if (fCorrelationR)
170   {
171     delete fCorrelationR;
172     fCorrelationR = 0;
173   }
174
175   if (fCorrelationLeading2Phi)
176   {
177     delete fCorrelationLeading2Phi;
178     fCorrelationLeading2Phi = 0;
179   }
180   
181   if (fCorrelationMultiplicity)
182   {
183     delete fCorrelationMultiplicity;
184     fCorrelationMultiplicity = 0;
185   }
186   
187   if (fEventCount)
188   {
189     delete fEventCount;
190     fEventCount = 0;
191   }
192
193   if (fEventCountDifferential)
194   {
195     delete fEventCountDifferential;
196     fEventCountDifferential = 0;
197   }
198   
199   if (fVertexContributors)
200   {
201     delete fVertexContributors;
202     fVertexContributors = 0;
203   }
204   
205   if (fCentralityDistribution)
206   {
207     delete fCentralityDistribution;
208     fCentralityDistribution = 0;
209   }
210 }
211
212 AliUEHist* AliUEHistograms::GetUEHist(Int_t id)
213 {
214   // returns AliUEHist object, useful for loops
215   
216   switch (id)
217   {
218     case 0: return fNumberDensitypT; break;
219     case 1: return fSumpT; break;
220     case 2: return fNumberDensityPhi; break;
221   }
222   
223   return 0;
224 }
225
226 //____________________________________________________________________
227 Int_t AliUEHistograms::CountParticles(TList* list, Float_t ptMin)
228 {
229   // counts the number of particles in the list with a pT above ptMin
230   // TODO eta cut needed here?
231   
232   Int_t count = 0;
233   for (Int_t j=0; j<list->GetEntries(); j++)
234     if (((AliVParticle*) list->At(j))->Pt() > ptMin)
235       count++;
236       
237   return count;
238 }
239   
240 //____________________________________________________________________
241 void AliUEHistograms::Fill(Int_t eventType, AliUEHist::CFStep step, AliVParticle* leading, TList* toward, TList* away, TList* min, TList* max)
242 {
243   // fills the UE event histograms
244   //
245   // this function needs the leading (track or jet or ...) and four lists of AliVParticles which contain the particles/tracks to be filled in the four regions
246   
247   // if leading is not set, just fill event statistics
248   if (leading)
249   {
250     Int_t multiplicity = 0;
251     
252     // TODO configurable?
253     Float_t ptMin = 0.15;
254     if (leading->Pt() > ptMin)
255       multiplicity++;
256     
257     multiplicity += CountParticles(toward, ptMin);
258     multiplicity += CountParticles(away, ptMin);
259     multiplicity += CountParticles(min, ptMin);
260     multiplicity += CountParticles(max, ptMin);
261      
262     FillRegion(AliUEHist::kToward, step, leading, toward, multiplicity);
263     FillRegion(AliUEHist::kAway,   step, leading, away, multiplicity);
264     FillRegion(AliUEHist::kMin,    step, leading, min, multiplicity);
265     FillRegion(AliUEHist::kMax,    step, leading, max, multiplicity);
266  
267     Double_t vars[2];
268     vars[0] = leading->Pt();
269     vars[1] = multiplicity;
270     for (Int_t i=0; i<fgkUEHists; i++)
271       if (GetUEHist(i))
272         GetUEHist(i)->GetEventHist()->Fill(vars, step);
273   
274     fEventCountDifferential->Fill(leading->Pt(), step, eventType);
275   }
276   
277   FillEvent(eventType, step);
278 }
279   
280 //____________________________________________________________________
281 void AliUEHistograms::FillRegion(AliUEHist::Region region, AliUEHist::CFStep step, AliVParticle* leading, TList* list, Int_t multiplicity)
282 {
283   // loops over AliVParticles in list and fills the given region at the given step
284   //
285   // See also Fill(...)
286
287   for (Int_t i=0; i<list->GetEntries(); i++)
288   {
289     AliVParticle* particle = (AliVParticle*) list->At(i);
290     
291     Double_t vars[5];
292     vars[0] = particle->Eta();
293     vars[1] = particle->Pt();
294     vars[2] = leading->Pt();
295     vars[3] = multiplicity;
296     vars[4] = leading->Phi() - particle->Phi();
297     if (vars[4] > 1.5 * TMath::Pi()) 
298       vars[4] -= TMath::TwoPi();
299     if (vars[4] < -0.5 * TMath::Pi())
300       vars[4] += TMath::TwoPi();
301
302     if (fNumberDensitypT)
303       fNumberDensitypT->GetTrackHist(region)->Fill(vars, step);
304       
305     if (fSumpT)
306       fSumpT->GetTrackHist(region)->Fill(vars, step, particle->Pt());
307     
308     // fill all in toward region (is anyway as function of delta phi!)
309     if (fNumberDensityPhi)
310       fNumberDensityPhi->GetTrackHist(AliUEHist::kToward)->Fill(vars, step);
311   }
312 }
313
314 //____________________________________________________________________
315 void AliUEHistograms::Fill(AliVParticle* leadingMC, AliVParticle* leadingReco)
316 {
317   // fills the correlation histograms
318   
319   if (leadingMC)
320   {
321     fCorrelationpT->Fill(leadingMC->Pt(), leadingReco->Pt());
322     if (leadingMC->Pt() > 0.5)
323     {
324       fCorrelationEta->Fill(leadingMC->Eta(), leadingReco->Eta());
325       fCorrelationPhi->Fill(leadingMC->Phi(), leadingReco->Phi());
326     }
327     
328     Float_t phiDiff = leadingMC->Phi() - leadingReco->Phi();
329     if (phiDiff > TMath::Pi())
330       phiDiff -= TMath::TwoPi();
331     if (phiDiff < -TMath::Pi())
332       phiDiff += TMath::TwoPi();
333       
334     Float_t etaDiff = leadingMC->Eta() - leadingReco->Eta();
335     
336     fCorrelationR->Fill(TMath::Sqrt(phiDiff * phiDiff + etaDiff * etaDiff), leadingMC->Pt());
337     fCorrelationLeading2Phi->Fill(phiDiff, leadingMC->Pt());
338   }
339   else
340   {
341     fCorrelationpT->Fill(1.0, leadingReco->Pt());
342     if (leadingReco->Pt() > 0.5)
343     {
344       fCorrelationEta->Fill(0.0, leadingReco->Eta());
345       fCorrelationPhi->Fill(0.0, leadingReco->Phi());
346     }
347   }
348 }
349
350 //____________________________________________________________________
351 void AliUEHistograms::FillCorrelations(Int_t eventType, Int_t centrality, AliUEHist::CFStep step, TSeqCollection* particles, TSeqCollection* mixed)
352 {
353   // fills the fNumberDensityPhi histogram
354   //
355   // this function need a list of AliVParticles which contain the particles/tracks to be filled
356   //
357   // if mixed is non-0, mixed events are filled, the trigger particle is from particles, the associated from mixed
358   
359   // if particles is not set, just fill event statistics
360   if (particles)
361   {
362     for (Int_t i=0; i<particles->GetEntries(); i++)
363     {
364       AliVParticle* triggerParticle = (AliVParticle*) particles->At(i);
365       Int_t jMax = particles->GetEntries();
366       if (mixed)
367         jMax = mixed->GetEntries();
368         
369       for (Int_t j=0; j<jMax; j++)
370       {
371         if (!mixed && i == j)
372           continue;
373       
374         AliVParticle* particle = 0;
375         if (!mixed)
376           particle = (AliVParticle*) particles->At(j);
377         else
378           particle = (AliVParticle*) mixed->At(j);
379         
380         Double_t vars[5];
381         vars[0] = triggerParticle->Eta() - particle->Eta();
382         vars[1] = particle->Pt();
383         vars[2] = triggerParticle->Pt();
384         vars[3] = centrality;
385         vars[4] = triggerParticle->Phi() - particle->Phi();
386         if (vars[4] > 1.5 * TMath::Pi()) 
387           vars[4] -= TMath::TwoPi();
388         if (vars[4] < -0.5 * TMath::Pi())
389           vars[4] += TMath::TwoPi();
390     
391         // fill all in toward region and no not use the other regions
392         fNumberDensityPhi->GetTrackHist(AliUEHist::kToward)->Fill(vars, step);
393       }    
394  
395       // once per trigger particle
396       Double_t vars[2];
397       vars[0] = triggerParticle->Pt();
398       vars[1] = centrality;
399       fNumberDensityPhi->GetEventHist()->Fill(vars, step);
400     
401       fEventCountDifferential->Fill(triggerParticle->Pt(), step, eventType);
402     }
403   }
404   
405   fCentralityDistribution->Fill(centrality);
406   FillEvent(eventType, step);
407 }
408   
409 //____________________________________________________________________
410 void AliUEHistograms::FillTrackingEfficiency(TObjArray* mc, TObjArray* recoPrim, TObjArray* recoAll, Int_t particleType)
411 {
412   // fills the tracking efficiency objects
413   //
414   // mc: all primary MC particles
415   // recoPrim: reconstructed primaries (again MC particles)
416   // recoAll: reconstructed (again MC particles)
417   // particleType is: 0 for pion, 1 for kaon, 2 for proton, 3 for others
418   
419   for (Int_t step=0; step<3; step++)
420   {
421     TObjArray* list = mc;
422     if (step == 1)
423       list = recoPrim;
424     else if (step == 2)
425       list = recoAll;
426       
427     for (Int_t i=0; i<list->GetEntries(); i++)
428     {
429       AliVParticle* particle = (AliVParticle*) list->At(i);
430       Double_t vars[3];
431       vars[0] = particle->Eta();
432       vars[1] = particle->Pt();
433       vars[2] = particleType;
434       
435       for (Int_t j=0; j<fgkUEHists; j++)
436         if (GetUEHist(j))
437           GetUEHist(j)->GetTrackHistEfficiency()->Fill(vars, step);
438     }
439   }
440 }
441
442 //____________________________________________________________________
443 void AliUEHistograms::FillEvent(Int_t eventType, Int_t step)
444 {
445   // fills the number of events at the given step and the given enty type
446   //
447   // WARNING: This function is called from Fill, so only call it for steps where Fill is not called
448   
449   fEventCount->Fill(step, eventType);
450 }
451
452 //____________________________________________________________________
453 void AliUEHistograms::SetEtaRange(Float_t etaMin, Float_t etaMax)
454 {
455   // sets eta min and max for all contained AliUEHist classes
456   
457   for (Int_t i=0; i<fgkUEHists; i++)
458     if (GetUEHist(i))
459       GetUEHist(i)->SetEtaRange(etaMin, etaMax);
460 }
461
462 //____________________________________________________________________
463 void AliUEHistograms::SetPtRange(Float_t ptMin, Float_t ptMax)
464 {
465   // sets pT min and max for all contained AliUEHist classes
466   
467   for (Int_t i=0; i<fgkUEHists; i++)
468     if (GetUEHist(i))
469       GetUEHist(i)->SetPtRange(ptMin, ptMax);
470 }
471
472 //____________________________________________________________________
473 void AliUEHistograms::SetContaminationEnhancement(TH1F* hist)
474 {
475   // sets the contamination enhancement histogram in all contained AliUEHist classes
476   
477   for (Int_t i=0; i<fgkUEHists; i++)
478     if (GetUEHist(i))
479       GetUEHist(i)->SetContaminationEnhancement(hist);
480 }  
481
482 //____________________________________________________________________
483 void AliUEHistograms::SetCombineMinMax(Bool_t flag)
484 {
485   // sets pT min and max for all contained AliUEHist classes
486   
487   for (Int_t i=0; i<fgkUEHists; i++)
488     if (GetUEHist(i))
489       GetUEHist(i)->SetCombineMinMax(flag);
490 }
491
492 //____________________________________________________________________
493 void AliUEHistograms::Correct(AliUEHistograms* corrections)
494 {
495   // corrects the contained histograms by calling AliUEHist::Correct
496   
497   for (Int_t i=0; i<fgkUEHists; i++)
498     if (GetUEHist(i))
499       GetUEHist(i)->Correct(corrections->GetUEHist(i));
500 }
501
502 //____________________________________________________________________
503 AliUEHistograms &AliUEHistograms::operator=(const AliUEHistograms &c)
504 {
505   // assigment operator
506
507   if (this != &c)
508     ((AliUEHistograms &) c).Copy(*this);
509
510   return *this;
511 }
512
513 //____________________________________________________________________
514 void AliUEHistograms::Copy(TObject& c) const
515 {
516   // copy function
517
518   AliUEHistograms& target = (AliUEHistograms &) c;
519
520   if (fNumberDensitypT)
521     target.fNumberDensitypT = dynamic_cast<AliUEHist*> (fNumberDensitypT->Clone());
522
523   if (fSumpT)
524     target.fSumpT = dynamic_cast<AliUEHist*> (fSumpT->Clone());
525
526   if (fNumberDensityPhi)
527     target.fNumberDensityPhi = dynamic_cast<AliUEHist*> (fNumberDensityPhi->Clone());
528
529   if (fCorrelationpT)
530     target.fCorrelationpT = dynamic_cast<TH2F*> (fCorrelationpT->Clone());
531
532   if (fCorrelationEta)
533     target.fCorrelationEta = dynamic_cast<TH2F*> (fCorrelationEta->Clone());
534
535   if (fCorrelationPhi)
536     target.fCorrelationPhi = dynamic_cast<TH2F*> (fCorrelationPhi->Clone());
537
538   if (fCorrelationR)
539     target.fCorrelationR = dynamic_cast<TH2F*> (fCorrelationR->Clone());
540
541   if (fCorrelationLeading2Phi)
542     target.fCorrelationLeading2Phi = dynamic_cast<TH2F*> (fCorrelationLeading2Phi->Clone());
543
544   if (fCorrelationMultiplicity)
545     target.fCorrelationMultiplicity = dynamic_cast<TH2F*> (fCorrelationMultiplicity->Clone());
546   
547   if (fEventCount)
548     target.fEventCount = dynamic_cast<TH2F*> (fEventCount->Clone());
549
550   if (fEventCountDifferential)
551     target.fEventCountDifferential = dynamic_cast<TH3F*> (fEventCountDifferential->Clone());
552     
553   if (fVertexContributors)
554     target.fVertexContributors = dynamic_cast<TH1F*> (fVertexContributors->Clone());
555
556   if (fCentralityDistribution)
557     target.fCentralityDistribution = dynamic_cast<TH1F*> (fCentralityDistribution->Clone());
558 }
559
560 //____________________________________________________________________
561 Long64_t AliUEHistograms::Merge(TCollection* list)
562 {
563   // Merge a list of AliUEHistograms objects with this (needed for
564   // PROOF). 
565   // Returns the number of merged objects (including this).
566
567   if (!list)
568     return 0;
569   
570   if (list->IsEmpty())
571     return 1;
572
573   TIterator* iter = list->MakeIterator();
574   TObject* obj;
575
576   // collections of objects
577   const Int_t kMaxLists = 13;
578   TList* lists[kMaxLists];
579   
580   for (Int_t i=0; i<kMaxLists; i++)
581     lists[i] = new TList;
582   
583   Int_t count = 0;
584   while ((obj = iter->Next())) {
585     
586     AliUEHistograms* entry = dynamic_cast<AliUEHistograms*> (obj);
587     if (entry == 0) 
588       continue;
589
590     if (entry->fNumberDensitypT)
591       lists[0]->Add(entry->fNumberDensitypT);
592     if (entry->fSumpT)
593       lists[1]->Add(entry->fSumpT);
594     if (entry->fNumberDensityPhi)
595       lists[2]->Add(entry->fNumberDensityPhi);
596     lists[3]->Add(entry->fCorrelationpT);
597     lists[4]->Add(entry->fCorrelationEta);
598     lists[5]->Add(entry->fCorrelationPhi);
599     lists[6]->Add(entry->fCorrelationR);
600     lists[7]->Add(entry->fCorrelationLeading2Phi);
601     lists[8]->Add(entry->fCorrelationMultiplicity);
602     lists[9]->Add(entry->fEventCount);
603     lists[10]->Add(entry->fEventCountDifferential);
604     lists[11]->Add(entry->fVertexContributors);
605     lists[12]->Add(entry->fCentralityDistribution);
606
607     count++;
608   }
609     
610   if (fNumberDensitypT)
611     fNumberDensitypT->Merge(lists[0]);
612   if (fSumpT)
613     fSumpT->Merge(lists[1]);
614   if (fNumberDensityPhi)
615     fNumberDensityPhi->Merge(lists[2]);
616   fCorrelationpT->Merge(lists[3]);
617   fCorrelationEta->Merge(lists[4]);
618   fCorrelationPhi->Merge(lists[5]);
619   fCorrelationR->Merge(lists[6]);
620   fCorrelationLeading2Phi->Merge(lists[7]);
621   fCorrelationMultiplicity->Merge(lists[8]);
622   fEventCount->Merge(lists[9]);
623   fEventCountDifferential->Merge(lists[10]);
624   fVertexContributors->Merge(lists[11]);
625   fCentralityDistribution->Merge(lists[12]);
626   
627   for (Int_t i=0; i<kMaxLists; i++)
628     delete lists[i];
629
630   return count+1;
631 }
632
633 void AliUEHistograms::CopyReconstructedData(AliUEHistograms* from)
634 {
635   // copies those histograms extracted from ESD to this object
636   
637   for (Int_t i=0; i<fgkUEHists; i++)
638     if (GetUEHist(i))
639       GetUEHist(i)->CopyReconstructedData(from->GetUEHist(i));
640 }
641
642 void AliUEHistograms::ExtendTrackingEfficiency()
643 {
644   // delegates to AliUEHists
645
646   for (Int_t i=0; i<fgkUEHists; i++)
647     if (GetUEHist(i))
648       GetUEHist(i)->ExtendTrackingEfficiency();
649 }
650