]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/FEMTOSCOPY/AliFemtoUser/AliFemtoCorrFctnPairFractions.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGCF / FEMTOSCOPY / AliFemtoUser / AliFemtoCorrFctnPairFractions.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 //                                                                            //
3 // AliFemtoCorrFctnPairFractions - A correlation function that analyzes            //
4 // two particle correlations with respect to the azimuthal angle (phi)        //
5 // and pseudorapidity (eta) difference                                        //
6 //                                                                            //
7 // Authors: Malgorzata Janik, majanik@cern.ch                                   //
8 //                                                                            //
9 ////////////////////////////////////////////////////////////////////////////////
10
11 #include "AliFemtoCorrFctnPairFractions.h"
12 #include "AliFemtoModelHiddenInfo.h"
13 //#include "AliFemtoHisto.hh"
14 #include <cstdio>
15 #include <TMath.h>
16
17 #ifdef __ROOT__ 
18 ClassImp(AliFemtoCorrFctnPairFractions)
19 #endif
20   
21 #define PIH 1.57079632679489656
22 #define PIT 6.28318530717958623
23
24 //____________________________
25 AliFemtoCorrFctnPairFractions::AliFemtoCorrFctnPairFractions(char* title):
26 AliFemtoCorrFctn(),
27   fPairFractions(0),
28   fPairFractionsDen(0),
29   fphiL(0),
30   fphiT(0)
31 {
32
33   //fphiL = (-(int)(aPhiBins/4)+0.5)*2.*TMath::Pi()/aPhiBins;
34   //fphiT = 2*TMath::Pi()+(-(int)(aPhiBins/4)+0.5)*2.*TMath::Pi()/aPhiBins;
35
36   TString  hname  = "hPairFraction"; hname+= title;
37   TString  htitle = "Pair Fraction "; htitle+= title;
38   fPairFractions = new TH1F(hname.Data(),htitle.Data(), 9, 0, 9);
39   fPairFractions->GetXaxis()->SetBinLabel(1,"#pi#pi, MC");
40   fPairFractions->GetXaxis()->SetBinLabel(2,"KK, MC");
41   fPairFractions->GetXaxis()->SetBinLabel(3,"pp, MC");
42   fPairFractions->GetXaxis()->SetBinLabel(4,"#pi K, MC");
43   fPairFractions->GetXaxis()->SetBinLabel(5,"#pi p, MC");
44   fPairFractions->GetXaxis()->SetBinLabel(6,"Kp, MC");
45   fPairFractions->GetXaxis()->SetBinLabel(7,"e+, MC");
46   fPairFractions->GetXaxis()->SetBinLabel(8,"#mu+, MC");
47   fPairFractions->GetXaxis()->SetBinLabel(9,"Other, MC");
48
49   hname  = "hPairFractionDen"; hname+= title;
50   htitle = "Pair Fraction in Mixing "; htitle+= title;
51   fPairFractionsDen = new TH1F(hname.Data(),htitle.Data(), 9, 0, 9);
52   fPairFractionsDen->GetXaxis()->SetBinLabel(1,"#pi#pi, MC");
53   fPairFractionsDen->GetXaxis()->SetBinLabel(2,"KK, MC");
54   fPairFractionsDen->GetXaxis()->SetBinLabel(3,"pp, MC");
55   fPairFractionsDen->GetXaxis()->SetBinLabel(4,"#pi K, MC");
56   fPairFractionsDen->GetXaxis()->SetBinLabel(5,"#pi p, MC");
57   fPairFractionsDen->GetXaxis()->SetBinLabel(6,"Kp, MC");
58   fPairFractionsDen->GetXaxis()->SetBinLabel(7,"e+, MC");
59   fPairFractionsDen->GetXaxis()->SetBinLabel(8,"#mu+, MC");
60   fPairFractionsDen->GetXaxis()->SetBinLabel(9,"Other, MC");
61
62   // to enable error bar calculation...
63   fPairFractions->Sumw2();
64   fPairFractionsDen->Sumw2();
65 }
66
67 //____________________________
68 AliFemtoCorrFctnPairFractions::AliFemtoCorrFctnPairFractions(const AliFemtoCorrFctnPairFractions& aCorrFctn) :
69   AliFemtoCorrFctn(),
70   fPairFractions(0),
71   fPairFractionsDen(0),
72   fphiL(0),
73   fphiT(0)
74 {
75   // copy constructor
76   if (aCorrFctn.fPairFractions)
77     fPairFractions = new TH1F(*aCorrFctn.fPairFractions);
78   else
79     fPairFractions = 0;
80
81  if (aCorrFctn.fPairFractions)
82     fPairFractions = new TH1F(*aCorrFctn.fPairFractions);
83   else
84     fPairFractions = 0;
85
86   fphiL = aCorrFctn.fphiL;
87   fphiT = aCorrFctn.fphiT;
88
89
90 }
91 //____________________________
92 AliFemtoCorrFctnPairFractions::~AliFemtoCorrFctnPairFractions(){
93   // destructor
94   if(fPairFractions)
95     delete fPairFractions;
96   if(fPairFractionsDen)
97     delete fPairFractionsDen;  
98 }
99 //_________________________
100 AliFemtoCorrFctnPairFractions& AliFemtoCorrFctnPairFractions::operator=(const AliFemtoCorrFctnPairFractions& aCorrFctn)
101 {
102   // assignment operator
103   if (this == &aCorrFctn)
104     return *this;
105
106   if (aCorrFctn.fPairFractions)
107     fPairFractions = new TH1F(*aCorrFctn.fPairFractions);
108   else
109     fPairFractions = 0;
110
111   if (aCorrFctn.fPairFractionsDen)
112     fPairFractionsDen = new TH1F(*aCorrFctn.fPairFractionsDen);
113   else
114     fPairFractionsDen = 0;
115   
116   fphiL = aCorrFctn.fphiL;
117   fphiT = aCorrFctn.fphiT;
118
119   return *this;
120 }
121 //_________________________
122 void AliFemtoCorrFctnPairFractions::Finish(){
123   // here is where we should normalize, fit, etc...
124   // we should NOT Draw() the histos (as I had done it below),
125   // since we want to insulate ourselves from root at this level
126   // of the code.  Do it instead at root command line with browser.
127   //  mShareNumerator->Draw();
128   // mShareDenominator->Draw();
129   // mRatio->Draw();
130
131 }
132
133 //____________________________
134 AliFemtoString AliFemtoCorrFctnPairFractions::Report(){
135   // create report
136   string stemp = "Pair Fractions Correlation Function Report:\n";
137   char ctemp[100];
138   snprintf(ctemp , 100, "Number of entries in numerator:\t%E\n",fPairFractions->GetEntries());
139   stemp += ctemp;
140   snprintf(ctemp , 100, "Number of entries in denominator:\t%E\n",fPairFractions->GetEntries());
141   stemp += ctemp;
142   //  stemp += mCoulombWeight->Report();
143   AliFemtoString returnThis = stemp;
144   return returnThis;
145 }
146 //____________________________
147 void AliFemtoCorrFctnPairFractions::AddRealPair( AliFemtoPair* pair){
148   // add real (effect) pair
149
150   //Applying pair cuts
151   if (fPairCut)
152     if (!fPairCut->Pass(pair)) return;
153
154
155
156   Int_t pdg1=0;
157   AliFemtoModelHiddenInfo *info1 = ( AliFemtoModelHiddenInfo *) pair->Track1()->GetHiddenInfo();
158   if(info1)pdg1 = info1->GetPDGPid();
159   Int_t pdg2=0;
160   AliFemtoModelHiddenInfo *info2 = ( AliFemtoModelHiddenInfo *) pair->Track2()->GetHiddenInfo();
161   if(info2)pdg2 = info2->GetPDGPid();
162
163   if(abs(pdg1)==211 && abs(pdg2)==211) //pi pi
164       fPairFractions->Fill(0.5);
165   else if(abs(pdg1)==321 && abs(pdg2)==321)// K K
166       fPairFractions->Fill(1.5);
167   else if(abs(pdg1)==2212 && abs(pdg2)==2212)// p p
168       fPairFractions->Fill(2.5);
169   else if((abs(pdg1)==211 && abs(pdg2)==321)||(abs(pdg1)==321 && abs(pdg2)==211))// pi K
170       fPairFractions->Fill(3.5);
171   else if((abs(pdg1)==211 && abs(pdg2)==2212)||(abs(pdg1)==2212 && abs(pdg2)==211))// pi p
172       fPairFractions->Fill(4.5);
173   else if((abs(pdg1)==321 && abs(pdg2)==2212)||(abs(pdg1)==2212 && abs(pdg2)==321))//K p
174       fPairFractions->Fill(5.5);
175   else if(abs(pdg1)==13 || abs(pdg2)==13)//one particle from the pair is electron
176       fPairFractions->Fill(6.5);
177  else if(abs(pdg1)==11 || abs(pdg2)==11)//one particle from the pair is muon
178       fPairFractions->Fill(7.5);
179   else //other
180     {
181       fPairFractions->Fill(8.5);
182     }
183   /*double phi1 = pair->Track1()->Track()->P().Phi();
184   double phi2 = pair->Track2()->Track()->P().Phi();
185   double eta1 = pair->Track1()->Track()->P().PseudoRapidity();
186   double eta2 = pair->Track2()->Track()->P().PseudoRapidity();
187
188   double phi1 = pair->Track1()->FourMomentum().Phi();
189   double phi2 = pair->Track2()->FourMomentum().Phi();
190   double eta1 = pair->Track1()->FourMomentum().PseudoRapidity();
191   double eta2 = pair->Track2()->FourMomentum().PseudoRapidity();
192
193   double dphi = phi1 - phi2;
194   while (dphi<fphiL) dphi+=PIT;
195   while (dphi>fphiT) dphi-=PIT;
196
197   double deta = eta1 - eta2;
198
199    double px1 = pair->Track1()->Track()->P().x();
200    double py1 = pair->Track1()->Track()->P().y();
201    double px2 = pair->Track2()->Track()->P().x();
202    double py2 = pair->Track2()->Track()->P().y();
203    double pt1 = TMath::Hypot(px1, py1);
204    double pt2 = TMath::Hypot(px2, py2);
205
206
207    double PionMass = 0.13956995;*/
208  
209 }
210 //____________________________
211 void AliFemtoCorrFctnPairFractions::AddMixedPair( AliFemtoPair* pair){
212   // add mixed (background) pair
213   if (fPairCut)
214     if (!fPairCut->Pass(pair)) return;
215
216
217
218   Int_t pdg1=0;
219   AliFemtoModelHiddenInfo *info1 = ( AliFemtoModelHiddenInfo *) pair->Track1()->GetHiddenInfo();
220   if(info1)pdg1 = info1->GetPDGPid();
221   Int_t pdg2=0;
222   AliFemtoModelHiddenInfo *info2 = ( AliFemtoModelHiddenInfo *) pair->Track2()->GetHiddenInfo();
223   if(info2)pdg2 = info2->GetPDGPid();
224
225   if(abs(pdg1)==211 && abs(pdg2)==211) //pi pi
226     fPairFractionsDen->Fill(0.5);
227   else if(abs(pdg1)==321 && abs(pdg2)==321)// K K
228     fPairFractionsDen->Fill(1.5);
229   else if(abs(pdg1)==2212 && abs(pdg2)==2212)// p p
230     fPairFractionsDen->Fill(2.5);
231   else if((abs(pdg1)==211 && abs(pdg2)==321)||(abs(pdg1)==321 && abs(pdg2)==211))// pi K
232     fPairFractionsDen->Fill(3.5);
233   else if((abs(pdg1)==211 && abs(pdg2)==2212)||(abs(pdg1)==2212 && abs(pdg2)==211))// pi p
234     fPairFractionsDen->Fill(4.5);
235   else if((abs(pdg1)==321 && abs(pdg2)==2212)||(abs(pdg1)==2212 && abs(pdg2)==321))//K p
236     fPairFractionsDen->Fill(5.5);
237   else if(abs(pdg1)==13 || abs(pdg2)==13)//one particle from the pair is electron
238     fPairFractionsDen->Fill(6.5);
239   else if(abs(pdg1)==11 || abs(pdg2)==11)//one particle from the pair is muon
240     fPairFractionsDen->Fill(7.5);
241   else //other
242     {
243       fPairFractionsDen->Fill(8.5);
244     }
245
246 }
247
248
249 void AliFemtoCorrFctnPairFractions::WriteHistos()
250 {
251   // Write out result histograms
252   fPairFractions->Write();
253   fPairFractionsDen->Write();
254 }
255
256 TList* AliFemtoCorrFctnPairFractions::GetOutputList()
257 {
258   // Prepare the list of objects to be written to the output
259   TList *tOutputList = new TList();
260
261   tOutputList->Add(fPairFractions);
262   tOutputList->Add(fPairFractionsDen);
263
264   return tOutputList;
265
266 }