]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FEMTOSCOPY/AliFemto/AliFemtoQinvCorrFctn.cxx
Lines getting the matched track moved to a method in AliCalorimeterUtils. Lines copie...
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / AliFemtoQinvCorrFctn.cxx
CommitLineData
d92ed900 1///////////////////////////////////////////////////////////////////////////
2// //
3// AliFemtoQinvCorrFctn: //
4// a simple Q-invariant correlation function //
5// //
6///////////////////////////////////////////////////////////////////////////
67427ff7 7
d0e92d9a 8#include "AliFemtoQinvCorrFctn.h"
9//#include "AliFemtoHisto.h"
67427ff7 10#include <cstdio>
11
12#ifdef __ROOT__
13ClassImp(AliFemtoQinvCorrFctn)
14#endif
15
16//____________________________
0215f606 17AliFemtoQinvCorrFctn::AliFemtoQinvCorrFctn(char* title, const int& nbins, const float& QinvLo, const float& QinvHi):
18 fNumerator(0),
19 fDenominator(0),
69c1c8ff 20 fRatio(0),
21 fkTMonitor(0)
0215f606 22{
67427ff7 23 // set up numerator
24 // title = "Num Qinv (MeV/c)";
9bd7eeb9 25 char tTitNum[101] = "Num";
3be563bf 26 strncat(tTitNum,title, 100);
d92ed900 27 fNumerator = new TH1D(tTitNum,title,nbins,QinvLo,QinvHi);
67427ff7 28 // set up denominator
29 //title = "Den Qinv (MeV/c)";
9bd7eeb9 30 char tTitDen[101] = "Den";
3be563bf 31 strncat(tTitDen,title, 100);
d92ed900 32 fDenominator = new TH1D(tTitDen,title,nbins,QinvLo,QinvHi);
67427ff7 33 // set up ratio
34 //title = "Ratio Qinv (MeV/c)";
9bd7eeb9 35 char tTitRat[101] = "Rat";
3be563bf 36 strncat(tTitRat,title, 100);
d92ed900 37 fRatio = new TH1D(tTitRat,title,nbins,QinvLo,QinvHi);
9bd7eeb9 38 char tTitkT[101] = "kTDep";
3be563bf 39 strncat(tTitkT,title, 100);
69c1c8ff 40 fkTMonitor = new TH1D(tTitkT,title,200,0.0,2.0);
67427ff7 41 // this next bit is unfortunately needed so that we can have many histos of same "title"
42 // it is neccessary if we typedef TH1D to TH1d (which we do)
43 //fNumerator->SetDirectory(0);
44 //fDenominator->SetDirectory(0);
45 //fRatio->SetDirectory(0);
46
47 // to enable error bar calculation...
48 fNumerator->Sumw2();
49 fDenominator->Sumw2();
50 fRatio->Sumw2();
69c1c8ff 51 fkTMonitor->Sumw2();
67427ff7 52}
53
0215f606 54//____________________________
55AliFemtoQinvCorrFctn::AliFemtoQinvCorrFctn(const AliFemtoQinvCorrFctn& aCorrFctn) :
fcda1d4e 56 AliFemtoCorrFctn(),
0215f606 57 fNumerator(0),
58 fDenominator(0),
69c1c8ff 59 fRatio(0),
60 fkTMonitor(0)
0215f606 61{
d92ed900 62 // copy constructor
0215f606 63 fNumerator = new TH1D(*aCorrFctn.fNumerator);
64 fDenominator = new TH1D(*aCorrFctn.fDenominator);
65 fRatio = new TH1D(*aCorrFctn.fRatio);
69c1c8ff 66 fkTMonitor = new TH1D(*aCorrFctn.fkTMonitor);
0215f606 67}
67427ff7 68//____________________________
69AliFemtoQinvCorrFctn::~AliFemtoQinvCorrFctn(){
d92ed900 70 // destructor
67427ff7 71 delete fNumerator;
72 delete fDenominator;
73 delete fRatio;
69c1c8ff 74 delete fkTMonitor;
67427ff7 75}
0215f606 76//_________________________
77AliFemtoQinvCorrFctn& AliFemtoQinvCorrFctn::operator=(const AliFemtoQinvCorrFctn& aCorrFctn)
78{
d92ed900 79 // assignment operator
0215f606 80 if (this == &aCorrFctn)
81 return *this;
82
83 if (fNumerator) delete fNumerator;
84 fNumerator = new TH1D(*aCorrFctn.fNumerator);
85 if (fDenominator) delete fDenominator;
86 fDenominator = new TH1D(*aCorrFctn.fDenominator);
87 if (fRatio) delete fRatio;
88 fRatio = new TH1D(*aCorrFctn.fRatio);
69c1c8ff 89 if (fkTMonitor) delete fkTMonitor;
90 fkTMonitor = new TH1D(*aCorrFctn.fkTMonitor);
0215f606 91
92 return *this;
93}
94
67427ff7 95//_________________________
96void AliFemtoQinvCorrFctn::Finish(){
97 // here is where we should normalize, fit, etc...
98 // we should NOT Draw() the histos (as I had done it below),
99 // since we want to insulate ourselves from root at this level
100 // of the code. Do it instead at root command line with browser.
101 // fNumerator->Draw();
102 //fDenominator->Draw();
103 //fRatio->Draw();
104 fRatio->Divide(fNumerator,fDenominator,1.0,1.0);
105
106}
107
108//____________________________
109AliFemtoString AliFemtoQinvCorrFctn::Report(){
d92ed900 110 // construct report
67427ff7 111 string stemp = "Qinv Correlation Function Report:\n";
112 char ctemp[100];
3be563bf 113 snprintf(ctemp , 100, "Number of entries in numerator:\t%E\n",fNumerator->GetEntries());
67427ff7 114 stemp += ctemp;
3be563bf 115 snprintf(ctemp , 100, "Number of entries in denominator:\t%E\n",fDenominator->GetEntries());
67427ff7 116 stemp += ctemp;
3be563bf 117 snprintf(ctemp , 100, "Number of entries in ratio:\t%E\n",fRatio->GetEntries());
67427ff7 118 stemp += ctemp;
119 // stemp += mCoulombWeight->Report();
120 AliFemtoString returnThis = stemp;
121 return returnThis;
122}
123//____________________________
d0e92d9a 124void AliFemtoQinvCorrFctn::AddRealPair(AliFemtoPair* pair){
d92ed900 125 // add true pair
ad795b76 126 if (fPairCut)
127 if (!fPairCut->Pass(pair)) return;
128
d92ed900 129 double tQinv = fabs(pair->QInv()); // note - qInv() will be negative for identical pairs...
130 fNumerator->Fill(tQinv);
69c1c8ff 131 fkTMonitor->Fill(pair->KT());
d92ed900 132 // cout << "AliFemtoQinvCorrFctn::AddRealPair : " << pair->qInv() << " " << tQinv <<
67427ff7 133 //" " << pair->track1().FourMomentum() << " " << pair->track2().FourMomentum() << endl;
134}
135//____________________________
d0e92d9a 136void AliFemtoQinvCorrFctn::AddMixedPair(AliFemtoPair* pair){
d92ed900 137 // add mixed (background) pair
ad795b76 138 if (fPairCut)
139 if (!fPairCut->Pass(pair)) return;
140
67427ff7 141 double weight = 1.0;
d92ed900 142 double tQinv = fabs(pair->QInv()); // note - qInv() will be negative for identical pairs...
143 fDenominator->Fill(tQinv,weight);
67427ff7 144}
0b3bd1ac 145//____________________________
146void AliFemtoQinvCorrFctn::Write(){
147 // Write out neccessary objects
148 fNumerator->Write();
149 fDenominator->Write();
69c1c8ff 150 fkTMonitor->Write();
0b3bd1ac 151}
152//______________________________
153TList* AliFemtoQinvCorrFctn::GetOutputList()
154{
155 // Prepare the list of objects to be written to the output
156 TList *tOutputList = new TList();
157
158 tOutputList->Add(fNumerator);
159 tOutputList->Add(fDenominator);
69c1c8ff 160 tOutputList->Add(fkTMonitor);
0b3bd1ac 161
162 return tOutputList;
163}
67427ff7 164
165