]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FEMTOSCOPY/AliFemto/CorrFctn/AliFemtoQinvCorrFctn.cxx
6b5eaea035a589cc29a5373cd5df56b21c4b99a8
[u/mrichter/AliRoot.git] / PWG2 / FEMTOSCOPY / AliFemto / CorrFctn / AliFemtoQinvCorrFctn.cxx
1 /***************************************************************************
2  *
3  * $Id$
4  *
5  * Author: Mike Lisa, Ohio State, lisa@mps.ohio-state.edu
6  ***************************************************************************
7  *
8  * Description: part of STAR HBT Framework: AliFemtoMaker package
9  *   a simple Q-invariant correlation function           
10  *
11  ***************************************************************************
12  *
13  * $Log$
14  * Revision 1.1.1.1  2007/03/07 10:14:49  mchojnacki
15  * First version on CVS
16  *
17  * Revision 1.4  2000/01/25 17:34:45  laue
18  * I. In order to run the stand alone version of the AliFemtoMaker the following
19  * changes have been done:
20  * a) all ClassDefs and ClassImps have been put into #ifdef __ROOT__ statements
21  * b) unnecessary includes of StMaker.h have been removed
22  * c) the subdirectory AliFemtoMaker/doc/Make has been created including everything
23  * needed for the stand alone version
24  *
25  * II. To reduce the amount of compiler warning
26  * a) some variables have been type casted
27  * b) some destructors have been declared as virtual
28  *
29  * Revision 1.3  1999/07/29 02:47:09  lisa
30  * 1) add OpeningAngle correlation function 2) add AliFemtoMcEventReader 3) make histos in CorrFctns do errors correctly
31  *
32  * Revision 1.2  1999/07/06 22:33:20  lisa
33  * Adjusted all to work in pro and new - dev itself is broken
34  *
35  * Revision 1.1.1.1  1999/06/29 16:02:57  lisa
36  * Installation of AliFemtoMaker
37  *
38  **************************************************************************/
39
40 #include "CorrFctn/AliFemtoQinvCorrFctn.h"
41 //#include "Infrastructure/AliFemtoHisto.h"
42 #include <cstdio>
43
44 #ifdef __ROOT__ 
45 ClassImp(AliFemtoQinvCorrFctn)
46 #endif
47
48 //____________________________
49 AliFemtoQinvCorrFctn::AliFemtoQinvCorrFctn(char* title, const int& nbins, const float& QinvLo, const float& QinvHi){
50   // set up numerator
51   //  title = "Num Qinv (MeV/c)";
52   char TitNum[100] = "Num";
53   strcat(TitNum,title);
54   fNumerator = new TH1D(TitNum,title,nbins,QinvLo,QinvHi);
55   // set up denominator
56   //title = "Den Qinv (MeV/c)";
57   char TitDen[100] = "Den";
58   strcat(TitDen,title);
59   fDenominator = new TH1D(TitDen,title,nbins,QinvLo,QinvHi);
60   // set up ratio
61   //title = "Ratio Qinv (MeV/c)";
62   char TitRat[100] = "Rat";
63   strcat(TitRat,title);
64   fRatio = new TH1D(TitRat,title,nbins,QinvLo,QinvHi);
65   // this next bit is unfortunately needed so that we can have many histos of same "title"
66   // it is neccessary if we typedef TH1D to TH1d (which we do)
67   //fNumerator->SetDirectory(0);
68   //fDenominator->SetDirectory(0);
69   //fRatio->SetDirectory(0);
70
71   // to enable error bar calculation...
72   fNumerator->Sumw2();
73   fDenominator->Sumw2();
74   fRatio->Sumw2();
75
76 }
77
78 //____________________________
79 AliFemtoQinvCorrFctn::~AliFemtoQinvCorrFctn(){
80   delete fNumerator;
81   delete fDenominator;
82   delete fRatio;
83 }
84 //_________________________
85 void AliFemtoQinvCorrFctn::Finish(){
86   // here is where we should normalize, fit, etc...
87   // we should NOT Draw() the histos (as I had done it below),
88   // since we want to insulate ourselves from root at this level
89   // of the code.  Do it instead at root command line with browser.
90   //  fNumerator->Draw();
91   //fDenominator->Draw();
92   //fRatio->Draw();
93   fRatio->Divide(fNumerator,fDenominator,1.0,1.0);
94
95 }
96
97 //____________________________
98 AliFemtoString AliFemtoQinvCorrFctn::Report(){
99   string stemp = "Qinv Correlation Function Report:\n";
100   char ctemp[100];
101   sprintf(ctemp,"Number of entries in numerator:\t%E\n",fNumerator->GetEntries());
102   stemp += ctemp;
103   sprintf(ctemp,"Number of entries in denominator:\t%E\n",fDenominator->GetEntries());
104   stemp += ctemp;
105   sprintf(ctemp,"Number of entries in ratio:\t%E\n",fRatio->GetEntries());
106   stemp += ctemp;
107   //  stemp += mCoulombWeight->Report();
108   AliFemtoString returnThis = stemp;
109   return returnThis;
110 }
111 //____________________________
112 void AliFemtoQinvCorrFctn::AddRealPair(const AliFemtoPair* pair){
113   double Qinv = fabs(pair->qInv());   // note - qInv() will be negative for identical pairs...
114   fNumerator->Fill(Qinv);
115   //  cout << "AliFemtoQinvCorrFctn::AddRealPair : " << pair->qInv() << " " << Qinv <<
116   //" " << pair->track1().FourMomentum() << " " << pair->track2().FourMomentum() << endl;
117 }
118 //____________________________
119 void AliFemtoQinvCorrFctn::AddMixedPair(const AliFemtoPair* pair){
120   double weight = 1.0;
121   double Qinv = fabs(pair->qInv());   // note - qInv() will be negative for identical pairs...
122   fDenominator->Fill(Qinv,weight);
123 }
124
125