StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
MinvCorrFctnM_vs_Phi.cxx
1 /***************************************************************************
2  *
3  * $Id: MinvCorrFctnM_vs_Phi.cxx,v 1.1 2000/05/25 22:09:15 laue Exp $
4  *
5  * Author: Frank Laue, Ohio State, laue@mps.ohio-state.edu
6  ***************************************************************************
7  *
8  * Description: part of STAR HBT Framework: StHbtMaker package
9  * A simple invariant-mass vs azimuth angle correlation function
10  *
11  ***************************************************************************
12  *
13  * $Log: MinvCorrFctnM_vs_Phi.cxx,v $
14  * Revision 1.1 2000/05/25 22:09:15 laue
15  * New correlation function. Minv vs azimuth angle
16  *
17  *
18  **************************************************************************/
19 
20 
21 #include "StHbtMaker/CorrFctn/MinvCorrFctnM_vs_Phi.h"
22 //#include "StHbtMaker/Infrastructure/StHbtHisto.hh"
23 #include <cstdio>
24 
25 #ifdef __ROOT__
26 ClassImp(MinvCorrFctnM_vs_Phi)
27 #endif
28 
29 //____________________________
30 MinvCorrFctnM_vs_Phi::MinvCorrFctnM_vs_Phi(char* title,
31  const int& nbins1, const float& MinvLo1, const float& MinvHi1,
32  const int& nbins2, const float& MinvLo2, const float& MinvHi2){
33  // set up numerator
34  char theTitle[100];
35  const char* TitNum = "MinvCorrFctnM_vs_Phi_Num";
36  sprintf(theTitle,"Num %s",title);
37  mNumerator = new StHbt2DHisto(TitNum,theTitle,nbins1,MinvLo1,MinvHi1,nbins2,MinvLo2,MinvHi2);
38  // set up denominator
39  const char* TitDen = "MinvCorrFctnM_vs_Phi_Den";
40  sprintf(theTitle,"Den %s",title);
41  mDenominator = new StHbt2DHisto(TitDen,theTitle,nbins1,MinvLo1,MinvHi1,nbins2,MinvLo2,MinvHi2);
42  // set up difference
43  const char* TitDif = "MinvCorrFctnM_vs_Phi_Dif";
44  sprintf(theTitle,"Dif %s",title);
45  mDifference = new StHbt2DHisto(TitDif,theTitle,nbins1,MinvLo1,MinvHi1,nbins2,MinvLo2,MinvHi2);
46  // this next bit is unfortunately needed so that we can have many histos of same "title"
47  // it is neccessary if we typedef StHbt1DHisto to TH1d (which we do)
48 
49  mNumerator->Sumw2();
50  mDenominator->Sumw2();
51  mDifference->Sumw2();
52 
53  mNumerator->SetDirectory(0);
54  mDenominator->SetDirectory(0);
55  mDifference->SetDirectory(0);
56 
57 }
58 
59 //____________________________
60 MinvCorrFctnM_vs_Phi::~MinvCorrFctnM_vs_Phi(){
61  delete mNumerator;
62  delete mDenominator;
63  delete mDifference;
64 }
65 //_________________________
66 void MinvCorrFctnM_vs_Phi::Finish(){
67  // here is where we should normalize, fit, etc...
68  // we should NOT Draw() the histos (as I had done it below),
69  // since we want to insulate ourselves from root at this level
70  // of the code. Do it instead at root command line with browser.
71  // mNumerator->Draw();
72  // mDenominator->Draw();
73  //mRatio->Draw();
74  double NumeratorInt = mNumerator->Integral();
75  double DenominatorInt = mDenominator->Integral();
76  mDifference->Add(mNumerator,mDenominator,1.0,-1*NumeratorInt/DenominatorInt);
77 
78 }
79 
80 //____________________________
81 StHbtString MinvCorrFctnM_vs_Phi::Report(){
82  string stemp = "Minv Correlation Function Report:\n";
83  char ctemp[100];
84  sprintf(ctemp,"Number of entries in numerator:\t%E\n",mNumerator->GetEntries());
85  stemp += ctemp;
86  sprintf(ctemp,"Number of entries in denominator:\t%E\n",mDenominator->GetEntries());
87  stemp += ctemp;
88  sprintf(ctemp,"Number of entries in difference:\t%E\n",mDifference->GetEntries());
89  stemp += ctemp;
90  StHbtString returnThis = stemp;
91  return returnThis;
92 }
93 //____________________________
94 void MinvCorrFctnM_vs_Phi::AddRealPair(const StHbtPair* pair){
95  mNumerator->Fill( pair->mInv(), pair->fourMomentumSum().vect().phi(), 1.);
96 }
97 //____________________________
98 void MinvCorrFctnM_vs_Phi::AddMixedPair(const StHbtPair* pair){
99  mDenominator->Fill( pair->mInv(), pair->fourMomentumSum().vect().phi(), 1.);
100 }
101 
102