Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:36

0001 #ifndef BackgroundHandler_h
0002 #define BackgroundHandler_h
0003 
0004 #include "MuonAnalysis/MomentumScaleCalibration/interface/Functions.h"
0005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0006 #include <vector>
0007 #include "MuonAnalysis/MomentumScaleCalibration/interface/MassWindow.h"
0008 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0009 #include "SimDataFormats/Track/interface/SimTrack.h"
0010 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0011 #include "DataFormats/Math/interface/LorentzVector.h"
0012 #include <CLHEP/Vector/LorentzVector.h>
0013 #include "DataFormats/MuonReco/interface/Muon.h"
0014 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0015 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0016 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
0017 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 // Unit test for testing BackgroundHandler
0021 class TestBackgroundHandler;
0022 
0023 /**
0024  * This class is used to handle the different background functions for the
0025  * different regions. <br>
0026  * It uses the backgroundFunctions defined in Functions.h and the
0027  * backgroundFunctionService defined in Functions.cc. <br>
0028  * More details are in the description of backgroundFunctionBase in Functions.h. <br>
0029  * <br>
0030  * A bool selects if to use the regions functions or the resonances functions.
0031  */
0032 class BackgroundHandler {
0033   // For tests
0034   friend class TestBackgroundHandler;
0035 
0036 public:
0037   BackgroundHandler(const std::vector<int>& identifiers,
0038                     const std::vector<double>& leftWindowBorders,
0039                     const std::vector<double>& rightWindowBorders,
0040                     const double* ResMass,
0041                     const double* massWindowHalfWidth);
0042   ~BackgroundHandler();
0043 
0044   /// Initialize the parNums to be used in the shifts of parval
0045   void initializeParNums();
0046 
0047   /// Returns the total number of parameters used for the regions
0048   inline int regionsParNum() { return parNumsResonances_[0]; }
0049 
0050   /// Check if the mass value is inside the given background region
0051   bool checkBackgroundWindow(const double& mass, const int iRegion) { return backgroundWindow_[iRegion].isIn(mass); }
0052 
0053   void countEventsInAllWindows(
0054       const std::vector<std::pair<reco::Particle::LorentzVector, reco::Particle::LorentzVector> >& muonPairs,
0055       const double& weight);
0056 
0057   /// Sets initial parameters for all the functions
0058   void setParameters(double* Start,
0059                      double* Step,
0060                      double* Mini,
0061                      double* Maxi,
0062                      int* ind,
0063                      TString* parname,
0064                      const std::vector<double>& parBgr,
0065                      const std::vector<int>& parBgrOrder,
0066                      const int muonType);
0067 
0068   /// returns true if the parameter is to be unlocked
0069   bool unlockParameter(const std::vector<int>& resfind, const unsigned int ipar);
0070 
0071   /// Returns the appropriate window borders depending on whether the background is being fitted and on the resonance
0072   std::pair<double, double> windowBorders(const bool doBackgroundFit, const int ires);
0073 
0074   /**
0075    * Returns the appropriate resMass value depending on whether the background is being fitted and on the resonance. <br>
0076    * 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,
0077    * for the Upsilons is the arithmetic mean of the Upsilon masses and the same for the J/Psi and Psi2S region.
0078    */
0079   double resMass(const bool doBackgroundFit, const int ires);
0080 
0081   /**
0082    * Computes the rescaled parameters from the regions functions to the
0083    * resonances functions. It takes into account the difference in intervals
0084    * and rescales the parameters so that the fraction of events is correctly accounted for. <br>
0085    * It uses the list of all muon pairs to compute the number of events in each resonance window.
0086    */
0087   void rescale(std::vector<double>& parBgr,
0088                const double* ResMass,
0089                const double* massWindowHalfWidth,
0090                const std::vector<std::pair<reco::Particle::LorentzVector, reco::Particle::LorentzVector> >& muonPairs,
0091                const double& weight = 1.);
0092 
0093   /**
0094    * Returns the background fraction parameter (parBgr[0], but shifted to the correct function) and
0095    * the value returned by the background function. <br>
0096    * Depending on the value of doBackgroundFit it returns the values for the regions or the resonances.
0097    */
0098   std::pair<double, double> backgroundFunction(const bool doBackgroundFit,
0099                                                const double* parval,
0100                                                const int resTotNum,
0101                                                const int ires,
0102                                                const bool* resConsidered,
0103                                                const double* ResMass,
0104                                                const double ResHalfWidth[],
0105                                                /* const int MuonType, const double & mass, const double & resEta ); */
0106                                                const int MuonType,
0107                                                const double& mass,
0108                                                const double& eta1,
0109                                                const double& eta2);
0110 
0111 private:
0112   /// Used to check the consistency of passed parameters
0113   void consistencyCheck(const std::vector<int>& identifiers,
0114                         const std::vector<double>& leftWindowBorders,
0115                         const std::vector<double>& rightWindowBorders) const noexcept(false);
0116 
0117   // Correspondence between regions and halfWidths used:
0118   // - for the Upsilons region we use the Upsilon
0119   // - for the J/Psi and Psi2S region we use the J/Psi
0120   int regToResHW_[3];
0121 
0122   // Correspondence between resonances and regions:
0123   // - Z -> region 0
0124   // - Uspilons -> region 1
0125   // - J/Psi and Psi2S -> region 2
0126   int resToReg_[6];
0127 
0128   // Used in the shifts of the parval
0129   // Contains 0 as the first number and Sum_0^(ires-1)(parNum(i)) for the rest,
0130   // so that calling parNums[ires] returns the sum of the number of parameters
0131   // of the previous functions (0 if none) and allows to shift the parval to the
0132   // parameters of the actual function.
0133   int parNumsRegions_[3];
0134   // These start from the correct parameters (take into account that the parRegions are
0135   // before the parResonances).
0136   int parNumsResonances_[6];
0137 
0138   // std::vector<double> leftWindowFactors_;
0139   // std::vector<double> rightWindowFactors_;
0140 
0141   std::vector<MassWindow> resonanceWindow_;
0142   std::vector<MassWindow> backgroundWindow_;
0143 };
0144 
0145 #endif  // BackgroundHandler_h