CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuScleFitUtils.h
Go to the documentation of this file.
1 #ifndef MuScleFitUtils_H
2 #define MuScleFitUtils_H
3 
11 #include <CLHEP/Vector/LorentzVector.h>
14 // #include "SimDataFormats/HepMCProduct/interface/HepMCProduct.h"
19 #include "TGraphErrors.h"
20 #include "TH2F.h"
21 #include "TMinuit.h"
22 
26 
27 #include <vector>
28 
29 // #include "Functions.h"
30 // class biasFunctionBase<std::vector<double> >;
31 // class scaleFunctionBase<double*>;
32 template <class T> class biasFunctionBase;
33 template <class T> class scaleFunctionBase;
34 class smearFunctionBase;
35 template <class T> class resolutionFunctionBase;
37 class BackgroundHandler;
38 
39 class SimTrack;
40 class TString;
41 class TTree;
42 
44 
46 public:
47 
48  // Constructor
49  // ----------
51 
52  // Destructor
53  // ----------
54  virtual ~MuScleFitUtils() {};
55 
56  // Operations
57  // ----------
58  static std::pair<SimTrack, SimTrack> findBestSimuRes( const std::vector<SimTrack>& simMuons );
59  static std::pair<lorentzVector, lorentzVector> findBestRecoRes( const std::vector<reco::LeafCandidate>& muons );
60  static std::pair <lorentzVector, lorentzVector> findGenMuFromRes( const reco::GenParticleCollection* genParticles);
61  static std::pair<lorentzVector, lorentzVector> findGenMuFromRes( const edm::HepMCProduct* evtMC );
62  static std::pair<lorentzVector, lorentzVector> findSimMuFromRes( const edm::Handle<edm::HepMCProduct> & evtMC,
63  const edm::Handle<edm::SimTrackContainer> & simTracks);
64 
65  static std::vector<TGraphErrors*> fitMass (TH2F* histo);
66  static std::vector<TGraphErrors*> fitReso (TH2F* histo);
67 
68  static lorentzVector applyScale( const lorentzVector & muon, const std::vector<double> & parval, const int charge );
69  static lorentzVector applyScale( const lorentzVector & muon, double* parval, const int charge );
70  static lorentzVector applyBias( const lorentzVector & muon, const int charge );
72  static lorentzVector fromPtEtaPhiToPxPyPz( const double* ptEtaPhiE );
73 
74  static void minimizeLikelihood();
75 
76  static double invDimuonMass( const lorentzVector & mu1, const lorentzVector & mu2 );
77  static double massResolution( const lorentzVector & mu1, const lorentzVector & mu2 );
78  static double massResolution( const lorentzVector & mu1, const lorentzVector & mu2, const std::vector<double> & parval );
79  static double massResolution( const lorentzVector & mu1, const lorentzVector & mu2, std::auto_ptr<double> parval );
80  static double massResolution( const lorentzVector & mu1, const lorentzVector & mu2, double* parval );
81  static double massResolution( const lorentzVector& mu1, const lorentzVector& mu2, const ResolutionFunction & resolFunc );
82 
83  static double massProb( const double & mass, const double & rapidity, const int ires, const double & massResol );
84  /* static double massProb( const double & mass, const double & resEta, const double & rapidity, const double & massResol, const std::vector<double> & parval, const bool doUseBkgrWindow = false ); */
85  /* static double massProb( const double & mass, const double & resEta, const double & rapidity, const double & massResol, double * parval, const bool doUseBkgrWindow = false ); */
86  static double massProb( const double & mass, const double & resEta, const double & rapidity, const double & massResol, const std::vector<double> & parval, const bool doUseBkgrWindow, const double & eta1, const double & eta2 );
87  static double massProb( const double & mass, const double & resEta, const double & rapidity, const double & massResol, double * parval, const bool doUseBkgrWindow, const double & eta1, const double & eta2 );
88  static double computeWeight( const double & mass, const int iev, const bool doUseBkgrWindow = false );
89 
90  static double deltaPhi( const double & phi1, const double & phi2 )
91  {
92  double deltaPhi = phi1 - phi2;
93  while(deltaPhi >= TMath::Pi()) deltaPhi -= 2*TMath::Pi();
94  while(deltaPhi < -TMath::Pi()) deltaPhi += 2*TMath::Pi();
95  return fabs(deltaPhi);
96  }
98  static double deltaPhiNoFabs( const double & phi1, const double & phi2 )
99  {
100  double deltaPhi = phi1 - phi2;
101  while(deltaPhi >= TMath::Pi()) deltaPhi -= 2*TMath::Pi();
102  while(deltaPhi < -TMath::Pi()) deltaPhi += 2*TMath::Pi();
103  return deltaPhi;
104  }
105  static double deltaR(const double & eta1, const double & eta2, const double & phi1, const double & phi2)
106  {
107  return sqrt( std::pow( eta1-eta2, 2 ) + std::pow( deltaPhi(phi1, phi2), 2 ) );
108  }
109 
110  static int debug; // debug option set by MuScleFit
111  static bool ResFound; // bool flag true if best resonance found (cuts on pt and eta)
112 
113  static const int totalResNum; // Total number of resonance: 6
114  static double massWindowHalfWidth[3][6]; // parameter set by MuScleFitUtils
115  static double ResGamma[6]; // parameter set by MuScleFitUtils
116  static double ResMass[6]; // parameter set by MuScleFitUtils
117  static double ResMinMass[6]; // parameter set by MuScleFitBase
118  static double crossSection[6];
119  static const double mMu2;
120  static const double muMass;
121 
122  // Array of the pdgId of resonances
123  static const unsigned int motherPdgIdArray[6];
124 
125  static unsigned int loopCounter; // parameter set by MuScleFit
126 
127  static int SmearType;
129  static int BiasType;
130  // No error, we take functions from the same group for scale and bias.
132  static int ResolFitType;
135  static int ScaleFitType;
138  static int BgrFitType;
139  // Three background regions:
140  // - one for the Z
141  // - one for the Upsilons
142  // - one for J/Psi and Psi2S
143  static const int backgroundFunctionsRegions;
144  // static backgroundFunctionBase * backgroundFunctionForRegion[];
145  // A background function for each resonance
146  // static backgroundFunctionBase * backgroundFunction[];
147 
148  // The Cross section handler takes care of computing the relative cross
149  // sections to be used depending on the resonances that are being fitted.
150  // This corresponds to a normalization of the signal pdf.
152 
153  // The background handler takes care of using the correct function in each
154  // window, use regions or resonance windows and rescale the fractions when needed
156 
157  // Parameters used to select whether to do a fit
158  static std::vector<int> doResolFit;
159  static std::vector<int> doScaleFit;
160  static std::vector<int> doCrossSectionFit;
161  static std::vector<int> doBackgroundFit;
162 
163  static int minuitLoop_;
164  static TH1D* likelihoodInLoop_;
165  static TH1D* signalProb_;
166  static TH1D* backgroundProb_;
167 
168  static bool duringMinos_;
169 
170  static std::vector<double> parSmear;
171  static std::vector<double> parBias;
172  static std::vector<double> parResol;
173  static std::vector<double> parResolStep;
174  static std::vector<double> parResolMin;
175  static std::vector<double> parResolMax;
176  static std::vector<double> parScale;
177  static std::vector<double> parScaleStep;
178  static std::vector<double> parScaleMin;
179  static std::vector<double> parScaleMax;
180  static std::vector<double> parCrossSection;
181  static std::vector<double> parBgr;
182  static std::vector<int> parResolFix;
183  static std::vector<int> parScaleFix;
184  static std::vector<int> parCrossSectionFix;
185  static std::vector<int> parBgrFix;
186  static std::vector<int> parResolOrder;
187  static std::vector<int> parScaleOrder;
188  static std::vector<int> parCrossSectionOrder;
189  static std::vector<int> parBgrOrder;
190  static std::vector<int> resfind;
191  static int FitStrategy;
192  static bool speedup; // parameter set by MuScleFit - whether to speedup processing
193  static double x[7][10000]; // smearing values set by MuScleFit constructor
194  static int goodmuon; // number of events with a usable resonance
195  static int counter_resprob;// number of times there are resolution problems
196  static double GLZValue[40][1001][1001]; // matrix with integral values of Lorentz * Gaussian
197  static double GLZNorm[40][1001]; // normalization values per each sigma
198  static double GLValue[6][1001][1001]; // matrix with integral values of Lorentz * Gaussian
199  static double GLNorm[6][1001]; // normalization values per each sigma
200  static double ResMaxSigma[6]; // max sigma of matrix
201  static double ResHalfWidth[6]; // halfwidth in matrix
202  static int nbins; // number of bins in matrix
203  static int MuonType; // 0, 1, 2 - 0 is GM, 1 is SM, 2 is track
204  static int MuonTypeForCheckMassWindow; // Reduced to be 0, 1 or 2. It is = MuonType when MuonType < 3, = 2 otherwise.
205 
206  static std::vector<std::vector<double> > parvalue;
207  // static std::map<unsigned int,std::vector<double> > parvalue;
208  static std::vector<int> parfix;
209  static std::vector<int> parorder;
210 
211  static std::vector<std::pair<lorentzVector,lorentzVector> > SavedPair;
212  static std::vector<std::pair<lorentzVector,lorentzVector> > ReducedSavedPair;
213  static std::vector<std::pair<lorentzVector,lorentzVector> > genPair;
214  static std::vector<std::pair<lorentzVector,lorentzVector> > simPair;
215 
216  static bool scaleFitNotDone_;
217 
219  // Pointer to the minuit object
220  static TMinuit * rminPtr_;
221  // Value stored to check whether to apply a new normalization to the likelihood
222  static double oldNormalization_;
223  static unsigned int normalizationChanged_;
224 
225  // This must be set to true if using events generated with Sherpa
226  static bool sherpa_;
227 
228  // Decide whether to use the rapidity bins for the Z
229  static bool rapidityBinsForZ_;
230 
231  static int iev_;
232 
233  static bool useProbsFile_;
234 
235  // Cuts on the muons to use in the fit
236  static bool separateRanges_;
237  static double minMuonPt_;
238  static double maxMuonPt_;
239  static double minMuonEtaFirstRange_;
240  static double maxMuonEtaFirstRange_;
241  static double minMuonEtaSecondRange_;
242  static double maxMuonEtaSecondRange_;
243  static double deltaPhiMinCut_;
244  static double deltaPhiMaxCut_;
245 
246  static bool debugMassResol_;
248  {
249  double dmdpt1;
250  double dmdpt2;
251  double dmdphi1;
252  double dmdphi2;
253  double dmdcotgth1;
254  double dmdcotgth2;
256 
257  // Fit accuracy and debug parameters
258  static bool startWithSimplex_;
259  static bool computeMinosErrors_;
260  static bool minimumShapePlots_;
261 
263  // static bool checkMassWindow( const double & mass, const int ires, const double & resMass, const double & leftFactor = 1., const double & rightFactor = 1. );
264  static bool checkMassWindow( const double & mass, const double & leftBorder, const double & rightBorder );
265 
267  static double probability( const double & mass, const double & massResol,
268  const double GLvalue[][1001][1001], const double GLnorm[][1001],
269  const int iRes, const int iY );
270 
271 protected:
272 
273 private:
274 
275  struct byPt {
276  bool operator() (const reco::Muon &a, const reco::Muon &b) const {
277  return a.pt() > b.pt();
278  }
279  };
280 
281 };
282 
283 extern "C" void likelihood (int& npar, double* grad, double& fval, double* xval, int flag);
284 
285 #endif
static double deltaPhiNoFabs(const double &phi1, const double &phi2)
Without fabs at the end, used to have a symmetric distribution for the resolution fits and variance c...
static std::vector< std::pair< lorentzVector, lorentzVector > > simPair
const double Pi
static std::vector< int > doScaleFit
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
static std::vector< int > doResolFit
static std::pair< lorentzVector, lorentzVector > findGenMuFromRes(const reco::GenParticleCollection *genParticles)
static double GLValue[6][1001][1001]
static std::vector< double > parBias
static std::vector< int > parCrossSectionOrder
static smearFunctionBase * smearFunction
static std::vector< int > parScaleOrder
static std::vector< int > doCrossSectionFit
static void minimizeLikelihood()
static double maxMuonEtaSecondRange_
static std::vector< int > parBgrOrder
static std::vector< int > parfix
static double deltaPhiMaxCut_
static unsigned int loopCounter
static std::vector< double > parResolMax
static std::pair< lorentzVector, lorentzVector > findSimMuFromRes(const edm::Handle< edm::HepMCProduct > &evtMC, const edm::Handle< edm::SimTrackContainer > &simTracks)
static int nbins
void likelihood(int &npar, double *grad, double &fval, double *xval, int flag)
static bool startWithSimplex_
static std::vector< int > doBackgroundFit
static std::vector< double > parResol
static bool debugMassResol_
int ires[2]
static double ResMinMass[6]
static TH1D * backgroundProb_
static BackgroundHandler * backgroundHandler
static double x[7][10000]
static int debug
static std::vector< TGraphErrors * > fitMass(TH2F *histo)
static double ResMass[6]
static double massWindowHalfWidth[3][6]
static bool speedup
static bool scaleFitNotDone_
static unsigned int normalizationChanged_
static bool ResFound
static scaleFunctionBase< std::vector< double > > * biasFunction
static bool checkMassWindow(const double &mass, const double &leftBorder, const double &rightBorder)
Method to check if the mass value is within the mass window of the i-th resonance.
static std::vector< int > parBgrFix
static std::vector< double > parResolMin
static bool minimumShapePlots_
virtual ~MuScleFitUtils()
static int MuonTypeForCheckMassWindow
static double minMuonEtaFirstRange_
static double massProb(const double &mass, const double &rapidity, const int ires, const double &massResol)
static struct MuScleFitUtils::massResolComponentsStruct massResolComponents
double charge(const std::vector< uint8_t > &Ampls)
static int BiasType
static std::vector< int > parScaleFix
static bool computeMinosErrors_
reco::Particle::LorentzVector lorentzVector
Definition: GenMuonPair.h:8
static scaleFunctionBase< std::vector< double > > * scaleFunctionForVec
static std::pair< SimTrack, SimTrack > findBestSimuRes(const std::vector< SimTrack > &simMuons)
static std::vector< std::vector< double > > parvalue
static double maxMuonPt_
static int ScaleFitType
static std::vector< std::pair< lorentzVector, lorentzVector > > genPair
static std::vector< std::pair< lorentzVector, lorentzVector > > ReducedSavedPair
static const int totalResNum
static const double muMass
static double massResolution(const lorentzVector &mu1, const lorentzVector &mu2)
static const unsigned int motherPdgIdArray[6]
static int minuitLoop_
static std::vector< double > parScaleMin
static lorentzVector applyBias(const lorentzVector &muon, const int charge)
static double GLZNorm[40][1001]
T sqrt(T t)
Definition: SSEVec.h:48
static std::vector< double > parBgr
static int SmearType
static lorentzVector fromPtEtaPhiToPxPyPz(const double *ptEtaPhiE)
static double ResMaxSigma[6]
static double computeWeight(const double &mass, const int iev, const bool doUseBkgrWindow=false)
static std::vector< double > parScaleStep
static std::vector< std::pair< lorentzVector, lorentzVector > > SavedPair
static double invDimuonMass(const lorentzVector &mu1, const lorentzVector &mu2)
static lorentzVector applyScale(const lorentzVector &muon, const std::vector< double > &parval, const int charge)
static TH1D * signalProb_
static std::vector< int > parResolFix
static double GLZValue[40][1001][1001]
static TMinuit * rminPtr_
static std::vector< double > parSmear
static std::vector< int > parResolOrder
static bool sherpa_
static int BgrFitType
static TH1D * likelihoodInLoop_
static const int backgroundFunctionsRegions
static std::vector< double > parScaleMax
static std::vector< double > parScale
static double oldNormalization_
static bool duringMinos_
static double ResGamma[6]
static double deltaPhiMinCut_
double b
Definition: hdecay.h:120
static bool rapidityBinsForZ_
static int ResolFitType
static int goodmuon
static bool separateRanges_
static double minMuonEtaSecondRange_
static int iev_
static std::vector< double > parCrossSection
static const double mMu2
static double GLNorm[6][1001]
static std::vector< TGraphErrors * > fitReso(TH2F *histo)
tuple muons
Definition: patZpeak.py:38
static lorentzVector applySmearing(const lorentzVector &muon)
static bool normalizeLikelihoodByEventNumber_
double a
Definition: hdecay.h:121
static int MuonType
static double deltaPhi(const double &phi1, const double &phi2)
static double probability(const double &mass, const double &massResol, const double GLvalue[][1001][1001], const double GLnorm[][1001], const int iRes, const int iY)
Computes the probability given the mass, mass resolution and the arrays with the probabilities and th...
static double minMuonPt_
static scaleFunctionBase< double * > * scaleFunction
static bool useProbsFile_
static std::vector< int > resfind
static double ResHalfWidth[6]
virtual float pt() const GCC11_FINAL
transverse momentum
static int counter_resprob
static resolutionFunctionBase< std::vector< double > > * resolutionFunctionForVec
static std::vector< int > parorder
static double crossSection[6]
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
static int FitStrategy
static double maxMuonEtaFirstRange_
static std::pair< lorentzVector, lorentzVector > findBestRecoRes(const std::vector< reco::LeafCandidate > &muons)
bool operator()(const reco::Muon &a, const reco::Muon &b) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
static double deltaR(const double &eta1, const double &eta2, const double &phi1, const double &phi2)
static CrossSectionHandler * crossSectionHandler
static std::vector< double > parResolStep
static std::vector< int > parCrossSectionFix
static resolutionFunctionBase< double * > * resolutionFunction