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