CMS 3D CMS Logo

MuScleFit.cc
Go to the documentation of this file.
1 
5 // \class MuScleFit
6 // Fitter of momentum scale and resolution from resonance decays to muon track pairs
7 //
8 // \author R. Bellan, C.Mariotti, S.Bolognesi - INFN Torino / T.Dorigo, M.De Mattia - INFN Padova
9 // revised S. Casasso, E. Migliore - UniTo & INFN Torino
10 //
11 // Recent additions:
12 // - several parameters allow a more flexible use, tests, and control handles
13 // for the likelihood. In particular, a set of integers controls the order
14 // with which parameters are released; another controls which parameters are
15 // fixed. A function allows to smear momenta and angles of muons from the
16 // resonance before any correction, using a set of random numbers generated
17 // once and for all in the constructor (this way the smearing remains the same
18 // for any given muon no matter how many times one loops and what corrections
19 // one applies).
20 // For a correct use of these flags, please see the function minimizeLikelihood() in
21 // MuScleFitUtils.cc
22 // - the fit now allows to extract resolution functions simultaneously with the
23 // momentum scale. So far a simple parametrization
24 // of muon momentum resolution and angle resolution has been implemented, but
25 // extensions are straightforward.
26 // - It is however advisable to fit separately resolution and scale. The suggested
27 // course of action is:
28 // 1) fit the scale with a simple parametrization
29 // 2) check results, fit with more complicated forms
30 // 3) verify which is a sufficiently accurate description of the data
31 // 4) fix scale parameters and fit for the resolution
32 // 5) go back to fitting the scale with resolution parameters fixed to fitted values
33 // - Also note that resolution fits may fail to converge due to instability of the
34 // probability distribution to the limit of large widths. Improvements here are
35 // advisable.
36 // - The treatment of signal windows in the Y region
37 // has to be refined because of overlaps. More work is needed here, assigning a different
38 // weight to different hypothesis of the resonance producing a given mass, if there are
39 // multiple candidates.
40 // - Also, larger windows are to be allowed for fits to SA muons.
41 // - File Probs_1000.root contains the probability distribution of lorentzians convoluted
42 // with gaussian smearing functions, for the six resonances. A 1000x1000 grid
43 // in mass,sigma has been computed (using root macro Probs.C).
44 // A wider interval of masses for each resonance should be computed, to be used for standalone muons
45 //
46 //
47 // Notes on additions, TD 31/3/08
48 //
49 // - background model: at least a couple of different models, with two parameters,
50 // should be included in the fitting procedure such that the function massprob(),
51 // which produces the probability in the likelihood computation, incorporates the
52 // probability that the event is from background. That is, if the fitting function
53 // knows the shape of the mass spectrum ( say, a falling exponential plus a gaussian
54 // signal) it becomes possible to fit the scale together with the background shape
55 // and normalization parameters. Of course, one should do one thing at a time: first
56 // a scale fit, then a shape fit with scale parameters fixed, and then a combined
57 // fit. Convergence issues should be handled case by case.
58 // - The correct implementation of the above idea requires a reorganization of pass
59 // parameters (in the cfg) and fit parameters. The user has to be able to smear,
60 // bias, fix parameters, choose scale fitting functions, resolution fitting functions,
61 // and background functions. It should be possible to separate the fit functions from
62 // the biasing ones, which would allow a more thorough testing.
63 // - all the above can be obtained by making the .cfg instructions heavier. Since this
64 // is a routine operated by experts only, it is a sensible choice.
65 // - One should thus envision the following:
66 // 1) a set of parameters controlling the biasing function: parBias()
67 // 2) a set of parameters controlling the smearing function: parSmear()
68 // 3) a set of parameters to define resolution modeling and initial values: parResol()
69 // 3b) parResol() gets fix and order bits by parResolFix() and parResolOrder()
70 // 4) a set of parameters to define scale modeling and initial values: parScale()
71 // 4b) parScale() gets fix and order bits by parScaleFix() and parScaleOrder()
72 // 5) a set of parameters controlling the background shape and normalization: parNorm()
73 // 5b) parNorm() gets fix and order bits by parNormFix() and parNormOrder()
74 // The likelihood parameters then become a vector which is dynamically composed of
75 // sets 3), 4), and 5): parval() = parResol()+parScale()+parNorm()
76 // - In order to study better the likelihood behavior it would be advisable to introduce
77 // some histogram filling on the last iteration of the likelihood. It is not clear
78 // how best to achieve that: probably the simplest way is to make a histogram filling
79 // function run just after the likelihood computation, such that the best value of the
80 // fit parameters is used.
81 // - The muon pair which we call our resonance must be chosen in a way which does not
82 // bias our likelihood: we cannot just choose the pair closest to a resonance.
83 //
84 //
85 // Notes on additions, T.Dorigo 22/12/2008
86 // ---------------------------------------
87 //
88 // - File Probs_new_1000_CTEQ.root now contains a set of 24 additional two-dim histograms,
89 // defining the probability distribution of Z boson decays as a function of measured mass
90 // and expected sigma in 24 different bins of Z rapidity, extracted from CTEQ 6 PDF (at
91 // Leading Order) from the convolution in the factorization integral. See programs CTEQ.cpp
92 // and Fits.C.
93 // - The probability for Z boson events now thus depends on the measured rapidity of the dimuon
94 // system. All functions in file MuScleFitUtils.cc have been suitably changed.
95 //
96 // ----------------------------------------------------------------------------------
97 // Modifications by M. De Mattia 13/3/2009
98 // ---------------------------------------
99 // - The histograms map was moved to a base class (MuScleFitBase) from which this one inherits.
100 //
101 // Modifications by M. De Mattia 20/7/2009
102 // ---------------------------------------
103 // - Reworked background fit based on ranges. See comments in the code for more details.
104 // ---------------------------------------------------------------------------------------------
105 // Base Class Headers
106 // ------------------
113 
114 #include <CLHEP/Vector/LorentzVector.h>
115 #include <vector>
116 
120 
121 #include "MuScleFitBase.h"
122 #include "Histograms.h"
123 #include "MuScleFitPlotter.h"
128 #include "MuScleFitMuonSelector.h"
129 
134 
137 
143 
145 
146 #include "HepPDT/defs.h"
147 #include "HepPDT/TableBuilder.hh"
148 #include "HepPDT/ParticleDataTable.hh"
149 
150 #include "HepMC/GenParticle.h"
151 #include "HepMC/GenEvent.h"
152 
157 
161 
162 #include "TFile.h"
163 #include "TTree.h"
164 #include "TMinuit.h"
165 
166 // To use callgrind for code profiling uncomment also the following define.
167 // #define USE_CALLGRIND
168 
169 #ifdef USE_CALLGRIND
170 #include "valgrind/callgrind.h"
171 #endif
172 
173 // To read likelihood distributions from the database.
174 //#include "CondFormats/RecoMuonObjects/interface/MuScleFitLikelihoodPdf.h"
175 //#include "CondFormats/DataRecord/interface/MuScleFitLikelihoodPdfRcd.h"
176 
177 namespace edm {
178  class ParameterSet;
179  class Event;
180  class EventSetup;
181 } // namespace edm
182 
184 public:
185  // Constructor
186  // -----------
188 
189  // Destructor
190  // ----------
191  ~MuScleFit() override;
192 
193  // Operations
194  // ----------
195  void beginOfJobInConstructor();
196  // void beginOfJob( const edm::EventSetup& eventSetup );
197  // virtual void beginOfJob();
198  void endOfJob() override;
199 
200  void startingNewLoop(unsigned int iLoop) override;
201 
202  edm::EDLooper::Status endOfLoop(const edm::EventSetup& eventSetup, unsigned int iLoop) override;
203  virtual void endOfFastLoop(const unsigned int iLoop);
204 
205  edm::EDLooper::Status duringLoop(const edm::Event& event, const edm::EventSetup& eventSetup) override;
210  virtual void duringFastLoop();
211 
212  template <typename T>
213  std::vector<MuScleFitMuon> fillMuonCollection(const std::vector<T>& tracks);
214 
215 private:
216 protected:
221  void selectMuons(const edm::Event& event);
227  void selectMuons(const int maxEvents, const TString& treeFileName);
228 
230  template <typename T>
231  void takeSelectedMuonType(const T& muon, std::vector<reco::Track>& tracks);
233  bool selGlobalMuon(const pat::Muon* aMuon);
234  bool selTrackerMuon(const pat::Muon* aMuon);
235 
237  bool checkDeltaR(reco::Particle::LorentzVector& genMu, reco::Particle::LorentzVector& recMu);
239  void fillComparisonHistograms(const reco::Particle::LorentzVector& genMu,
240  const reco::Particle::LorentzVector& recoMu,
241  const std::string& inputName,
242  const int charge);
243 
245  void applySmearing(reco::Particle::LorentzVector& mu);
247  void applyBias(reco::Particle::LorentzVector& mu, const int charge);
248 
253  void checkParameters();
254 
256 
257  // Counters
258  // --------
263 
264  bool ifHepMC;
265  bool ifGenPart;
266 
267  // Constants
268  // ---------
269  double minResMass_hwindow[6];
270  double maxResMass_hwindow[6];
271 
272  // Total number of loops
273  // ---------------------
274  unsigned int maxLoopNumber;
275  unsigned int loopCounter;
276 
277  bool fastLoop;
278 
280 
281  // The reconstructed muon 4-momenta to be put in the tree
282  // ------------------------------------------------------
285  int iev;
287 
290  bool PATmuons_;
292 
293  // Input Root Tree file name. If empty events are read from the edm root file.
295  // Output Root Tree file name. If not empty events are dumped to this file at the end of the last iteration.
297  // Maximum number of events from root tree. It works in the same way as the maxEvents to configure a input source.
299 
302  std::vector<std::string> triggerPath_;
305 
306  // input collections for PU related infos
309 
310  std::unique_ptr<MuScleFitMuonSelector> muonSelector_;
311 };
312 
313 template <typename T>
314 std::vector<MuScleFitMuon> MuScleFit::fillMuonCollection(const std::vector<T>& tracks) {
315  std::vector<MuScleFitMuon> muons;
316  typename std::vector<T>::const_iterator track;
317  for (track = tracks.begin(); track != tracks.end(); ++track) {
320  track->px(), track->py(), track->pz(), sqrt(track->p() * track->p() + MuScleFitUtils::mMu2));
321  // Apply smearing if needed, and then bias
322  // ---------------------------------------
324  if (debug_ > 0)
325  std::cout << std::setprecision(9) << "Muon #" << MuScleFitUtils::goodmuon << ": initial value Pt = " << mu.Pt()
326  << std::endl;
327 
328  applySmearing(mu);
329  applyBias(mu, track->charge());
330  if (debug_ > 0)
331  std::cout << "track charge: " << track->charge() << std::endl;
332 
333  Double_t hitsTk = track->innerTrack()->hitPattern().numberOfValidTrackerHits();
334  Double_t hitsMuon = track->innerTrack()->hitPattern().numberOfValidMuonHits();
335  Double_t ptError = track->innerTrack()->ptError();
336  MuScleFitMuon muon(mu, track->charge(), ptError, hitsTk, hitsMuon, false);
337  if (debug_ > 0) {
338  std::cout << "[MuScleFit::fillMuonCollection]" << std::endl;
339  std::cout << " muon = " << muon << std::endl;
340  }
341 
342  // Store modified muon
343  // -------------------
344  muons.push_back(muon);
345  }
346  return muons;
347 }
348 
349 template <typename T>
350 void MuScleFit::takeSelectedMuonType(const T& muon, std::vector<reco::Track>& tracks) {
351  // std::cout<<"muon "<<muon->isGlobalMuon()<<muon->isStandAloneMuon()<<muon->isTrackerMuon()<<std::endl;
352  //NNBB: one muon can be of many kinds at once but with the theMuonType_ we are sure
353  // to avoid double counting of the same muon
354  if (muon->isGlobalMuon() && theMuonType_ == 1)
355  tracks.push_back(*(muon->globalTrack()));
356  else if (muon->isStandAloneMuon() && theMuonType_ == 2)
357  tracks.push_back(*(muon->outerTrack()));
358  else if (muon->isTrackerMuon() && theMuonType_ == 3)
359  tracks.push_back(*(muon->innerTrack()));
360 
361  else if (theMuonType_ == 10 && !(muon->isStandAloneMuon())) //particular case!!
362  tracks.push_back(*(muon->innerTrack()));
363  else if (theMuonType_ == 11 && muon->isGlobalMuon())
364  tracks.push_back(*(muon->innerTrack()));
365  else if (theMuonType_ == 13 && muon->isTrackerMuon())
366  tracks.push_back(*(muon->innerTrack()));
367 }
368 
369 // Constructor
370 // -----------
371 MuScleFit::MuScleFit(const edm::ParameterSet& pset) : MuScleFitBase(pset), totalEvents_(0) {
373  if (debug_ > 0)
374  std::cout << "[MuScleFit]: Constructor" << std::endl;
375 
376  if ((theMuonType_ < -4 || theMuonType_ > 5) && theMuonType_ < 10) {
377  std::cout << "[MuScleFit]: Unknown muon type! Aborting." << std::endl;
378  abort();
379  }
380 
381  loopCounter = 0;
382 
383  // Boundaries for h-function computation (to be improved!)
384  // -------------------------------------------------------
385  minResMass_hwindow[0] = 71.1876; // 76.;
386  maxResMass_hwindow[0] = 111.188; // 106.;
387  minResMass_hwindow[1] = 10.15;
388  maxResMass_hwindow[1] = 10.55;
389  minResMass_hwindow[2] = 9.8;
390  maxResMass_hwindow[2] = 10.2;
391  minResMass_hwindow[3] = 9.25;
392  maxResMass_hwindow[3] = 9.65;
393  minResMass_hwindow[4] = 3.58;
394  maxResMass_hwindow[4] = 3.78;
395  minResMass_hwindow[5] = 3.0;
396  maxResMass_hwindow[5] = 3.2;
397 
398  // Max number of loops (if > 2 then try to minimize likelihood more than once)
399  // ---------------------------------------------------------------------------
400  maxLoopNumber = pset.getUntrackedParameter<int>("maxLoopNumber", 2);
401  fastLoop = pset.getUntrackedParameter<bool>("FastLoop", true);
402 
403  // Selection of fits according to loop
404  MuScleFitUtils::doResolFit = pset.getParameter<std::vector<int> >("doResolFit");
405  MuScleFitUtils::doScaleFit = pset.getParameter<std::vector<int> >("doScaleFit");
406  MuScleFitUtils::doCrossSectionFit = pset.getParameter<std::vector<int> >("doCrossSectionFit");
407  MuScleFitUtils::doBackgroundFit = pset.getParameter<std::vector<int> >("doBackgroundFit");
408 
409  // Bias and smear types
410  // --------------------
411  int biasType = pset.getParameter<int>("BiasType");
412  MuScleFitUtils::BiasType = biasType;
413  // No error, the scale functions are used also for the bias
415  int smearType = pset.getParameter<int>("SmearType");
416  MuScleFitUtils::SmearType = smearType;
418 
419  // Fit types
420  // ---------
421  int resolFitType = pset.getParameter<int>("ResolFitType");
422  MuScleFitUtils::ResolFitType = resolFitType;
425  int scaleType = pset.getParameter<int>("ScaleFitType");
426  MuScleFitUtils::ScaleFitType = scaleType;
429 
430  // Initial parameters values
431  // -------------------------
432  MuScleFitUtils::parBias = pset.getParameter<std::vector<double> >("parBias");
433  MuScleFitUtils::parSmear = pset.getParameter<std::vector<double> >("parSmear");
434  MuScleFitUtils::parResol = pset.getParameter<std::vector<double> >("parResol");
436  pset.getUntrackedParameter<std::vector<double> >("parResolStep", std::vector<double>());
437  MuScleFitUtils::parResolMin = pset.getUntrackedParameter<std::vector<double> >("parResolMin", std::vector<double>());
438  MuScleFitUtils::parResolMax = pset.getUntrackedParameter<std::vector<double> >("parResolMax", std::vector<double>());
439  MuScleFitUtils::parScale = pset.getParameter<std::vector<double> >("parScale");
441  pset.getUntrackedParameter<std::vector<double> >("parScaleStep", std::vector<double>());
442  MuScleFitUtils::parScaleMin = pset.getUntrackedParameter<std::vector<double> >("parScaleMin", std::vector<double>());
443  MuScleFitUtils::parScaleMax = pset.getUntrackedParameter<std::vector<double> >("parScaleMax", std::vector<double>());
444  MuScleFitUtils::parCrossSection = pset.getParameter<std::vector<double> >("parCrossSection");
445  MuScleFitUtils::parBgr = pset.getParameter<std::vector<double> >("parBgr");
446  MuScleFitUtils::parResolFix = pset.getParameter<std::vector<int> >("parResolFix");
447  MuScleFitUtils::parScaleFix = pset.getParameter<std::vector<int> >("parScaleFix");
448  MuScleFitUtils::parCrossSectionFix = pset.getParameter<std::vector<int> >("parCrossSectionFix");
449  MuScleFitUtils::parBgrFix = pset.getParameter<std::vector<int> >("parBgrFix");
450  MuScleFitUtils::parResolOrder = pset.getParameter<std::vector<int> >("parResolOrder");
451  MuScleFitUtils::parScaleOrder = pset.getParameter<std::vector<int> >("parScaleOrder");
452  MuScleFitUtils::parCrossSectionOrder = pset.getParameter<std::vector<int> >("parCrossSectionOrder");
453  MuScleFitUtils::parBgrOrder = pset.getParameter<std::vector<int> >("parBgrOrder");
454 
455  MuScleFitUtils::resfind = pset.getParameter<std::vector<int> >("resfind");
456  MuScleFitUtils::FitStrategy = pset.getParameter<int>("FitStrategy");
457 
458  // Option to skip unnecessary stuff
459  // --------------------------------
460  MuScleFitUtils::speedup = pset.getParameter<bool>("speedup");
461 
462  // Option to skip simTracks comparison
463  compareToSimTracks_ = pset.getParameter<bool>("compareToSimTracks");
464  simTracksCollection_ = pset.getUntrackedParameter<edm::InputTag>("SimTracksCollection", edm::InputTag("g4SimHits"));
465 
466  triggerResultsLabel_ = pset.getUntrackedParameter<std::string>("TriggerResultsLabel");
467  triggerResultsProcess_ = pset.getUntrackedParameter<std::string>("TriggerResultsProcess");
468  triggerPath_ = pset.getUntrackedParameter<std::vector<std::string> >("TriggerPath");
469  negateTrigger_ = pset.getUntrackedParameter<bool>("NegateTrigger", false);
470  saveAllToTree_ = pset.getUntrackedParameter<bool>("SaveAllToTree", false);
471 
472  // input collections for PU related infos
473  puInfoSrc_ = pset.getUntrackedParameter<edm::InputTag>("PileUpSummaryInfo");
474  vertexSrc_ = pset.getUntrackedParameter<edm::InputTag>("PrimaryVertexCollection");
475 
476  PATmuons_ = pset.getUntrackedParameter<bool>("PATmuons", false);
477  genParticlesName_ = pset.getUntrackedParameter<std::string>("GenParticlesName", "genParticles");
478 
479  // Use the probability file or not. If not it will perform a simpler selection taking the muon pair with
480  // invariant mass closer to the pdf value and will crash if some fit is attempted.
481  MuScleFitUtils::useProbsFile_ = pset.getUntrackedParameter<bool>("UseProbsFile", true);
482 
483  // This must be set to true if using events generated with Sherpa
484  MuScleFitUtils::sherpa_ = pset.getUntrackedParameter<bool>("Sherpa", false);
485 
486  MuScleFitUtils::rapidityBinsForZ_ = pset.getUntrackedParameter<bool>("RapidityBinsForZ", true);
487 
488  // Set the cuts on muons to be used in the fit
489  MuScleFitUtils::separateRanges_ = pset.getUntrackedParameter<bool>("SeparateRanges", true);
490  MuScleFitUtils::maxMuonPt_ = pset.getUntrackedParameter<double>("MaxMuonPt", 100000000.);
491  MuScleFitUtils::minMuonPt_ = pset.getUntrackedParameter<double>("MinMuonPt", 0.);
492  MuScleFitUtils::minMuonEtaFirstRange_ = pset.getUntrackedParameter<double>("MinMuonEtaFirstRange", -6.);
493  MuScleFitUtils::maxMuonEtaFirstRange_ = pset.getUntrackedParameter<double>("MaxMuonEtaFirstRange", 6.);
494  MuScleFitUtils::minMuonEtaSecondRange_ = pset.getUntrackedParameter<double>("MinMuonEtaSecondRange", -100.);
495  MuScleFitUtils::maxMuonEtaSecondRange_ = pset.getUntrackedParameter<double>("MaxMuonEtaSecondRange", 100.);
496  MuScleFitUtils::deltaPhiMinCut_ = pset.getUntrackedParameter<double>("DeltaPhiMinCut", -100.);
497  MuScleFitUtils::deltaPhiMaxCut_ = pset.getUntrackedParameter<double>("DeltaPhiMaxCut", 100.);
498 
499  MuScleFitUtils::debugMassResol_ = pset.getUntrackedParameter<bool>("DebugMassResol", false);
500  // MuScleFitUtils::massResolComponentsStruct MuScleFitUtils::massResolComponents;
501 
502  // Check for parameters consistency
503  // it will abort in case of errors.
504  checkParameters();
505 
506  // Generate array of gaussian-distributed numbers for smearing
507  // -----------------------------------------------------------
508  if (MuScleFitUtils::SmearType > 0) {
509  std::cout << "[MuScleFit-Constructor]: Generating random values for smearing" << std::endl;
510  TF1 G("G", "[0]*exp(-0.5*pow(x,2))", -5., 5.);
511  double norm = 1 / sqrt(2 * TMath::Pi());
512  G.SetParameter(0, norm);
513  for (int i = 0; i < 10000; i++) {
514  for (int j = 0; j < 7; j++) {
515  MuScleFitUtils::x[j][i] = G.GetRandom();
516  }
517  }
518  }
520 
521  if (theMuonType_ > 0 && theMuonType_ < 4) {
524  } else if (theMuonType_ == 0 || theMuonType_ == 4 || theMuonType_ == 5 || theMuonType_ >= 10 || theMuonType_ == -1 ||
525  theMuonType_ == -2 || theMuonType_ == -3 || theMuonType_ == -4) {
528  } else {
529  std::cout << "Wrong muon type " << theMuonType_ << std::endl;
530  exit(1);
531  }
532 
533  // When using standalone muons switch to the single Z pdf
534  if (theMuonType_ == 2) {
535  MuScleFitUtils::rapidityBinsForZ_ = false;
536  }
537 
538  // Initialize ResMaxSigma And ResHalfWidth - 0 = global, 1 = SM, 2 = tracker
539  // -------------------------------------------------------------------------
558 
560  theMuonType_,
561  PATmuons_,
562  MuScleFitUtils::resfind,
563  MuScleFitUtils::speedup,
564  genParticlesName_,
565  compareToSimTracks_,
566  simTracksCollection_,
567  MuScleFitUtils::sherpa_,
568  debug_));
569 
571  new BackgroundHandler(pset.getParameter<std::vector<int> >("BgrFitType"),
572  pset.getParameter<std::vector<double> >("LeftWindowBorder"),
573  pset.getParameter<std::vector<double> >("RightWindowBorder"),
576 
578  new CrossSectionHandler(MuScleFitUtils::parCrossSection, MuScleFitUtils::resfind);
579 
580  // Build cross section scale factors
581  // MuScleFitUtils::resfind
582 
584  pset.getUntrackedParameter<bool>("NormalizeLikelihoodByEventNumber", true);
585  if (debug_ > 0)
586  std::cout << "End of MuScleFit constructor" << std::endl;
587 
588  inputRootTreeFileName_ = pset.getParameter<std::string>("InputRootTreeFileName");
589  outputRootTreeFileName_ = pset.getParameter<std::string>("OutputRootTreeFileName");
590  maxEventsFromRootTree_ = pset.getParameter<int>("MaxEventsFromRootTree");
591 
592  MuScleFitUtils::startWithSimplex_ = pset.getParameter<bool>("StartWithSimplex");
593  MuScleFitUtils::computeMinosErrors_ = pset.getParameter<bool>("ComputeMinosErrors");
594  MuScleFitUtils::minimumShapePlots_ = pset.getParameter<bool>("MinimumShapePlots");
595 
597 }
598 
599 // Destructor
600 // ----------
602  if (debug_ > 0)
603  std::cout << "[MuScleFit]: Destructor" << std::endl;
604  std::cout << "Total number of analyzed events = " << totalEvents_ << std::endl;
605 
606  if (!(outputRootTreeFileName_.empty())) {
607  // Save the events to a root tree unless we are reading from the edm root file and the SavedPair size is different from the totalEvents_
609  std::cout << "Saving muon pairs to root tree" << std::endl;
610  RootTreeHandler rootTreeHandler;
612  // rootTreeHandler.writeTree(outputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, 0, saveAllToTree_);
613  if (debug_ > 0) {
614  std::vector<MuonPair>::const_iterator it = muonPairs_.begin();
615  std::cout << "[MuScleFit::~MuScleFit] (Destructor)" << std::endl;
616  for (; it < muonPairs_.end(); ++it) {
617  std::cout << " Debugging pairs that are going to be written to file" << std::endl;
618  std::cout << " muon1 = " << it->mu1 << std::endl;
619  std::cout << " muon2 = " << it->mu2 << std::endl;
620  }
621  }
623  } else {
624  // rootTreeHandler.writeTree(outputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, &(MuScleFitUtils::genPair), saveAllToTree_ );
625  rootTreeHandler.writeTree(
627  }
628  } else {
629  std::cout << "ERROR: events in the vector = " << MuScleFitUtils::SavedPair.size()
630  << " != totalEvents = " << totalEvents_ << std::endl;
631  }
632  }
633 }
634 
635 // Begin job
636 // ---------
638 // void MuScleFit::beginOfJob ()
639 // void MuScleFit::beginOfJob( const edm::EventSetup& eventSetup )
640 {
641  if (debug_ > 0)
642  std::cout << "[MuScleFit]: beginOfJob" << std::endl;
643  //if(maxLoopNumber>1)
646  }
647 
648  if (debug_ > 0)
649  std::cout << "[MuScleFit]: beginOfJob" << std::endl;
650 
651  // Create the root file
652  // --------------------
653  for (unsigned int i = 0; i < (maxLoopNumber); i++) {
654  std::stringstream ss;
655  ss << i;
656  std::string rootFileName = ss.str() + "_" + theRootFileName_;
657  theFiles_.push_back(new TFile(rootFileName.c_str(), "RECREATE"));
658  }
659  if (debug_ > 0)
660  std::cout << "[MuScleFit]: Root file created" << std::endl;
661 
662  std::cout << "creating plotter" << std::endl;
664  plotter->debug = debug_;
665 }
666 
667 // End of job method
668 // -----------------
670  if (debug_ > 0)
671  std::cout << "[MuScleFit]: endOfJob" << std::endl;
672 }
673 
674 // New loop
675 // --------
676 void MuScleFit::startingNewLoop(unsigned int iLoop) {
677  if (debug_ > 0)
678  std::cout << "[MuScleFit]: Starting loop # " << iLoop << std::endl;
679 
680  // Number of muons used
681  // --------------------
683 
684  // Counters for problem std::cout-ing
685  // -----------------------------
687 
688  // Create the root file
689  // --------------------
690  fillHistoMap(theFiles_[iLoop], iLoop);
691 
692  loopCounter = iLoop;
694 
695  iev = 0;
697 
699 }
700 
701 // End of loop routine
702 // -------------------
703 edm::EDLooper::Status MuScleFit::endOfLoop(const edm::EventSetup& eventSetup, unsigned int iLoop) {
704  unsigned int iFastLoop = 1;
705 
706  // Read the events from the root tree if requested
707  if (!(inputRootTreeFileName_.empty())) {
709  // When reading from local file all the loops are done here
711  iFastLoop = 0;
712  } else {
713  endOfFastLoop(iLoop);
714  }
715 
716  // If a fastLoop is required we do all the remaining iterations here
717  if (fastLoop == true) {
718  for (; iFastLoop < maxLoopNumber; ++iFastLoop) {
719  std::cout << "Starting fast loop number " << iFastLoop << std::endl;
720 
721  // In the first loop is called by the framework
722  // if( iFastLoop > 0 ) {
723  startingNewLoop(iFastLoop);
724  // }
725 
726  // std::vector<std::pair<lorentzVector,lorentzVector> >::const_iterator it = MuScleFitUtils::SavedPair.begin();
727  // for( ; it != SavedPair.end(); ++it ) {
728  while (iev < totalEvents_) {
729  if (iev % 50000 == 0) {
730  std::cout << "Fast looping on event number " << iev << std::endl;
731  }
732  // This reads muons from SavedPair using iev to keep track of the event
733  duringFastLoop();
734  }
735  std::cout << "End of fast loop number " << iFastLoop << ". Ran on " << iev << " events" << std::endl;
736  endOfFastLoop(iFastLoop);
737  }
738  }
739 
740  if (iFastLoop >= maxLoopNumber - 1) {
741  return kStop;
742  } else {
743  return kContinue;
744  }
745 }
746 
747 void MuScleFit::endOfFastLoop(const unsigned int iLoop) {
748  // std::cout<< "Inside endOfFastLoop, iLoop = " << iLoop << " and loopCounter = " << loopCounter << std::endl;
749 
750  if (loopCounter == 0) {
751  // plotter->writeHistoMap();
752  // The destructor will call the writeHistoMap after the cd to the output file
753  delete plotter;
754  }
755 
756  std::cout << "Ending loop # " << iLoop << std::endl;
757 
758  // Write the histos to file
759  // ------------------------
760  // theFiles_[iLoop]->cd();
761  writeHistoMap(iLoop);
762 
763  // Likelihood minimization to compute corrections
764  // ----------------------------------------------
765  // theFiles_[iLoop]->cd();
766  TDirectory* likelihoodDir = theFiles_[iLoop]->mkdir("likelihood");
767  likelihoodDir->cd();
769 
770  // ATTENTION, this was put BEFORE the minimizeLikelihood. Check for problems.
771  theFiles_[iLoop]->Close();
772  // ATTENTION: Check that this delete does not give any problem
773  delete theFiles_[iLoop];
774 
775  // Clear the histos
776  // ----------------
777  clearHistoMap();
778 }
779 
780 // Stuff to do during loop
781 // -----------------------
784  event.getByLabel(edm::InputTag(triggerResultsLabel_.c_str(), "", triggerResultsProcess_.c_str()), triggerResults);
785  //event.getByLabel(InputTag(triggerResultsLabel_),triggerResults);
786  bool isFired = false;
787 
788  if (triggerPath_[0].empty())
789  isFired = true;
790  else if (triggerPath_[0] == "All") {
791  isFired = triggerResults->accept();
792  if (debug_ > 0)
793  std::cout << "Trigger " << isFired << std::endl;
794  } else {
795  bool changed;
797  hltConfig.init(event.getRun(), eventSetup, triggerResultsProcess_, changed);
798 
799  const edm::TriggerNames& triggerNames = event.triggerNames(*triggerResults);
800 
801  for (unsigned i = 0; i < triggerNames.size(); i++) {
802  const std::string& hltName = triggerNames.triggerName(i);
803 
804  // match the path in the pset with the true name of the trigger
805  for (unsigned int ipath = 0; ipath < triggerPath_.size(); ipath++) {
806  if (hltName.find(triggerPath_[ipath]) != std::string::npos) {
807  unsigned int triggerIndex(hltConfig.triggerIndex(hltName));
808 
809  // triggerIndex must be less than the size of HLTR or you get a CMSException: _M_range_check
810  if (triggerIndex < triggerResults->size()) {
811  isFired = triggerResults->accept(triggerIndex);
812  if (debug_ > 0)
813  std::cout << triggerPath_[ipath] << " " << hltName << " " << isFired << std::endl;
814  }
815  } // end if (matching the path in the pset with the true trigger name
816  }
817  }
818  }
819 
820  if (negateTrigger_ && isFired)
821  return kContinue;
822  else if (!(negateTrigger_) && !isFired)
823  return kContinue;
824 
825 #ifdef USE_CALLGRIND
826  CALLGRIND_START_INSTRUMENTATION;
827 #endif
828 
829  if (debug_ > 0) {
830  std::cout << "[MuScleFit-duringLoop]: loopCounter = " << loopCounter << " Run: " << event.id().run()
831  << " Event: " << event.id().event() << std::endl;
832  }
833 
834  // On the first iteration we read the bank, otherwise we fetch the information from the muon tree
835  // ------------------------------------ Important Note --------------------------------------- //
836  // The fillMuonCollection method applies any smearing or bias to the muons, so we NEVER use
837  // unbiased muons.
838  // ----------------------------------------------------------------------------------------------
839  if (loopCounter == 0) {
840  if (!fastLoop || inputRootTreeFileName_.empty()) {
841  if (debug_ > 0)
842  std::cout << "Reading from edm event" << std::endl;
843  selectMuons(event);
844  duringFastLoop();
845  ++totalEvents_;
846  }
847  }
848 
849  return kContinue;
850 
851 #ifdef USE_CALLGRIND
852  CALLGRIND_STOP_INSTRUMENTATION;
853  CALLGRIND_DUMP_STATS;
854 #endif
855 }
856 
860 
861  std::vector<MuScleFitMuon> muons;
862  muonSelector_->selectMuons(event, muons, genMuonPairs_, MuScleFitUtils::simPair, plotter);
863  // plotter->fillRec(muons); // @EM method already invoked inside MuScleFitMuonSelector::selectMuons()
864 
865  if (debug_ > 0) {
866  std::cout << "[MuScleFit::selectMuons] Debugging muons collections after call to muonSelector_->selectMuons"
867  << std::endl;
868  int iMu = 0;
869  for (std::vector<MuScleFitMuon>::const_iterator it = muons.begin(); it < muons.end(); ++it) {
870  std::cout << " - muon n. " << iMu << " = " << (*it) << std::endl;
871  ++iMu;
872  }
873  }
874 
875  // Find the two muons from the resonance, and set ResFound bool
876  // ------------------------------------------------------------
877  std::pair<MuScleFitMuon, MuScleFitMuon> recMuFromBestRes = MuScleFitUtils::findBestRecoRes(muons);
878 
880  if (debug_ > 0) {
881  std::cout << std::setprecision(9) << "Pt after findbestrecores: " << (recMuFromBestRes.first).Pt() << " "
882  << (recMuFromBestRes.second).Pt() << std::endl;
883  std::cout << "recMu1 = " << recMu1 << std::endl;
884  std::cout << "recMu2 = " << recMu2 << std::endl;
885  }
886  recMu1 = recMuFromBestRes.first.p4();
887  recMu2 = recMuFromBestRes.second.p4();
888  recMuScleMu1 = recMuFromBestRes.first;
889  recMuScleMu2 = recMuFromBestRes.second;
890 
891  if (debug_ > 0) {
892  std::cout << "after recMu1 = " << recMu1 << std::endl;
893  std::cout << "after recMu2 = " << recMu2 << std::endl;
894  std::cout << "mu1.pt = " << recMu1.Pt() << std::endl;
895  std::cout << "mu2.pt = " << recMu2.Pt() << std::endl;
896  std::cout << "after recMuScleMu1 = " << recMuScleMu1 << std::endl;
897  std::cout << "after recMuScleMu2 = " << recMuScleMu2 << std::endl;
898  }
899  MuScleFitUtils::SavedPair.push_back(std::make_pair(recMu1, recMu2));
901  } else {
902  MuScleFitUtils::SavedPair.push_back(std::make_pair(lorentzVector(0., 0., 0., 0.), lorentzVector(0., 0., 0., 0.)));
904  }
905  // Save the events also in the external tree so that it can be saved late
906 
907  // Fetch extra information (per event)
908  UInt_t the_NVtx(0);
909  Int_t the_numPUvtx(0);
910  Float_t the_TrueNumInteractions(0);
911 
912  // Fill pile-up related informations
913  // --------------------------------
915  event.getByLabel(puInfoSrc_, puInfo);
916  if (puInfo.isValid()) {
917  std::vector<PileupSummaryInfo>::const_iterator PVI;
918  for (PVI = puInfo->begin(); PVI != puInfo->end(); ++PVI) {
919  int BX = PVI->getBunchCrossing();
920  if (BX == 0) { // "0" is the in-time crossing, negative values are the early crossings, positive are late
921  the_TrueNumInteractions = PVI->getTrueNumInteractions();
922  the_numPUvtx = PVI->getPU_NumInteractions();
923  }
924  }
925  }
926 
928  event.getByLabel(vertexSrc_, vertices);
929  if (vertices.isValid()) {
930  std::vector<reco::Vertex>::const_iterator itv;
931  // now, count vertices
932  for (itv = vertices->begin(); itv != vertices->end(); ++itv) {
933  // require that the vertex meets certain criteria
934  if (itv->ndof() < 5)
935  continue;
936  if (fabs(itv->z()) > 50.0)
937  continue;
938  if (fabs(itv->position().rho()) > 2.0)
939  continue;
940  ++the_NVtx;
941  }
942  }
943 
944  // get the MC event weight
946  event.getByLabel("generator", genEvtInfo);
947  double the_genEvtweight = 1.;
948  if (genEvtInfo.isValid()) {
949  the_genEvtweight = genEvtInfo->weight();
950  }
951 
952  muonPairs_.push_back(MuonPair(
956  event.run(), event.id().event(), the_genEvtweight, the_numPUvtx, the_TrueNumInteractions, the_NVtx)));
957  // Fill the internal genPair tree from the external one
958  if (MuScleFitUtils::speedup == false) {
959  MuScleFitUtils::genPair.push_back(std::make_pair(genMuonPairs_.back().mu1.p4(), genMuonPairs_.back().mu2.p4()));
960  }
961 }
962 
963 void MuScleFit::selectMuons(const int maxEvents, const TString& treeFileName) {
964  std::cout << "Reading muon pairs from Root Tree in " << treeFileName << std::endl;
965  RootTreeHandler rootTreeHandler;
966  std::vector<std::pair<unsigned int, unsigned long long> > evtRun;
968  rootTreeHandler.readTree(
970  } else {
971  rootTreeHandler.readTree(maxEvents,
974  theMuonType_,
975  &evtRun,
977  }
978  // Now loop on all the pairs and apply any smearing and bias if needed
979  std::vector<std::pair<unsigned int, unsigned long long> >::iterator evtRunIt = evtRun.begin();
980  std::vector<std::pair<MuScleFitMuon, MuScleFitMuon> >::iterator it = MuScleFitUtils::SavedPairMuScleFitMuons.begin();
981  std::vector<std::pair<MuScleFitMuon, MuScleFitMuon> >::iterator genIt;
982  if (MuScleFitUtils::speedup == false)
983  genIt = MuScleFitUtils::genMuscleFitPair.begin();
984  for (; it != MuScleFitUtils::SavedPairMuScleFitMuons.end(); ++it, ++evtRunIt) {
985  // Apply any cut if requested
986  // Note that cuts here are only applied to already selected muons. They should not be used unless
987  // you are sure that the difference is negligible (e.g. the number of events with > 2 muons is negligible).
988  double pt1 = it->first.pt();
989  //std::cout << "pt1 = " << pt1 << std::endl;
990  double pt2 = it->second.pt();
991  //std::cout << "pt2 = " << pt2 << std::endl;
992  double eta1 = it->first.eta();
993  //std::cout << "eta1 = " << eta1 << std::endl;
994  double eta2 = it->second.eta();
995  //std::cout << "eta2 = " << eta2 << std::endl;
996  // If they don't pass the cuts set to null vectors
997  bool dontPass = false;
998  bool eta1InFirstRange;
999  bool eta2InFirstRange;
1000  bool eta1InSecondRange;
1001  bool eta2InSecondRange;
1002 
1003  int ch1 = it->first.charge();
1004  int ch2 = it->second.charge();
1005 
1009  eta1InSecondRange =
1011  eta2InSecondRange =
1013 
1014  // This is my logic, which should be erroneous, but certainly simpler...
1017  ((eta1InFirstRange && eta2InSecondRange && ch1 >= ch2) ||
1018  (eta1InSecondRange && eta2InFirstRange && ch1 < ch2))))
1019  dontPass = true;
1020  } else {
1023  eta1InSecondRange =
1025  eta2InSecondRange =
1029  (((eta1InFirstRange && !eta2InFirstRange) && (eta2InSecondRange && !eta1InSecondRange) && ch1 >= ch2) ||
1030  ((eta2InFirstRange && !eta1InFirstRange) && (eta1InSecondRange && !eta2InSecondRange) && ch1 < ch2))))
1031  dontPass = true;
1032  }
1033 
1034  // Additional check on deltaPhi
1035  double deltaPhi = MuScleFitUtils::deltaPhi(it->first.phi(), it->second.phi());
1036  if ((deltaPhi <= MuScleFitUtils::deltaPhiMinCut_) || (deltaPhi >= MuScleFitUtils::deltaPhiMaxCut_))
1037  dontPass = true;
1038 
1039  lorentzVector vec1 = it->first.p4();
1040  lorentzVector vec2 = it->second.p4();
1041  if (ch1 >= ch2) {
1042  lorentzVector vectemp = vec1;
1043  vec1 = vec2;
1044  vec2 = vectemp;
1045  }
1046 
1047  if (!dontPass) {
1048  // First is always mu-, second mu+
1049  if ((MuScleFitUtils::SmearType != 0) || (MuScleFitUtils::BiasType != 0)) {
1050  applySmearing(vec1);
1051  applyBias(vec1, -1);
1052  applySmearing(vec2);
1053  applyBias(vec2, 1);
1054  }
1055 
1056  MuScleFitUtils::SavedPair.push_back(std::make_pair(vec1, vec2));
1057  }
1058 
1059  //FIXME: we loose the additional information besides the 4-momenta
1060  muonPairs_.push_back(MuonPair(
1061  MuScleFitMuon(vec1, -1),
1062  MuScleFitMuon(vec2, +1),
1064  (*evtRunIt).first, (*evtRunIt).second, 0, 0, 0, 0)) // FIXME: order of event and run number mixed up!
1065  );
1066 
1067  // Fill the internal genPair tree from the external one
1068  if (!MuScleFitUtils::speedup) {
1069  MuScleFitUtils::genPair.push_back(std::make_pair(genIt->first.p4(), genIt->second.p4()));
1070  genMuonPairs_.push_back(GenMuonPair(genIt->first.p4(), genIt->second.p4(), 0));
1071  ++genIt;
1072  }
1073  }
1075  if (!(MuScleFitUtils::speedup)) {
1077  }
1078 }
1079 
1081  // On loops>0 the two muons are directly obtained from the SavedMuon array
1082  // -----------------------------------------------------------------------
1083  MuScleFitUtils::ResFound = false;
1085  recMu2 = (MuScleFitUtils::SavedPair[iev].second);
1086 
1087  //std::cout << "iev = " << iev << ", recMu1 pt = " << recMu1.Pt() << ", recMu2 pt = " << recMu2.Pt() << std::endl;
1088 
1089  if (recMu1.Pt() > 0 && recMu2.Pt() > 0) {
1090  MuScleFitUtils::ResFound = true;
1091  if (debug_ > 0)
1092  std::cout << "Ev = " << iev << ": found muons in tree with Pt = " << recMu1.Pt() << " " << recMu2.Pt()
1093  << std::endl;
1094  }
1095 
1096  if (debug_ > 0)
1097  std::cout << "About to start lik par correction and histo filling; ResFound is " << MuScleFitUtils::ResFound
1098  << std::endl;
1099  // If resonance found, do the hard work
1100  // ------------------------------------
1102  // Find weight and reference mass for this muon pair
1103  // -------------------------------------------------
1104  // The last parameter = true means that we want to use always the background window to compute the weight,
1105  // otherwise the probability will be filled only for the resonance region.
1106  double weight = MuScleFitUtils::computeWeight((recMu1 + recMu2).mass(), iev, true);
1107  if (debug_ > 0) {
1108  std::cout << "Loop #" << loopCounter << "Event #" << iev << ": before correction Pt1 = " << recMu1.Pt()
1109  << " Pt2 = " << recMu2.Pt() << std::endl;
1110  }
1111  // For successive iterations, correct the muons only if the previous iteration was a scale fit.
1112  // --------------------------------------------------------------------------------------------
1113  if (loopCounter > 0) {
1117  }
1118  }
1119  if (debug_ > 0) {
1120  std::cout << "Loop #" << loopCounter << "Event #" << iev << ": after correction Pt1 = " << recMu1.Pt()
1121  << " Pt2 = " << recMu2.Pt() << std::endl;
1122  }
1123 
1125 
1126  //Fill histograms
1127  //------------------
1128 
1129  mapHisto_["hRecBestMu"]->Fill(recMu1, -1, weight);
1130  mapHisto_["hRecBestMuVSEta"]->Fill(recMu1);
1131  mapHisto_["hRecBestMu"]->Fill(recMu2, +1, weight);
1132  mapHisto_["hRecBestMuVSEta"]->Fill(recMu2);
1133  mapHisto_["hDeltaRecBestMu"]->Fill(recMu1, recMu2);
1134  // Reconstructed resonance
1135  mapHisto_["hRecBestRes"]->Fill(bestRecRes, +1, weight);
1136  mapHisto_["hRecBestResAllEvents"]->Fill(bestRecRes, +1, 1.);
1137  // // Fill histogram of Res mass vs muon variables
1138  // mapHisto_["hRecBestResVSMu"]->Fill (recMu1, bestRecRes, -1);
1139  // mapHisto_["hRecBestResVSMu"]->Fill (recMu2, bestRecRes, +1);
1140  // // Fill also the mass mu+/mu- comparisons
1141  // mapHisto_["hRecBestResVSMu"]->Fill(recMu1, recMu2, bestRecRes);
1142 
1143  mapHisto_["hRecBestResVSMu"]->Fill(recMu1, bestRecRes, -1, weight);
1144  mapHisto_["hRecBestResVSMu"]->Fill(recMu2, bestRecRes, +1, weight);
1145  // Fill also the mass mu+/mu- comparisons
1146  mapHisto_["hRecBestResVSMu"]->Fill(recMu1, recMu2, bestRecRes, weight);
1147 
1148  //-- rc 2010 filling histograms for mu+ /mu- ------
1149  // mapHisto_["hRecBestResVSMuMinus"]->Fill (recMu1, bestRecRes, -1);
1150  // mapHisto_["hRecBestResVSMuPlus"]->Fill (recMu2, bestRecRes, +1);
1151 
1152  //-- rc 2010 filling histograms MassVsMuEtaPhi------
1153  // mapHisto_["hRecBestResVSMuEtaPhi"]->Fill (recMu1, bestRecRes,-1);
1154  // mapHisto_["hRecBestResVSMuEtaPhi"]->Fill (recMu2, bestRecRes,+1);
1155 
1156  // Fill histogram of Res mass vs Res variables
1157  // mapHisto_["hRecBestResVSRes"]->Fill (bestRecRes, bestRecRes, +1);
1158  mapHisto_["hRecBestResVSRes"]->Fill(bestRecRes, bestRecRes, +1, weight);
1159 
1160  std::vector<double>* parval;
1161  std::vector<double> initpar;
1162  // Store a pointer to the vector of parameters of the last iteration, or the initial
1163  // parameters if this is the first iteration
1164  if (loopCounter == 0) {
1165  initpar = MuScleFitUtils::parResol;
1166  initpar.insert(initpar.end(), MuScleFitUtils::parScale.begin(), MuScleFitUtils::parScale.end());
1167  initpar.insert(initpar.end(), MuScleFitUtils::parCrossSection.begin(), MuScleFitUtils::parCrossSection.end());
1168  initpar.insert(initpar.end(), MuScleFitUtils::parBgr.begin(), MuScleFitUtils::parBgr.end());
1169  parval = &initpar;
1170  } else {
1171  parval = &(MuScleFitUtils::parvalue[loopCounter - 1]);
1172  }
1173 
1174  //Compute pt resolution w.r.t generated and simulated muons
1175  //--------------------------------------------------------
1176  if (!MuScleFitUtils::speedup) {
1177  //first is always mu-, second is always mu+
1180  }
1183  }
1184  if (compareToSimTracks_) {
1185  //first is always mu-, second is always mu+
1188  }
1191  }
1192  }
1193  }
1194 
1195  // ATTENTION: this was done only when a matching was found. Moved it outside because, genInfo or not, we still want to see the resolution function
1196  // Fill also the resolution histogramsm using the resolution functions:
1197  // the parameters are those from the last iteration, as the muons up to this point have also the corrections from the same iteration.
1198  // Need to use a different array (ForVec), containing functors able to operate on std::vector<double>
1199  mapHisto_["hFunctionResolPt"]->Fill(
1200  recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu1.Pt(), recMu1.Eta(), *parval), -1);
1201  mapHisto_["hFunctionResolCotgTheta"]->Fill(
1202  recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu1.Pt(), recMu1.Eta(), *parval), -1);
1203  mapHisto_["hFunctionResolPhi"]->Fill(
1204  recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu1.Pt(), recMu1.Eta(), *parval), -1);
1205  mapHisto_["hFunctionResolPt"]->Fill(
1206  recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu2.Pt(), recMu2.Eta(), *parval), +1);
1207  mapHisto_["hFunctionResolCotgTheta"]->Fill(
1208  recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu2.Pt(), recMu2.Eta(), *parval), +1);
1209  mapHisto_["hFunctionResolPhi"]->Fill(
1210  recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu2.Pt(), recMu2.Eta(), *parval), +1);
1211 
1212  // Compute likelihood histograms
1213  // -----------------------------
1214  if (debug_ > 0)
1215  std::cout << "mass = " << bestRecRes.mass() << std::endl;
1216  if (weight != 0.) {
1217  double massResol;
1218  double prob;
1219  double deltalike;
1220  if (loopCounter == 0) {
1221  std::vector<double> initpar;
1222  for (int i = 0; i < (int)(MuScleFitUtils::parResol.size()); i++) {
1223  initpar.push_back(MuScleFitUtils::parResol[i]);
1224  }
1225  for (int i = 0; i < (int)(MuScleFitUtils::parScale.size()); i++) {
1226  initpar.push_back(MuScleFitUtils::parScale[i]);
1227  }
1228  // for (int i=0; i<(int)(MuScleFitUtils::parCrossSection.size()); i++) {
1229  // initpar.push_back(MuScleFitUtils::parCrossSection[i]);
1230  // }
1232 
1233  for (int i = 0; i < (int)(MuScleFitUtils::parBgr.size()); i++) {
1234  initpar.push_back(MuScleFitUtils::parBgr[i]);
1235  }
1236  massResol = MuScleFitUtils::massResolution(recMu1, recMu2, initpar);
1237  // prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(), massResol, initpar, true );
1238  prob = MuScleFitUtils::massProb(bestRecRes.mass(),
1239  bestRecRes.Eta(),
1240  bestRecRes.Rapidity(),
1241  massResol,
1242  initpar,
1243  true,
1244  recMu1.eta(),
1245  recMu2.eta());
1246  } else {
1248  // prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(),
1249  // massResol, MuScleFitUtils::parvalue[loopCounter-1], true );
1250  prob = MuScleFitUtils::massProb(bestRecRes.mass(),
1251  bestRecRes.Eta(),
1252  bestRecRes.Rapidity(),
1253  massResol,
1255  true,
1256  recMu1.eta(),
1257  recMu2.eta());
1258  }
1259  if (debug_ > 0)
1260  std::cout << "inside weight: mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
1261  if (prob > 0) {
1262  if (debug_ > 0)
1263  std::cout << "inside prob: mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
1264 
1265  deltalike = log(prob) * weight; // NB maximum likelihood --> deltalike is maximized
1266  mapHisto_["hLikeVSMu"]->Fill(recMu1, deltalike);
1267  mapHisto_["hLikeVSMu"]->Fill(recMu2, deltalike);
1268  mapHisto_["hLikeVSMuMinus"]->Fill(recMu1, deltalike);
1269  mapHisto_["hLikeVSMuPlus"]->Fill(recMu2, deltalike);
1270 
1271  double recoMass = (recMu1 + recMu2).mass();
1272  if (recoMass != 0) {
1273  // IMPORTANT: massResol is not a relative resolution
1274  mapHisto_["hResolMassVSMu"]->Fill(recMu1, massResol, -1);
1275  mapHisto_["hResolMassVSMu"]->Fill(recMu2, massResol, +1);
1276  mapHisto_["hFunctionResolMassVSMu"]->Fill(recMu1, massResol / recoMass, -1);
1277  mapHisto_["hFunctionResolMassVSMu"]->Fill(recMu2, massResol / recoMass, +1);
1278  }
1279 
1281  mapHisto_["hdMdPt1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdpt1, -1);
1282  mapHisto_["hdMdPt2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdpt2, +1);
1283  mapHisto_["hdMdPhi1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdphi1, -1);
1284  mapHisto_["hdMdPhi2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdphi2, +1);
1285  mapHisto_["hdMdCotgTh1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdcotgth1, -1);
1286  mapHisto_["hdMdCotgTh2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdcotgth2, +1);
1287  }
1288 
1289  if (!MuScleFitUtils::speedup) {
1290  double genMass = (MuScleFitUtils::genPair[iev].first + MuScleFitUtils::genPair[iev].second).mass();
1291  // Fill the mass resolution (computed from MC), we use the covariance class to compute the variance
1292  if (genMass != 0) {
1293  mapHisto_["hGenResVSMu"]->Fill((MuScleFitUtils::genPair[iev].first),
1295  -1);
1296  mapHisto_["hGenResVSMu"]->Fill((MuScleFitUtils::genPair[iev].second),
1298  +1);
1299  double diffMass = (recoMass - genMass) / genMass;
1300  // double diffMass = recoMass - genMass;
1301  // Fill if for both muons
1302  double pt1 = recMu1.pt();
1303  double eta1 = recMu1.eta();
1304  double pt2 = recMu2.pt();
1305  double eta2 = recMu2.eta();
1306  // This is to avoid nan
1307  if (diffMass == diffMass) {
1308  // Mass relative difference vs Pt and Eta. To be used to extract the true mass resolution
1309  mapHisto_["hDeltaMassOverGenMassVsPt"]->Fill(pt1, diffMass);
1310  mapHisto_["hDeltaMassOverGenMassVsPt"]->Fill(pt2, diffMass);
1311  mapHisto_["hDeltaMassOverGenMassVsEta"]->Fill(eta1, diffMass);
1312  mapHisto_["hDeltaMassOverGenMassVsEta"]->Fill(eta2, diffMass);
1313  // This is used for the covariance comparison
1314  mapHisto_["hMassResolutionVsPtEta"]->Fill(pt1, eta1, diffMass, diffMass);
1315  mapHisto_["hMassResolutionVsPtEta"]->Fill(pt2, eta2, diffMass, diffMass);
1316  } else {
1317  std::cout << "Error, there is a nan: recoMass = " << recoMass << ", genMass = " << genMass << std::endl;
1318  }
1319  }
1320  // Fill with mass resolution from resolution function
1322  mapHisto_["hFunctionResolMass"]->Fill(recMu1, std::pow(massRes, 2), -1);
1323  mapHisto_["hFunctionResolMass"]->Fill(recMu2, std::pow(massRes, 2), +1);
1324  }
1325 
1326  mapHisto_["hMass_P"]->Fill(bestRecRes.mass(), prob);
1327  if (debug_ > 0)
1328  std::cout << "mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
1329  mapHisto_["hMass_fine_P"]->Fill(bestRecRes.mass(), prob);
1330 
1331  mapHisto_["hMassProbVsRes"]->Fill(bestRecRes, bestRecRes, +1, prob);
1332  mapHisto_["hMassProbVsMu"]->Fill(recMu1, bestRecRes, -1, prob);
1333  mapHisto_["hMassProbVsMu"]->Fill(recMu2, bestRecRes, +1, prob);
1334  mapHisto_["hMassProbVsRes_fine"]->Fill(bestRecRes, bestRecRes, +1, prob);
1335  mapHisto_["hMassProbVsMu_fine"]->Fill(recMu1, bestRecRes, -1, prob);
1336  mapHisto_["hMassProbVsMu_fine"]->Fill(recMu2, bestRecRes, +1, prob);
1337  }
1338  }
1339  } // end if ResFound
1340 
1341  // Fill the pair
1342  // -------------
1343  if (loopCounter > 0) {
1344  if (debug_ > 0)
1345  std::cout << "[MuScleFit]: filling the pair" << std::endl;
1346  MuScleFitUtils::SavedPair[iev] = std::make_pair(recMu1, recMu2);
1347  }
1348 
1349  iev++;
1351 
1352  // return kContinue;
1353 }
1354 
1356  //first is always mu-, second is always mu+
1357  double deltaR =
1358  sqrt(MuScleFitUtils::deltaPhi(recMu.Phi(), genMu.Phi()) * MuScleFitUtils::deltaPhi(recMu.Phi(), genMu.Phi()) +
1359  ((recMu.Eta() - genMu.Eta()) * (recMu.Eta() - genMu.Eta())));
1360  if (deltaR < 0.01)
1361  return true;
1362  else if (debug_ > 0) {
1363  std::cout << "Reco muon " << recMu << " with eta " << recMu.Eta() << " and phi " << recMu.Phi() << std::endl
1364  << " DOES NOT MATCH with generated muon from resonance: " << std::endl
1365  << genMu << " with eta " << genMu.Eta() << " and phi " << genMu.Phi() << std::endl;
1366  }
1367  return false;
1368 }
1369 
1371  const reco::Particle::LorentzVector& recMu,
1372  const std::string& inputName,
1373  const int charge) {
1374  std::string name(inputName + "VSMu");
1375  mapHisto_["hResolPt" + name]->Fill(recMu, (-genMu.Pt() + recMu.Pt()) / genMu.Pt(), charge);
1376  mapHisto_["hResolTheta" + name]->Fill(recMu, (-genMu.Theta() + recMu.Theta()), charge);
1377  mapHisto_["hResolCotgTheta" + name]->Fill(
1378  recMu, (-cos(genMu.Theta()) / sin(genMu.Theta()) + cos(recMu.Theta()) / sin(recMu.Theta())), charge);
1379  mapHisto_["hResolEta" + name]->Fill(recMu, (-genMu.Eta() + recMu.Eta()), charge);
1380  mapHisto_["hResolPhi" + name]->Fill(recMu, MuScleFitUtils::deltaPhiNoFabs(recMu.Phi(), genMu.Phi()), charge);
1381 
1382  // Fill only if it was matched to a genMu and this muon is valid
1383  if ((genMu.Pt() != 0) && (recMu.Pt() != 0)) {
1384  mapHisto_["hPtRecoVsPt" + inputName]->Fill(genMu.Pt(), recMu.Pt());
1385  }
1386 }
1387 
1389  if (MuScleFitUtils::SmearType > 0) {
1391  if (debug_ > 0)
1392  std::cout << "Muon #" << MuScleFitUtils::goodmuon << ": after smearing Pt = " << mu.Pt() << std::endl;
1393  }
1394 }
1395 
1397  if (MuScleFitUtils::BiasType > 0) {
1398  mu = MuScleFitUtils::applyBias(mu, charge);
1399  if (debug_ > 0)
1400  std::cout << "Muon #" << MuScleFitUtils::goodmuon << ": after bias Pt = " << mu.Pt() << std::endl;
1401  }
1402 }
1403 
1404 // Simple method to check parameters consistency. It aborts the job if the parameters
1405 // are not consistent.
1407  // Fits selection dimension check
1409  std::cout << "[MuScleFit-Constructor]: wrong size of resolution fits selector = "
1410  << MuScleFitUtils::doResolFit.size() << std::endl;
1411  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1412  abort();
1413  }
1415  std::cout << "[MuScleFit-Constructor]: wrong size of scale fits selector = " << MuScleFitUtils::doScaleFit.size()
1416  << std::endl;
1417  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1418  abort();
1419  }
1421  std::cout << "[MuScleFit-Constructor]: wrong size of cross section fits selector = "
1422  << MuScleFitUtils::doCrossSectionFit.size() << std::endl;
1423  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1424  abort();
1425  }
1427  std::cout << "[MuScleFit-Constructor]: wrong size of background fits selector = "
1428  << MuScleFitUtils::doBackgroundFit.size() << std::endl;
1429  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1430  abort();
1431  }
1432 
1433  // Bias parameters: dimension check
1434  // --------------------------------
1435  if ((MuScleFitUtils::BiasType == 1 && MuScleFitUtils::parBias.size() != 2) || // linear in pt
1436  (MuScleFitUtils::BiasType == 2 && MuScleFitUtils::parBias.size() != 2) || // linear in |eta|
1437  (MuScleFitUtils::BiasType == 3 && MuScleFitUtils::parBias.size() != 4) || // sinusoidal in phi
1438  (MuScleFitUtils::BiasType == 4 && MuScleFitUtils::parBias.size() != 3) || // linear in pt and |eta|
1439  (MuScleFitUtils::BiasType == 5 && MuScleFitUtils::parBias.size() != 3) || // linear in pt and sinusoidal in phi
1440  (MuScleFitUtils::BiasType == 6 &&
1441  MuScleFitUtils::parBias.size() != 3) || // linear in |eta| and sinusoidal in phi
1442  (MuScleFitUtils::BiasType == 7 &&
1443  MuScleFitUtils::parBias.size() != 4) || // linear in pt and |eta| and sinusoidal in phi
1444  (MuScleFitUtils::BiasType == 8 && MuScleFitUtils::parBias.size() != 4) || // linear in pt and parabolic in |eta|
1445  (MuScleFitUtils::BiasType == 9 && MuScleFitUtils::parBias.size() != 2) || // exponential in pt
1446  (MuScleFitUtils::BiasType == 10 && MuScleFitUtils::parBias.size() != 3) || // parabolic in pt
1447  (MuScleFitUtils::BiasType == 11 &&
1448  MuScleFitUtils::parBias.size() != 4) || // linear in pt and sin in phi with chg
1449  (MuScleFitUtils::BiasType == 12 &&
1450  MuScleFitUtils::parBias.size() != 6) || // linear in pt and para in plus sin in phi with chg
1451  (MuScleFitUtils::BiasType == 13 &&
1452  MuScleFitUtils::parBias.size() != 8) || // linear in pt and para in plus sin in phi with chg
1454  MuScleFitUtils::BiasType > 13) {
1455  std::cout << "[MuScleFit-Constructor]: Wrong bias type or number of parameters: aborting!" << std::endl;
1456  abort();
1457  }
1458  // Smear parameters: dimension check
1459  // ---------------------------------
1468  std::cout << "[MuScleFit-Constructor]: Wrong smear type or number of parameters: aborting!" << std::endl;
1469  abort();
1470  }
1471  // Protect against bad size of parameters
1472  // --------------------------------------
1475  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Resol: aborting!" << std::endl;
1476  abort();
1477  }
1480  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Scale: aborting!" << std::endl;
1481  abort();
1482  }
1485  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Bgr: aborting!" << std::endl;
1486  abort();
1487  }
1490  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Bgr: aborting!" << std::endl;
1491  abort();
1492  }
1493 
1494  // Protect against an incorrect number of resonances
1495  // -------------------------------------------------
1496  if (MuScleFitUtils::resfind.size() != 6) {
1497  std::cout << "[MuScleFit-Constructor]: resfind must have 6 elements (1 Z, 3 Y, 2 Psi): aborting!" << std::endl;
1498  abort();
1499  }
1500 }
1501 
1503  reco::TrackRef iTrack = aMuon->innerTrack();
1504  const reco::HitPattern& p = iTrack->hitPattern();
1505 
1506  reco::TrackRef gTrack = aMuon->globalTrack();
1507  const reco::HitPattern& q = gTrack->hitPattern();
1508 
1509  return ( //isMuonInAccept(aMuon) &&// no acceptance cuts!
1510  iTrack->found() > 11 && gTrack->chi2() / gTrack->ndof() < 20.0 && q.numberOfValidMuonHits() > 0 &&
1511  iTrack->chi2() / iTrack->ndof() < 4.0 && aMuon->muonID("TrackerMuonArbitrated") &&
1512  aMuon->muonID("TMLastStationAngTight") && p.pixelLayersWithMeasurement() > 1 &&
1513  fabs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
1514  fabs(iTrack->dz()) < 15.0); //should be done w.r.t. PV!
1515 }
1516 
1518  reco::TrackRef iTrack = aMuon->innerTrack();
1519  const reco::HitPattern& p = iTrack->hitPattern();
1520 
1521  return ( //isMuonInAccept(aMuon) // no acceptance cuts!
1522  iTrack->found() > 11 && iTrack->chi2() / iTrack->ndof() < 4.0 && aMuon->muonID("TrackerMuonArbitrated") &&
1523  aMuon->muonID("TMLastStationAngTight") && p.pixelLayersWithMeasurement() > 1 &&
1524  fabs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
1525  fabs(iTrack->dz()) < 15.0); //should be done w.r.t. PV!
1526 }
1527 
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...
size
Write out results.
static std::vector< std::pair< lorentzVector, lorentzVector > > simPair
const double Pi
static std::vector< int > doScaleFit
static std::vector< int > doResolFit
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
static std::vector< double > parBias
void selectMuons(const edm::Event &event)
Definition: MuScleFit.cc:857
static std::vector< int > parCrossSectionOrder
std::vector< GenMuonPair > genMuonPairs_
Stores the genMuon pairs and the motherId prior to the creation of the internal tree.
Definition: MuScleFitBase.h:82
static smearFunctionBase * smearFunction
static std::vector< int > parScaleOrder
edm::InputTag theMuonLabel_
Definition: MuScleFitBase.h:47
static std::vector< int > doCrossSectionFit
static void minimizeLikelihood()
static double maxMuonEtaSecondRange_
std::vector< MuScleFitMuon > fillMuonCollection(const std::vector< T > &tracks)
Definition: MuScleFit.cc:314
static std::vector< int > parBgrOrder
bool selGlobalMuon(const pat::Muon *aMuon)
Function for onia selections.
Definition: MuScleFit.cc:1502
static unsigned int loopCounter
static double deltaPhiMaxCut_
static std::vector< double > parResolMax
static std::vector< std::pair< MuScleFitMuon, MuScleFitMuon > > SavedPairMuScleFitMuons
unsigned int loopCounter
Definition: MuScleFit.cc:275
static std::vector< std::pair< MuScleFitMuon, MuScleFitMuon > > genMuscleFitPair
int totalEvents_
Definition: MuScleFit.cc:286
~MuScleFit() override
Definition: MuScleFit.cc:601
static bool startWithSimplex_
std::string triggerResultsProcess_
Definition: MuScleFit.cc:301
static std::vector< int > doBackgroundFit
static std::vector< double > parResol
static bool debugMassResol_
bool saveAllToTree_
Definition: MuScleFit.cc:304
virtual void duringFastLoop()
Definition: MuScleFit.cc:1080
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void checkParameters()
Definition: MuScleFit.cc:1406
bool accept() const
Has at least one path accepted the event?
static BackgroundHandler * backgroundHandler
static double x[7][10000]
static int debug
static double ResMass[6]
edm::InputTag vertexSrc_
Definition: MuScleFit.cc:308
static double massWindowHalfWidth[3][6]
static bool speedup
#define DEFINE_FWK_LOOPER(type)
bool muonID(const std::string &name) const
Definition: weight.py:1
Run const & getRun() const
Definition: Event.cc:108
std::map< std::string, Histograms * > mapHisto_
The map of histograms.
Definition: MuScleFitBase.h:77
static bool ResFound
void readTree(const int maxEvents, const TString &fileName, MuonPairVector *savedPair, const int muonType, std::vector< std::pair< unsigned int, unsigned long long > > *evtRun, MuonPairVector *genPair=nullptr)
void fillTreeRec(const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &savedPairs)
Used when running on the root tree containing preselected muon pairs.
static scaleFunctionBase< std::vector< double > > * biasFunction
bool selTrackerMuon(const pat::Muon *aMuon)
Definition: MuScleFit.cc:1517
bool ifHepMC
Definition: MuScleFit.cc:264
void startingNewLoop(unsigned int iLoop) override
Definition: MuScleFit.cc:676
Strings::size_type size() const
Definition: TriggerNames.cc:31
static std::vector< int > parBgrFix
static std::vector< double > parResolMin
static bool minimumShapePlots_
double weight() const
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:492
edm::InputTag simTracksCollection_
Definition: MuScleFit.cc:289
static int MuonTypeForCheckMassWindow
static double minMuonEtaFirstRange_
static double massProb(const double &mass, const double &rapidity, const int ires, const double &massResol)
void fillComparisonHistograms(const reco::Particle::LorentzVector &genMu, const reco::Particle::LorentzVector &recoMu, const std::string &inputName, const int charge)
Fill the reco vs gen and reco vs sim comparison histograms.
Definition: MuScleFit.cc:1370
static struct MuScleFitUtils::massResolComponentsStruct massResolComponents
static int BiasType
std::vector< std::string > triggerPath_
Definition: MuScleFit.cc:302
static std::vector< int > parScaleFix
MuScleFitMuon recMuScleMu1
Definition: MuScleFit.cc:284
static bool computeMinosErrors_
reco::Particle::LorentzVector lorentzVector
Definition: GenMuonPair.h:9
static scaleFunctionBase< std::vector< double > > * scaleFunctionForVec
static std::vector< std::vector< double > > parvalue
unsigned int triggerIndex(const std::string &triggerName) const
slot position of trigger path in trigger table (0 to size-1)
resolutionFunctionBase< double * > * resolutionFunctionService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier.
Definition: Functions.cc:70
U second(std::pair< T, U > const &p)
void writeTree(const TString &fileName, const std::vector< MuonPair > *savedPair, const int muonType=0, const std::vector< GenMuonPair > *genPair=nullptr, const bool saveAll=false)
static double maxMuonPt_
static std::pair< MuScleFitMuon, MuScleFitMuon > findBestRecoRes(const std::vector< MuScleFitMuon > &muons)
static int ScaleFitType
static std::vector< std::pair< lorentzVector, lorentzVector > > genPair
std::string outputRootTreeFileName_
Definition: MuScleFit.cc:296
std::string theGenInfoRootFileName_
Definition: MuScleFitBase.h:49
int numberOfSimVertices
Definition: MuScleFit.cc:261
reco::TrackRef innerTrack() const override
reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
Definition: Muon.h:72
resolutionFunctionBase< std::vector< double > > * resolutionFunctionVecService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier when receiving a std::...
Definition: Functions.cc:92
static double massResolution(const lorentzVector &mu1, const lorentzVector &mu2)
reco::Particle::LorentzVector recMu2
Definition: MuScleFit.cc:283
bool PATmuons_
Definition: MuScleFit.cc:290
static std::vector< double > parScaleMin
std::vector< double > vec1
Definition: HCALResponse.h:15
reco::Particle::LorentzVector recMu1
Definition: MuScleFit.cc:283
static lorentzVector applyBias(const lorentzVector &muon, const int charge)
T sqrt(T t)
Definition: SSEVec.h:19
static std::vector< double > parBgr
unsigned int maxLoopNumber
Definition: MuScleFit.cc:274
static int SmearType
MuonServiceProxy * theService
Definition: MuScleFit.cc:255
std::string genParticlesName_
Definition: MuScleFit.cc:291
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
scaleFunctionBase< double * > * scaleFunctionService(const int identifier)
Service to build the scale functor corresponding to the passed identifier.
Definition: Functions.cc:3
RunNumber_t run() const
Definition: Event.h:107
void clearHistoMap()
Clean the histograms map.
MuScleFitMuon recMuScleMu2
Definition: MuScleFit.cc:284
bool compareToSimTracks_
Definition: MuScleFit.cc:288
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
MuScleFit(const edm::ParameterSet &pset)
Definition: MuScleFit.cc:371
static std::string const triggerResults
Definition: EdmProvDump.cc:45
static lorentzVector applyScale(const lorentzVector &muon, const std::vector< double > &parval, const int charge)
static std::vector< int > parResolFix
bool isValid() const
Definition: HandleBase.h:70
static std::vector< double > parSmear
void beginOfJobInConstructor()
Definition: MuScleFit.cc:637
smearFunctionBase * smearFunctionService(const int identifier)
Service to build the smearing functor corresponding to the passed identifier.
Definition: Functions.cc:37
static std::vector< int > parResolOrder
static bool sherpa_
void endOfJob() override
Definition: MuScleFit.cc:669
bool fastLoop
Definition: MuScleFit.cc:277
edm::EDLooper::Status duringLoop(const edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: MuScleFit.cc:782
int numberOfSimMuons
Definition: MuScleFit.cc:260
static std::vector< double > parScaleMax
scaleFunctionBase< std::vector< double > > * scaleFunctionVecService(const int identifier)
Service to build the scale functor corresponding to the passed identifier when receiving a std::vecto...
Definition: Functions.cc:20
void writeHistoMap(const unsigned int iLoop)
Save the histograms map to file.
static std::vector< double > parScale
static double oldNormalization_
void applySmearing(reco::Particle::LorentzVector &mu)
Apply the smearing if needed using the function in MuScleFitUtils.
Definition: MuScleFit.cc:1388
int numberOfSimTracks
Definition: MuScleFit.cc:259
reco::TrackRef globalTrack() const override
reference to Track reconstructed in both tracked and muon detector (reimplemented from reco::Muon) ...
Definition: Muon.h:80
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:22
std::string theRootFileName_
Definition: MuScleFitBase.h:48
void fillTreeGen(const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &genPairs)
edm::InputTag puInfoSrc_
Definition: MuScleFit.cc:307
bool init(const edm::Run &iRun, const edm::EventSetup &iSetup, const std::string &processName, bool &changed)
d&#39;tor
static double deltaPhiMinCut_
std::string triggerResultsLabel_
Definition: MuScleFit.cc:300
static bool rapidityBinsForZ_
MuScleFitPlotter * plotter
Definition: MuScleFit.cc:279
void fillHistoMap(TFile *outputFile, unsigned int iLoop)
Create the histograms map.
Definition: MuScleFitBase.cc:9
std::vector< TFile * > theFiles_
The files were the histograms are saved.
Definition: MuScleFitBase.h:74
double minResMass_hwindow[6]
Definition: MuScleFit.cc:269
static int ResolFitType
static int goodmuon
static bool separateRanges_
static double minMuonEtaSecondRange_
static int iev_
static std::vector< double > parCrossSection
static const double mMu2
HLT enums.
static lorentzVector applySmearing(const lorentzVector &muon)
static bool normalizeLikelihoodByEventNumber_
std::string inputRootTreeFileName_
Definition: MuScleFit.cc:294
bool checkDeltaR(reco::Particle::LorentzVector &genMu, reco::Particle::LorentzVector &recMu)
Check if two lorentzVector are near in deltaR.
Definition: MuScleFit.cc:1355
int maxEventsFromRootTree_
Definition: MuScleFit.cc:298
double maxResMass_hwindow[6]
Definition: MuScleFit.cc:270
static int MuonType
std::vector< MuonPair > muonPairs_
Used to store the muon pairs plus run and event number prior to the creation of the internal tree...
Definition: MuScleFitBase.h:80
static double deltaPhi(const double &phi1, const double &phi2)
static double minMuonPt_
static scaleFunctionBase< double * > * scaleFunction
int numberOfValidMuonHits() const
Definition: HitPattern.h:793
static bool useProbsFile_
static std::vector< int > resfind
void readProbabilityDistributionsFromFile()
Read probability distributions from a local root file.
void takeSelectedMuonType(const T &muon, std::vector< reco::Track > &tracks)
Template method used to fill the track collection starting from reco::muons or pat::muons.
Definition: MuScleFit.cc:350
std::vector< vec1 > vec2
Definition: HCALResponse.h:16
long double T
static int counter_resprob
static resolutionFunctionBase< std::vector< double > > * resolutionFunctionForVec
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
static int FitStrategy
Analysis-level muon class.
Definition: Muon.h:51
static double maxMuonEtaFirstRange_
std::unique_ptr< MuScleFitMuonSelector > muonSelector_
Definition: MuScleFit.cc:310
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
void applyBias(reco::Particle::LorentzVector &mu, const int charge)
Apply the bias if needed using the function in MuScleFitUtils.
Definition: MuScleFit.cc:1396
bool negateTrigger_
Definition: MuScleFit.cc:303
bool ifGenPart
Definition: MuScleFit.cc:265
static CrossSectionHandler * crossSectionHandler
edm::EDLooper::Status endOfLoop(const edm::EventSetup &eventSetup, unsigned int iLoop) override
Definition: MuScleFit.cc:703
static std::vector< double > parResolStep
Definition: event.py:1
static std::vector< int > parCrossSectionFix
static resolutionFunctionBase< double * > * resolutionFunction
void addParameters(std::vector< double > &initpar)
Inputs the vars in a vector.
def exit(msg="")
virtual void endOfFastLoop(const unsigned int iLoop)
Definition: MuScleFit.cc:747
int numberOfEwkZ
Definition: MuScleFit.cc:262