]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemtoUser/AliFemtoTPCInnerCorrFctn.cxx
Bring AliFemto up to date with latest code developements
[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   fDTPCNumerator(0),
42   fDTPCDenominator(0)
43 {
44   // copy constructor
45   if (aCorrFctn.fDTPCNumerator)
46     fDTPCNumerator = new TH2D(*aCorrFctn.fDTPCNumerator);
47   if (aCorrFctn.fDTPCDenominator)
48     fDTPCDenominator = new TH2D(*aCorrFctn.fDTPCDenominator);
49 }
50 //____________________________
51 AliFemtoTPCInnerCorrFctn::~AliFemtoTPCInnerCorrFctn(){
52   // destructor
53   delete fDTPCNumerator;
54   delete fDTPCDenominator;
55 }
56 //_________________________
57 AliFemtoTPCInnerCorrFctn& AliFemtoTPCInnerCorrFctn::operator=(const AliFemtoTPCInnerCorrFctn& aCorrFctn)
58 {
59   // assignment operator
60   if (this == &aCorrFctn)
61     return *this;
62
63   if (aCorrFctn.fDTPCNumerator)
64     fDTPCNumerator = new TH2D(*aCorrFctn.fDTPCNumerator);
65   else
66     fDTPCNumerator = 0;
67   if (aCorrFctn.fDTPCDenominator)
68     fDTPCDenominator = new TH2D(*aCorrFctn.fDTPCDenominator);
69   else
70     fDTPCDenominator = 0;
71
72   return *this;
73 }
74 //_________________________
75 void AliFemtoTPCInnerCorrFctn::Finish(){
76   // here is where we should normalize, fit, etc...
77   // we should NOT Draw() the histos (as I had done it below),
78   // since we want to insulate ourselves from root at this level
79   // of the code.  Do it instead at root command line with browser.
80   //  mShareNumerator->Draw();
81   //mShareDenominator->Draw();
82   //mRatio->Draw();
83
84 }
85
86 //____________________________
87 AliFemtoString AliFemtoTPCInnerCorrFctn::Report(){
88   // create report
89   string stemp = "Entrace TPC distance Correlation Function Report:\n";
90   char ctemp[100];
91   sprintf(ctemp,"Number of entries in numerator:\t%E\n",fDTPCNumerator->GetEntries());
92   stemp += ctemp;
93   sprintf(ctemp,"Number of entries in denominator:\t%E\n",fDTPCDenominator->GetEntries());
94   stemp += ctemp;
95   //  stemp += mCoulombWeight->Report();
96   AliFemtoString returnThis = stemp;
97   return returnThis;
98 }
99 //____________________________
100 void AliFemtoTPCInnerCorrFctn::AddRealPair( AliFemtoPair* pair){
101   // add real (effect) pair
102   double tQinv = fabs(pair->QInv());   // note - qInv() will be negative for identical pairs...
103   double distx = pair->Track1()->Track()->NominalTpcEntrancePoint().x() - pair->Track2()->Track()->NominalTpcEntrancePoint().x();
104   double disty = pair->Track1()->Track()->NominalTpcEntrancePoint().y() - pair->Track2()->Track()->NominalTpcEntrancePoint().y();
105   double distz = pair->Track1()->Track()->NominalTpcEntrancePoint().z() - pair->Track2()->Track()->NominalTpcEntrancePoint().z();
106   double dist = sqrt(distx*distx + disty*disty + distz*distz);
107
108   fDTPCNumerator->Fill(tQinv, dist);
109 //   cout << "AliFemtoTPCInnerCorrFctn::AddRealPair : " << tQinv << " " << dist << endl;
110 //   cout << distx << " " << disty << " " << distz << endl;
111 }
112 //____________________________
113 void AliFemtoTPCInnerCorrFctn::AddMixedPair( AliFemtoPair* pair){
114   // add mixed (background) pair
115   double tQinv = fabs(pair->QInv());   // note - qInv() will be negative for identical pairs...
116   double distx = pair->Track1()->Track()->NominalTpcEntrancePoint().x() - pair->Track2()->Track()->NominalTpcEntrancePoint().x();
117   double disty = pair->Track1()->Track()->NominalTpcEntrancePoint().y() - pair->Track2()->Track()->NominalTpcEntrancePoint().y();
118   double distz = pair->Track1()->Track()->NominalTpcEntrancePoint().z() - pair->Track2()->Track()->NominalTpcEntrancePoint().z();
119   double dist = sqrt(distx*distx + disty*disty + distz*distz);
120
121   fDTPCDenominator->Fill(tQinv,dist);
122 }
123
124
125 void AliFemtoTPCInnerCorrFctn::WriteHistos()
126 {
127   // Write out result histograms
128   fDTPCNumerator->Write();
129   fDTPCDenominator->Write();
130 }
131 //______________________________
132 TList* AliFemtoTPCInnerCorrFctn::GetOutputList()
133 {
134   // Prepare the list of objects to be written to the output
135   TList *tOutputList = new TList();
136
137   tOutputList->Add(fDTPCNumerator); 
138   tOutputList->Add(fDTPCDenominator);  
139
140   return tOutputList;
141 }