1 ////////////////////////////////////////////////////////////////////////////////
3 /// AliFemtoTPCInnerCorrFctn - A correlation function that saves the ///
4 /// distance at the entrance to the TPC between two tracks as a function ///
6 /// Authors: Adam Kisiel kisiel@mps.ohio-state.edu ///
8 ////////////////////////////////////////////////////////////////////////////////
10 #include "AliFemtoTPCInnerCorrFctn.h"
11 //#include "AliFemtoHisto.hh"
15 ClassImp(AliFemtoTPCInnerCorrFctn)
18 //____________________________
19 AliFemtoTPCInnerCorrFctn::AliFemtoTPCInnerCorrFctn(char* title, const int& nbins, const float& QinvLo, const float& QinvHi):
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);
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);
34 // to enable error bar calculation...
35 fDTPCNumerator->Sumw2();
36 fDTPCDenominator->Sumw2();
39 //____________________________
40 AliFemtoTPCInnerCorrFctn::AliFemtoTPCInnerCorrFctn(const AliFemtoTPCInnerCorrFctn& aCorrFctn) :
45 if (aCorrFctn.fDTPCNumerator)
46 fDTPCNumerator = new TH2D(*aCorrFctn.fDTPCNumerator);
47 if (aCorrFctn.fDTPCDenominator)
48 fDTPCDenominator = new TH2D(*aCorrFctn.fDTPCDenominator);
50 //____________________________
51 AliFemtoTPCInnerCorrFctn::~AliFemtoTPCInnerCorrFctn(){
53 delete fDTPCNumerator;
54 delete fDTPCDenominator;
56 //_________________________
57 AliFemtoTPCInnerCorrFctn& AliFemtoTPCInnerCorrFctn::operator=(const AliFemtoTPCInnerCorrFctn& aCorrFctn)
59 // assignment operator
60 if (this == &aCorrFctn)
63 if (aCorrFctn.fDTPCNumerator)
64 fDTPCNumerator = new TH2D(*aCorrFctn.fDTPCNumerator);
67 if (aCorrFctn.fDTPCDenominator)
68 fDTPCDenominator = new TH2D(*aCorrFctn.fDTPCDenominator);
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();
86 //____________________________
87 AliFemtoString AliFemtoTPCInnerCorrFctn::Report(){
89 string stemp = "Entrace TPC distance Correlation Function Report:\n";
91 sprintf(ctemp,"Number of entries in numerator:\t%E\n",fDTPCNumerator->GetEntries());
93 sprintf(ctemp,"Number of entries in denominator:\t%E\n",fDTPCDenominator->GetEntries());
95 // stemp += mCoulombWeight->Report();
96 AliFemtoString returnThis = stemp;
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);
108 fDTPCNumerator->Fill(tQinv, dist);
109 // cout << "AliFemtoTPCInnerCorrFctn::AddRealPair : " << tQinv << " " << dist << endl;
110 // cout << distx << " " << disty << " " << distz << endl;
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);
121 fDTPCDenominator->Fill(tQinv,dist);
125 void AliFemtoTPCInnerCorrFctn::WriteHistos()
127 // Write out result histograms
128 fDTPCNumerator->Write();
129 fDTPCDenominator->Write();
131 //______________________________
132 TList* AliFemtoTPCInnerCorrFctn::GetOutputList()
134 // Prepare the list of objects to be written to the output
135 TList *tOutputList = new TList();
137 tOutputList->Add(fDTPCNumerator);
138 tOutputList->Add(fDTPCDenominator);