]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemtoUser/AliFemtoTPCInnerCorrFctn.cxx
Add DEta DPhi correlation code
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemtoUser / AliFemtoTPCInnerCorrFctn.cxx
1 ////////////////////////////////////////////////////////////////////////////////
2 ///                                                                          ///
3 /// AliFemtoTPCInnerCorrFctn - A correlation function that saves the         ///
4 /// distance at the entrance to the TPC between two tracks as a function     ///
5 /// of qinv                                                                  ///
6 /// Authors: Adam Kisiel kisiel@mps.ohio-state.edu                           ///
7 ///                                                                          ///
8 ////////////////////////////////////////////////////////////////////////////////
9
10 #include "AliFemtoTPCInnerCorrFctn.h"
11 //#include "AliFemtoHisto.hh"
12 #include <cstdio>
13
14 #ifdef __ROOT__ 
15 ClassImp(AliFemtoTPCInnerCorrFctn)
16 #endif
17
18 //____________________________
19 AliFemtoTPCInnerCorrFctn::AliFemtoTPCInnerCorrFctn(char* title, const int& nbins, const float& QinvLo, const float& QinvHi):
20   fDTPCNumerator(0),
21   fDTPCDenominator(0)
22 {
23   // set up numerator
24   //  title = "Num Qinv (MeV/c)";
25   char tTitNum[100] = "NumDTPC";
26   strcat(tTitNum,title);
27   fDTPCNumerator = new TH2D(tTitNum,title,nbins,QinvLo,QinvHi,100,0.0,20.0);
28   // set up denominator
29   //title = "Den Qinv (MeV/c)";
30   char tTitDen[100] = "DenDTPC";
31   strcat(tTitDen,title);
32   fDTPCDenominator = new TH2D(tTitDen,title,nbins,QinvLo,QinvHi,100,0.0,20.0);
33
34   // to enable error bar calculation...
35   fDTPCNumerator->Sumw2();
36   fDTPCDenominator->Sumw2();
37 }
38
39 //____________________________
40 AliFemtoTPCInnerCorrFctn::AliFemtoTPCInnerCorrFctn(const AliFemtoTPCInnerCorrFctn& aCorrFctn) :
41   AliFemtoCorrFctn(),
42   fDTPCNumerator(0),
43   fDTPCDenominator(0)
44 {
45   // copy constructor
46   if (aCorrFctn.fDTPCNumerator)
47     fDTPCNumerator = new TH2D(*aCorrFctn.fDTPCNumerator);
48   if (aCorrFctn.fDTPCDenominator)
49     fDTPCDenominator = new TH2D(*aCorrFctn.fDTPCDenominator);
50 }
51 //____________________________
52 AliFemtoTPCInnerCorrFctn::~AliFemtoTPCInnerCorrFctn(){
53   // destructor
54   delete fDTPCNumerator;
55   delete fDTPCDenominator;
56 }
57 //_________________________
58 AliFemtoTPCInnerCorrFctn& AliFemtoTPCInnerCorrFctn::operator=(const AliFemtoTPCInnerCorrFctn& aCorrFctn)
59 {
60   // assignment operator
61   if (this == &aCorrFctn)
62     return *this;
63
64   if (aCorrFctn.fDTPCNumerator)
65     fDTPCNumerator = new TH2D(*aCorrFctn.fDTPCNumerator);
66   else
67     fDTPCNumerator = 0;
68   if (aCorrFctn.fDTPCDenominator)
69     fDTPCDenominator = new TH2D(*aCorrFctn.fDTPCDenominator);
70   else
71     fDTPCDenominator = 0;
72
73   return *this;
74 }
75 //_________________________
76 void AliFemtoTPCInnerCorrFctn::Finish(){
77   // here is where we should normalize, fit, etc...
78   // we should NOT Draw() the histos (as I had done it below),
79   // since we want to insulate ourselves from root at this level
80   // of the code.  Do it instead at root command line with browser.
81   //  mShareNumerator->Draw();
82   //mShareDenominator->Draw();
83   //mRatio->Draw();
84
85 }
86
87 //____________________________
88 AliFemtoString AliFemtoTPCInnerCorrFctn::Report(){
89   // create report
90   string stemp = "Entrace TPC distance Correlation Function Report:\n";
91   char ctemp[100];
92   sprintf(ctemp,"Number of entries in numerator:\t%E\n",fDTPCNumerator->GetEntries());
93   stemp += ctemp;
94   sprintf(ctemp,"Number of entries in denominator:\t%E\n",fDTPCDenominator->GetEntries());
95   stemp += ctemp;
96   //  stemp += mCoulombWeight->Report();
97   AliFemtoString returnThis = stemp;
98   return returnThis;
99 }
100 //____________________________
101 void AliFemtoTPCInnerCorrFctn::AddRealPair( AliFemtoPair* pair){
102   // add real (effect) pair
103   double tQinv = fabs(pair->QInv());   // note - qInv() will be negative for identical pairs...
104   double distx = pair->Track1()->Track()->NominalTpcEntrancePoint().x() - pair->Track2()->Track()->NominalTpcEntrancePoint().x();
105   double disty = pair->Track1()->Track()->NominalTpcEntrancePoint().y() - pair->Track2()->Track()->NominalTpcEntrancePoint().y();
106   double distz = pair->Track1()->Track()->NominalTpcEntrancePoint().z() - pair->Track2()->Track()->NominalTpcEntrancePoint().z();
107   double dist = sqrt(distx*distx + disty*disty + distz*distz);
108
109   fDTPCNumerator->Fill(tQinv, dist);
110 }
111 //____________________________
112 void AliFemtoTPCInnerCorrFctn::AddMixedPair( AliFemtoPair* pair){
113   // add mixed (background) pair
114   double tQinv = fabs(pair->QInv());   // note - qInv() will be negative for identical pairs...
115   double distx = pair->Track1()->Track()->NominalTpcEntrancePoint().x() - pair->Track2()->Track()->NominalTpcEntrancePoint().x();
116   double disty = pair->Track1()->Track()->NominalTpcEntrancePoint().y() - pair->Track2()->Track()->NominalTpcEntrancePoint().y();
117   double distz = pair->Track1()->Track()->NominalTpcEntrancePoint().z() - pair->Track2()->Track()->NominalTpcEntrancePoint().z();
118   double dist = sqrt(distx*distx + disty*disty + distz*distz);
119
120   fDTPCDenominator->Fill(tQinv,dist);
121 }
122
123
124 void AliFemtoTPCInnerCorrFctn::WriteHistos()
125 {
126   // Write out result histograms
127   fDTPCNumerator->Write();
128   fDTPCDenominator->Write();
129 }
130 //______________________________
131 TList* AliFemtoTPCInnerCorrFctn::GetOutputList()
132 {
133   // Prepare the list of objects to be written to the output
134   TList *tOutputList = new TList();
135
136   tOutputList->Add(fDTPCNumerator); 
137   tOutputList->Add(fDTPCDenominator);  
138
139   return tOutputList;
140 }