1 ////////////////////////////////////////////////////////////////////////////////
3 // AliFemtoGammaMonitor - A correlation function that analyzes //
4 // two particle correlations with respect to the azimuthal angle (phi) //
5 // and pseudorapidity (eta) difference //
7 // Authors: Adam Kisiel Adam.Kisiel@cern.ch //
9 ////////////////////////////////////////////////////////////////////////////////
11 #include "AliFemtoCorrFctnGammaMonitor.h"
12 //#include "AliFemtoHisto.hh"
17 ClassImp(AliFemtoCorrFctnGammaMonitor)
20 //____________________________
21 AliFemtoCorrFctnGammaMonitor::AliFemtoCorrFctnGammaMonitor(char* title, const int& aMinvBins=20, const int& aDThetaBins=20):
29 char tTitNumD[101] = "NumPMinvTheta";
30 strncat(tTitNumD,title, 100);
31 fNumPMinvDTheta = new TH2D(tTitNumD,title,aMinvBins,0.0,0.2,aDThetaBins,0.0,0.2);
33 char tTitDenD[101] = "DenPMinvTheta";
34 strncat(tTitDenD,title, 100);
35 fDenPMinvDTheta = new TH2D(tTitDenD,title,aMinvBins,0.0,0.2,aDThetaBins,0.0,0.2);
38 char tTitNumR[101] = "NumNMinvTheta";
39 strncat(tTitNumR,title, 100);
40 fNumNMinvDTheta = new TH2D(tTitNumR,title,aMinvBins,0.0,0.2,aDThetaBins,0.0,0.2);
42 char tTitDenR[101] = "DenNMinvTheta";
43 strncat(tTitDenR,title, 100);
44 fDenNMinvDTheta = new TH2D(tTitDenR,title,aMinvBins,0.0,0.2,aDThetaBins,0.0,0.2);
46 // to enable error bar calculation...
47 fNumPMinvDTheta->Sumw2();
48 fDenPMinvDTheta->Sumw2();
49 fNumNMinvDTheta->Sumw2();
50 fDenNMinvDTheta->Sumw2();
53 //____________________________
54 AliFemtoCorrFctnGammaMonitor::AliFemtoCorrFctnGammaMonitor(const AliFemtoCorrFctnGammaMonitor& aCorrFctn) :
62 if (aCorrFctn.fNumPMinvDTheta)
63 fNumPMinvDTheta = new TH2D(*aCorrFctn.fNumPMinvDTheta);
66 if (aCorrFctn.fDenPMinvDTheta)
67 fDenPMinvDTheta = new TH2D(*aCorrFctn.fDenPMinvDTheta);
71 if (aCorrFctn.fNumNMinvDTheta)
72 fNumNMinvDTheta = new TH2D(*aCorrFctn.fNumNMinvDTheta);
75 if (aCorrFctn.fDenNMinvDTheta)
76 fDenNMinvDTheta = new TH2D(*aCorrFctn.fDenNMinvDTheta);
81 //____________________________
82 AliFemtoCorrFctnGammaMonitor::~AliFemtoCorrFctnGammaMonitor(){
84 delete fNumPMinvDTheta;
85 delete fDenPMinvDTheta;
86 delete fNumNMinvDTheta;
87 delete fDenNMinvDTheta;
89 //_________________________
90 AliFemtoCorrFctnGammaMonitor& AliFemtoCorrFctnGammaMonitor::operator=(const AliFemtoCorrFctnGammaMonitor& aCorrFctn)
92 // assignment operator
93 if (this == &aCorrFctn)
96 if (aCorrFctn.fNumPMinvDTheta)
97 fNumPMinvDTheta = new TH2D(*aCorrFctn.fNumPMinvDTheta);
100 if (aCorrFctn.fDenPMinvDTheta)
101 fDenPMinvDTheta = new TH2D(*aCorrFctn.fDenPMinvDTheta);
105 if (aCorrFctn.fNumNMinvDTheta)
106 fNumNMinvDTheta = new TH2D(*aCorrFctn.fNumNMinvDTheta);
109 if (aCorrFctn.fDenNMinvDTheta)
110 fDenNMinvDTheta = new TH2D(*aCorrFctn.fDenNMinvDTheta);
117 //_________________________
118 void AliFemtoCorrFctnGammaMonitor::Finish(){
119 // here is where we should normalize, fit, etc...
120 // we should NOT Draw() the histos (as I had done it below),
121 // since we want to insulate ourselves from root at this level
122 // of the code. Do it instead at root command line with browser.
123 // mShareNumerator->Draw();
124 //mShareDenominator->Draw();
129 //____________________________
130 AliFemtoString AliFemtoCorrFctnGammaMonitor::Report(){
132 string stemp = "Gamma Monitor Function Report:\n";
134 snprintf(ctemp , 100, "Number of entries in numerator:\t%E\n",fNumPMinvDTheta->GetEntries());
136 snprintf(ctemp , 100, "Number of entries in denominator:\t%E\n",fDenPMinvDTheta->GetEntries());
138 // stemp += mCoulombWeight->Report();
139 AliFemtoString returnThis = stemp;
142 //____________________________
143 void AliFemtoCorrFctnGammaMonitor::AddRealPair( AliFemtoPair* pair){
144 // add real (effect) pair
145 double me = 0.000511;
147 double theta1 = pair->Track1()->Track()->P().Theta();
148 double theta2 = pair->Track2()->Track()->P().Theta();
149 double dtheta = TMath::Abs(theta1 - theta2);
151 double e1 = TMath::Sqrt(me*me + pair->Track1()->Track()->P().Mag2());
152 double e2 = TMath::Sqrt(me*me + pair->Track2()->Track()->P().Mag2());
154 double minv = 2*me*me + 2*(e1*e2 -
155 pair->Track1()->Track()->P().x()*pair->Track2()->Track()->P().x() -
156 pair->Track1()->Track()->P().y()*pair->Track2()->Track()->P().y() -
157 pair->Track1()->Track()->P().z()*pair->Track2()->Track()->P().z());
159 if (pair->KSide()>0.0)
160 fNumPMinvDTheta->Fill(minv, dtheta);
162 fNumNMinvDTheta->Fill(minv, dtheta);
164 //____________________________
165 void AliFemtoCorrFctnGammaMonitor::AddMixedPair( AliFemtoPair* pair){
166 // add mixed (background) pair
167 double me = 0.000511;
169 double theta1 = pair->Track1()->Track()->P().Theta();
170 double theta2 = pair->Track2()->Track()->P().Theta();
171 double dtheta = TMath::Abs(theta1 - theta2);
173 double e1 = TMath::Sqrt(me*me + pair->Track1()->Track()->P().Mag2());
174 double e2 = TMath::Sqrt(me*me + pair->Track2()->Track()->P().Mag2());
176 double minv = 2*me*me + 2*(e1*e2 -
177 pair->Track1()->Track()->P().x()*pair->Track2()->Track()->P().x() -
178 pair->Track1()->Track()->P().y()*pair->Track2()->Track()->P().y() -
179 pair->Track1()->Track()->P().z()*pair->Track2()->Track()->P().z());
181 if (pair->KSide()>0.0)
182 fDenPMinvDTheta->Fill(minv, dtheta);
184 fDenNMinvDTheta->Fill(minv, dtheta);
188 void AliFemtoCorrFctnGammaMonitor::WriteHistos()
190 // Write out result histograms
191 fNumPMinvDTheta->Write();
192 fDenPMinvDTheta->Write();
193 fNumNMinvDTheta->Write();
194 fDenNMinvDTheta->Write();
197 TList* AliFemtoCorrFctnGammaMonitor::GetOutputList()
199 // Prepare the list of objects to be written to the output
200 TList *tOutputList = new TList();
202 tOutputList->Add(fNumPMinvDTheta);
203 tOutputList->Add(fDenPMinvDTheta);
204 tOutputList->Add(fNumNMinvDTheta);
205 tOutputList->Add(fDenNMinvDTheta);