CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
167 // To use callgrind for code profiling uncomment also the following define.
168 // #define USE_CALLGRIND
169 
170 #ifdef USE_CALLGRIND
171 #include "valgrind/callgrind.h"
172 #endif
173 
174 
175 // To read likelihood distributions from the database.
176  //#include "CondFormats/RecoMuonObjects/interface/MuScleFitLikelihoodPdf.h"
177  //#include "CondFormats/DataRecord/interface/MuScleFitLikelihoodPdfRcd.h"
178 
179 namespace edm {
180  class ParameterSet;
181  class Event;
182  class EventSetup;
183 }
184 
185 
187 {
188  public:
189  // Constructor
190  // -----------
191  MuScleFit( const edm::ParameterSet& pset );
192 
193  // Destructor
194  // ----------
195  virtual ~MuScleFit();
196 
197  // Operations
198  // ----------
200  // void beginOfJob( const edm::EventSetup& eventSetup );
201  // virtual void beginOfJob();
202  virtual void endOfJob() override;
203 
204  virtual void startingNewLoop( unsigned int iLoop ) override;
205 
206  virtual edm::EDLooper::Status endOfLoop( const edm::EventSetup& eventSetup, unsigned int iLoop ) override;
207  virtual void endOfFastLoop( const unsigned int iLoop );
208 
209  virtual edm::EDLooper::Status duringLoop( const edm::Event & event, const edm::EventSetup& eventSetup ) override;
214  virtual void duringFastLoop();
215 
216  template<typename T>
217  std::vector<MuScleFitMuon> fillMuonCollection( const std::vector<T>& tracks );
218 
219  private:
220 
221  protected:
226  void selectMuons(const edm::Event & event);
232  void selectMuons(const int maxEvents, const TString & treeFileName);
233 
235  template<typename T>
236  void takeSelectedMuonType(const T & muon, std::vector<reco::Track> & tracks);
238  bool selGlobalMuon(const pat::Muon* aMuon);
239  bool selTrackerMuon(const pat::Muon* aMuon);
240 
244  void fillComparisonHistograms( const reco::Particle::LorentzVector & genMu, const reco::Particle::LorentzVector & recoMu, const std::string & inputName, const int charge );
245 
249  void applyBias( reco::Particle::LorentzVector & mu, const int charge );
250 
255  void checkParameters();
256 
258 
259  // Counters
260  // --------
265 
266  bool ifHepMC;
267  bool ifGenPart;
268 
269  // Constants
270  // ---------
273 
274  // Total number of loops
275  // ---------------------
276  unsigned int maxLoopNumber;
277  unsigned int loopCounter;
278 
279  bool fastLoop;
280 
282 
283  // The reconstructed muon 4-momenta to be put in the tree
284  // ------------------------------------------------------
287  int iev;
289 
292  bool PATmuons_;
294 
295  // Input Root Tree file name. If empty events are read from the edm root file.
297  // Output Root Tree file name. If not empty events are dumped to this file at the end of the last iteration.
299  // Maximum number of events from root tree. It works in the same way as the maxEvents to configure a input source.
301 
304  std::vector<std::string> triggerPath_;
307 
308  // input collections for PU related infos
311 
312  std::auto_ptr<MuScleFitMuonSelector> muonSelector_;
313 };
314 
315 template<typename T>
316 std::vector<MuScleFitMuon> MuScleFit::fillMuonCollection( const std::vector<T>& tracks )
317 {
318  std::vector<MuScleFitMuon> muons;
319  typename std::vector<T>::const_iterator track;
320  for( track = tracks.begin(); track != tracks.end(); ++track ) {
322  mu = reco::Particle::LorentzVector(track->px(),track->py(),track->pz(),
323  sqrt(track->p()*track->p() + MuScleFitUtils::mMu2));
324  // Apply smearing if needed, and then bias
325  // ---------------------------------------
327  if (debug_>0)
328  std::cout <<std::setprecision(9)<< "Muon #" << MuScleFitUtils::goodmuon
329  << ": initial value Pt = " << mu.Pt() << std::endl;
330 
331  applySmearing(mu);
332  applyBias(mu, track->charge());
333  if (debug_>0) std::cout<<"track charge: "<<track->charge()<<std::endl;
334 
335  Double_t hitsTk = track->innerTrack()->hitPattern().numberOfValidTrackerHits();
336  Double_t hitsMuon = track->innerTrack()->hitPattern().numberOfValidMuonHits();
337  Double_t ptError = track->innerTrack()->ptError();
338  MuScleFitMuon muon(mu,track->charge(),ptError,hitsTk,hitsMuon,false);
339  if (debug_>0) {
340  std::cout<<"[MuScleFit::fillMuonCollection]"<<std::endl;
341  std::cout<<" muon = "<<muon<<std::endl;
342  }
343 
344  // Store modified muon
345  // -------------------
346  muons.push_back (muon);
347  }
348  return muons;
349 }
350 
351 
352 template<typename T>
353 void MuScleFit::takeSelectedMuonType(const T & muon, std::vector<reco::Track> & tracks)
354 {
355  // std::cout<<"muon "<<muon->isGlobalMuon()<<muon->isStandAloneMuon()<<muon->isTrackerMuon()<<std::endl;
356  //NNBB: one muon can be of many kinds at once but with the theMuonType_ we are sure
357  // to avoid double counting of the same muon
358  if(muon->isGlobalMuon() && theMuonType_==1)
359  tracks.push_back(*(muon->globalTrack()));
360  else if(muon->isStandAloneMuon() && theMuonType_==2)
361  tracks.push_back(*(muon->outerTrack()));
362  else if(muon->isTrackerMuon() && theMuonType_==3)
363  tracks.push_back(*(muon->innerTrack()));
364 
365  else if( theMuonType_ == 10 && !(muon->isStandAloneMuon()) ) //particular case!!
366  tracks.push_back(*(muon->innerTrack()));
367  else if( theMuonType_ == 11 && muon->isGlobalMuon() )
368  tracks.push_back(*(muon->innerTrack()));
369  else if( theMuonType_ == 13 && muon->isTrackerMuon() )
370  tracks.push_back(*(muon->innerTrack()));
371 }
372 
373 
374 // Constructor
375 // -----------
377  MuScleFitBase( pset ),
378  totalEvents_(0)
379 {
381  if (debug_>0) std::cout << "[MuScleFit]: Constructor" << std::endl;
382 
383  if ((theMuonType_<-4 || theMuonType_>5) && theMuonType_<10) {
384  std::cout << "[MuScleFit]: Unknown muon type! Aborting." << std::endl;
385  abort();
386  }
387 
388  loopCounter = 0;
389 
390  // Boundaries for h-function computation (to be improved!)
391  // -------------------------------------------------------
392  minResMass_hwindow[0] = 71.1876; // 76.;
393  maxResMass_hwindow[0] = 111.188; // 106.;
394  minResMass_hwindow[1] = 10.15;
395  maxResMass_hwindow[1] = 10.55;
396  minResMass_hwindow[2] = 9.8;
397  maxResMass_hwindow[2] = 10.2;
398  minResMass_hwindow[3] = 9.25;
399  maxResMass_hwindow[3] = 9.65;
400  minResMass_hwindow[4] = 3.58;
401  maxResMass_hwindow[4] = 3.78;
402  minResMass_hwindow[5] = 3.0;
403  maxResMass_hwindow[5] = 3.2;
404 
405  // Max number of loops (if > 2 then try to minimize likelihood more than once)
406  // ---------------------------------------------------------------------------
407  maxLoopNumber = pset.getUntrackedParameter<int>("maxLoopNumber", 2);
408  fastLoop = pset.getUntrackedParameter<bool>("FastLoop", true);
409 
410  // Selection of fits according to loop
411  MuScleFitUtils::doResolFit = pset.getParameter<std::vector<int> >("doResolFit");
412  MuScleFitUtils::doScaleFit = pset.getParameter<std::vector<int> >("doScaleFit");
413  MuScleFitUtils::doCrossSectionFit = pset.getParameter<std::vector<int> >("doCrossSectionFit");
414  MuScleFitUtils::doBackgroundFit = pset.getParameter<std::vector<int> >("doBackgroundFit");
415 
416  // Bias and smear types
417  // --------------------
418  int biasType = pset.getParameter<int>("BiasType");
419  MuScleFitUtils::BiasType = biasType;
420  // No error, the scale functions are used also for the bias
422  int smearType = pset.getParameter<int>("SmearType");
423  MuScleFitUtils::SmearType = smearType;
425 
426  // Fit types
427  // ---------
428  int resolFitType = pset.getParameter<int>("ResolFitType");
429  MuScleFitUtils::ResolFitType = resolFitType;
432  int scaleType = pset.getParameter<int>("ScaleFitType");
433  MuScleFitUtils::ScaleFitType = scaleType;
436 
437  // Initial parameters values
438  // -------------------------
439  MuScleFitUtils::parBias = pset.getParameter<std::vector<double> >("parBias");
440  MuScleFitUtils::parSmear = pset.getParameter<std::vector<double> >("parSmear");
441  MuScleFitUtils::parResol = pset.getParameter<std::vector<double> >("parResol");
442  MuScleFitUtils::parResolStep = pset.getUntrackedParameter<std::vector<double> >("parResolStep", std::vector<double>());
443  MuScleFitUtils::parResolMin = pset.getUntrackedParameter<std::vector<double> >("parResolMin", std::vector<double>());
444  MuScleFitUtils::parResolMax = pset.getUntrackedParameter<std::vector<double> >("parResolMax", std::vector<double>());
445  MuScleFitUtils::parScale = pset.getParameter<std::vector<double> >("parScale");
446  MuScleFitUtils::parScaleStep = pset.getUntrackedParameter<std::vector<double> >("parScaleStep", std::vector<double>());
447  MuScleFitUtils::parScaleMin = pset.getUntrackedParameter<std::vector<double> >("parScaleMin", std::vector<double>());
448  MuScleFitUtils::parScaleMax = pset.getUntrackedParameter<std::vector<double> >("parScaleMax", std::vector<double>());
449  MuScleFitUtils::parCrossSection = pset.getParameter<std::vector<double> >("parCrossSection");
450  MuScleFitUtils::parBgr = pset.getParameter<std::vector<double> >("parBgr");
451  MuScleFitUtils::parResolFix = pset.getParameter<std::vector<int> >("parResolFix");
452  MuScleFitUtils::parScaleFix = pset.getParameter<std::vector<int> >("parScaleFix");
453  MuScleFitUtils::parCrossSectionFix = pset.getParameter<std::vector<int> >("parCrossSectionFix");
454  MuScleFitUtils::parBgrFix = pset.getParameter<std::vector<int> >("parBgrFix");
455  MuScleFitUtils::parResolOrder = pset.getParameter<std::vector<int> >("parResolOrder");
456  MuScleFitUtils::parScaleOrder = pset.getParameter<std::vector<int> >("parScaleOrder");
457  MuScleFitUtils::parCrossSectionOrder = pset.getParameter<std::vector<int> >("parCrossSectionOrder");
458  MuScleFitUtils::parBgrOrder = pset.getParameter<std::vector<int> >("parBgrOrder");
459 
460  MuScleFitUtils::resfind = pset.getParameter<std::vector<int> >("resfind");
461  MuScleFitUtils::FitStrategy = pset.getParameter<int>("FitStrategy");
462 
463  // Option to skip unnecessary stuff
464  // --------------------------------
465  MuScleFitUtils::speedup = pset.getParameter<bool>("speedup");
466 
467  // Option to skip simTracks comparison
468  compareToSimTracks_ = pset.getParameter<bool>("compareToSimTracks");
469  simTracksCollection_ = pset.getUntrackedParameter<edm::InputTag>("SimTracksCollection", edm::InputTag("g4SimHits"));
470 
471  triggerResultsLabel_ = pset.getUntrackedParameter<std::string>("TriggerResultsLabel");
472  triggerResultsProcess_ = pset.getUntrackedParameter<std::string>("TriggerResultsProcess");
473  triggerPath_ = pset.getUntrackedParameter<std::vector<std::string> >("TriggerPath");
474  negateTrigger_ = pset.getUntrackedParameter<bool>("NegateTrigger", false);
475  saveAllToTree_ = pset.getUntrackedParameter<bool>("SaveAllToTree", false);
476 
477  // input collections for PU related infos
478  puInfoSrc_ = pset.getUntrackedParameter<edm::InputTag>("PileUpSummaryInfo");
479  vertexSrc_ = pset.getUntrackedParameter<edm::InputTag>("PrimaryVertexCollection");
480 
481 
482  PATmuons_ = pset.getUntrackedParameter<bool>("PATmuons", false);
483  genParticlesName_ = pset.getUntrackedParameter<std::string>("GenParticlesName", "genParticles");
484 
485  // Use the probability file or not. If not it will perform a simpler selection taking the muon pair with
486  // invariant mass closer to the pdf value and will crash if some fit is attempted.
487  MuScleFitUtils::useProbsFile_ = pset.getUntrackedParameter<bool>("UseProbsFile", true);
488 
489  // This must be set to true if using events generated with Sherpa
490  MuScleFitUtils::sherpa_ = pset.getUntrackedParameter<bool>("Sherpa", false);
491 
492  MuScleFitUtils::rapidityBinsForZ_ = pset.getUntrackedParameter<bool>("RapidityBinsForZ", true);
493 
494  // Set the cuts on muons to be used in the fit
495  MuScleFitUtils::separateRanges_ = pset.getUntrackedParameter<bool>("SeparateRanges", true);
496  MuScleFitUtils::maxMuonPt_ = pset.getUntrackedParameter<double>("MaxMuonPt", 100000000.);
497  MuScleFitUtils::minMuonPt_ = pset.getUntrackedParameter<double>("MinMuonPt", 0.);
498  MuScleFitUtils::minMuonEtaFirstRange_ = pset.getUntrackedParameter<double>("MinMuonEtaFirstRange", -6.);
499  MuScleFitUtils::maxMuonEtaFirstRange_ = pset.getUntrackedParameter<double>("MaxMuonEtaFirstRange", 6.);
500  MuScleFitUtils::minMuonEtaSecondRange_ = pset.getUntrackedParameter<double>("MinMuonEtaSecondRange", -100.);
501  MuScleFitUtils::maxMuonEtaSecondRange_ = pset.getUntrackedParameter<double>("MaxMuonEtaSecondRange", 100.);
502  MuScleFitUtils::deltaPhiMinCut_ = pset.getUntrackedParameter<double>("DeltaPhiMinCut", -100.);
503  MuScleFitUtils::deltaPhiMaxCut_ = pset.getUntrackedParameter<double>("DeltaPhiMaxCut", 100.);
504 
505  MuScleFitUtils::debugMassResol_ = pset.getUntrackedParameter<bool>("DebugMassResol", false);
506  // MuScleFitUtils::massResolComponentsStruct MuScleFitUtils::massResolComponents;
507 
508  // Check for parameters consistency
509  // it will abort in case of errors.
510  checkParameters();
511 
512  // Generate array of gaussian-distributed numbers for smearing
513  // -----------------------------------------------------------
514  if (MuScleFitUtils::SmearType>0) {
515  std::cout << "[MuScleFit-Constructor]: Generating random values for smearing" << std::endl;
516  TF1 G("G", "[0]*exp(-0.5*pow(x,2))", -5., 5.);
517  double norm = 1/sqrt(2*TMath::Pi());
518  G.SetParameter (0,norm);
519  for (int i=0; i<10000; i++) {
520  for (int j=0; j<7; j++) {
521  MuScleFitUtils::x[j][i] = G.GetRandom();
522  }
523  }
524  }
526 
527  if(theMuonType_ > 0 && theMuonType_ < 4) {
530  }
531  else if(theMuonType_ == 0 || theMuonType_ == 4 || theMuonType_ == 5 || theMuonType_ >= 10 || theMuonType_==-1 || theMuonType_==-2 || theMuonType_==-3 || theMuonType_==-4) {
534  }
535  else{
536  std::cout<<"Wrong muon type "<<theMuonType_<<std::endl;
537  exit(1);
538  }
539 
540  // When using standalone muons switch to the single Z pdf
541  if( theMuonType_ == 2 ) {
542  MuScleFitUtils::rapidityBinsForZ_ = false;
543  }
544 
545  // Initialize ResMaxSigma And ResHalfWidth - 0 = global, 1 = SM, 2 = tracker
546  // -------------------------------------------------------------------------
565 
567  MuScleFitUtils::resfind,
568  MuScleFitUtils::speedup, genParticlesName_,
569  compareToSimTracks_, simTracksCollection_,
570  MuScleFitUtils::sherpa_, debug_));
571 
572  MuScleFitUtils::backgroundHandler = new BackgroundHandler( pset.getParameter<std::vector<int> >("BgrFitType"),
573  pset.getParameter<std::vector<double> >("LeftWindowBorder"),
574  pset.getParameter<std::vector<double> >("RightWindowBorder"),
577 
579 
580  // Build cross section scale factors
581  // MuScleFitUtils::resfind
582 
583  MuScleFitUtils::normalizeLikelihoodByEventNumber_ = pset.getUntrackedParameter<bool>("NormalizeLikelihoodByEventNumber", true);
584  if(debug_>0) std::cout << "End of MuScleFit constructor" << std::endl;
585 
586  inputRootTreeFileName_ = pset.getParameter<std::string>("InputRootTreeFileName");
587  outputRootTreeFileName_ = pset.getParameter<std::string>("OutputRootTreeFileName");
588  maxEventsFromRootTree_ = pset.getParameter<int>("MaxEventsFromRootTree");
589 
590  MuScleFitUtils::startWithSimplex_ = pset.getParameter<bool>("StartWithSimplex");
591  MuScleFitUtils::computeMinosErrors_ = pset.getParameter<bool>("ComputeMinosErrors");
592  MuScleFitUtils::minimumShapePlots_ = pset.getParameter<bool>("MinimumShapePlots");
593 
595 }
596 
597 // Destructor
598 // ----------
600  if (debug_>0) std::cout << "[MuScleFit]: Destructor" << std::endl;
601  std::cout << "Total number of analyzed events = " << totalEvents_ << std::endl;
602 
603  if( !(outputRootTreeFileName_.empty()) ) {
604  // 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_
605  if( !(inputRootTreeFileName_.empty() && (int(MuScleFitUtils::SavedPair.size()) != totalEvents_)) ) {
606  std::cout << "Saving muon pairs to root tree" << std::endl;
607  RootTreeHandler rootTreeHandler;
609  // rootTreeHandler.writeTree(outputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, 0, saveAllToTree_);
610  if (debug_>0) {
611  std::vector<MuonPair>::const_iterator it = muonPairs_.begin();
612  std::cout<<"[MuScleFit::~MuScleFit] (Destructor)"<<std::endl;
613  for (; it<muonPairs_.end();++it){
614  std::cout<<" Debugging pairs that are going to be written to file"<<std::endl;
615  std::cout<<" muon1 = "<<it->mu1<<std::endl;
616  std::cout<<" muon2 = "<<it->mu2<<std::endl;
617  }
618  }
620  }
621  else {
622  // rootTreeHandler.writeTree(outputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, &(MuScleFitUtils::genPair), saveAllToTree_ );
624  }
625  }
626  else {
627  std::cout << "ERROR: events in the vector = " << MuScleFitUtils::SavedPair.size() << " != totalEvents = " << totalEvents_ << std::endl;
628  }
629  }
630 }
631 
632 // Begin job
633 // ---------
635 // void MuScleFit::beginOfJob ()
636 // void MuScleFit::beginOfJob( const edm::EventSetup& eventSetup )
637 {
638  if (debug_>0) std::cout << "[MuScleFit]: beginOfJob" << std::endl;
639  //if(maxLoopNumber>1)
642  }
643 
644  if (debug_>0) std::cout << "[MuScleFit]: beginOfJob" << std::endl;
645 
646  // Create the root file
647  // --------------------
648  for (unsigned int i=0; i<(maxLoopNumber); i++) {
649  std::stringstream ss;
650  ss << i;
651  std::string rootFileName = ss.str() + "_" + theRootFileName_;
652  theFiles_.push_back (new TFile(rootFileName.c_str(), "RECREATE"));
653  }
654  if (debug_>0) std::cout << "[MuScleFit]: Root file created" << std::endl;
655 
656  std::cout << "creating plotter" << std::endl;
658  plotter->debug = debug_;
659 }
660 
661 // End of job method
662 // -----------------
664  if (debug_>0) std::cout << "[MuScleFit]: endOfJob" << std::endl;
665 }
666 
667 // New loop
668 // --------
669 void MuScleFit::startingNewLoop( unsigned int iLoop )
670 {
671  if (debug_>0) std::cout << "[MuScleFit]: Starting loop # " << iLoop << std::endl;
672 
673  // Number of muons used
674  // --------------------
676 
677  // Counters for problem std::cout-ing
678  // -----------------------------
680 
681  // Create the root file
682  // --------------------
683  fillHistoMap(theFiles_[iLoop], iLoop);
684 
685  loopCounter = iLoop;
687 
688  iev = 0;
690 
692 }
693 
694 // End of loop routine
695 // -------------------
696 edm::EDLooper::Status MuScleFit::endOfLoop( const edm::EventSetup& eventSetup, unsigned int iLoop )
697 {
698  unsigned int iFastLoop = 1;
699 
700  // Read the events from the root tree if requested
701  if( !(inputRootTreeFileName_.empty()) ) {
703  // When reading from local file all the loops are done here
705  iFastLoop = 0;
706  }
707  else {
708  endOfFastLoop(iLoop);
709  }
710 
711  // If a fastLoop is required we do all the remaining iterations here
712  if( fastLoop == true ) {
713  for( ; iFastLoop<maxLoopNumber; ++iFastLoop ) {
714 
715  std::cout << "Starting fast loop number " << iFastLoop << std::endl;
716 
717  // In the first loop is called by the framework
718  // if( iFastLoop > 0 ) {
719  startingNewLoop(iFastLoop);
720  // }
721 
722  // std::vector<std::pair<lorentzVector,lorentzVector> >::const_iterator it = MuScleFitUtils::SavedPair.begin();
723  // for( ; it != SavedPair.end(); ++it ) {
724  while( iev<totalEvents_ ) {
725  if( iev%50000 == 0 ) {
726  std::cout << "Fast looping on event number " << iev << std::endl;
727  }
728  // This reads muons from SavedPair using iev to keep track of the event
729  duringFastLoop();
730  }
731  std::cout << "End of fast loop number " << iFastLoop << ". Ran on " << iev << " events" << std::endl;
732  endOfFastLoop(iFastLoop);
733  }
734  }
735 
736  if (iFastLoop>=maxLoopNumber-1) {
737  return kStop;
738  } else {
739  return kContinue;
740  }
741 }
742 
743 void MuScleFit::endOfFastLoop( const unsigned int iLoop )
744 {
745  // std::cout<< "Inside endOfFastLoop, iLoop = " << iLoop << " and loopCounter = " << loopCounter << std::endl;
746 
747  if( loopCounter == 0 ) {
748  // plotter->writeHistoMap();
749  // The destructor will call the writeHistoMap after the cd to the output file
750  delete plotter;
751  }
752 
753  std::cout << "Ending loop # " << iLoop << std::endl;
754 
755  // Write the histos to file
756  // ------------------------
757  // theFiles_[iLoop]->cd();
758  writeHistoMap(iLoop);
759 
760  // Likelihood minimization to compute corrections
761  // ----------------------------------------------
762  // theFiles_[iLoop]->cd();
763  TDirectory * likelihoodDir = theFiles_[iLoop]->mkdir("likelihood");
764  likelihoodDir->cd();
766 
767  // ATTENTION, this was put BEFORE the minimizeLikelihood. Check for problems.
768  theFiles_[iLoop]->Close();
769  // ATTENTION: Check that this delete does not give any problem
770  delete theFiles_[iLoop];
771 
772  // Clear the histos
773  // ----------------
774  clearHistoMap();
775 }
776 
777 // Stuff to do during loop
778 // -----------------------
780 {
782  event.getByLabel(edm::InputTag(triggerResultsLabel_.c_str(), "", triggerResultsProcess_.c_str()), triggerResults);
783  //event.getByLabel(InputTag(triggerResultsLabel_),triggerResults);
784  bool isFired = false;
785 
786  if(triggerPath_[0] == "")
787  isFired = true;
788  else if(triggerPath_[0] == "All"){
789  isFired =triggerResults->accept();
790  if(debug_>0)
791  std::cout<<"Trigger "<<isFired<<std::endl;
792  }
793  else{
794  bool changed;
796  hltConfig.init(event.getRun(), eventSetup, triggerResultsProcess_, changed);
797 
798 
799  const edm::TriggerNames triggerNames = event.triggerNames(*triggerResults);
800 
801  for (unsigned i=0; i<triggerNames.size(); i++) {
802  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 
821  if( negateTrigger_ && isFired ) return kContinue;
822  else if( !(negateTrigger_) && !isFired ) return kContinue;
823 
824 #ifdef USE_CALLGRIND
825  CALLGRIND_START_INSTRUMENTATION;
826 #endif
827 
828  if (debug_>0) {
829  std::cout << "[MuScleFit-duringLoop]: loopCounter = " << loopCounter
830  << " Run: " << event.id().run() << " Event: " << event.id().event() << std::endl;
831  }
832 
833  // On the first iteration we read the bank, otherwise we fetch the information from the muon tree
834  // ------------------------------------ Important Note --------------------------------------- //
835  // The fillMuonCollection method applies any smearing or bias to the muons, so we NEVER use
836  // unbiased muons.
837  // ----------------------------------------------------------------------------------------------
838  if( loopCounter == 0 ) {
839 
840  if( !fastLoop || inputRootTreeFileName_.empty() ) {
841  if( debug_ > 0 ) std::cout << "Reading from edm event" << std::endl;
842  selectMuons(event);
843  duringFastLoop();
844  ++totalEvents_;
845  }
846  }
847 
848  return kContinue;
849 
850 #ifdef USE_CALLGRIND
851  CALLGRIND_STOP_INSTRUMENTATION;
852  CALLGRIND_DUMP_STATS;
853 #endif
854 }
855 
857 {
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"<<std::endl;
867  int iMu=0;
868  for (std::vector<MuScleFitMuon>::const_iterator it = muons.begin(); it < muons.end(); ++it) {
869  std::cout<<" - muon n. "<<iMu<<" = "<<(*it)<<std::endl;
870  ++iMu;
871  }
872  }
873 
874  // Find the two muons from the resonance, and set ResFound bool
875  // ------------------------------------------------------------
876  std::pair<MuScleFitMuon, MuScleFitMuon> recMuFromBestRes =
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 ) );
900  MuScleFitUtils::SavedPairMuScleFitMuons.push_back( std::make_pair( recMuScleMu1, recMuScleMu2 ) );
901  } else {
902  MuScleFitUtils::SavedPair.push_back( std::make_pair( lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.) ) );
903  MuScleFitUtils::SavedPairMuScleFitMuons.push_back( std::make_pair( MuScleFitMuon() , MuScleFitMuon() ) );
904  }
905  // Save the events also in the external tree so that it can be saved late
906 
907 
908 
909  // Fetch extra information (per event)
910  UInt_t the_NVtx(0);
911  Int_t the_numPUvtx(0);
912  Float_t the_TrueNumInteractions(0);
913 
914  // Fill pile-up related informations
915  // --------------------------------
917  event.getByLabel(puInfoSrc_, puInfo);
918  if ( puInfo.isValid() ) {
919  std::vector<PileupSummaryInfo>::const_iterator PVI;
920  for(PVI = puInfo->begin(); PVI != puInfo->end(); ++PVI) {
921  int BX = PVI->getBunchCrossing();
922  if(BX == 0) { // "0" is the in-time crossing, negative values are the early crossings, positive are late
923  the_TrueNumInteractions = PVI->getTrueNumInteractions();
924  the_numPUvtx = PVI->getPU_NumInteractions();
925  }
926  }
927  }
928 
930  event.getByLabel(vertexSrc_, vertices);
931  if ( vertices.isValid() ) {
932  std::vector<reco::Vertex>::const_iterator itv;
933  // now, count vertices
934  for (itv = vertices->begin(); itv != vertices->end(); ++itv) {
935  // require that the vertex meets certain criteria
936  if(itv->ndof()<5) continue;
937  if(fabs(itv->z())>50.0) continue;
938  if(fabs(itv->position().rho())>2.0) continue;
939  ++the_NVtx;
940  }
941  }
942 
943  // get the MC event weight
945  event.getByLabel("generator", genEvtInfo);
946  double the_genEvtweight = 1.;
947  if ( genEvtInfo.isValid() ) {
948  the_genEvtweight = genEvtInfo->weight();
949  }
950 
952  MuScleFitEvent(event.run(), event.id().event(), the_genEvtweight, the_numPUvtx, the_TrueNumInteractions, the_NVtx)
953  ));
954  // Fill the internal genPair tree from the external one
955  if( MuScleFitUtils::speedup == false ) {
956  MuScleFitUtils::genPair.push_back(std::make_pair( genMuonPairs_.back().mu1.p4(), genMuonPairs_.back().mu2.p4() ));
957  }
958 }
959 
960 void MuScleFit::selectMuons(const int maxEvents, const TString & treeFileName)
961 {
962  std::cout << "Reading muon pairs from Root Tree in " << treeFileName << std::endl;
963  RootTreeHandler rootTreeHandler;
964  std::vector<std::pair<unsigned int, unsigned long long> > evtRun;
966  rootTreeHandler.readTree(maxEvents, inputRootTreeFileName_, &(MuScleFitUtils::SavedPair), theMuonType_, &evtRun);
967  }
968  else {
970  }
971  // Now loop on all the pairs and apply any smearing and bias if needed
972  std::vector<std::pair<unsigned int, unsigned long long> >::iterator evtRunIt = evtRun.begin();
973  std::vector<std::pair<lorentzVector,lorentzVector> >::iterator it = MuScleFitUtils::SavedPair.begin();
974  std::vector<std::pair<lorentzVector,lorentzVector> >::iterator genIt;
975  if(MuScleFitUtils::speedup == false) genIt = MuScleFitUtils::genPair.begin();
976  for( ; it != MuScleFitUtils::SavedPair.end(); ++it, ++evtRunIt ) {
977 
978  // Apply any cut if requested
979  // Note that cuts here are only applied to already selected muons. They should not be used unless
980  // you are sure that the difference is negligible (e.g. the number of events with > 2 muons is negligible).
981  double pt1 = it->first.pt();
982  //std::cout << "pt1 = " << pt1 << std::endl;
983  double pt2 = it->second.pt();
984  //std::cout << "pt2 = " << pt2 << std::endl;
985  double eta1 = it->first.eta();
986  //std::cout << "eta1 = " << eta1 << std::endl;
987  double eta2 = it->second.eta();
988  //std::cout << "eta2 = " << eta2 << std::endl;
989  // If they don't pass the cuts set to null vectors
990  bool dontPass = false;
991  bool eta1InFirstRange;
992  bool eta2InFirstRange;
993  bool eta1InSecondRange;
994  bool eta2InSecondRange;
995 
1001 
1002  // This is my logic, which should be erroneous, but certainly simpler...
1005  eta1InFirstRange && eta2InSecondRange ) ) {
1006  dontPass = true;
1007  }
1008  }
1009  else {
1010  eta1 = fabs(eta1);
1011  eta2 = fabs(eta2);
1018  ( ((eta1InFirstRange && !eta2InFirstRange) && (eta2InSecondRange && !eta1InSecondRange)) ||
1019  ((eta2InFirstRange && !eta1InFirstRange) && (eta1InSecondRange && !eta2InSecondRange)) )) ) {
1020  dontPass = true;
1021  }
1022  }
1023 
1024  // Additional check on deltaPhi
1025  double deltaPhi = MuScleFitUtils::deltaPhi(it->first.phi(), it->second.phi());
1026  if( (deltaPhi <= MuScleFitUtils::deltaPhiMinCut_) || (deltaPhi >= MuScleFitUtils::deltaPhiMaxCut_) ) dontPass = true;
1027 
1028 
1029  if( dontPass ) {
1030  // std::cout << "removing muons not passing cuts" << std::endl;
1031  it->first = reco::Particle::LorentzVector(0,0,0,0);
1032  it->second = reco::Particle::LorentzVector(0,0,0,0);
1033  }
1034 
1035 
1036  // First is always mu-, second mu+
1037  if( (MuScleFitUtils::SmearType != 0) || (MuScleFitUtils::BiasType != 0) ) {
1038  applySmearing(it->first);
1039  applyBias(it->first, -1);
1040  applySmearing(it->second);
1041  applyBias(it->second, 1);
1042  }
1043 
1044  //FIXME: we loose the additional information besides the 4-momenta
1045  muonPairs_.push_back(MuonPair(MuScleFitMuon(it->first,-1), MuScleFitMuon(it->second,+1),
1046  MuScleFitEvent((*evtRunIt).first, (*evtRunIt).second, 0, 0, 0, 0) ) // FIXME: order of event and run number mixed up!
1047  );
1048 
1049 
1050  // Fill the internal genPair tree from the external one
1051  if( MuScleFitUtils::speedup == false ) {
1052  genMuonPairs_.push_back(GenMuonPair(genIt->first, genIt->second, 0));
1053  ++genIt;
1054  }
1055  }
1057  if( !(MuScleFitUtils::speedup) ) {
1059  }
1060 }
1061 
1063 {
1064  // On loops>0 the two muons are directly obtained from the SavedMuon array
1065  // -----------------------------------------------------------------------
1066  MuScleFitUtils::ResFound = false;
1068  recMu2 = (MuScleFitUtils::SavedPair[iev].second);
1069 
1070  //std::cout << "iev = " << iev << ", recMu1 pt = " << recMu1.Pt() << ", recMu2 pt = " << recMu2.Pt() << std::endl;
1071 
1072  if (recMu1.Pt()>0 && recMu2.Pt()>0) {
1073  MuScleFitUtils::ResFound = true;
1074  if (debug_>0) std::cout << "Ev = " << iev << ": found muons in tree with Pt = "
1075  << recMu1.Pt() << " " << recMu2.Pt() << std::endl;
1076  }
1077 
1078  if( debug_>0 ) std::cout << "About to start lik par correction and histo filling; ResFound is "
1079  << MuScleFitUtils::ResFound << std::endl;
1080  // If resonance found, do the hard work
1081  // ------------------------------------
1082  if( MuScleFitUtils::ResFound ) {
1083 
1084  // Find weight and reference mass for this muon pair
1085  // -------------------------------------------------
1086  // The last parameter = true means that we want to use always the background window to compute the weight,
1087  // otherwise the probability will be filled only for the resonance region.
1088  double weight = MuScleFitUtils::computeWeight( (recMu1+recMu2).mass(), iev, true );
1089  if (debug_>0) {
1090  std::cout << "Loop #" << loopCounter << "Event #" << iev << ": before correction Pt1 = "
1091  << recMu1.Pt() << " Pt2 = " << recMu2.Pt() << std::endl;
1092  }
1093  // For successive iterations, correct the muons only if the previous iteration was a scale fit.
1094  // --------------------------------------------------------------------------------------------
1095  if ( loopCounter>0 ) {
1099  }
1100  }
1101  if (debug_>0) {
1102  std::cout << "Loop #" << loopCounter << "Event #" << iev << ": after correction Pt1 = "
1103  << recMu1.Pt() << " Pt2 = " << recMu2.Pt() << std::endl;
1104  }
1105 
1107 
1108  //Fill histograms
1109  //------------------
1110 
1111  mapHisto_["hRecBestMu"]->Fill(recMu1, -1,weight);
1112  mapHisto_["hRecBestMuVSEta"]->Fill(recMu1);
1113  mapHisto_["hRecBestMu"]->Fill(recMu2, +1,weight);
1114  mapHisto_["hRecBestMuVSEta"]->Fill(recMu2);
1115  mapHisto_["hDeltaRecBestMu"]->Fill(recMu1, recMu2);
1116  // Reconstructed resonance
1117  mapHisto_["hRecBestRes"]->Fill(bestRecRes,+1, weight);
1118  mapHisto_["hRecBestResAllEvents"]->Fill(bestRecRes,+1, 1.);
1119 // // Fill histogram of Res mass vs muon variables
1120 // mapHisto_["hRecBestResVSMu"]->Fill (recMu1, bestRecRes, -1);
1121 // mapHisto_["hRecBestResVSMu"]->Fill (recMu2, bestRecRes, +1);
1122 // // Fill also the mass mu+/mu- comparisons
1123 // mapHisto_["hRecBestResVSMu"]->Fill(recMu1, recMu2, bestRecRes);
1124 
1125  mapHisto_["hRecBestResVSMu"]->Fill (recMu1, bestRecRes, -1, weight);
1126  mapHisto_["hRecBestResVSMu"]->Fill (recMu2, bestRecRes, +1, weight);
1127  // Fill also the mass mu+/mu- comparisons
1128  mapHisto_["hRecBestResVSMu"]->Fill(recMu1, recMu2, bestRecRes, weight);
1129 
1130  //-- rc 2010 filling histograms for mu+ /mu- ------
1131  // mapHisto_["hRecBestResVSMuMinus"]->Fill (recMu1, bestRecRes, -1);
1132  // mapHisto_["hRecBestResVSMuPlus"]->Fill (recMu2, bestRecRes, +1);
1133 
1134  //-- rc 2010 filling histograms MassVsMuEtaPhi------
1135  // mapHisto_["hRecBestResVSMuEtaPhi"]->Fill (recMu1, bestRecRes,-1);
1136  // mapHisto_["hRecBestResVSMuEtaPhi"]->Fill (recMu2, bestRecRes,+1);
1137 
1138  // Fill histogram of Res mass vs Res variables
1139  // mapHisto_["hRecBestResVSRes"]->Fill (bestRecRes, bestRecRes, +1);
1140  mapHisto_["hRecBestResVSRes"]->Fill (bestRecRes, bestRecRes, +1, weight);
1141 
1142 
1143 
1144 
1145 
1146 
1147  std::vector<double> * parval;
1148  std::vector<double> initpar;
1149  // Store a pointer to the vector of parameters of the last iteration, or the initial
1150  // parameters if this is the first iteration
1151  if (loopCounter==0) {
1152  initpar = MuScleFitUtils::parResol;
1153  initpar.insert( initpar.end(), MuScleFitUtils::parScale.begin(), MuScleFitUtils::parScale.end() );
1154  initpar.insert( initpar.end(), MuScleFitUtils::parCrossSection.begin(), MuScleFitUtils::parCrossSection.end() );
1155  initpar.insert( initpar.end(), MuScleFitUtils::parBgr.begin(), MuScleFitUtils::parBgr.end() );
1156  parval = &initpar;
1157  } else {
1158  parval = &(MuScleFitUtils::parvalue[loopCounter-1]);
1159  }
1160 
1161  //Compute pt resolution w.r.t generated and simulated muons
1162  //--------------------------------------------------------
1163  if( !MuScleFitUtils::speedup ) {
1164 
1165  //first is always mu-, second is always mu+
1168  }
1171  }
1172  if( compareToSimTracks_ ) {
1173  //first is always mu-, second is always mu+
1176  }
1179  }
1180  }
1181  }
1182 
1183  // 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
1184  // Fill also the resolution histogramsm using the resolution functions:
1185  // the parameters are those from the last iteration, as the muons up to this point have also the corrections from the same iteration.
1186  // Need to use a different array (ForVec), containing functors able to operate on std::vector<double>
1187  mapHisto_["hFunctionResolPt"]->Fill( recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu1.Pt(), recMu1.Eta(), *parval ), -1 );
1188  mapHisto_["hFunctionResolCotgTheta"]->Fill( recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu1.Pt(), recMu1.Eta(), *parval ), -1 );
1189  mapHisto_["hFunctionResolPhi"]->Fill( recMu1, MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu1.Pt(), recMu1.Eta(), *parval ), -1 );
1190  mapHisto_["hFunctionResolPt"]->Fill( recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaPt(recMu2.Pt(), recMu2.Eta(), *parval ), +1 );
1191  mapHisto_["hFunctionResolCotgTheta"]->Fill( recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaCotgTh(recMu2.Pt(), recMu2.Eta(), *parval ), +1 );
1192  mapHisto_["hFunctionResolPhi"]->Fill( recMu2, MuScleFitUtils::resolutionFunctionForVec->sigmaPhi(recMu2.Pt(), recMu2.Eta(), *parval ), +1 );
1193 
1194  // Compute likelihood histograms
1195  // -----------------------------
1196  if( debug_ > 0 ) std::cout << "mass = " << bestRecRes.mass() << std::endl;
1197  if (weight!=0.) {
1198  double massResol;
1199  double prob;
1200  double deltalike;
1201  if (loopCounter==0) {
1202  std::vector<double> initpar;
1203  for (int i=0; i<(int)(MuScleFitUtils::parResol.size()); i++) {
1204  initpar.push_back(MuScleFitUtils::parResol[i]);
1205  }
1206  for (int i=0; i<(int)(MuScleFitUtils::parScale.size()); i++) {
1207  initpar.push_back(MuScleFitUtils::parScale[i]);
1208  }
1209 // for (int i=0; i<(int)(MuScleFitUtils::parCrossSection.size()); i++) {
1210 // initpar.push_back(MuScleFitUtils::parCrossSection[i]);
1211 // }
1213 
1214  for (int i=0; i<(int)(MuScleFitUtils::parBgr.size()); i++) {
1215  initpar.push_back(MuScleFitUtils::parBgr[i]);
1216  }
1217  massResol = MuScleFitUtils::massResolution( recMu1, recMu2, initpar );
1218  // prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(), massResol, initpar, true );
1219  prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(), massResol,
1220  initpar, true, recMu1.eta(), recMu2.eta() );
1221  } else {
1224  // prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(),
1225  // massResol, MuScleFitUtils::parvalue[loopCounter-1], true );
1226  prob = MuScleFitUtils::massProb( bestRecRes.mass(), bestRecRes.Eta(), bestRecRes.Rapidity(),
1227  massResol, MuScleFitUtils::parvalue[loopCounter-1], true,
1228  recMu1.eta(), recMu2.eta() );
1229  }
1230  if( debug_ > 0 ) std::cout << "inside weight: mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
1231  if (prob>0) {
1232  if( debug_ > 0 ) std::cout << "inside prob: mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
1233 
1234  deltalike = log(prob)*weight; // NB maximum likelihood --> deltalike is maximized
1235  mapHisto_["hLikeVSMu"]->Fill(recMu1, deltalike);
1236  mapHisto_["hLikeVSMu"]->Fill(recMu2, deltalike);
1237  mapHisto_["hLikeVSMuMinus"]->Fill(recMu1, deltalike);
1238  mapHisto_["hLikeVSMuPlus"]->Fill(recMu2, deltalike);
1239 
1240  double recoMass = (recMu1+recMu2).mass();
1241  if( recoMass != 0 ) {
1242  // IMPORTANT: massResol is not a relative resolution
1243  mapHisto_["hResolMassVSMu"]->Fill(recMu1, massResol, -1);
1244  mapHisto_["hResolMassVSMu"]->Fill(recMu2, massResol, +1);
1245  mapHisto_["hFunctionResolMassVSMu"]->Fill(recMu1, massResol/recoMass, -1);
1246  mapHisto_["hFunctionResolMassVSMu"]->Fill(recMu2, massResol/recoMass, +1);
1247  }
1248 
1250  mapHisto_["hdMdPt1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdpt1, -1);
1251  mapHisto_["hdMdPt2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdpt2, +1);
1252  mapHisto_["hdMdPhi1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdphi1, -1);
1253  mapHisto_["hdMdPhi2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdphi2, +1);
1254  mapHisto_["hdMdCotgTh1"]->Fill(recMu1, MuScleFitUtils::massResolComponents.dmdcotgth1, -1);
1255  mapHisto_["hdMdCotgTh2"]->Fill(recMu2, MuScleFitUtils::massResolComponents.dmdcotgth2, +1);
1256  }
1257 
1258  if( !MuScleFitUtils::speedup ) {
1259  double genMass = (MuScleFitUtils::genPair[iev].first + MuScleFitUtils::genPair[iev].second).mass();
1260  // Fill the mass resolution (computed from MC), we use the covariance class to compute the variance
1261  if( genMass != 0 ) {
1263  mapHisto_["hGenResVSMu"]->Fill((MuScleFitUtils::genPair[iev].second), (MuScleFitUtils::genPair[iev].first + MuScleFitUtils::genPair[iev].second), +1);
1264  double diffMass = (recoMass - genMass)/genMass;
1265  // double diffMass = recoMass - genMass;
1266  // Fill if for both muons
1267  double pt1 = recMu1.pt();
1268  double eta1 = recMu1.eta();
1269  double pt2 = recMu2.pt();
1270  double eta2 = recMu2.eta();
1271  // This is to avoid nan
1272  if( diffMass == diffMass ) {
1273  // Mass relative difference vs Pt and Eta. To be used to extract the true mass resolution
1274  mapHisto_["hDeltaMassOverGenMassVsPt"]->Fill(pt1, diffMass);
1275  mapHisto_["hDeltaMassOverGenMassVsPt"]->Fill(pt2, diffMass);
1276  mapHisto_["hDeltaMassOverGenMassVsEta"]->Fill(eta1, diffMass);
1277  mapHisto_["hDeltaMassOverGenMassVsEta"]->Fill(eta2, diffMass);
1278  // This is used for the covariance comparison
1279  mapHisto_["hMassResolutionVsPtEta"]->Fill(pt1, eta1, diffMass, diffMass);
1280  mapHisto_["hMassResolutionVsPtEta"]->Fill(pt2, eta2, diffMass, diffMass);
1281  }
1282  else {
1283  std::cout << "Error, there is a nan: recoMass = " << recoMass << ", genMass = " << genMass << std::endl;
1284  }
1285  }
1286  // Fill with mass resolution from resolution function
1288  mapHisto_["hFunctionResolMass"]->Fill( recMu1, std::pow(massRes,2), -1 );
1289  mapHisto_["hFunctionResolMass"]->Fill( recMu2, std::pow(massRes,2), +1 );
1290  }
1291 
1292  mapHisto_["hMass_P"]->Fill(bestRecRes.mass(), prob);
1293  if( debug_ > 0 ) std::cout << "mass = " << bestRecRes.mass() << ", prob = " << prob << std::endl;
1294  mapHisto_["hMass_fine_P"]->Fill(bestRecRes.mass(), prob);
1295 
1296  mapHisto_["hMassProbVsRes"]->Fill(bestRecRes, bestRecRes, +1, prob);
1297  mapHisto_["hMassProbVsMu"]->Fill(recMu1, bestRecRes, -1, prob);
1298  mapHisto_["hMassProbVsMu"]->Fill(recMu2, bestRecRes, +1, prob);
1299  mapHisto_["hMassProbVsRes_fine"]->Fill(bestRecRes, bestRecRes, +1, prob);
1300  mapHisto_["hMassProbVsMu_fine"]->Fill(recMu1, bestRecRes, -1, prob);
1301  mapHisto_["hMassProbVsMu_fine"]->Fill(recMu2, bestRecRes, +1, prob);
1302  }
1303  }
1304  } // end if ResFound
1305 
1306  // Fill the pair
1307  // -------------
1308  if (loopCounter>0) {
1309  if (debug_>0) std::cout << "[MuScleFit]: filling the pair" << std::endl;
1310  MuScleFitUtils::SavedPair[iev] = std::make_pair( recMu1, recMu2 );
1311  }
1312 
1313  iev++;
1315 
1316  // return kContinue;
1317 }
1318 
1320  //first is always mu-, second is always mu+
1321  double deltaR = sqrt(MuScleFitUtils::deltaPhi(recMu.Phi(),genMu.Phi()) * MuScleFitUtils::deltaPhi(recMu.Phi(),genMu.Phi()) +
1322  ((recMu.Eta()-genMu.Eta()) * (recMu.Eta()-genMu.Eta())));
1323  if(deltaR<0.01)
1324  return true;
1325  else if( debug_ > 0 ) {
1326  std::cout<<"Reco muon "<<recMu<<" with eta "<<recMu.Eta()<<" and phi "<<recMu.Phi()<<std::endl
1327  <<" DOES NOT MATCH with generated muon from resonance: "<<std::endl
1328  <<genMu<<" with eta "<<genMu.Eta()<<" and phi "<<genMu.Phi()<<std::endl;
1329  }
1330  return false;
1331 }
1332 
1334  const std::string & inputName, const int charge )
1335 {
1336  std::string name(inputName + "VSMu");
1337  mapHisto_["hResolPt"+name]->Fill(recMu, (-genMu.Pt()+recMu.Pt())/genMu.Pt(), charge);
1338  mapHisto_["hResolTheta"+name]->Fill(recMu, (-genMu.Theta()+recMu.Theta()), charge);
1339  mapHisto_["hResolCotgTheta"+name]->Fill(recMu,(-cos(genMu.Theta())/sin(genMu.Theta())
1340  +cos(recMu.Theta())/sin(recMu.Theta())), charge);
1341  mapHisto_["hResolEta"+name]->Fill(recMu, (-genMu.Eta()+recMu.Eta()),charge);
1342  mapHisto_["hResolPhi"+name]->Fill(recMu, MuScleFitUtils::deltaPhiNoFabs(recMu.Phi(), genMu.Phi()), charge);
1343 
1344  // Fill only if it was matched to a genMu and this muon is valid
1345  if( (genMu.Pt() != 0) && (recMu.Pt() != 0) ) {
1346  mapHisto_["hPtRecoVsPt"+inputName]->Fill(genMu.Pt(), recMu.Pt());
1347  }
1348 }
1349 
1351 {
1352  if( MuScleFitUtils::SmearType>0 ) {
1353  mu = MuScleFitUtils::applySmearing( mu );
1354  if (debug_>0) std::cout << "Muon #" << MuScleFitUtils::goodmuon
1355  << ": after smearing Pt = " << mu.Pt() << std::endl;
1356  }
1357 }
1358 
1360 {
1361  if( MuScleFitUtils::BiasType>0 ) {
1362  mu = MuScleFitUtils::applyBias( mu, charge );
1363  if (debug_>0) std::cout << "Muon #" << MuScleFitUtils::goodmuon
1364  << ": after bias Pt = " << mu.Pt() << std::endl;
1365  }
1366 }
1367 
1368 // Simple method to check parameters consistency. It aborts the job if the parameters
1369 // are not consistent.
1371 
1372  // Fits selection dimension check
1374  std::cout << "[MuScleFit-Constructor]: wrong size of resolution fits selector = " << MuScleFitUtils::doResolFit.size() << std::endl;
1375  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1376  abort();
1377  }
1379  std::cout << "[MuScleFit-Constructor]: wrong size of scale fits selector = " << MuScleFitUtils::doScaleFit.size() << std::endl;
1380  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1381  abort();
1382  }
1384  std::cout << "[MuScleFit-Constructor]: wrong size of cross section fits selector = " << MuScleFitUtils::doCrossSectionFit.size() << std::endl;
1385  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1386  abort();
1387  }
1389  std::cout << "[MuScleFit-Constructor]: wrong size of background fits selector = " << MuScleFitUtils::doBackgroundFit.size() << std::endl;
1390  std::cout << "it must have as many values as the number of loops, which is = " << maxLoopNumber << std::endl;
1391  abort();
1392  }
1393 
1394  // Bias parameters: dimension check
1395  // --------------------------------
1396  if ((MuScleFitUtils::BiasType==1 && MuScleFitUtils::parBias.size()!=2) || // linear in pt
1397  (MuScleFitUtils::BiasType==2 && MuScleFitUtils::parBias.size()!=2) || // linear in |eta|
1398  (MuScleFitUtils::BiasType==3 && MuScleFitUtils::parBias.size()!=4) || // sinusoidal in phi
1399  (MuScleFitUtils::BiasType==4 && MuScleFitUtils::parBias.size()!=3) || // linear in pt and |eta|
1400  (MuScleFitUtils::BiasType==5 && MuScleFitUtils::parBias.size()!=3) || // linear in pt and sinusoidal in phi
1401  (MuScleFitUtils::BiasType==6 && MuScleFitUtils::parBias.size()!=3) || // linear in |eta| and sinusoidal in phi
1402  (MuScleFitUtils::BiasType==7 && MuScleFitUtils::parBias.size()!=4) || // linear in pt and |eta| and sinusoidal in phi
1403  (MuScleFitUtils::BiasType==8 && MuScleFitUtils::parBias.size()!=4) || // linear in pt and parabolic in |eta|
1404  (MuScleFitUtils::BiasType==9 && MuScleFitUtils::parBias.size()!=2) || // exponential in pt
1405  (MuScleFitUtils::BiasType==10 && MuScleFitUtils::parBias.size()!=3) || // parabolic in pt
1406  (MuScleFitUtils::BiasType==11 && MuScleFitUtils::parBias.size()!=4) || // linear in pt and sin in phi with chg
1407  (MuScleFitUtils::BiasType==12 && MuScleFitUtils::parBias.size()!=6) || // linear in pt and para in plus sin in phi with chg
1408  (MuScleFitUtils::BiasType==13 && MuScleFitUtils::parBias.size()!=8) || // linear in pt and para in plus sin in phi with chg
1409  MuScleFitUtils::BiasType<0 || MuScleFitUtils::BiasType>13) {
1410  std::cout << "[MuScleFit-Constructor]: Wrong bias type or number of parameters: aborting!" << std::endl;
1411  abort();
1412  }
1413  // Smear parameters: dimension check
1414  // ---------------------------------
1422  MuScleFitUtils::SmearType<0 || MuScleFitUtils::SmearType>7) {
1423  std::cout << "[MuScleFit-Constructor]: Wrong smear type or number of parameters: aborting!" << std::endl;
1424  abort();
1425  }
1426  // Protect against bad size of parameters
1427  // --------------------------------------
1430  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Resol: aborting!" << std::endl;
1431  abort();
1432  }
1435  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Scale: aborting!" << std::endl;
1436  abort();
1437  }
1440  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Bgr: aborting!" << std::endl;
1441  abort();
1442  }
1445  std::cout << "[MuScleFit-Constructor]: Mismatch in number of parameters for Bgr: aborting!" << std::endl;
1446  abort();
1447  }
1448 
1449  // Protect against an incorrect number of resonances
1450  // -------------------------------------------------
1451  if (MuScleFitUtils::resfind.size()!=6) {
1452  std::cout << "[MuScleFit-Constructor]: resfind must have 6 elements (1 Z, 3 Y, 2 Psi): aborting!" << std::endl;
1453  abort();
1454  }
1455 }
1456 
1458 
1459  reco::TrackRef iTrack = aMuon->innerTrack();
1460  const reco::HitPattern& p = iTrack->hitPattern();
1461 
1462  reco::TrackRef gTrack = aMuon->globalTrack();
1463  const reco::HitPattern& q = gTrack->hitPattern();
1464 
1465  return (//isMuonInAccept(aMuon) &&// no acceptance cuts!
1466  iTrack->found() > 11 &&
1467  gTrack->chi2()/gTrack->ndof() < 20.0 &&
1468  q.numberOfValidMuonHits() > 0 &&
1469  iTrack->chi2()/iTrack->ndof() < 4.0 &&
1470  aMuon->muonID("TrackerMuonArbitrated") &&
1471  aMuon->muonID("TMLastStationAngTight") &&
1472  p.pixelLayersWithMeasurement() > 1 &&
1473  fabs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
1474  fabs(iTrack->dz()) < 15.0 );//should be done w.r.t. PV!
1475 }
1476 
1477 
1479 
1480  reco::TrackRef iTrack = aMuon->innerTrack();
1481  const reco::HitPattern& p = iTrack->hitPattern();
1482 
1483  return (//isMuonInAccept(aMuon) // no acceptance cuts!
1484  iTrack->found() > 11 &&
1485  iTrack->chi2()/iTrack->ndof() < 4.0 &&
1486  aMuon->muonID("TrackerMuonArbitrated") &&
1487  aMuon->muonID("TMLastStationAngTight") &&
1488  p.pixelLayersWithMeasurement() > 1 &&
1489  fabs(iTrack->dxy()) < 3.0 && //should be done w.r.t. PV!
1490  fabs(iTrack->dz()) < 15.0 );//should be done w.r.t. PV!
1491 
1492 }
1493 
1494 
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
static std::vector< int > doResolFit
T getParameter(std::string const &) const
scaleFunctionBase< double * > * scaleFunctionService(const int identifier)
Service to build the scale functor corresponding to the passed identifier.
Definition: Functions.cc:3
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
static std::vector< double > parBias
void selectMuons(const edm::Event &event)
Definition: MuScleFit.cc:856
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
reco::TrackRef innerTrack() const
reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
Definition: Muon.h:72
static smearFunctionBase * smearFunction
static std::vector< int > parScaleOrder
edm::InputTag theMuonLabel_
Definition: MuScleFitBase.h:46
static std::vector< int > doCrossSectionFit
static void minimizeLikelihood()
static double maxMuonEtaSecondRange_
std::vector< MuScleFitMuon > fillMuonCollection(const std::vector< T > &tracks)
Definition: MuScleFit.cc:316
static std::vector< int > parBgrOrder
bool selGlobalMuon(const pat::Muon *aMuon)
Function for onia selections.
Definition: MuScleFit.cc:1457
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:277
int totalEvents_
Definition: MuScleFit.cc:288
static bool startWithSimplex_
std::string triggerResultsProcess_
Definition: MuScleFit.cc:303
static std::vector< int > doBackgroundFit
static std::vector< double > parResol
static bool debugMassResol_
bool saveAllToTree_
Definition: MuScleFit.cc:306
virtual void duringFastLoop()
Definition: MuScleFit.cc:1062
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
void checkParameters()
Definition: MuScleFit.cc:1370
static BackgroundHandler * backgroundHandler
static double x[7][10000]
static int debug
static double ResMass[6]
edm::InputTag vertexSrc_
Definition: MuScleFit.cc:310
static double massWindowHalfWidth[3][6]
static bool speedup
#define DEFINE_FWK_LOOPER(type)
bool muonID(const std::string &name) const
Run const & getRun() const
Definition: Event.cc:63
std::map< std::string, Histograms * > mapHisto_
The map of histograms.
Definition: MuScleFitBase.h:77
static bool ResFound
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:1478
bool ifHepMC
Definition: MuScleFit.cc:266
virtual void startingNewLoop(unsigned int iLoop) override
Definition: MuScleFit.cc:669
Strings::size_type size() const
Definition: TriggerNames.cc:39
static std::vector< int > parBgrFix
static std::vector< double > parResolMin
static bool minimumShapePlots_
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:458
edm::InputTag simTracksCollection_
Definition: MuScleFit.cc:291
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:1333
static struct MuScleFitUtils::massResolComponentsStruct massResolComponents
Strings const & triggerNames() const
Definition: TriggerNames.cc:24
static int BiasType
std::vector< std::string > triggerPath_
Definition: MuScleFit.cc:304
static std::vector< int > parScaleFix
static bool computeMinosErrors_
MuScleFitMuon recMuScleMu1
Definition: MuScleFit.cc:286
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)
U second(std::pair< T, U > const &p)
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:298
std::string theGenInfoRootFileName_
Definition: MuScleFitBase.h:48
int numberOfSimVertices
Definition: MuScleFit.cc:263
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:50
static double massResolution(const lorentzVector &mu1, const lorentzVector &mu2)
reco::Particle::LorentzVector recMu2
Definition: MuScleFit.cc:285
bool PATmuons_
Definition: MuScleFit.cc:292
static std::vector< double > parScaleMin
reco::Particle::LorentzVector recMu1
Definition: MuScleFit.cc:285
static lorentzVector applyBias(const lorentzVector &muon, const int charge)
T sqrt(T t)
Definition: SSEVec.h:48
static std::vector< double > parBgr
unsigned int maxLoopNumber
Definition: MuScleFit.cc:276
static int SmearType
MuonServiceProxy * theService
Definition: MuScleFit.cc:257
std::string genParticlesName_
Definition: MuScleFit.cc:293
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
RunNumber_t run() const
Definition: Event.h:87
void clearHistoMap()
Clean the histograms map.
int j
Definition: DBlmapReader.cc:9
MuScleFitMuon recMuScleMu2
Definition: MuScleFit.cc:286
bool compareToSimTracks_
Definition: MuScleFit.cc:290
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
const int mu
Definition: Constants.h:22
MuScleFit(const edm::ParameterSet &pset)
Definition: MuScleFit.cc:376
static std::string const triggerResults
Definition: EdmProvDump.cc:40
static lorentzVector applyScale(const lorentzVector &muon, const std::vector< double > &parval, const int charge)
static std::vector< int > parResolFix
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:75
static std::vector< double > parSmear
void beginOfJobInConstructor()
Definition: MuScleFit.cc:634
smearFunctionBase * smearFunctionService(const int identifier)
Service to build the smearing functor corresponding to the passed identifier.
Definition: Functions.cc:23
reco::TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector (reimplemented from reco::Muon) ...
Definition: Muon.h:80
static std::vector< int > parResolOrder
static bool sherpa_
std::auto_ptr< MuScleFitMuonSelector > muonSelector_
Definition: MuScleFit.cc:312
virtual void endOfJob() override
Definition: MuScleFit.cc:663
bool fastLoop
Definition: MuScleFit.cc:279
virtual edm::EDLooper::Status duringLoop(const edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: MuScleFit.cc:779
tuple G
Definition: callgraph.py:12
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
int numberOfSimMuons
Definition: MuScleFit.cc:262
void writeTree(const TString &fileName, const std::vector< MuonPair > *savedPair, const int muonType=0, const std::vector< GenMuonPair > *genPair=0, const bool saveAll=false)
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:13
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:1350
virtual ~MuScleFit()
Definition: MuScleFit.cc:599
int numberOfSimTracks
Definition: MuScleFit.cc:261
tuple tracks
Definition: testEve_cfg.py:39
std::string const & triggerName(unsigned int index) const
Definition: TriggerNames.cc:27
std::string theRootFileName_
Definition: MuScleFitBase.h:47
void fillTreeGen(const std::vector< std::pair< reco::Particle::LorentzVector, reco::Particle::LorentzVector > > &genPairs)
edm::InputTag puInfoSrc_
Definition: MuScleFit.cc:309
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:302
static bool rapidityBinsForZ_
MuScleFitPlotter * plotter
Definition: MuScleFit.cc:281
void fillHistoMap(TFile *outputFile, unsigned int iLoop)
Create the histograms map.
Definition: MuScleFitBase.cc:7
std::vector< TFile * > theFiles_
The files were the histograms are saved.
Definition: MuScleFitBase.h:74
double minResMass_hwindow[6]
Definition: MuScleFit.cc:271
static int ResolFitType
static int goodmuon
static bool separateRanges_
static double minMuonEtaSecondRange_
static int iev_
static std::vector< double > parCrossSection
static const double mMu2
tuple muons
Definition: patZpeak.py:38
static lorentzVector applySmearing(const lorentzVector &muon)
static bool normalizeLikelihoodByEventNumber_
std::string inputRootTreeFileName_
Definition: MuScleFit.cc:296
bool checkDeltaR(reco::Particle::LorentzVector &genMu, reco::Particle::LorentzVector &recMu)
Check if two lorentzVector are near in deltaR.
Definition: MuScleFit.cc:1319
resolutionFunctionBase< double * > * resolutionFunctionService(const int identifier)
Service to build the resolution functor corresponding to the passed identifier.
Definition: Functions.cc:38
int maxEventsFromRootTree_
Definition: MuScleFit.cc:300
double maxResMass_hwindow[6]
Definition: MuScleFit.cc:272
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)
tuple cout
Definition: gather_cfg.py:121
static double minMuonPt_
static scaleFunctionBase< double * > * scaleFunction
int numberOfValidMuonHits() const
Definition: HitPattern.h:744
static bool useProbsFile_
int weight
Definition: histoStyle.py:50
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:353
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:49
static double maxMuonEtaFirstRange_
tuple size
Write out results.
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void applyBias(reco::Particle::LorentzVector &mu, const int charge)
Apply the bias if needed using the function in MuScleFitUtils.
Definition: MuScleFit.cc:1359
bool negateTrigger_
Definition: MuScleFit.cc:305
bool ifGenPart
Definition: MuScleFit.cc:267
static CrossSectionHandler * crossSectionHandler
virtual edm::EDLooper::Status endOfLoop(const edm::EventSetup &eventSetup, unsigned int iLoop) override
Definition: MuScleFit.cc:696
static std::vector< double > parResolStep
static std::vector< int > parCrossSectionFix
static resolutionFunctionBase< double * > * resolutionFunction
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=0)
void addParameters(std::vector< double > &initpar)
Inputs the vars in a vector.
virtual void endOfFastLoop(const unsigned int iLoop)
Definition: MuScleFit.cc:743
int numberOfEwkZ
Definition: MuScleFit.cc:264