CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes | Friends
BackgroundHandler Class Reference

#include <BackgroundHandler.h>

Public Member Functions

std::pair< double, double > backgroundFunction (const bool doBackgroundFit, const double *parval, const int resTotNum, const int ires, const bool *resConsidered, const double *ResMass, const double ResHalfWidth[], const int MuonType, const double &mass, const double &eta1, const double &eta2)
 
 BackgroundHandler (const std::vector< int > &identifiers, const std::vector< double > &leftWindowBorders, const std::vector< double > &rightWindowBorders, const double *ResMass, const double *massWindowHalfWidth)
 
bool checkBackgroundWindow (const double &mass, const int iRegion)
 Check if the mass value is inside the given background region. More...
 
void countEventsInAllWindows (const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &muonPairs, const double &weight)
 
void initializeParNums ()
 Initialize the parNums to be used in the shifts of parval. More...
 
int regionsParNum ()
 Returns the total number of parameters used for the regions. More...
 
void rescale (std::vector< double > &parBgr, const double *ResMass, const double *massWindowHalfWidth, const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &muonPairs, const double &weight=1.)
 
double resMass (const bool doBackgroundFit, const int ires)
 
void setParameters (double *Start, double *Step, double *Mini, double *Maxi, int *ind, TString *parname, const std::vector< double > &parBgr, const std::vector< int > &parBgrOrder, const int muonType)
 Sets initial parameters for all the functions. More...
 
bool unlockParameter (const std::vector< int > &resfind, const unsigned int ipar)
 returns true if the parameter is to be unlocked More...
 
std::pair< double, double > windowBorders (const bool doBackgroundFit, const int ires)
 Returns the appropriate window borders depending on whether the background is being fitted and on the resonance. More...
 
 ~BackgroundHandler ()
 

Private Member Functions

void consistencyCheck (const std::vector< int > &identifiers, const std::vector< double > &leftWindowBorders, const std::vector< double > &rightWindowBorders) const (false)
 Used to check the consistency of passed parameters. More...
 

Private Attributes

std::vector< MassWindowbackgroundWindow_
 
int parNumsRegions_ [3]
 
int parNumsResonances_ [6]
 
int regToResHW_ [3]
 
std::vector< MassWindowresonanceWindow_
 
int resToReg_ [6]
 

Friends

class TestBackgroundHandler
 

Detailed Description

This class is used to handle the different background functions for the different regions.
It uses the backgroundFunctions defined in Functions.h and the backgroundFunctionService defined in Functions.cc.
More details are in the description of backgroundFunctionBase in Functions.h.

A bool selects if to use the regions functions or the resonances functions.

Definition at line 32 of file BackgroundHandler.h.

Constructor & Destructor Documentation

BackgroundHandler::BackgroundHandler ( const std::vector< int > &  identifiers,
const std::vector< double > &  leftWindowBorders,
const std::vector< double > &  rightWindowBorders,
const double *  ResMass,
const double *  massWindowHalfWidth 
)

Definition at line 11 of file BackgroundHandler.cc.

References backgroundFunctionService(), backgroundWindow_, consistencyCheck(), mps_fire::i, initializeParNums(), ResonanceBuilder::mass, regToResHW_, resonanceWindow_, and resToReg_.

16 {
17  // : leftWindowFactors_(leftWindowFactors),
18  // rightWindowFactors_(rightWindowFactors)
19 
20  // Define the correspondence between regions and halfWidth to use
21  // Defines also the function type to use (but they are checked to be consistent over a region)
22  regToResHW_[0] = 0; // Region 0 use the one from Z
23  regToResHW_[1] = 3; // Region 1 use the one from Upsilon1S
24  regToResHW_[2] = 5; // Region 2 use the one from J/Psi
25 
26  // Define the correspondence between resonances and regions
27  resToReg_[0] = 0; // Z
28  resToReg_[1] = 1; // Upsilon3S
29  resToReg_[2] = 1; // Upsilon2S
30  resToReg_[3] = 1; // Upsilon1S
31  resToReg_[4] = 2; // Psi2S
32  resToReg_[5] = 2; // J/Psi
33 
34  // Throws cms::Exception("Configuration") in case the parameters are not what is expected
35  consistencyCheck(identifiers, leftWindowBorders, rightWindowBorders);
36 
37  // Build the resonance windows
38  for( unsigned int i=0; i<6; ++i ) {
39  double mass = ResMass[i];
40  double lowerLimit = mass - massWindowHalfWidth[i];
41  double upperLimit = mass + massWindowHalfWidth[i];
42  resonanceWindow_.push_back( MassWindow( mass, lowerLimit, upperLimit,
43  std::vector<unsigned int>(1,i),
45  lowerLimit,
46  upperLimit) ) );
47  }
48 
49  // Build the background windows
50  // ----------------------------
51  // Compute the mass center of each region
52  double resMassForRegion[3];
53  resMassForRegion[0] = ResMass[0];
54  resMassForRegion[1] = (ResMass[1]+ResMass[2]+ResMass[3])/3;
55  resMassForRegion[2] = (ResMass[4]+ResMass[5])/2;
56 
57  // Define which resonance is in which background window
58  std::vector<std::vector<unsigned int> > indexes;
59  indexes.push_back(std::vector<unsigned int>(1,0));
60  indexes.push_back(std::vector<unsigned int>());
61  for( int i=1; i<=3; ++i ) { indexes[1].push_back(i); }
62  indexes.push_back(std::vector<unsigned int>());
63  for( int i=4; i<=5; ++i ) { indexes[2].push_back(i); }
64 
65  unsigned int i=0;
66  typedef std::vector<unsigned int> indexVec;
67  for(auto const& index : indexes) {
68  // double lowerLimit = resMassForRegion[i] - massWindowHalfWidth[regToResHW_[i]]*leftWindowFactors[i];
69  // double upperLimit = resMassForRegion[i] + massWindowHalfWidth[regToResHW_[i]]*rightWindowFactors[i];
70  // backgroundWindow_.push_back( MassWindow( resMassForRegion[i], lowerLimit, upperLimit, index,
71  // backgroundFunctionService(identifiers[i], lowerLimit, upperLimit ) ) );
72  backgroundWindow_.push_back( MassWindow( resMassForRegion[i], leftWindowBorders[i], rightWindowBorders[i], index,
73  backgroundFunctionService(identifiers[i], leftWindowBorders[i], rightWindowBorders[i] ) ) );
74  ++i;
75  }
76  // Initialize the parNums to be used in the shifts of parval
78 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
void consistencyCheck(const std::vector< int > &identifiers, const std::vector< double > &leftWindowBorders, const std::vector< double > &rightWindowBorders) const (false)
Used to check the consistency of passed parameters.
backgroundFunctionBase * backgroundFunctionService(const int identifier, const double &lowerLimit, const double &upperLimit)
Service to build the background functor corresponding to the passed identifier.
Definition: Functions.cc:62
void initializeParNums()
Initialize the parNums to be used in the shifts of parval.
BackgroundHandler::~BackgroundHandler ( )

Definition at line 80 of file BackgroundHandler.cc.

81 {
82 }

Member Function Documentation

std::pair< double, double > BackgroundHandler::backgroundFunction ( const bool  doBackgroundFit,
const double *  parval,
const int  resTotNum,
const int  ires,
const bool *  resConsidered,
const double *  ResMass,
const double  ResHalfWidth[],
const int  MuonType,
const double &  mass,
const double &  eta1,
const double &  eta2 
)

Returns the background fraction parameter (parBgr[0], but shifted to the correct function) and the value returned by the background function.
Depending on the value of doBackgroundFit it returns the values for the regions or the resonances.

Definition at line 225 of file BackgroundHandler.cc.

References backgroundWindow_, ires, parNumsRegions_, parNumsResonances_, resonanceWindow_, and resToReg_.

Referenced by checkBackgroundWindow(), and MuScleFitUtils::massProb().

230 {
231  if( doBackgroundFit ) {
232  // Return the values for the region
233  int iReg = resToReg_[ires];
234  // return std::make_pair( parval[parNumsRegions_[iReg]] * backgroundWindow_[iReg].backgroundFunction()->fracVsEta(&(parval[parNumsRegions_[iReg]]), resEta),
235  return std::make_pair( parval[parNumsRegions_[iReg]] * backgroundWindow_[iReg].backgroundFunction()->fracVsEta(&(parval[parNumsRegions_[iReg]]), eta1, eta2 ),
236  (*(backgroundWindow_[iReg].backgroundFunction()))( &(parval[parNumsRegions_[iReg]]), mass, eta1, eta2 ) );
237  // return std::make_pair( backgroundWindow_[iReg].backgroundFunction()->fracVsEta(&(parval[parNumsRegions_[iReg]]), eta1, eta2),
238  // (*(backgroundWindow_[iReg].backgroundFunction()))( &(parval[parNumsRegions_[iReg]]), mass, eta1, eta2 ) );
239  }
240  // Return the values for the resonance
241  // return std::make_pair( parval[parNumsResonances_[ires]] * resonanceWindow_[ires].backgroundFunction()->fracVsEta(&(parval[parNumsResonances_[ires]]), resEta),
242  // (*(resonanceWindow_[ires].backgroundFunction()))( &(parval[parNumsResonances_[ires]]), mass, resEta ) );
243  return std::make_pair( parval[parNumsResonances_[ires]] * resonanceWindow_[ires].backgroundFunction()->fracVsEta(&(parval[parNumsResonances_[ires]]), eta1, eta2),
244  (*(resonanceWindow_[ires].backgroundFunction()))( &(parval[parNumsResonances_[ires]]), mass, eta1, eta2 ) );
245 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
int ires[2]
std::pair< double, double > backgroundFunction(const bool doBackgroundFit, const double *parval, const int resTotNum, const int ires, const bool *resConsidered, const double *ResMass, const double ResHalfWidth[], const int MuonType, const double &mass, const double &eta1, const double &eta2)
bool BackgroundHandler::checkBackgroundWindow ( const double &  mass,
const int  iRegion 
)
inline

Check if the mass value is inside the given background region.

Definition at line 55 of file BackgroundHandler.h.

References backgroundFunction(), backgroundWindow_, consistencyCheck(), countEventsInAllWindows(), ires, ResonanceBuilder::mass, noexcept, rescale(), resMass(), setParameters(), cmsPerfCommons::Step, unlockParameter(), and windowBorders().

56  {
57  return backgroundWindow_[iRegion].isIn(mass);
58  }
std::vector< MassWindow > backgroundWindow_
void BackgroundHandler::consistencyCheck ( const std::vector< int > &  identifiers,
const std::vector< double > &  leftWindowBorders,
const std::vector< double > &  rightWindowBorders 
) const
private

Used to check the consistency of passed parameters.

Definition at line 273 of file BackgroundHandler.cc.

References Exception.

Referenced by BackgroundHandler(), and checkBackgroundWindow().

276 {
277  if( leftWindowBorders.size() != rightWindowBorders.size() ) {
278  throw cms::Exception("Configuration") << "BackgroundHandler::BackgroundHandler: leftWindowBorders.size() = " << leftWindowBorders.size()
279  << " != rightWindowBorders.size() = " << rightWindowBorders.size() << std::endl;
280  }
281  if( leftWindowBorders.size() != 3 ) {
282  throw cms::Exception("Configuration") << "BackgroundHandler::BackgroundHandler: leftWindowBorders.size() = rightWindowBorders.size() = "
283  << leftWindowBorders.size() << " != 3" << std::endl;
284  }
285  if( identifiers.size() != 3 ) {
286  throw cms::Exception("Configuration") << "BackgroundHandler::BackgroundHandler: identifiers must match the number of regions = 3" << std::endl;
287  }
288 }
void BackgroundHandler::countEventsInAllWindows ( const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &  muonPairs,
const double &  weight 
)

Definition at line 247 of file BackgroundHandler.cc.

References backgroundWindow_, resonanceWindow_, and mps_merge::weight.

Referenced by checkBackgroundWindow(), and rescale().

249 {
250  // First reset all the counters
251  for(auto& resonanceWindow : resonanceWindow_) {
252  resonanceWindow.resetCounter();
253  }
254  // Count events in background windows
255  for(auto& backgroundWindow : backgroundWindow_) {
256  backgroundWindow.resetCounter();
257  }
258 
259  // Now count the events in each window
260  std::pair<lorentzVector,lorentzVector> muonPair;
261  for(auto const& muonPair : muonPairs) {
262  // Count events in resonance windows
263  for(auto& resonanceWindow : resonanceWindow_) {
264  resonanceWindow.count((muonPair.first + muonPair.second).mass(), weight);
265  }
266  // Count events in background windows
267  for(auto& backgroundWindow : backgroundWindow_) {
268  backgroundWindow.count((muonPair.first + muonPair.second).mass(), weight);
269  }
270  }
271 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
void BackgroundHandler::initializeParNums ( )

Initialize the parNums to be used in the shifts of parval.

Definition at line 84 of file BackgroundHandler.cc.

References backgroundWindow_, mps_fire::i, parNumsRegions_, parNumsResonances_, and resonanceWindow_.

Referenced by BackgroundHandler().

85 {
86  // Initialize the parNums to be used in the shifts of parval
87  parNumsRegions_[0] = 0;
88  for( unsigned int i=1; i<backgroundWindow_.size() ; ++i ) {
89  parNumsRegions_[i] = parNumsRegions_[i-1] + backgroundWindow_[i-1].backgroundFunction()->parNum();
90  }
91  parNumsResonances_[0] = parNumsRegions_[2]+backgroundWindow_[2].backgroundFunction()->parNum();
92  for( unsigned int i=1; i<resonanceWindow_.size() ; ++i ) {
93  parNumsResonances_[i] = parNumsResonances_[i-1] + resonanceWindow_[i-1].backgroundFunction()->parNum();
94  }
95 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
int BackgroundHandler::regionsParNum ( )
inline

Returns the total number of parameters used for the regions.

Definition at line 49 of file BackgroundHandler.h.

References parNumsResonances_.

Referenced by MuScleFitUtils::minimizeLikelihood().

50  {
51  return parNumsResonances_[0];
52  }
void BackgroundHandler::rescale ( std::vector< double > &  parBgr,
const double *  ResMass,
const double *  massWindowHalfWidth,
const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &  muonPairs,
const double &  weight = 1. 
)

Computes the rescaled parameters from the regions functions to the resonances functions. It takes into account the difference in intervals and rescales the parameters so that the fraction of events is correctly accounted for.
It uses the list of all muon pairs to compute the number of events in each resonance window.

Definition at line 170 of file BackgroundHandler.cc.

References backgroundWindow_, countEventsInAllWindows(), gather_cfg::cout, parNumsRegions_, parNumsResonances_, resonanceWindow_, and resToReg_.

Referenced by checkBackgroundWindow(), and MuScleFitUtils::minimizeLikelihood().

173 {
174  countEventsInAllWindows(muonPairs, weight);
175 
176  // Loop on all regions and on all the resonances of each region and compute the background fraction
177  // for each resonance window.
178  unsigned int iRegion = 0;
179  for(auto const& backgroundWindow : backgroundWindow_)
180  {
181  // Iterator pointing to the first parameter of this background function in the full set of parameters
182  std::vector<double>::const_iterator parBgrIt = (parBgr.begin()+parNumsRegions_[iRegion]);
183  TF1 * backgroundFunctionForIntegral = backgroundWindow.backgroundFunction()->functionForIntegral(parBgrIt);
184  // WARNING: this expects the background fraction parameter to be parBgr[0] for all the background functions.
185  double kOld = *parBgrIt;
186  double Nbw = backgroundWindow.events();
187  double Ibw = backgroundFunctionForIntegral->Integral(backgroundWindow.lowerBound(),
188  backgroundWindow.upperBound());
189 
190  // index is the index of the resonance in the background window
191  for( unsigned int index : *(backgroundWindow.indexes()) )
192  {
193  // First set all parameters of the resonance window background function to those of the corresponding region
194  for( int iPar = 0; iPar < resonanceWindow_[index].backgroundFunction()->parNum(); ++iPar ) {
195  parBgr[parNumsResonances_[index]+iPar] = parBgr[parNumsRegions_[resToReg_[index]]+iPar];
196  }
197  // Estimated fraction of events in the resonance window
198  double Irw = backgroundFunctionForIntegral->Integral(resonanceWindow_[index].lowerBound(),
199  resonanceWindow_[index].upperBound());
200  double Nrw = resonanceWindow_[index].events();
201 
202  // Ibw is 1 (to avoid effects from approximation errors we set it to 1 and do not include it in the computation).
203  if( Nrw != 0 ) parBgr[parNumsResonances_[index]] = kOld*Nbw/Nrw*Irw;
204  else parBgr[parNumsResonances_[index]] = 0.;
205 
206  // Protect against fluctuations of number of events which could cause the fraction to go above 1.
207  if( parBgr[parNumsResonances_[index]] > 1. ) parBgr[parNumsResonances_[index]] = 1.;
208 
209  double kNew = parBgr[parNumsResonances_[index]];
210  std::cout << "For resonance = " << index << std::endl;
211  std::cout << "backgroundWindow.lowerBound = " << backgroundWindow.lowerBound() << std::endl;
212  std::cout << "backgroundWindow.upperBound = " << backgroundWindow.upperBound() << std::endl;
213  std::cout << "parNumsResonances_["<<index<<"] = " << parNumsResonances_[index] << std::endl;
214  std::cout << "Nbw = " << Nbw << ", Ibw = " << Ibw << std::endl;
215  std::cout << "Nrw = " << Nrw << ", Irw = " << Irw << std::endl;
216  std::cout << "k = " << kOld << ", k' = " << parBgr[parNumsResonances_[index]] << std::endl;
217  std::cout << "background fraction in background window = Nbw*k = " << Nbw*kOld << std::endl;
218  std::cout << "background fraction in resonance window = Nrw*k' = " << Nrw*kNew << std::endl;
219  }
220  ++iRegion;
221  delete backgroundFunctionForIntegral;
222  }
223 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
Definition: weight.py:1
void countEventsInAllWindows(const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &muonPairs, const double &weight)
double BackgroundHandler::resMass ( const bool  doBackgroundFit,
const int  ires 
)

Returns the appropriate resMass value depending on whether the background is being fitted and on the resonance.
The resMass used for the region is the mean of the mass of the corresponding resonances, so for the Z is the same Z mass, for the Upsilons is the arithmetic mean of the Upsilon masses and the same for the J/Psi and Psi2S region.

Definition at line 158 of file BackgroundHandler.cc.

References backgroundWindow_, ires, resonanceWindow_, and resToReg_.

Referenced by checkBackgroundWindow(), and MuScleFitUtils::massProb().

159 {
160  if( doBackgroundFit ) {
161  // Fitting the background: use the regions
162  return backgroundWindow_[resToReg_[ires]].mass();
163  }
164  else {
165  // Not fitting the background: use the resonances
166  return resonanceWindow_[ires].mass();
167  }
168 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
int ires[2]
void BackgroundHandler::setParameters ( double *  Start,
double *  Step,
double *  Mini,
double *  Maxi,
int *  ind,
TString *  parname,
const std::vector< double > &  parBgr,
const std::vector< int > &  parBgrOrder,
const int  muonType 
)

Sets initial parameters for all the functions.

Definition at line 97 of file BackgroundHandler.cc.

References backgroundWindow_, parNumsRegions_, parNumsResonances_, resonanceWindow_, and edm::shift.

Referenced by checkBackgroundWindow(), and MuScleFitUtils::minimizeLikelihood().

98 {
99  std::vector<double>::const_iterator parBgrIt = parBgr.begin();
100  std::vector<int>::const_iterator parBgrOrderIt = parBgrOrder.begin();
101  // Set the parameters for the regions only if this is not a rescaling
102  for( int iReg = 0; iReg < 3; ++iReg ) {
103  int shift = parNumsRegions_[iReg];
104  backgroundWindow_[iReg].backgroundFunction()->setParameters( &(Start[shift]), &(Step[shift]), &(Mini[shift]),
105  &(Maxi[shift]), &(ind[shift]), &(parname[shift]),
106  parBgrIt+shift, parBgrOrderIt+shift, muonType );
107  }
108  for( int iRes = 0; iRes < 6; ++iRes ) {
109  // parNumsResonances is already shifted for the regions parameters
110  int shift = parNumsResonances_[iRes];
111  resonanceWindow_[iRes].backgroundFunction()->setParameters( &(Start[shift]), &(Step[shift]), &(Mini[shift]),
112  &(Maxi[shift]), &(ind[shift]), &(parname[shift]),
113  parBgrIt+shift, parBgrOrderIt+shift, muonType );
114  }
115 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
static unsigned int const shift
bool BackgroundHandler::unlockParameter ( const std::vector< int > &  resfind,
const unsigned int  ipar 
)

returns true if the parameter is to be unlocked

Definition at line 117 of file BackgroundHandler.cc.

References parNumsRegions_, and parNumsResonances_.

Referenced by checkBackgroundWindow(), and MuScleFitUtils::minimizeLikelihood().

118 {
119  // parNumsRegions_ are shifted: [1] contains the number of parameters for 0 and so on.
120  if( ipar < unsigned(parNumsRegions_[1]) && resfind[0] > 0 ) {
121  return true;
122  }
123  if( ipar >= unsigned(parNumsRegions_[1]) && ipar < unsigned(parNumsRegions_[2]) && ( resfind[1] > 0 || resfind[2] > 0 || resfind[3] > 0 ) ) {
124  return true;
125  }
126  // The first of parNumsResonances_ has the sum of parNums of the regions.
127  if( ipar >= unsigned(parNumsRegions_[2]) && ipar < unsigned(parNumsResonances_[0]) && ( resfind[4] > 0 || resfind[5] > 0 ) ) {
128  return true;
129  }
130  return false;
131 }
std::pair< double, double > BackgroundHandler::windowBorders ( const bool  doBackgroundFit,
const int  ires 
)

Returns the appropriate window borders depending on whether the background is being fitted and on the resonance.

Definition at line 145 of file BackgroundHandler.cc.

References backgroundWindow_, resonanceWindow_, and resToReg_.

Referenced by checkBackgroundWindow(), MuScleFitUtils::computeWeight(), MuScleFitUtils::massProb(), and MuScleFitUtils::minimizeLikelihood().

146 {
147  if( doBackgroundFit ) {
148  // Fitting the background: use the regions
149  return std::make_pair(backgroundWindow_[resToReg_[ires]].lowerBound(), backgroundWindow_[resToReg_[ires]].upperBound());
150  }
151  else {
152  // Not fitting the background: use the resonances
153  // return std::make_pair(1.,1.);
154  return std::make_pair(resonanceWindow_[ires].lowerBound(), resonanceWindow_[ires].upperBound());
155  }
156 }
std::vector< MassWindow > backgroundWindow_
std::vector< MassWindow > resonanceWindow_
int ires[2]

Friends And Related Function Documentation

friend class TestBackgroundHandler
friend

Definition at line 35 of file BackgroundHandler.h.

Member Data Documentation

std::vector<MassWindow> BackgroundHandler::backgroundWindow_
private
int BackgroundHandler::parNumsRegions_[3]
private
int BackgroundHandler::parNumsResonances_[6]
private
int BackgroundHandler::regToResHW_[3]
private

Definition at line 108 of file BackgroundHandler.h.

Referenced by BackgroundHandler().

std::vector<MassWindow> BackgroundHandler::resonanceWindow_
private
int BackgroundHandler::resToReg_[6]
private