]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctnNonIdDR.cxx
Lines getting the matched track moved to a method in AliCalorimeterUtils. Lines copie...
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / AliFemtoCorrFctnNonIdDR.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoCorrFctnNonIdDR - correlation function for non-identical particles //
4 // uses k* as a function variable. Stores the correlation function separately //
5 // for positive and negative signs of k* projections into out, side and long  //
6 // directions, enabling the calculations of double ratios                     //
7 //                                                                            //
8 ////////////////////////////////////////////////////////////////////////////////
9
10 #include "AliFemtoCorrFctnNonIdDR.h"
11 //#include "AliFemtoHisto.h"
12 #include <cstdio>
13
14 #ifdef __ROOT__ 
15 ClassImp(AliFemtoCorrFctnNonIdDR)
16 #endif
17
18 //____________________________
19 AliFemtoCorrFctnNonIdDR::AliFemtoCorrFctnNonIdDR(char* title, const int& nbins, const float& QinvLo, const float& QinvHi):
20   fNumOutP(0), 
21   fNumOutN(0),  
22   fNumSideP(0), 
23   fNumSideN(0), 
24   fNumLongP(0), 
25   fNumLongN(0), 
26   fDenOutP(0),  
27   fDenOutN(0),  
28   fDenSideP(0), 
29   fDenSideN(0), 
30   fDenLongP(0), 
31   fDenLongN(0)
32 {
33   // Default constructor
34   // set up numerators
35   char bufname[200];
36   snprintf(bufname, 200, "NumOutP%s", title);
37   fNumOutP = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
38   snprintf(bufname, 200, "NumOutN%s", title);
39   fNumOutN = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
40   snprintf(bufname, 200, "NumSideP%s", title);
41   fNumSideP = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
42   snprintf(bufname, 200, "NumSideN%s", title);
43   fNumSideN = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
44   snprintf(bufname, 200, "NumLongP%s", title);
45   fNumLongP = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
46   snprintf(bufname, 200, "NumLongN%s", title);
47   fNumLongN = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
48
49   // set up denominators
50   snprintf(bufname, 200, "DenOutP%s", title);
51   fDenOutP = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
52   snprintf(bufname, 200, "DenOutN%s", title);
53   fDenOutN = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
54   snprintf(bufname, 200, "DenSideP%s", title);
55   fDenSideP = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
56   snprintf(bufname, 200, "DenSideN%s", title);
57   fDenSideN = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
58   snprintf(bufname, 200, "DenLongP%s", title);
59   fDenLongP = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
60   snprintf(bufname, 200, "DenLongN%s", title);
61   fDenLongN = new TH1D(bufname,title,nbins,QinvLo,QinvHi);
62
63   // to enable error bar calculation...
64   fNumOutP->Sumw2(); 
65   fNumOutN->Sumw2();  
66   fNumSideP->Sumw2(); 
67   fNumSideN->Sumw2(); 
68   fNumLongP->Sumw2(); 
69   fNumLongN->Sumw2(); 
70   fDenOutP->Sumw2();  
71   fDenOutN->Sumw2();  
72   fDenSideP->Sumw2(); 
73   fDenSideN->Sumw2(); 
74   fDenLongP->Sumw2(); 
75   fDenLongN->Sumw2();
76 }
77
78 //____________________________
79 AliFemtoCorrFctnNonIdDR::AliFemtoCorrFctnNonIdDR(const AliFemtoCorrFctnNonIdDR& aCorrFctn) :
80   AliFemtoCorrFctn(),
81   fNumOutP(0), 
82   fNumOutN(0),  
83   fNumSideP(0), 
84   fNumSideN(0), 
85   fNumLongP(0), 
86   fNumLongN(0), 
87   fDenOutP(0),  
88   fDenOutN(0),  
89   fDenSideP(0), 
90   fDenSideN(0), 
91   fDenLongP(0), 
92   fDenLongN(0)
93 {
94   // copy constructor
95   if (aCorrFctn.fNumOutP)
96     fNumOutP = new TH1D(*aCorrFctn.fNumOutP);
97   if (aCorrFctn.fNumOutN)
98     fNumOutN = new TH1D(*aCorrFctn.fNumOutN);
99   if (aCorrFctn.fNumSideP)
100     fNumSideP = new TH1D(*aCorrFctn.fNumSideP);
101   if (aCorrFctn.fNumSideN)
102     fNumSideN = new TH1D(*aCorrFctn.fNumSideN);
103   if (aCorrFctn.fNumLongP)
104     fNumLongP = new TH1D(*aCorrFctn.fNumLongP);
105   if (aCorrFctn.fNumLongN)
106     fNumLongN = new TH1D(*aCorrFctn.fNumLongN);
107
108   if (aCorrFctn.fDenOutP)
109     fDenOutP = new TH1D(*aCorrFctn.fDenOutP);
110   if (aCorrFctn.fDenOutN)
111     fDenOutN = new TH1D(*aCorrFctn.fDenOutN);
112   if (aCorrFctn.fDenSideP)
113     fDenSideP = new TH1D(*aCorrFctn.fDenSideP);
114   if (aCorrFctn.fDenSideN)
115     fDenSideN = new TH1D(*aCorrFctn.fDenSideN);
116   if (aCorrFctn.fDenLongP)
117     fDenLongP = new TH1D(*aCorrFctn.fDenLongP);
118   if (aCorrFctn.fDenLongN)
119     fDenLongN = new TH1D(*aCorrFctn.fDenLongN);
120 }
121 //____________________________
122 AliFemtoCorrFctnNonIdDR::~AliFemtoCorrFctnNonIdDR(){
123   delete fNumOutP; 
124   delete fNumOutN;  
125   delete fNumSideP; 
126   delete fNumSideN; 
127   delete fNumLongP; 
128   delete fNumLongN; 
129   delete fDenOutP;  
130   delete fDenOutN;  
131   delete fDenSideP; 
132   delete fDenSideN; 
133   delete fDenLongP; 
134   delete fDenLongN;
135 }
136 //_________________________
137 AliFemtoCorrFctnNonIdDR& AliFemtoCorrFctnNonIdDR::operator=(const AliFemtoCorrFctnNonIdDR& aCorrFctn)
138 {
139   // assignment operator
140   if (this == &aCorrFctn)
141     return *this;
142
143   if (aCorrFctn.fNumOutP)
144     fNumOutP = new TH1D(*aCorrFctn.fNumOutP);
145   if (aCorrFctn.fNumOutN)
146     fNumOutN = new TH1D(*aCorrFctn.fNumOutN);
147   if (aCorrFctn.fNumSideP)
148     fNumSideP = new TH1D(*aCorrFctn.fNumSideP);
149   if (aCorrFctn.fNumSideN)
150     fNumSideN = new TH1D(*aCorrFctn.fNumSideN);
151   if (aCorrFctn.fNumLongP)
152     fNumLongP = new TH1D(*aCorrFctn.fNumLongP);
153   if (aCorrFctn.fNumLongN)
154     fNumLongN = new TH1D(*aCorrFctn.fNumLongN);
155
156   if (aCorrFctn.fDenOutP)
157     fDenOutP = new TH1D(*aCorrFctn.fDenOutP);
158   if (aCorrFctn.fDenOutN)
159     fDenOutN = new TH1D(*aCorrFctn.fDenOutN);
160   if (aCorrFctn.fDenSideP)
161     fDenSideP = new TH1D(*aCorrFctn.fDenSideP);
162   if (aCorrFctn.fDenSideN)
163     fDenSideN = new TH1D(*aCorrFctn.fDenSideN);
164   if (aCorrFctn.fDenLongP)
165     fDenLongP = new TH1D(*aCorrFctn.fDenLongP);
166   if (aCorrFctn.fDenLongN)
167     fDenLongN = new TH1D(*aCorrFctn.fDenLongN);
168
169   return *this;
170 }
171
172 //_________________________
173 void AliFemtoCorrFctnNonIdDR::Finish(){
174   // here is where we should normalize, fit, etc...
175   // we should NOT Draw() the histos (as I had done it below),
176   // since we want to insulate ourselves from root at this level
177   // of the code.  Do it instead at root command line with browser.
178   //  fNumerator->Draw();
179   //fDenominator->Draw();
180   //fRatio->Draw();
181   //  fRatio->Divide(fNumerator,fDenominator,1.0,1.0);
182
183 }
184
185 //____________________________
186 AliFemtoString AliFemtoCorrFctnNonIdDR::Report(){
187   // construct report
188   string stemp = "Non-identical particles Correlation Function Report:\n";
189   char ctemp[100];
190   snprintf(ctemp , 100, "Number of entries in numerators:\t%E\n",fNumOutP->GetEntries()+fNumOutN->GetEntries());
191   stemp += ctemp;
192   snprintf(ctemp , 100, "Number of entries in denominators:\t%E\n",fDenOutP->GetEntries()+fDenOutN->GetEntries());
193   stemp += ctemp;
194   //  stemp += mCoulombWeight->Report();
195   AliFemtoString returnThis = stemp;
196   return returnThis;
197 }
198 //____________________________
199 void AliFemtoCorrFctnNonIdDR::AddRealPair(AliFemtoPair* pair){
200   // add true pair
201   double tKStar = pair->KStar();
202   if (pair->KOut()>0.0)
203     fNumOutP->Fill(tKStar);
204   else
205     fNumOutN->Fill(tKStar);
206
207   if (pair->KSide()>0.0)
208     fNumSideP->Fill(tKStar);
209   else
210     fNumSideN->Fill(tKStar);
211
212   if (pair->KLong()>0.0)
213     fNumLongP->Fill(tKStar);
214   else
215     fNumLongN->Fill(tKStar);
216
217 }
218 //____________________________
219 void AliFemtoCorrFctnNonIdDR::AddMixedPair(AliFemtoPair* pair){
220   // add mixed (background) pair
221   double tKStar = pair->KStar();
222   if (pair->KOut()>0.0)
223     fDenOutP->Fill(tKStar);
224   else
225     fDenOutN->Fill(tKStar);
226
227   if (pair->KSide()>0.0)
228     fDenSideP->Fill(tKStar);
229   else
230     fDenSideN->Fill(tKStar);
231
232   if (pair->KLong()>0.0)
233     fDenLongP->Fill(tKStar);
234   else
235     fDenLongN->Fill(tKStar);
236 }
237 //____________________________
238 void AliFemtoCorrFctnNonIdDR::Write(){
239   fNumOutP->Write(); 
240   fNumOutN->Write();  
241   fNumSideP->Write(); 
242   fNumSideN->Write(); 
243   fNumLongP->Write(); 
244   fNumLongN->Write(); 
245   fDenOutP->Write();  
246   fDenOutN->Write();  
247   fDenSideP->Write(); 
248   fDenSideN->Write(); 
249   fDenLongP->Write(); 
250   fDenLongN->Write();
251 }
252
253 TList* AliFemtoCorrFctnNonIdDR::GetOutputList()
254 {
255   // Prepare the list of objects to be written to the output
256   TList *tOutputList = new TList();
257
258   tOutputList->Add(fNumOutP); 
259   tOutputList->Add(fNumOutN);  
260   tOutputList->Add(fNumSideP); 
261   tOutputList->Add(fNumSideN); 
262   tOutputList->Add(fNumLongP); 
263   tOutputList->Add(fNumLongN); 
264   tOutputList->Add(fDenOutP);  
265   tOutputList->Add(fDenOutN);  
266   tOutputList->Add(fDenSideP); 
267   tOutputList->Add(fDenSideN); 
268   tOutputList->Add(fDenLongP); 
269   tOutputList->Add(fDenLongN);
270
271   return tOutputList;
272 }