laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
Lau2DAbsHistDP.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2013 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
29 #include "Lau2DAbsHistDP.hh"
30 
31 #include "LauBifurcatedGaussPdf.hh"
32 #include "LauDaughters.hh"
33 #include "LauKinematics.hh"
34 #include "LauRandom.hh"
35 
36 #include "TAxis.h"
37 #include "TH2.h"
38 #include "TRandom.h"
39 #include "TSystem.h"
40 
41 #include <iostream>
42 
44  Bool_t useUpperHalfOnly,
45  Bool_t squareDP ) :
46  kinematics_( ( daughters != 0 ) ? daughters->getKinematics() : 0 ),
47  upperHalf_( useUpperHalfOnly ),
48  squareDP_( squareDP )
49 {
50  if ( squareDP && ! daughters->squareDP() ) {
51  // The histogram provided is defined in the square DP but the
52  // kinematics object has calculation of the square DP
53  // co-ordinates disabled, so need to enable it,
54  // which requires a bit of a unpleasant const_cast
55  std::cerr << "WARNING in Lau2DAbsHistDP constructor : forcing kinematics to calculate the required square DP co-ordinates"
56  << std::endl;
57  LauKinematics* kine = const_cast<LauKinematics*>( kinematics_ );
58  kine->squareDP( kTRUE );
59  }
60 }
61 
63 {
64 }
65 
67 {
68  TRandom* random = LauRandom::randomFun();
69 
70  Int_t nBinsX = static_cast<Int_t>( hist->GetNbinsX() );
71  Int_t nBinsY = static_cast<Int_t>( hist->GetNbinsY() );
72 
73  for ( Int_t i( 0 ); i < nBinsX; ++i ) {
74  for ( Int_t j( 0 ); j < nBinsY; ++j ) {
75  Double_t currentContent = hist->GetBinContent( i + 1, j + 1 );
76  Double_t currentError = hist->GetBinError( i + 1, j + 1 );
77  Double_t newContent = random->Gaus( currentContent, currentError );
78  if ( newContent < 0.0 ) {
79  hist->SetBinContent( i + 1, j + 1, 0.0 );
80  } else {
81  hist->SetBinContent( i + 1, j + 1, newContent );
82  }
83  }
84  }
85 }
86 
87 void Lau2DAbsHistDP::doBinFluctuation( TH2* hist, const TH2* errorHi, const TH2* errorLo )
88 {
89  LauParameter* mean = new LauParameter( "mean", 0.5, 0.0, 1.0, kFALSE );
90  LauParameter* sigL = new LauParameter( "sigmaL", 0.5, 0.0, 1.0, kFALSE );
91  LauParameter* sigR = new LauParameter( "sigmaR", 0.5, 0.0, 1.0, kFALSE );
92 
93  std::vector<LauAbsRValue*> pars( 3 );
94  pars[0] = mean;
95  pars[1] = sigL;
96  pars[2] = sigR;
97 
98  const TString varName( "tmp" );
99 
100  Int_t nBinsX = static_cast<Int_t>( hist->GetNbinsX() );
101  Int_t nBinsY = static_cast<Int_t>( hist->GetNbinsY() );
102 
103  LauFitData genData;
104 
105  for ( Int_t i( 0 ); i < nBinsX; ++i ) {
106  for ( Int_t j( 0 ); j < nBinsY; ++j ) {
107  const Double_t currentContent = hist->GetBinContent( i + 1, j + 1 );
108  const Double_t currentErrorLo = errorLo->GetBinContent( i + 1, j + 1 );
109  const Double_t currentErrorHi = errorHi->GetBinContent( i + 1, j + 1 );
110 
111  mean->value( currentContent );
112  sigL->value( currentErrorLo );
113  sigR->value( currentErrorHi );
114 
115  const Double_t minVal = TMath::Max( 0.0, currentContent - 5.0 * currentErrorLo );
116  const Double_t maxVal = TMath::Min( 1.0, currentContent + 5.0 * currentErrorHi );
117 
118  LauBifurcatedGaussPdf bfgaus( varName, pars, minVal, maxVal );
119  bfgaus.heightUpToDate( kFALSE );
120  genData = bfgaus.generate( 0 );
121 
122  const Double_t newContent = genData[varName];
123 
124  hist->SetBinContent( i + 1, j + 1, newContent );
125  }
126  }
127 
128  delete pars[0];
129  delete pars[1];
130  delete pars[2];
131  pars.clear();
132 }
133 
134 void Lau2DAbsHistDP::raiseOrLowerBins( TH2* hist, const Double_t avEff, const Double_t avEffError )
135 {
136  TRandom* random = LauRandom::randomFun();
137 
138  Double_t curAvg = this->computeAverageContents( hist );
139  Double_t newAvg = random->Gaus( avEff, avEffError );
140 
141  hist->Scale( newAvg / curAvg );
142 }
143 
144 Double_t Lau2DAbsHistDP::computeAverageContents( const TH2* hist ) const
145 {
146  Double_t totalContent( 0.0 );
147  Int_t binsWithinDPBoundary( 0 );
148 
149  Int_t nBinsX = static_cast<Int_t>( hist->GetNbinsX() );
150  Int_t nBinsY = static_cast<Int_t>( hist->GetNbinsY() );
151 
152  // Loop through the bins and include any that have their centre or any
153  // of the four corners within the kinematic boundary
154  for ( Int_t i( 0 ); i < nBinsX; ++i ) {
155  Double_t binXCentre = hist->GetXaxis()->GetBinCenter( i + 1 );
156  Double_t binXLowerEdge = hist->GetXaxis()->GetBinLowEdge( i + 1 );
157  Double_t binXUpperEdge = hist->GetXaxis()->GetBinUpEdge( i + 1 );
158 
159  for ( Int_t j( 0 ); j < nBinsY; ++j ) {
160  Double_t binYCentre = hist->GetYaxis()->GetBinCenter( i + 1 );
161  Double_t binYLowerEdge = hist->GetYaxis()->GetBinLowEdge( i + 1 );
162  Double_t binYUpperEdge = hist->GetYaxis()->GetBinUpEdge( i + 1 );
163 
164  if ( this->withinDPBoundaries( binXCentre, binYCentre ) ||
165  this->withinDPBoundaries( binXLowerEdge, binYLowerEdge ) ||
166  this->withinDPBoundaries( binXUpperEdge, binYUpperEdge ) ||
167  this->withinDPBoundaries( binXLowerEdge, binYUpperEdge ) ||
168  this->withinDPBoundaries( binXUpperEdge, binYLowerEdge ) ) {
169 
170  totalContent += hist->GetBinContent( i + 1, j + 1 );
171  ++binsWithinDPBoundary;
172  }
173  }
174  }
175 
176  return totalContent / binsWithinDPBoundary;
177 }
178 
179 Bool_t Lau2DAbsHistDP::withinDPBoundaries( Double_t x, Double_t y ) const
180 {
181  return squareDP_ ? kinematics_->withinSqDPLimits( x, y ) : kinematics_->withinDPLimits( x, y );
182 }
183 
184 void Lau2DAbsHistDP::getUpperHalf( Double_t& x, Double_t& y ) const
185 {
186  if ( upperHalf_ == kTRUE ) {
187  if ( squareDP_ == kFALSE && x > y ) {
188  Double_t temp = y;
189  y = x;
190  x = temp;
191  } else if ( squareDP_ == kTRUE && y > 0.5 ) {
192  y = 1.0 - y;
193  }
194  }
195 }
File containing declaration of LauBifurcatedGaussPdf class.
Bool_t squareDP_
Boolean for using square DP variables.
File containing LauRandom namespace.
virtual LauFitData generate(const LauKinematics *kinematics)
Generate an event from the PDF.
Definition: LauAbsPdf.cc:338
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
void doBinFluctuation(TH2 *hist)
Fluctuate the contents of each histogram bin independently, in accordance with their errors.
Bool_t withinDPBoundaries(Double_t x, Double_t y) const
Check whether the given co-ordinates are within the kinematic boundary.
Double_t value() const
The value of the parameter.
Bool_t squareDP() const
Determine to use or not the square Dalitz plot.
Double_t computeAverageContents(const TH2 *hist) const
Compute the average bin content for bins within the kinematic boundary.
File containing declaration of Lau2DAbsHistDP class.
Bool_t withinDPLimits(const Double_t m13Sq, const Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot.
virtual ~Lau2DAbsHistDP()
Copy constructor.
void raiseOrLowerBins(TH2 *hist, const Double_t avEff, const Double_t avEffError)
Rescale the histogram bin contents based on the desired average efficiency and its uncertainty.
void squareDP(const Bool_t calcSquareDPCoords)
Enable/disable the calculation of square Dalitz plot co-ordinates.
Bool_t upperHalf_
Boolean for using the upper half of DP.
Class for defining a bifurcated Gaussian PDF.
File containing declaration of LauDaughters class.
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:33
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
Lau2DAbsHistDP(const LauDaughters *daughters, Bool_t useUpperHalfOnly=kFALSE, Bool_t squareDP=kFALSE)
Constructor.
const LauKinematics * kinematics_
Kinematics used to check events are within DP boundary.
void getUpperHalf(Double_t &x, Double_t &y) const
If only using the upper half of the (symmetric) DP then transform into the correct half.
std::map< TString, Double_t > LauFitData
Type for holding event data.
Class for calculating 3-body kinematic quantities.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
virtual Bool_t heightUpToDate() const
Check if the maximum height of the PDF is up to date.
Definition: LauAbsPdf.hh:287
File containing declaration of LauKinematics class.
Bool_t withinSqDPLimits(const Double_t mPrime, const Double_t thetaPrime) const
Check whether a given (m',theta') point is within the kinematic limits of the Dalitz plot.