]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemtoUser/AliFemtoAvgSepCorrFctn.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemtoUser / AliFemtoAvgSepCorrFctn.cxx
1 ///////////////////////////////////////////////////////////////////////////
2 //                                                                       //
3 // AliFemtoAvgSepCorrFctn:                                                 //
4 // a simple Q-invariant correlation function                             // 
5 //                                                                       //
6 ///////////////////////////////////////////////////////////////////////////
7
8 #include "AliFemtoAvgSepCorrFctn.h"
9 //#include "AliFemtoHisto.h"
10 #include <cstdio>
11
12 #ifdef __ROOT__ 
13 ClassImp(AliFemtoAvgSepCorrFctn)
14 #endif
15
16 //____________________________
17 AliFemtoAvgSepCorrFctn::AliFemtoAvgSepCorrFctn(char* title, const int& nbins, const float& Low, const float& High):
18 fNumerator(0), //2 tracks
19   fDenominator(0),
20   fNumeratorPos(0), //track + V0
21   fDenominatorPos(0),
22   fNumeratorNeg(0),
23   fDenominatorNeg(0),
24   fNumeratorPosPos(0), //2 V0s
25   fDenominatorPosPos(0),
26   fNumeratorPosNeg(0),
27   fDenominatorPosNeg(0),
28   fNumeratorNegPos(0),
29   fDenominatorNegPos(0),
30   fNumeratorNegNeg(0),
31   fDenominatorNegNeg(0),
32   fRatio(0),
33   fPairType(kTracks)
34 {
35   // set up numerator
36   //  title = "Num Qinv (MeV/c)";
37   char tTitNum[101] = "Num";
38   strncat(tTitNum,title, 100);
39   fNumerator = new TH1D(tTitNum,title,nbins,Low,High);
40   // set up denominator
41   //title = "Den Qinv (MeV/c)";
42   char tTitDen[101] = "Den";
43   strncat(tTitDen,title, 100);
44   fDenominator = new TH1D(tTitDen,title,nbins,Low,High);
45   // set up ratio
46   //title = "Ratio Qinv (MeV/c)";
47   char tTitRat[101] = "Rat";
48   strncat(tTitRat,title, 100);
49   fRatio = new TH1D(tTitRat,title,nbins,Low,High);
50
51
52   char tTitNumPos[101] = "NumV0TrackPos";
53   strncat(tTitNumPos,title, 100);
54   fNumeratorPos = new TH1D(tTitNumPos,title,nbins,Low,High);
55   char tTitDenPos[101] = "DenV0TrackPos";
56   strncat(tTitDenPos,title, 100);
57   fDenominatorPos = new TH1D(tTitDenPos,title,nbins,Low,High);
58   char tTitNumNeg[101] = "NumV0TrackNeg";
59   strncat(tTitNumNeg,title, 100);
60   fNumeratorNeg = new TH1D(tTitNumNeg,title,nbins,Low,High);
61   char tTitDenNeg[101] = "DenV0TrackNeg";
62   strncat(tTitDenNeg,title, 100);
63   fDenominatorNeg = new TH1D(tTitDenNeg,title,nbins,Low,High);
64
65
66   char tTitNumPosPos[101] = "NumV0sPosPos";
67   strncat(tTitNumPosPos,title, 100);
68   fNumeratorPosPos = new TH1D(tTitNumPosPos,title,nbins,Low,High);
69   char tTitDenPosPos[101] = "DenV0sPosPos";
70   strncat(tTitDenPosPos,title, 100);
71   fDenominatorPosPos = new TH1D(tTitDenPosPos,title,nbins,Low,High);
72   char tTitNumPosNeg[101] = "NumV0sPosNeg";
73   strncat(tTitNumPosNeg,title, 100);
74   fNumeratorPosNeg = new TH1D(tTitNumPosNeg,title,nbins,Low,High);
75   char tTitDenPosNeg[101] = "DenV0sPosNeg";
76   strncat(tTitDenPosNeg,title, 100);
77   fDenominatorPosNeg = new TH1D(tTitDenPosNeg,title,nbins,Low,High);
78   char tTitNumNegPos[101] = "NumV0sNegPos";
79   strncat(tTitNumNegPos,title, 100);
80   fNumeratorNegPos = new TH1D(tTitNumNegPos,title,nbins,Low,High);
81   char tTitDenNegPos[101] = "DenV0sNegPos";
82   strncat(tTitDenNegPos,title, 100);
83   fDenominatorNegPos = new TH1D(tTitDenNegPos,title,nbins,Low,High);
84   char tTitNumNegNeg[101] = "NumV0sNegNeg";
85   strncat(tTitNumNegNeg,title, 100);
86   fNumeratorNegNeg = new TH1D(tTitNumNegNeg,title,nbins,Low,High);
87   char tTitDenNegNeg[101] = "DenV0sNegNeg";
88   strncat(tTitDenNegNeg,title, 100);
89   fDenominatorNegNeg = new TH1D(tTitDenNegNeg,title,nbins,Low,High);
90
91
92   // to enable error bar calculation...
93   fNumerator->Sumw2();
94   fDenominator->Sumw2();
95
96   fNumeratorPos->Sumw2();
97   fDenominatorPos->Sumw2();
98   fNumeratorNeg->Sumw2();
99   fDenominatorNeg->Sumw2();
100
101   fNumeratorPosPos->Sumw2();
102   fDenominatorPosPos->Sumw2();
103   fNumeratorPosNeg->Sumw2();
104   fDenominatorPosNeg->Sumw2();
105   fNumeratorNegPos->Sumw2();
106   fDenominatorNegPos->Sumw2();
107   fNumeratorNegNeg->Sumw2();
108   fDenominatorNegNeg->Sumw2();
109
110   fRatio->Sumw2();
111 }
112
113 //____________________________
114 AliFemtoAvgSepCorrFctn::AliFemtoAvgSepCorrFctn(const AliFemtoAvgSepCorrFctn& aCorrFctn) :
115   AliFemtoCorrFctn(),
116   fNumerator(0),
117   fDenominator(0),
118   fNumeratorPos(0),
119   fDenominatorPos(0),
120   fNumeratorNeg(0),
121   fDenominatorNeg(0),
122   fNumeratorPosPos(0), 
123   fDenominatorPosPos(0),
124   fNumeratorPosNeg(0),
125   fDenominatorPosNeg(0),
126   fNumeratorNegPos(0),
127   fDenominatorNegPos(0),
128   fNumeratorNegNeg(0),
129   fDenominatorNegNeg(0),
130   fRatio(0),
131   fPairType(kTracks)
132 {
133   // copy constructor
134   fNumerator = new TH1D(*aCorrFctn.fNumerator);
135   fDenominator = new TH1D(*aCorrFctn.fDenominator);
136
137   fNumeratorPos = new TH1D(*aCorrFctn.fNumerator);
138   fDenominatorPos = new TH1D(*aCorrFctn.fDenominator);
139   fNumeratorNeg = new TH1D(*aCorrFctn.fNumerator);
140   fDenominatorNeg = new TH1D(*aCorrFctn.fDenominator);
141
142   fNumeratorPosPos = new TH1D(*aCorrFctn.fNumerator);
143   fDenominatorPosPos = new TH1D(*aCorrFctn.fDenominator);
144   fNumeratorPosNeg = new TH1D(*aCorrFctn.fNumerator);
145   fDenominatorPosNeg = new TH1D(*aCorrFctn.fDenominator);
146   fNumeratorNegPos = new TH1D(*aCorrFctn.fNumerator);
147   fDenominatorNegPos = new TH1D(*aCorrFctn.fDenominator);
148   fNumeratorNegNeg = new TH1D(*aCorrFctn.fNumerator);
149   fDenominatorNegNeg = new TH1D(*aCorrFctn.fDenominator);
150
151   fRatio = new TH1D(*aCorrFctn.fRatio);
152 }
153 //____________________________
154 AliFemtoAvgSepCorrFctn::~AliFemtoAvgSepCorrFctn(){
155   // destructor
156   delete fNumerator;
157   delete fDenominator;
158
159   delete fNumeratorPos;
160   delete fDenominatorPos;
161   delete fNumeratorNeg;
162   delete fDenominatorNeg;
163
164   delete fNumeratorPosPos;
165   delete fDenominatorPosPos;
166   delete fNumeratorPosNeg;
167   delete fDenominatorPosNeg;
168   delete fNumeratorNegPos;
169   delete fDenominatorNegPos;
170   delete fNumeratorNegNeg;
171   delete fDenominatorNegNeg;
172
173   delete fRatio;
174 }
175 //_________________________
176 AliFemtoAvgSepCorrFctn& AliFemtoAvgSepCorrFctn::operator=(const AliFemtoAvgSepCorrFctn& aCorrFctn)
177 {
178   // assignment operator
179   if (this == &aCorrFctn)
180     return *this;
181
182   if (fNumerator) delete fNumerator;
183   fNumerator = new TH1D(*aCorrFctn.fNumerator);
184   if (fDenominator) delete fDenominator;
185   fDenominator = new TH1D(*aCorrFctn.fDenominator);
186
187   if (fNumeratorPos) delete fNumeratorPos;
188   fNumeratorPos = new TH1D(*aCorrFctn.fNumeratorPos);
189   if (fDenominatorPos) delete fDenominatorPos;
190   fDenominatorPos = new TH1D(*aCorrFctn.fDenominatorPos);
191   if (fNumeratorNeg) delete fNumeratorNeg;
192   fNumeratorNeg = new TH1D(*aCorrFctn.fNumeratorNeg);
193   if (fDenominatorNeg) delete fDenominatorNeg;
194   fDenominatorNeg = new TH1D(*aCorrFctn.fDenominatorNeg);
195
196   if (fNumeratorPosPos) delete fNumeratorPosPos;
197   fNumeratorPosPos = new TH1D(*aCorrFctn.fNumeratorPosPos);
198   if (fDenominatorPosPos) delete fDenominatorPosPos;
199   fDenominatorPosPos = new TH1D(*aCorrFctn.fDenominatorPosPos);
200   if (fNumeratorPosNeg) delete fNumeratorPosNeg;
201   fNumeratorPosNeg = new TH1D(*aCorrFctn.fNumeratorPosNeg);
202   if (fDenominatorPosNeg) delete fDenominatorPosNeg;
203   fDenominatorPosNeg = new TH1D(*aCorrFctn.fDenominatorPosNeg);
204   if (fNumeratorNegPos) delete fNumeratorNegPos;
205   fNumeratorNegPos = new TH1D(*aCorrFctn.fNumeratorNegPos);
206   if (fDenominatorNegPos) delete fDenominatorNegPos;
207   fDenominatorNegPos = new TH1D(*aCorrFctn.fDenominatorNegPos);
208   if (fNumeratorNegNeg) delete fNumeratorNegNeg;
209   fNumeratorNegNeg = new TH1D(*aCorrFctn.fNumeratorNegNeg);
210   if (fDenominatorNegNeg) delete fDenominatorNegNeg;
211   fDenominatorNegNeg = new TH1D(*aCorrFctn.fDenominatorNegNeg);
212
213   if (fRatio) delete fRatio;
214   fRatio = new TH1D(*aCorrFctn.fRatio);
215
216   return *this;
217 }
218
219 //_________________________
220 void AliFemtoAvgSepCorrFctn::Finish(){
221   // here is where we should normalize, fit, etc...
222   // we should NOT Draw() the histos (as I had done it below),
223   // since we want to insulate ourselves from root at this level
224   // of the code.  Do it instead at root command line with browser.
225   //  fNumerator->Draw();
226   //fDenominator->Draw();
227   //fRatio->Draw();
228   fRatio->Divide(fNumerator,fDenominator,1.0,1.0);
229
230 }
231
232 //____________________________
233 AliFemtoString AliFemtoAvgSepCorrFctn::Report(){
234   // construct report
235   string stemp = "Qinv Correlation Function Report:\n";
236   char ctemp[100];
237   snprintf(ctemp , 100, "Number of entries in numerator:\t%E\n",fNumerator->GetEntries());
238   stemp += ctemp;
239   snprintf(ctemp , 100, "Number of entries in denominator:\t%E\n",fDenominator->GetEntries());
240   stemp += ctemp;
241   snprintf(ctemp , 100, "Number of entries in ratio:\t%E\n",fRatio->GetEntries());
242   stemp += ctemp;
243   //  stemp += mCoulombWeight->Report();
244   AliFemtoString returnThis = stemp;
245   return returnThis;
246 }
247 //____________________________
248 void AliFemtoAvgSepCorrFctn::AddRealPair(AliFemtoPair* pair){
249   // add true pair
250   if (fPairCut)
251     if (!fPairCut->Pass(pair)) return;
252   
253   double avgSep=0;
254   AliFemtoThreeVector first, second, tmp;
255
256
257   if(fPairType == 0) //2 tracks
258     {
259       for(int i=0; i<8 ;i++)
260         {
261           tmp = pair->Track1()->Track()->NominalTpcPoint(i);
262           first.SetX((double)(tmp.x()));
263           first.SetY((double)tmp.y());
264           first.SetZ((double)tmp.z());
265       
266           tmp = pair->Track2()->Track()->NominalTpcPoint(i);
267           second.SetX((double)tmp.x());
268           second.SetY((double)tmp.y());
269           second.SetZ((double)tmp.z()); 
270
271           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
272         }
273       avgSep /= 8;
274       fNumerator->Fill(avgSep);
275     }
276   else if(fPairType == 1) // track + V0
277     {
278       for(int i=0; i<8 ;i++)
279         {
280           tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
281           first.SetX((double)(tmp.x()));
282           first.SetY((double)tmp.y());
283           first.SetZ((double)tmp.z());
284           //cout<<"V0pos x: "<<tmp.x()<<", y: "<<tmp.y()<<", z: "<<tmp.z()<<endl;
285           //cout<<"V0pos x: "<<first.x()<<", y: "<<first.y()<<", z: "<<first.z()<<endl;
286       
287           tmp = pair->Track2()->Track()->NominalTpcPoint(i);
288           second.SetX((double)tmp.x());
289           second.SetY((double)tmp.y());
290           second.SetZ((double)tmp.z());
291           //cout<<"Track x: "<<tmp.x()<<", y: "<<tmp.y()<<", z: "<<tmp.z()<<endl;
292           //cout<<"Track x: "<<second.x()<<", y: "<<second.y()<<", z: "<<second.z()<<endl;
293
294           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
295         }
296       //cout<<"****************************************************"<<endl;
297       avgSep /= 8;
298       fNumeratorPos->Fill(avgSep);
299       //cout<<"Track + Pos V0 Avg Sep: "<<avgSep<<endl;
300       avgSep = 0;
301
302       for(int i=0; i<8 ;i++)
303         {
304           tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
305           first.SetX((double)(tmp.x()));
306           first.SetY((double)tmp.y());
307           first.SetZ((double)tmp.z());
308           //cout<<"V0 X: "<<tmp.x()<<endl;
309       
310           tmp = pair->Track2()->Track()->NominalTpcPoint(i);
311           second.SetX((double)tmp.x());
312           second.SetY((double)tmp.y());
313           second.SetZ((double)tmp.z()); 
314           //cout<<"Track X: "<<tmp.x()<<endl;
315
316           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
317         }
318       fNumeratorNeg->Fill(avgSep);
319
320     }
321   else if(fPairType == 2) 
322     {
323       for(int i=0; i<8 ;i++)
324         {
325           tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
326           //cout<<"X pos: "<<tmp.x()<<endl;
327           first.SetX((double)(tmp.x()));
328           first.SetY((double)tmp.y());
329           first.SetZ((double)tmp.z());
330       
331           tmp = pair->Track2()->V0()->NominalTpcPointPos(i);
332           second.SetX((double)tmp.x());
333           second.SetY((double)tmp.y());
334           second.SetZ((double)tmp.z()); 
335
336           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
337         }
338       avgSep /= 8;
339       fNumeratorPosPos->Fill(avgSep);
340       //cout<<"PovV0 + PosV0 Avg Sep: "<<avgSep<<endl;
341
342       avgSep = 0;
343
344       for(int i=0; i<8 ;i++)
345         {
346           tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
347           first.SetX((double)(tmp.x()));
348           first.SetY((double)tmp.y());
349           first.SetZ((double)tmp.z());
350       
351           tmp = pair->Track2()->V0()->NominalTpcPointNeg(i);
352           //cout<<"X neg: "<<tmp.x()<<endl;
353           second.SetX((double)tmp.x());
354           second.SetY((double)tmp.y());
355           second.SetZ((double)tmp.z()); 
356
357           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
358         }
359       avgSep /= 8;
360       fNumeratorPosNeg->Fill(avgSep);
361
362       avgSep = 0;
363
364       for(int i=0; i<8 ;i++)
365         {
366           tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
367           first.SetX((double)(tmp.x()));
368           first.SetY((double)tmp.y());
369           first.SetZ((double)tmp.z());
370       
371           tmp = pair->Track2()->V0()->NominalTpcPointPos(i);
372           second.SetX((double)tmp.x());
373           second.SetY((double)tmp.y());
374           second.SetZ((double)tmp.z()); 
375
376           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
377         }
378       avgSep /= 8;
379       fNumeratorNegPos->Fill(avgSep);
380       
381
382       avgSep = 0;
383  
384      for(int i=0; i<8 ;i++)
385         {
386           tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
387           first.SetX((double)(tmp.x()));
388           first.SetY((double)tmp.y());
389           first.SetZ((double)tmp.z());
390       
391           tmp = pair->Track2()->V0()->NominalTpcPointNeg(i);
392           second.SetX((double)tmp.x());
393           second.SetY((double)tmp.y());
394           second.SetZ((double)tmp.z()); 
395
396           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
397         }
398       avgSep /= 8;
399       fNumeratorNegNeg->Fill(avgSep);
400
401
402     } 
403   
404
405   //2 daughters of first V0
406
407
408   
409  
410
411 }
412 //____________________________
413 void AliFemtoAvgSepCorrFctn::AddMixedPair(AliFemtoPair* pair){
414   // add mixed (background) pair
415   if (fPairCut)
416     if (!fPairCut->Pass(pair)) return;
417   
418   double avgSep=0;
419   AliFemtoThreeVector first, second, tmp;
420
421
422   if(fPairType == 0) //2 tracks
423     {
424       for(int i=0; i<8 ;i++)
425         {
426           tmp = pair->Track1()->Track()->NominalTpcPoint(i);
427           first.SetX((double)(tmp.x()));
428           first.SetY((double)tmp.y());
429           first.SetZ((double)tmp.z());
430       
431           tmp = pair->Track2()->Track()->NominalTpcPoint(i);
432           second.SetX((double)tmp.x());
433           second.SetY((double)tmp.y());
434           second.SetZ((double)tmp.z()); 
435
436           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
437         }
438       avgSep /= 8;
439       fDenominator->Fill(avgSep);
440     }
441   else if(fPairType == 1) // track + V0
442     {
443       for(int i=0; i<8 ;i++)
444         {
445           tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
446           first.SetX((double)(tmp.x()));
447           first.SetY((double)tmp.y());
448           first.SetZ((double)tmp.z());
449       
450           tmp = pair->Track2()->Track()->NominalTpcPoint(i);
451           second.SetX((double)tmp.x());
452           second.SetY((double)tmp.y());
453           second.SetZ((double)tmp.z()); 
454
455           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
456         }
457       avgSep /= 8;
458       fDenominatorPos->Fill(avgSep);
459
460       avgSep = 0;
461
462       for(int i=0; i<8 ;i++)
463         {
464           tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
465           first.SetX((double)(tmp.x()));
466           first.SetY((double)tmp.y());
467           first.SetZ((double)tmp.z());
468       
469           tmp = pair->Track2()->Track()->NominalTpcPoint(i);
470           second.SetX((double)tmp.x());
471           second.SetY((double)tmp.y());
472           second.SetZ((double)tmp.z()); 
473
474           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
475         }
476       fDenominatorNeg->Fill(avgSep);   
477
478     }
479   else if(fPairType == 2) 
480     {
481       for(int i=0; i<8 ;i++)
482         {
483           
484           tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
485           first.SetX((double)(tmp.x()));
486           first.SetY((double)tmp.y());
487           first.SetZ((double)tmp.z());
488       
489           tmp = pair->Track2()->V0()->NominalTpcPointPos(i);
490           second.SetX((double)tmp.x());
491           second.SetY((double)tmp.y());
492           second.SetZ((double)tmp.z()); 
493
494           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
495         }
496       avgSep /= 8;
497       fDenominatorPosPos->Fill(avgSep);
498
499       avgSep = 0;
500
501       for(int i=0; i<8 ;i++)
502         {
503           tmp = pair->Track1()->V0()->NominalTpcPointPos(i);
504           first.SetX((double)(tmp.x()));
505           first.SetY((double)tmp.y());
506           first.SetZ((double)tmp.z());
507       
508           tmp = pair->Track2()->V0()->NominalTpcPointNeg(i);
509           second.SetX((double)tmp.x());
510           second.SetY((double)tmp.y());
511           second.SetZ((double)tmp.z()); 
512
513           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
514         }
515       avgSep /= 8;
516       fDenominatorPosNeg->Fill(avgSep);
517
518       avgSep = 0;
519
520       for(int i=0; i<8 ;i++)
521         {
522           tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
523           first.SetX((double)(tmp.x()));
524           first.SetY((double)tmp.y());
525           first.SetZ((double)tmp.z());
526       
527           tmp = pair->Track2()->V0()->NominalTpcPointPos(i);
528           second.SetX((double)tmp.x());
529           second.SetY((double)tmp.y());
530           second.SetZ((double)tmp.z()); 
531
532           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
533         }
534       avgSep /= 8;
535       fDenominatorNegPos->Fill(avgSep);
536
537       avgSep = 0;
538  
539      for(int i=0; i<8 ;i++)
540         {
541           tmp = pair->Track1()->V0()->NominalTpcPointNeg(i);
542           first.SetX((double)(tmp.x()));
543           first.SetY((double)tmp.y());
544           first.SetZ((double)tmp.z());
545       
546           tmp = pair->Track2()->V0()->NominalTpcPointNeg(i);
547           second.SetX((double)tmp.x());
548           second.SetY((double)tmp.y());
549           second.SetZ((double)tmp.z()); 
550
551           avgSep += TMath::Sqrt(((double)first.x()-(double)second.x())*((double)first.x()-(double)second.x())+((double)first.y()-(double)second.y())*((double)first.y()-second.y())+((double)first.z()-(double)second.z())*((double)first.z()-(double)second.z()));
552         }
553       avgSep /= 8;
554       fDenominatorNegNeg->Fill(avgSep);
555
556
557     } 
558
559 }
560 //____________________________
561 void AliFemtoAvgSepCorrFctn::Write(){
562   // Write out neccessary objects
563   if(fPairType == kTracks)
564     {
565       fNumerator->Write(); 
566       fDenominator->Write();  
567       fRatio->Write();
568     }
569   else if(fPairType == kTrackV0)
570     {
571       fNumeratorPos->Write(); 
572       fDenominatorPos->Write();  
573       fNumeratorNeg->Write(); 
574       fDenominatorNeg->Write();  
575     }
576   else if(fPairType == kV0s)
577     {
578       fNumeratorPosPos->Write(); 
579       fDenominatorPosPos->Write();  
580       fNumeratorPosNeg->Write(); 
581       fDenominatorPosNeg->Write();  
582       fNumeratorNegPos->Write(); 
583       fDenominatorNegPos->Write();  
584       fNumeratorNegNeg->Write(); 
585       fDenominatorNegNeg->Write();  
586     }
587
588  
589 }
590 //______________________________
591 TList* AliFemtoAvgSepCorrFctn::GetOutputList()
592 {
593   // Prepare the list of objects to be written to the output
594   TList *tOutputList = new TList();
595
596   if(fPairType == kTracks)
597     {
598       tOutputList->Add(fNumerator);
599       tOutputList->Add(fDenominator);
600       tOutputList->Add(fRatio);
601     }
602   else if(fPairType == kTrackV0)
603     {
604       tOutputList->Add(fNumeratorPos); 
605       tOutputList->Add(fDenominatorPos);  
606       tOutputList->Add(fNumeratorNeg); 
607       tOutputList->Add(fDenominatorNeg);  
608     }
609   else if(fPairType == kV0s)
610     {
611       tOutputList->Add(fNumeratorPosPos); 
612       tOutputList->Add(fDenominatorPosPos);  
613       tOutputList->Add(fNumeratorPosNeg); 
614       tOutputList->Add(fDenominatorPosNeg);  
615       tOutputList->Add(fNumeratorNegPos); 
616       tOutputList->Add(fDenominatorNegPos);  
617       tOutputList->Add(fNumeratorNegNeg); 
618       tOutputList->Add(fDenominatorNegNeg);  
619     }
620   return tOutputList;
621 }
622
623
624 void AliFemtoAvgSepCorrFctn::SetPairType(AliFemtoPairType pairtype)
625 {
626   fPairType = pairtype;
627 }