CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BackgroundHandler.h
Go to the documentation of this file.
1 #ifndef BackgroundHandler_h
2 #define BackgroundHandler_h
3 
6 #include <vector>
12 #include <CLHEP/Vector/LorentzVector.h>
19 
20 // Unit test for testing BackgroundHandler
21 class TestBackgroundHandler;
22 
33 {
34  // For tests
35  friend class TestBackgroundHandler;
36 
37 public:
38  BackgroundHandler( const std::vector<int> & identifiers,
39  const std::vector<double> & leftWindowBorders,
40  const std::vector<double> & rightWindowBorders,
41  const double * ResMass,
42  const double * massWindowHalfWidth );
44 
46  void initializeParNums();
47 
49  inline int regionsParNum()
50  {
51  return parNumsResonances_[0];
52  }
53 
55  bool checkBackgroundWindow(const double & mass, const int iRegion)
56  {
57  return backgroundWindow_[iRegion].isIn(mass);
58  }
59 
60  void countEventsInAllWindows(const std::vector<std::pair<reco::Particle::LorentzVector,reco::Particle::LorentzVector> > & muonPairs,
61  const double & weight);
62 
64  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);
65 
67  bool unlockParameter(const std::vector<int> & resfind, const unsigned int ipar);
68 
70  std::pair<double, double> windowBorders( const bool doBackgroundFit, const int ires );
71 
77  double resMass( const bool doBackgroundFit, const int ires );
78 
85  void rescale( std::vector<double> & parBgr, const double * ResMass, const double * massWindowHalfWidth,
86  const std::vector<std::pair<reco::Particle::LorentzVector,reco::Particle::LorentzVector> > & muonPairs,
87  const double & weight = 1. );
88 
94  std::pair<double, double> backgroundFunction( const bool doBackgroundFit,
95  const double * parval, const int resTotNum, const int ires,
96  const bool * resConsidered, const double * ResMass, const double ResHalfWidth[],
97  /* const int MuonType, const double & mass, const double & resEta ); */
98  const int MuonType, const double & mass, const double & eta1, const double & eta2 );
99 private:
101  void consistencyCheck( const std::vector<int> & identifiers,
102  const std::vector<double> & leftWindowBorders,
103  const std::vector<double> & rightWindowBorders ) const throw(cms::Exception);
104 
105  // Correspondence between regions and halfWidths used:
106  // - for the Upsilons region we use the Upsilon
107  // - for the J/Psi and Psi2S region we use the J/Psi
108  int regToResHW_[3];
109 
110  // Correspondence between resonances and regions:
111  // - Z -> region 0
112  // - Uspilons -> region 1
113  // - J/Psi and Psi2S -> region 2
114  int resToReg_[6];
115 
116  // Used in the shifts of the parval
117  // Contains 0 as the first number and Sum_0^(ires-1)(parNum(i)) for the rest,
118  // so that calling parNums[ires] returns the sum of the number of parameters
119  // of the previous functions (0 if none) and allows to shift the parval to the
120  // parameters of the actual function.
122  // These start from the correct parameters (take into account that the parRegions are
123  // before the parResonances).
125 
126  // std::vector<double> leftWindowFactors_;
127  // std::vector<double> rightWindowFactors_;
128 
129  std::vector<MassWindow> resonanceWindow_;
130  std::vector<MassWindow> backgroundWindow_;
131 };
132 
133 #endif // BackgroundHandler_h
std::vector< MassWindow > backgroundWindow_
int regionsParNum()
Returns the total number of parameters used for the regions.
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)
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 consistencyCheck(const std::vector< int > &identifiers, const std::vector< double > &leftWindowBorders, const std::vector< double > &rightWindowBorders) const
Used to check the consistency of passed parameters.
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.
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...
bool unlockParameter(const std::vector< int > &resfind, const unsigned int ipar)
returns true if the parameter is to be unlocked
void countEventsInAllWindows(const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &muonPairs, const double &weight)
bool checkBackgroundWindow(const double &mass, const int iRegion)
Check if the mass value is inside the given background region.
friend class TestBackgroundHandler
BackgroundHandler(const std::vector< int > &identifiers, const std::vector< double > &leftWindowBorders, const std::vector< double > &rightWindowBorders, const double *ResMass, const double *massWindowHalfWidth)
void initializeParNums()
Initialize the parNums to be used in the shifts of parval.