]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/AliFemtoCorrFctn3DSpherical.cxx
Fixing pair cut
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / AliFemtoCorrFctn3DSpherical.cxx
1 ///////////////////////////////////////////////////////////////////////////
2 //                                                                       //
3 // AliFemtoCorrFctn3DSpherical: a class to calculate 3D correlation      //
4 // for pairs of identical particles, binned in spherical coordinates.    //
5 // In analysis the function should be first created in a macro, then     //
6 // added to the analysis, and at the end of the macro the procedure to   //
7 // write out histograms should be called.                                //
8 //                                                                       //
9 ///////////////////////////////////////////////////////////////////////////
10
11 #include "AliFemtoCorrFctn3DSpherical.h"
12 #include <TMath.h>
13 #include <cstdio>
14
15 #ifdef __ROOT__ 
16 ClassImp(AliFemtoCorrFctn3DSpherical)
17 #endif
18
19 //____________________________
20   AliFemtoCorrFctn3DSpherical::AliFemtoCorrFctn3DSpherical(char* title, const int& nqbins, const float& QLo, const float& QHi, const int& nphibins, const int& ncthetabins):
21   fNumerator(0),
22   fDenominator(0) //,
23                                                           //  fPairCut(0x0)
24 {
25   // set up numerator
26   char tTitNum[100] = "Num";
27   strcat(tTitNum,title);
28   fNumerator = new TH3D(tTitNum,title,nqbins,QLo,QHi,nphibins,-TMath::Pi(),TMath::Pi(),ncthetabins,-1.0,1.0);
29   // set up denominator
30   char tTitDen[100] = "Den";
31   strcat(tTitDen,title);
32   fDenominator = new TH3D(tTitDen,title,nqbins,QLo,QHi,nphibins,-TMath::Pi(),TMath::Pi(),ncthetabins,-1.0,1.0);
33
34   // to enable error bar calculation...
35   fNumerator->Sumw2();
36   fDenominator->Sumw2();
37 }
38
39 AliFemtoCorrFctn3DSpherical::AliFemtoCorrFctn3DSpherical(const AliFemtoCorrFctn3DSpherical& aCorrFctn) :
40   AliFemtoCorrFctn(aCorrFctn),
41   fNumerator(0),
42   fDenominator(0) //,
43                                                         //  fPairCut(0x0)
44 {
45   // Copy constructor
46   fNumerator = new TH3D(*aCorrFctn.fNumerator);
47   fDenominator = new TH3D(*aCorrFctn.fDenominator);
48   //  fPairCut = aCorrFctn.fPairCut;
49 }
50 //____________________________
51 AliFemtoCorrFctn3DSpherical::~AliFemtoCorrFctn3DSpherical(){
52   // Destructor
53   delete fNumerator;
54   delete fDenominator;
55 }
56 //_________________________
57 AliFemtoCorrFctn3DSpherical& AliFemtoCorrFctn3DSpherical::operator=(const AliFemtoCorrFctn3DSpherical& aCorrFctn)
58 {
59   // assignment operator
60   if (this == &aCorrFctn)
61     return *this;
62
63   if (fNumerator) delete fNumerator;
64   fNumerator = new TH3D(*aCorrFctn.fNumerator);
65   if (fDenominator) delete fDenominator;
66   fDenominator = new TH3D(*aCorrFctn.fDenominator);
67   
68   //  fPairCut = aCorrFctn.fPairCut;
69   
70   return *this;
71 }
72
73 //_________________________
74 void AliFemtoCorrFctn3DSpherical::WriteOutHistos(){
75   // Write out all histograms to file
76   fNumerator->Write();
77   fDenominator->Write();
78 }
79 //______________________________
80 TList* AliFemtoCorrFctn3DSpherical::GetOutputList()
81 {
82   // Prepare the list of objects to be written to the output
83   TList *tOutputList = new TList();
84
85   tOutputList->Add(fNumerator); 
86   tOutputList->Add(fDenominator);  
87
88   return tOutputList;
89 }
90
91 //_________________________
92 void AliFemtoCorrFctn3DSpherical::Finish(){
93   // here is where we should normalize, fit, etc...
94 }
95
96 //____________________________
97 AliFemtoString AliFemtoCorrFctn3DSpherical::Report(){
98   // Construct the report
99   string stemp = "PRF Frame Spherical 3D Correlation Function Report:\n";
100   char ctemp[100];
101   sprintf(ctemp,"Number of entries in numerator:\t%E\n",fNumerator->GetEntries());
102   stemp += ctemp;
103   sprintf(ctemp,"Number of entries in denominator:\t%E\n",fDenominator->GetEntries());
104   stemp += ctemp;
105
106   if (fPairCut){
107     sprintf(ctemp,"Here is the PairCut specific to this CorrFctn\n");
108     stemp += ctemp;
109     stemp += fPairCut->Report();
110   }
111   else{
112     sprintf(ctemp,"No PairCut specific to this CorrFctn\n");
113     stemp += ctemp;
114   }
115
116   //  
117   AliFemtoString returnThis = stemp;
118   return returnThis;
119 }
120 //____________________________
121 void AliFemtoCorrFctn3DSpherical::AddRealPair( AliFemtoPair* pair){
122   // perform operations on real pairs
123   if (fPairCut){
124     if (!(fPairCut->Pass(pair))) return;
125   }
126
127   double tKO = pair->KOut();
128   double tKS = pair->KSide();
129   double tKL = pair->KLong();
130
131   double tKR = sqrt(tKO*tKO + tKS*tKS + tKL*tKL);
132   double tKC;
133   if ( fabs(tKR) < 1e-10 ) tKC = 0.0;
134   else tKC=tKL/tKR;
135   double tKP=atan2(tKS,tKO);
136
137   fNumerator->Fill(tKR,tKP,tKC);
138 }
139 //____________________________
140 void AliFemtoCorrFctn3DSpherical::AddMixedPair( AliFemtoPair* pair){
141   // perform operations on mixed pairs
142   if (fPairCut){
143     if (!(fPairCut->Pass(pair))) return;
144   }
145
146   double tKO = pair->KOut();
147   double tKS = pair->KSide();
148   double tKL = pair->KLong();
149
150   double tKR = sqrt(tKO*tKO + tKS*tKS + tKL*tKL);
151   double tKC;
152   if ( fabs(tKR) < 1e-10 ) tKC = 0.0;
153   else tKC=tKL/tKR;
154   double tKP=atan2(tKS,tKO);
155
156   fDenominator->Fill(tKR,tKP,tKC);
157 }
158
159