]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemtoUser/AliFemtoModelCorrFctnDEtaDPhi.cxx
Migration of PWG2/FEMTOSCOPY to PWGCF/FEMTOSCOPY
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemtoUser / AliFemtoModelCorrFctnDEtaDPhi.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoModelCorrFctnDEtaDPhi - A correlation function that analyzes            //
4 // two particle correlations with respect to the azimuthal angle (phi)        //
5 // and pseudorapidity (eta) difference                                        //
6 //                                                                            //
7 // Authors: Adam Kisiel Adam.Kisiel@cern.ch                                   //
8 //                                                                            //
9 ////////////////////////////////////////////////////////////////////////////////
10
11 #include "AliFemtoModelCorrFctnDEtaDPhi.h"
12 //#include "AliFemtoHisto.hh"
13 #include <cstdio>
14 #include <TMath.h>
15
16 #ifdef __ROOT__ 
17 ClassImp(AliFemtoModelCorrFctnDEtaDPhi)
18 #endif
19
20 /*And some Model libraries..*/
21 //1
22 //#include "AliFemtoModelGausLCMSFreezeOutGenerator.h"
23 //#include "AliFemtoModelHiddenInfo.h"
24 //2
25 #include "AliFemtoModelGausLCMSFreezeOutGenerator.h"
26 #include "AliFemtoModelHiddenInfo.h"
27 #include "AliFemtoPair.h"
28 #include "AliFemtoModelManager.h"
29 //
30
31 //____________________________
32 AliFemtoModelCorrFctnDEtaDPhi::AliFemtoModelCorrFctnDEtaDPhi(char* title, const int& aPhiBins=20, const int& aEtaBins=20):
33   AliFemtoModelCorrFctn(),
34   fDPhiDEtaNumeratorTrue(0),
35   fDPhiDEtaNumeratorFake(0),
36   fDPhiDEtaDenominator(0),
37   fDPhiDEtaColNumerator(0),
38   fDPhiDEtaColDenominator(0),
39   fDPhiNumeratorTrue(0),
40   fDPhiNumeratorFake(0),
41   fDPhiDenominator(0),
42   fDCosNumeratorTrue(0),
43   fDCosNumeratorFake(0),
44   fDCosDenominator(0),
45   fDPhiPtNumerator(0),
46   fDPhiPtDenominator(0),
47   fDCosPtNumerator(0),
48   fDCosPtDenominator(0)
49 {
50   // set up numerator
51   char tTitNumDT[101] = "NumDPhiDEtaTrue";
52   strncat(tTitNumDT,title, 100);
53   fDPhiDEtaNumeratorTrue = new TH2D(tTitNumDT,title,aPhiBins,-0.5*TMath::Pi(),1.5*TMath::Pi(),aEtaBins,-2.0,2.0);
54
55   char tTitNumDF[101] = "NumDPhiDEtaFake";
56   strncat(tTitNumDF,title, 100);
57   fDPhiDEtaNumeratorFake = new TH2D(tTitNumDF,title,aPhiBins,-0.5*TMath::Pi(),1.5*TMath::Pi(),aEtaBins,-2.0,2.0);
58
59
60   // set up denominator
61   char tTitDenD[101] = "DenDPhiDEta";
62   strncat(tTitDenD,title, 100);
63   fDPhiDEtaDenominator = new TH2D(tTitDenD,title,aPhiBins,-0.5*TMath::Pi(),1.5*TMath::Pi(),aEtaBins,-2.0,2.0);
64
65   // set up numerator
66   char tTitNumR[101] = "NumDPhiDEtaCol";
67   strncat(tTitNumR,title, 100);
68   fDPhiDEtaColNumerator = new TH2D(tTitNumR,title,aPhiBins,-0.5*TMath::Pi(),1.5*TMath::Pi(),aEtaBins,-2.0,2.0);
69   // set up denominator
70   char tTitDenR[101] = "DenDPhiDEtaCol";
71   strncat(tTitDenR,title, 100);
72   fDPhiDEtaColDenominator = new TH2D(tTitDenR,title,aPhiBins,-0.5*TMath::Pi(),1.5*TMath::Pi(),aEtaBins,-2.0,2.0);
73
74   // set up numerator
75   char tTitNumDPhiT[101] = "NumDPhiTrue";
76   strncat(tTitNumDPhiT,title, 100);
77   fDPhiNumeratorTrue = new TH1D(tTitNumDPhiT,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi());
78   // set up numerator
79   char tTitNumDPhiF[101] = "NumDPhiFake";
80   strncat(tTitNumDPhiF,title, 100);
81   fDPhiNumeratorFake = new TH1D(tTitNumDPhiF,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi());
82
83   // set up denominator
84   char tTitDenDPhi[101] = "DenDPhi";
85   strncat(tTitDenDPhi,title, 100);
86   fDPhiDenominator = new TH1D(tTitDenDPhi,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi());
87
88   // set up numerator
89   char tTitNumDCosT[101] = "NumDCosTrue";
90   strncat(tTitNumDCosT,title, 100);
91   fDCosNumeratorTrue = new TH1D(tTitNumDCosT,title,aPhiBins*2,-1.0,1.0);
92   // set up numerator
93   char tTitNumDCosF[101] = "NumDCosFake";
94   strncat(tTitNumDCosF,title, 100);
95   fDCosNumeratorFake = new TH1D(tTitNumDCosF,title,aPhiBins*2,-1.0,1.0);
96
97   // set up denominator
98   char tTitDenDCos[101] = "DenDCos";
99   strncat(tTitDenDCos,title, 100);
100   fDCosDenominator = new TH1D(tTitDenDCos,title,aPhiBins*2,-1.0,1.0);
101
102   // set up numerator
103   char tTitNumDPhiPt[101] = "NumDPhiPt";
104   strncat(tTitNumDPhiPt,title, 100);
105   fDPhiPtNumerator = new TH2D(tTitNumDPhiPt,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi(), 30, 0.0, 3.0);
106   // set up denominator
107   char tTitDenDPhiPt[101] = "DenDPhiPt";
108   strncat(tTitDenDPhiPt,title, 100);
109   fDPhiPtDenominator = new TH2D(tTitDenDPhiPt,title,aPhiBins*2,-0.5*TMath::Pi(),1.5*TMath::Pi(), 30, 0.0, 3.0);
110
111   // set up numerator
112   char tTitNumDCosPt[101] = "NumDCosPt";
113   strncat(tTitNumDCosPt,title, 100);
114   fDCosPtNumerator = new TH2D(tTitNumDCosPt,title,aPhiBins*2,-1.0,1.0, 30, 0.0, 3.0);
115   // set up denominator
116   char tTitDenDCosPt[101] = "DenDCosPt";
117   strncat(tTitDenDCosPt,title, 100);
118   fDCosPtDenominator = new TH2D(tTitDenDCosPt,title,aPhiBins*2,-1.0,1.0, 30, 0.0, 3.0);
119
120   // to enable error bar calculation...
121   fDPhiDEtaNumeratorTrue->Sumw2();
122   fDPhiDEtaNumeratorFake->Sumw2();
123   fDPhiDEtaDenominator->Sumw2();
124   fDPhiDEtaColNumerator->Sumw2();
125   fDPhiDEtaColDenominator->Sumw2();
126   fDPhiNumeratorTrue->Sumw2();
127   fDPhiNumeratorFake->Sumw2();
128   fDPhiDenominator->Sumw2();
129   fDCosNumeratorTrue->Sumw2();
130   fDCosNumeratorFake->Sumw2();
131   fDCosDenominator->Sumw2();
132   fDPhiPtNumerator->Sumw2();
133   fDPhiPtDenominator->Sumw2();
134   fDCosPtNumerator->Sumw2();
135   fDCosPtDenominator->Sumw2();
136
137 }
138
139 //____________________________
140 AliFemtoModelCorrFctnDEtaDPhi::AliFemtoModelCorrFctnDEtaDPhi(const AliFemtoModelCorrFctnDEtaDPhi& aCorrFctn) :
141   AliFemtoModelCorrFctn(),
142   fDPhiDEtaNumeratorTrue(0),
143   fDPhiDEtaNumeratorFake(0),
144   fDPhiDEtaDenominator(0),
145   fDPhiDEtaColNumerator(0),
146   fDPhiDEtaColDenominator(0),
147   fDPhiNumeratorTrue(0),
148   fDPhiNumeratorFake(0),
149   fDPhiDenominator(0),
150   fDCosNumeratorTrue(0),
151   fDCosNumeratorFake(0),
152   fDCosDenominator(0),
153   fDPhiPtNumerator(0),
154   fDPhiPtDenominator(0),
155   fDCosPtNumerator(0),
156   fDCosPtDenominator(0)
157 {
158   // copy constructor
159   if (aCorrFctn.fDPhiDEtaNumeratorTrue)
160     fDPhiDEtaNumeratorTrue = new TH2D(*aCorrFctn.fDPhiDEtaNumeratorTrue);
161   else
162     fDPhiDEtaNumeratorTrue = 0;
163   if (aCorrFctn.fDPhiDEtaNumeratorFake)
164     fDPhiDEtaNumeratorFake = new TH2D(*aCorrFctn.fDPhiDEtaNumeratorFake);
165   else
166     fDPhiDEtaNumeratorFake = 0;
167
168   if (aCorrFctn.fDPhiDEtaDenominator)
169     fDPhiDEtaDenominator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
170   else
171     fDPhiDEtaDenominator = 0;
172
173   if (aCorrFctn.fDPhiDEtaColNumerator)
174     fDPhiDEtaColNumerator = new TH2D(*aCorrFctn.fDPhiDEtaColNumerator);
175   else
176     fDPhiDEtaColNumerator = 0;
177   if (aCorrFctn.fDPhiDEtaColDenominator)
178     fDPhiDEtaColDenominator = new TH2D(*aCorrFctn.fDPhiDEtaColDenominator);
179   else
180     fDPhiDEtaColDenominator = 0;
181
182   if (aCorrFctn.fDPhiNumeratorTrue)
183     fDPhiNumeratorTrue = new TH1D(*aCorrFctn.fDPhiNumeratorTrue);
184   else
185     fDPhiNumeratorTrue = 0;
186
187   if (aCorrFctn.fDPhiNumeratorFake)
188     fDPhiNumeratorFake = new TH1D(*aCorrFctn.fDPhiNumeratorFake);
189   else
190     fDPhiNumeratorFake = 0;
191
192   if (aCorrFctn.fDPhiDenominator)
193     fDPhiDenominator = new TH1D(*aCorrFctn.fDPhiDenominator);
194   else
195     fDPhiDenominator = 0;
196
197   if (aCorrFctn.fDCosNumeratorTrue)
198     fDCosNumeratorTrue = new TH1D(*aCorrFctn.fDCosNumeratorTrue);
199   else
200     fDCosNumeratorTrue = 0;
201   if (aCorrFctn.fDCosNumeratorFake)
202     fDCosNumeratorFake = new TH1D(*aCorrFctn.fDCosNumeratorFake);
203   else
204     fDCosNumeratorFake = 0;
205
206   if (aCorrFctn.fDCosDenominator)
207     fDCosDenominator = new TH1D(*aCorrFctn.fDCosDenominator);
208   else
209     fDCosDenominator = 0;
210
211   if (aCorrFctn.fDPhiPtNumerator)
212     fDPhiPtNumerator = new TH2D(*aCorrFctn.fDPhiPtNumerator);
213   else
214     fDPhiPtNumerator = 0;
215   if (aCorrFctn.fDPhiPtDenominator)
216     fDPhiPtDenominator = new TH2D(*aCorrFctn.fDPhiPtDenominator);
217   else
218     fDPhiPtDenominator = 0;
219
220   if (aCorrFctn.fDCosPtNumerator)
221     fDCosPtNumerator = new TH2D(*aCorrFctn.fDCosPtNumerator);
222   else
223     fDCosPtNumerator = 0;
224   if (aCorrFctn.fDCosPtDenominator)
225     fDCosPtDenominator = new TH2D(*aCorrFctn.fDCosPtDenominator);
226   else
227     fDCosPtDenominator = 0;
228 }
229 //____________________________
230 AliFemtoModelCorrFctnDEtaDPhi::~AliFemtoModelCorrFctnDEtaDPhi(){
231   // destructor
232   if(fDPhiDEtaNumeratorTrue)      delete fDPhiDEtaNumeratorTrue;
233   if(fDPhiDEtaNumeratorFake)      delete fDPhiDEtaNumeratorFake;
234   if(fDPhiDEtaDenominator)    delete fDPhiDEtaDenominator;
235   if(fDPhiDEtaColNumerator)   delete fDPhiDEtaColNumerator;
236   if(fDPhiDEtaColDenominator) delete fDPhiDEtaColDenominator;
237   if(fDPhiNumeratorTrue)       delete fDPhiNumeratorTrue;
238   if(fDPhiNumeratorFake)       delete fDPhiNumeratorFake;
239   if(fDPhiDenominator)     delete fDPhiDenominator;
240   if(fDCosNumeratorTrue)       delete fDCosNumeratorTrue;
241   if(fDCosNumeratorFake)       delete fDCosNumeratorFake;
242   if(fDCosDenominator)     delete fDCosDenominator;
243   if(fDPhiPtNumerator)     delete fDPhiPtNumerator;
244   if(fDPhiPtDenominator)   delete fDPhiPtDenominator;
245   if(fDCosPtNumerator)     delete fDCosPtNumerator;
246   if(fDCosPtDenominator)   delete fDCosPtDenominator;
247 }
248 //_________________________
249 AliFemtoModelCorrFctnDEtaDPhi& AliFemtoModelCorrFctnDEtaDPhi::operator=(const AliFemtoModelCorrFctnDEtaDPhi& aCorrFctn)
250 {
251   // assignment operator
252   if (this == &aCorrFctn)
253     return *this;
254
255   if (aCorrFctn.fDPhiDEtaNumeratorTrue)
256     fDPhiDEtaNumeratorTrue = new TH2D(*aCorrFctn.fDPhiDEtaNumeratorTrue);
257   else
258     fDPhiDEtaNumeratorTrue = 0;
259
260   if (aCorrFctn.fDPhiDEtaNumeratorFake)
261     fDPhiDEtaNumeratorFake = new TH2D(*aCorrFctn.fDPhiDEtaNumeratorFake);
262   else
263     fDPhiDEtaNumeratorFake = 0;
264
265   if (aCorrFctn.fDPhiDEtaDenominator)
266     fDPhiDEtaDenominator = new TH2D(*aCorrFctn.fDPhiDEtaDenominator);
267   else
268     fDPhiDEtaDenominator = 0;
269
270   if (aCorrFctn.fDPhiDEtaColNumerator)
271     fDPhiDEtaColNumerator = new TH2D(*aCorrFctn.fDPhiDEtaColNumerator);
272   else
273     fDPhiDEtaColNumerator = 0;
274   if (aCorrFctn.fDPhiDEtaColDenominator)
275     fDPhiDEtaColDenominator = new TH2D(*aCorrFctn.fDPhiDEtaColDenominator);
276   else
277     fDPhiDEtaColDenominator = 0;
278
279   if (aCorrFctn.fDPhiNumeratorTrue)
280     fDPhiNumeratorTrue = new TH1D(*aCorrFctn.fDPhiNumeratorTrue);
281   else
282     fDPhiNumeratorTrue = 0;
283   if (aCorrFctn.fDPhiNumeratorFake)
284     fDPhiNumeratorFake = new TH1D(*aCorrFctn.fDPhiNumeratorFake);
285   else
286     fDPhiNumeratorFake = 0;
287
288   if (aCorrFctn.fDPhiDenominator)
289     fDPhiDenominator = new TH1D(*aCorrFctn.fDPhiDenominator);
290   else
291     fDPhiDenominator = 0;
292
293   if (aCorrFctn.fDCosNumeratorTrue)
294     fDCosNumeratorTrue = new TH1D(*aCorrFctn.fDCosNumeratorTrue);
295   else
296     fDCosNumeratorTrue = 0;
297   if (aCorrFctn.fDCosNumeratorFake)
298     fDCosNumeratorFake = new TH1D(*aCorrFctn.fDCosNumeratorFake);
299   else
300     fDCosNumeratorFake = 0;
301
302   if (aCorrFctn.fDCosDenominator)
303     fDCosDenominator = new TH1D(*aCorrFctn.fDCosDenominator);
304   else
305     fDCosDenominator = 0;
306
307   if (aCorrFctn.fDPhiPtNumerator)
308     fDPhiPtNumerator = new TH2D(*aCorrFctn.fDPhiPtNumerator);
309   else
310     fDPhiPtNumerator = 0;
311   if (aCorrFctn.fDPhiPtDenominator)
312     fDPhiPtDenominator = new TH2D(*aCorrFctn.fDPhiPtDenominator);
313   else
314     fDPhiPtDenominator = 0;
315
316   if (aCorrFctn.fDCosPtNumerator)
317     fDCosPtNumerator = new TH2D(*aCorrFctn.fDCosPtNumerator);
318   else
319     fDCosPtNumerator = 0;
320   if (aCorrFctn.fDCosPtDenominator)
321     fDCosPtDenominator = new TH2D(*aCorrFctn.fDCosPtDenominator);
322   else
323     fDCosPtDenominator = 0;
324
325   return *this;
326 }
327 //_________________________
328 void AliFemtoModelCorrFctnDEtaDPhi::Finish(){
329   // here is where we should normalize, fit, etc...
330   // we should NOT Draw() the histos (as I had done it below),
331   // since we want to insulate ourselves from root at this level
332   // of the code.  Do it instead at root command line with browser.
333   //  mShareNumerator->Draw();
334   //mShareDenominator->Draw();
335   //mRatio->Draw();
336
337 }
338
339 //____________________________
340 AliFemtoString AliFemtoModelCorrFctnDEtaDPhi::Report(){
341   // create report
342   string stemp = "TPC Ncls Correlation Function Report:\n";
343   char ctemp[100];
344   snprintf(ctemp,100,"Number of entries in numerator true:\t%E\n",fDPhiDEtaNumeratorTrue->GetEntries());
345   snprintf(ctemp,100,"Number of entries in numerator fake:\t%E\n",fDPhiDEtaNumeratorFake->GetEntries());
346   stemp += ctemp;
347   snprintf(ctemp,100,"Number of entries in denominator:\t%E\n",fDPhiDEtaDenominator->GetEntries());
348   stemp += ctemp;
349   //  stemp += mCoulombWeight->Report();
350   AliFemtoString returnThis = stemp;
351   return returnThis;
352 }
353 //____________________________
354 void AliFemtoModelCorrFctnDEtaDPhi::AddRealPair( AliFemtoPair* pair){
355   // add real (effect) pair
356   double phi1 = pair->Track1()->Track()->P().Phi();
357   double phi2 = pair->Track2()->Track()->P().Phi();
358   double eta1 = pair->Track1()->Track()->P().PseudoRapidity();
359   double eta2 = pair->Track2()->Track()->P().PseudoRapidity();
360
361   double dphi = phi1 - phi2;
362   while (dphi<-TMath::Pi()/2) dphi+=TMath::Pi()*2;
363   while (dphi>3*TMath::Pi()/2) dphi-=TMath::Pi()*2;
364
365   double deta = eta1 - eta2;
366
367   double px1 = pair->Track1()->Track()->P().x();
368   double py1 = pair->Track1()->Track()->P().y();
369   double pz1 = pair->Track1()->Track()->P().z();
370
371   double px2 = pair->Track2()->Track()->P().x();
372   double py2 = pair->Track2()->Track()->P().y();
373   double pz2 = pair->Track2()->Track()->P().z();
374
375   double pt1 = TMath::Hypot(px1, py1);
376   double pt2 = TMath::Hypot(px2, py2);
377   double ptmin = pt1>pt2 ? pt2 : pt1;
378
379   double cosphi = (px1*px2 + py1*py2 + pz1*pz2)/
380     sqrt((px1*px1 + py1*py1 + pz1*pz1)*(px2*px2 + py2*py2 + pz2*pz2));
381
382 /*weights*/
383   Double_t weight = fManager->GetWeight(pair);
384   fDPhiDEtaNumeratorTrue->Fill(dphi, deta,weight);
385
386   if (cosphi > 0) {
387     fDPhiDEtaColNumerator->Fill(dphi, deta);
388   }
389   else {
390     fDPhiDEtaColNumerator->Fill(dphi, -eta1-eta2);
391   }
392
393   fDPhiNumeratorTrue->Fill(dphi,weight);
394   fDCosNumeratorTrue->Fill(cosphi,weight);
395
396   fDPhiPtNumerator->Fill(dphi, ptmin);
397   fDCosPtNumerator->Fill(cosphi, ptmin);
398
399 }
400 //____________________________
401 void AliFemtoModelCorrFctnDEtaDPhi::AddMixedPair( AliFemtoPair* pair){
402   // add mixed (background) pair
403   double phi1 = pair->Track1()->Track()->P().Phi();
404   double phi2 = pair->Track2()->Track()->P().Phi();
405   double eta1 = pair->Track1()->Track()->P().PseudoRapidity();
406   double eta2 = pair->Track2()->Track()->P().PseudoRapidity();
407
408   double dphi = phi1 - phi2;
409   while (dphi<-TMath::Pi()/2) dphi+=TMath::Pi()*2;
410   while (dphi>3*TMath::Pi()/2) dphi-=TMath::Pi()*2;
411
412   double deta = eta1 - eta2;
413
414   double px1 = pair->Track1()->Track()->P().x();
415   double py1 = pair->Track1()->Track()->P().y();
416   double pz1 = pair->Track1()->Track()->P().z();
417
418   double px2 = pair->Track2()->Track()->P().x();
419   double py2 = pair->Track2()->Track()->P().y();
420   double pz2 = pair->Track2()->Track()->P().z();
421
422   double pt1 = TMath::Hypot(px1, py1);
423   double pt2 = TMath::Hypot(px2, py2);
424   double ptmin = pt1>pt2 ? pt2 : pt1;
425
426   double cosphi = (px1*px2 + py1*py2 + pz1*pz2)/
427     sqrt((px1*px1 + py1*py1 + pz1*pz1)*(px2*px2 + py2*py2 + pz2*pz2));
428
429
430   Double_t weight = fManager->GetWeight(pair);
431   fDPhiDEtaNumeratorFake->Fill(dphi, deta,weight);
432
433   fDPhiDEtaDenominator->Fill(dphi, deta,1.0);
434
435   if (cosphi > 0) {
436     fDPhiDEtaColDenominator->Fill(dphi, deta);
437   }
438   else {
439     fDPhiDEtaColDenominator->Fill(dphi, -eta1-eta2);
440   }
441
442   fDPhiNumeratorFake->Fill(dphi,weight);
443   fDCosNumeratorFake->Fill(cosphi,weight);
444
445   fDPhiDenominator->Fill(dphi,1.0);
446   fDCosDenominator->Fill(cosphi,1.0);
447
448   fDPhiPtDenominator->Fill(dphi, ptmin);
449   fDCosPtDenominator->Fill(cosphi, ptmin);
450 }
451
452
453 void AliFemtoModelCorrFctnDEtaDPhi::WriteHistos()
454 {
455   // Write out result histograms
456   fDPhiDEtaNumeratorTrue->Write();
457   fDPhiDEtaNumeratorFake->Write();
458   fDPhiDEtaDenominator->Write();
459   fDPhiDEtaColNumerator->Write();
460   fDPhiDEtaColDenominator->Write();
461   fDPhiNumeratorTrue->Write();
462   fDPhiNumeratorFake->Write();
463   fDPhiDenominator->Write();
464   fDCosNumeratorTrue->Write();
465   fDCosNumeratorFake->Write();
466   fDCosDenominator->Write();
467   fDPhiPtNumerator->Write();
468   fDPhiPtDenominator->Write();
469   fDCosPtNumerator->Write();
470   fDCosPtDenominator->Write();
471 }
472
473 TList* AliFemtoModelCorrFctnDEtaDPhi::GetOutputList()
474 {
475   // Prepare the list of objects to be written to the output
476   TList *tOutputList = new TList();
477
478   tOutputList->Add(fDPhiDEtaNumeratorTrue);
479   tOutputList->Add(fDPhiDEtaNumeratorFake);
480   tOutputList->Add(fDPhiDEtaDenominator);
481   tOutputList->Add(fDPhiDEtaColNumerator);
482   tOutputList->Add(fDPhiDEtaColDenominator);
483   tOutputList->Add(fDPhiNumeratorTrue);
484   tOutputList->Add(fDPhiNumeratorFake);
485   tOutputList->Add(fDPhiDenominator);
486   tOutputList->Add(fDCosNumeratorTrue);
487   tOutputList->Add(fDCosNumeratorFake);
488   tOutputList->Add(fDCosDenominator);
489   tOutputList->Add(fDPhiPtNumerator);
490   tOutputList->Add(fDPhiPtDenominator);
491   tOutputList->Add(fDCosPtNumerator);
492   tOutputList->Add(fDCosPtDenominator);
493
494   return tOutputList;
495
496 }