CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
TestCorrection Class Reference

#include <MuonAnalysis/MomentumScaleCalibration/plugins/TestCorrection.cc>

Inheritance diagram for TestCorrection:
edm::EDAnalyzer MuScleFitBase edm::EDAnalyzer MuScleFitBase edm::EDConsumerBase edm::EDConsumerBase

Public Member Functions

 TestCorrection (const edm::ParameterSet &)
 
 TestCorrection (const edm::ParameterSet &)
 
 ~TestCorrection ()
 
 ~TestCorrection ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
lorentzVector correctMuon (const lorentzVector &muon)
 
lorentzVector correctMuon (const lorentzVector &muon)
 
virtual void endJob ()
 
virtual void endJob () override
 
template<typename T >
std::vector< MuScleFitMuonfillMuonCollection (const std::vector< T > &tracks)
 
template<typename T >
std::vector< MuScleFitMuonfillMuonCollection (const std::vector< T > &tracks)
 
virtual void initialize (const edm::EventSetup &)
 
virtual void initialize (const edm::EventSetup &)
 
- Private Member Functions inherited from MuScleFitBase
 MuScleFitBase (const edm::ParameterSet &iConfig)
 
virtual ~MuScleFitBase ()(false)
 
void clearHistoMap ()
 Clean the histograms map. More...
 
void fillHistoMap (TFile *outputFile, unsigned int iLoop)
 Create the histograms map. More...
 
void readProbabilityDistributionsFromFile ()
 Read probability distributions from a local root file. More...
 
void writeHistoMap (const unsigned int iLoop)
 Save the histograms map to file. More...
 

Private Attributes

std::auto_ptr< BackgroundFunctionbackground_
 
TH1F * correctedPt_
 
TProfile * correctedPtVsEta_
 
std::auto_ptr< MomentumScaleCorrectorcorrector_
 
int eventCounter_
 
edm::EDGetTokenT< reco::MuonCollectionglbMuonsToken_
 
std::auto_ptr< ResolutionFunctionresolution_
 
edm::EDGetTokenT< reco::TrackCollectionsaMuonsToken_
 
edm::EDGetTokenT< reco::TrackCollectiontracksToken_
 
TH1F * uncorrectedPt_
 
TProfile * uncorrectedPtVsEta_
 
- Private Attributes inherited from MuScleFitBase
int debug_
 
std::vector< GenMuonPairgenMuonPairs_
 Stores the genMuon pairs and the motherId prior to the creation of the internal tree. More...
 
std::map< std::string, Histograms * > mapHisto_
 The map of histograms. More...
 
std::vector< MuonPairmuonPairs_
 Used to store the muon pairs plus run and event number prior to the creation of the internal tree. More...
 
std::string probabilitiesFile_
 
std::string probabilitiesFileInPath_
 
std::vector< TFile * > theFiles_
 The files were the histograms are saved. More...
 
std::string theGenInfoRootFileName_
 
edm::InputTag theMuonLabel_
 
int theMuonType_
 
std::string theRootFileName_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 63 of file TestCorrection.cc.

Constructor & Destructor Documentation

TestCorrection::TestCorrection ( const edm::ParameterSet iConfig)
explicit

Definition at line 148 of file TestCorrection.cc.

References background_, correctedPt_, correctedPtVsEta_, corrector_, gather_cfg::cout, eventCounter_, MuScleFitBase::fillHistoMap(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), download_sqlite_cfg::outputFile, MuScleFitUtils::parResol, MuScleFitUtils::resfind, resolution_, MuScleFitUtils::resolutionFunction, MuScleFitUtils::resolutionFunctionForVec, resolutionFunctionVecService(), AlCaHLTBitMon_QueryRunRegistry::string, MuScleFitBase::theFiles_, MuScleFitBase::theRootFileName_, uncorrectedPt_, and uncorrectedPtVsEta_.

148  :
149  MuScleFitBase( iConfig ),
150  glbMuonsToken_(mayConsume<reco::MuonCollection>(theMuonLabel_)),
151  saMuonsToken_(mayConsume<reco::TrackCollection>(theMuonLabel_)),
152  tracksToken_(mayConsume<reco::TrackCollection>(theMuonLabel_))
153 {
154  //now do what ever initialization is needed
155  TFile * outputFile = new TFile(theRootFileName_.c_str(), "RECREATE");
156  theFiles_.push_back(outputFile);
157  // outputFile_ = new TFile(theRootFileName_.c_str(), "RECREATE");
158  // outputFile_->cd();
159  outputFile->cd();
160  MuScleFitUtils::resfind = iConfig.getParameter<std::vector<int> >("resfind");
161  fillHistoMap(outputFile, 0);
162  uncorrectedPt_ = new TH1F("uncorrectedPt", "uncorrected pt", 1000, 0, 100);
163  uncorrectedPtVsEta_ = new TProfile("uncorrectedPtVsEta", "uncorrected pt vs eta", 1000, 0, 100, -3., 3.);
164  correctedPt_ = new TH1F("correctedPt", "corrected pt", 1000, 0, 100);
165  correctedPtVsEta_ = new TProfile("correctedPtVsEta", "corrected pt vs eta", 1000, 0, 100, -3., 3.);
166  eventCounter_ = 0;
167  // Create the corrector and set the parameters
168  corrector_.reset(new MomentumScaleCorrector( iConfig.getUntrackedParameter<std::string>("CorrectionsIdentifier") ) );
169  std::cout << "corrector_ = " << &*corrector_ << std::endl;
170  resolution_.reset(new ResolutionFunction(iConfig.getUntrackedParameter<std::string>("ResolutionsIdentifier") ) );
171  std::cout << "resolution_ = " << &*resolution_ << std::endl;
172  background_.reset(new BackgroundFunction(iConfig.getUntrackedParameter<std::string>("BackgroundIdentifier") ) );
173 
174  // Initialize the parameters of MuScleFitUtils from those saved in the functions.
175  // MuScleFitUtils::parScale = corrector_.getFunction(0)->parameters();
178 
179  MuScleFitUtils::parResol = resolution_->parameters();
180 }
T getParameter(std::string const &) const
MuScleFitBase(const edm::ParameterSet &iConfig)
Definition: MuScleFitBase.h:21
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag theMuonLabel_
Definition: MuScleFitBase.h:46
TProfile * correctedPtVsEta_
static std::vector< double > parResol
std::auto_ptr< BackgroundFunction > background_
edm::EDGetTokenT< reco::TrackCollection > saMuonsToken_
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
TProfile * uncorrectedPtVsEta_
std::string theRootFileName_
Definition: MuScleFitBase.h:47
std::auto_ptr< MomentumScaleCorrector > corrector_
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
edm::EDGetTokenT< reco::MuonCollection > glbMuonsToken_
std::auto_ptr< ResolutionFunction > resolution_
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
static std::vector< int > resfind
static resolutionFunctionBase< std::vector< double > > * resolutionFunctionForVec
static resolutionFunctionBase< double * > * resolutionFunction
TestCorrection::~TestCorrection ( )

Definition at line 182 of file TestCorrection.cc.

References svgfig::canvas(), correctedPt_, correctedPtVsEta_, gather_cfg::cout, eventCounter_, listHistos::legend, MuScleFitBase::theFiles_, uncorrectedPt_, uncorrectedPtVsEta_, and MuScleFitBase::writeHistoMap().

183 {
184  theFiles_[0]->cd();
185  TCanvas canvas("ptComparison","pt comparison", 1000, 800);
186  canvas.cd();
187  uncorrectedPt_->GetXaxis()->SetTitle("Pt(GeV)");
188  correctedPt_->SetLineColor(kRed);
189  TLegend * legend = new TLegend(0.7,0.71,0.98,1.);
190  legend->SetTextSize(0.02);
191  legend->SetFillColor(0); // Have a white background
192  legend->AddEntry(uncorrectedPt_, "original pt");
193  legend->AddEntry(correctedPt_, "corrected pt");
194  uncorrectedPt_->Draw();
195  correctedPt_->Draw("same");
196  legend->Draw("same");
197 
198  canvas.Write();
199  uncorrectedPt_->Write();
200  uncorrectedPtVsEta_->Write();
201  correctedPt_->Write();
202  correctedPtVsEta_->Write();
203 
204  writeHistoMap(0);
205  theFiles_[0]->Close();
206 
207  std::cout << "Total analyzed events = " << eventCounter_ << std::endl;
208 }
TProfile * correctedPtVsEta_
TProfile * uncorrectedPtVsEta_
void writeHistoMap(const unsigned int iLoop)
Save the histograms map to file.
std::vector< TFile * > theFiles_
The files were the histograms are saved.
Definition: MuScleFitBase.h:74
def canvas(sub, attr)
Definition: svgfig.py:481
TestCorrection::TestCorrection ( const edm::ParameterSet )
explicit
TestCorrection::~TestCorrection ( )

Member Function Documentation

virtual void TestCorrection::analyze ( const edm::Event ,
const edm::EventSetup  
)
privatevirtual
void TestCorrection::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Definition at line 215 of file TestCorrection.cc.

References funct::abs(), correctedPt_, correctedPtVsEta_, correctMuon(), gather_cfg::cout, eventCounter_, fillMuonCollection(), MuScleFitUtils::findBestRecoRes(), edm::Event::getByToken(), glbMuonsToken_, initialize(), MuScleFitBase::mapHisto_, electronCleaner_cfi::muons, MuScleFitUtils::ResFound, saMuonsToken_, MuScleFitUtils::SavedPair, MuScleFitBase::theMuonType_, l1t::tracks, tracksToken_, uncorrectedPt_, and uncorrectedPtVsEta_.

215  {
216  using namespace edm;
217 
218  initialize(iSetup);
219 
220  ++eventCounter_;
221  if ( eventCounter_%100 == 0 ) {
222  std::cout << "Event number " << eventCounter_ << std::endl;
223  }
224 
225  // Take the reco-muons, depending on the type selected in the cfg
226  // --------------------------------------------------------------
227 
228  std::vector<MuScleFitMuon> muons;
229 
230  if (theMuonType_==1) { // GlobalMuons
232  iEvent.getByToken (glbMuonsToken_, glbMuons);
233  muons = fillMuonCollection(*glbMuons);
234  }
235  else if (theMuonType_==2) { // StandaloneMuons
237  iEvent.getByToken (saMuonsToken_, saMuons);
238  muons = fillMuonCollection(*saMuons);
239  }
240  else if (theMuonType_==3) { // Tracker tracks
242  iEvent.getByToken (tracksToken_, tracks);
243  muons = fillMuonCollection(*tracks);
244  }
245 
246  // Find the two muons from the resonance, and set ResFound bool
247  // ------------------------------------------------------------
248  std::pair <MuScleFitMuon, MuScleFitMuon> recMuFromBestRes =
251  MuScleFitUtils::SavedPair.push_back( std::make_pair (recMuFromBestRes.first.p4(), recMuFromBestRes.second.p4()) );
252  } else {
253  MuScleFitUtils::SavedPair.push_back( std::make_pair (lorentzVector(0.,0.,0.,0.), lorentzVector(0.,0.,0.,0.)) );
254  }
255 
256  // If resonance found, do the hard work
257  // ------------------------------------
259 
260  // Find weight and reference mass for this muon pair
261  // -------------------------------------------------
262  // double weight = MuScleFitUtils::computeWeight ((recMu1+recMu2).mass());
263 
264  // Use the correction function to correct the pt scale of the muons. Note that this takes into
265  // account the corrections from all iterations.
266  lorentzVector recMu1;
267  recMu1 = correctMuon(recMu1);
268  lorentzVector recMu2;
269  recMu2 = correctMuon(recMu2);
270 
271  reco::Particle::LorentzVector bestRecRes (recMu1+recMu2);
272 
273  //Fill histograms
274  //------------------
275  mapHisto_["hRecBestMu"]->Fill(recMu1);
276  if ((std::abs(recMu1.eta())<2.5) && (recMu1.pt()>2.5)) {
277  mapHisto_["hRecBestMu_Acc"]->Fill(recMu1);
278  }
279  mapHisto_["hRecBestMu"]->Fill(recMu2);
280  if ((std::abs(recMu2.eta())<2.5) && (recMu2.pt()>2.5)) {
281  mapHisto_["hRecBestMu_Acc"]->Fill(recMu2);
282  }
283  mapHisto_["hDeltaRecBestMu"]->Fill(recMu1, recMu2);
284 
285  mapHisto_["hRecBestRes"]->Fill(bestRecRes);
286  if ((std::abs(recMu1.eta())<2.5) && (recMu1.pt()>2.5) && (std::abs(recMu2.eta())<2.5) && (recMu2.pt()>2.5)){
287  mapHisto_["hRecBestRes_Acc"]->Fill(bestRecRes);
288  // Fill histogram of Res mass vs muon variable
289  mapHisto_["hRecBestResVSMu"]->Fill (recMu1, bestRecRes, -1);
290  mapHisto_["hRecBestResVSMu"]->Fill (recMu2, bestRecRes, +1);
291  }
292  }
293 
294  // Loop on the recMuons
295  std::vector<MuScleFitMuon>::const_iterator recMuon = muons.begin();
296  int muonCount = 0;
297  for ( ; recMuon!=muons.end(); ++recMuon, ++muonCount ) {
298 
299  // Fill the histogram with uncorrected pt values
300  uncorrectedPt_->Fill(recMuon->pt());
301  uncorrectedPtVsEta_->Fill(recMuon->pt(), recMuon->eta());
302 
303  // Fill the histogram with corrected pt values
304  std::cout << "correcting muon["<<muonCount<<"] with pt = " << recMuon->pt() << std::endl;
305  double corrPt = (*corrector_)(*recMuon);
306  std::cout << "to pt = " << corrPt << std::endl;
307  correctedPt_->Fill(corrPt);
308  correctedPtVsEta_->Fill(corrPt, recMuon->eta());
309  // correctedPt_->Fill(recMuon->pt());
310  }
311 }
lorentzVector correctMuon(const lorentzVector &muon)
TProfile * correctedPtVsEta_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
virtual void initialize(const edm::EventSetup &)
edm::EDGetTokenT< reco::TrackCollection > saMuonsToken_
std::map< std::string, Histograms * > mapHisto_
The map of histograms.
Definition: MuScleFitBase.h:77
static bool ResFound
reco::Particle::LorentzVector lorentzVector
Definition: GenMuonPair.h:9
static std::pair< MuScleFitMuon, MuScleFitMuon > findBestRecoRes(const std::vector< MuScleFitMuon > &muons)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static std::vector< std::pair< lorentzVector, lorentzVector > > SavedPair
TProfile * uncorrectedPtVsEta_
std::vector< MuScleFitMuon > fillMuonCollection(const std::vector< T > &tracks)
edm::EDGetTokenT< reco::MuonCollection > glbMuonsToken_
HLT enums.
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
lorentzVector TestCorrection::correctMuon ( const lorentzVector muon)
private
lorentzVector TestCorrection::correctMuon ( const lorentzVector muon)
private

Definition at line 313 of file TestCorrection.cc.

References corrector_, and MuScleFitUtils::fromPtEtaPhiToPxPyPz().

Referenced by analyze(), and fillMuonCollection().

313  {
314  double corrPt = corrector_->correct(muon);
315  double ptEtaPhiE[4] = {corrPt, muon.Eta(), muon.Phi(), muon.E()};
316  return MuScleFitUtils::fromPtEtaPhiToPxPyPz(ptEtaPhiE);
317 }
static lorentzVector fromPtEtaPhiToPxPyPz(const double *ptEtaPhiE)
std::auto_ptr< MomentumScaleCorrector > corrector_
virtual void TestCorrection::endJob ( void  )
inlineprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 69 of file TestCorrection.h.

69 {};
virtual void TestCorrection::endJob ( void  )
inlineoverrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 71 of file TestCorrection.cc.

71 {};
template<typename T >
std::vector<MuScleFitMuon> TestCorrection::fillMuonCollection ( const std::vector< T > &  tracks)
inlineprivate

Definition at line 71 of file TestCorrection.h.

References background_, correctedPt_, correctedPtVsEta_, correctMuon(), corrector_, gather_cfg::cout, MuScleFitBase::debug_, eventCounter_, RPCpg::mu, metsig::muon, electronCleaner_cfi::muons, resolution_, mathSSE::sqrt(), HiIsolationCommonParameters_cff::track, uncorrectedPt_, and uncorrectedPtVsEta_.

72  {
73  std::vector<MuScleFitMuon> muons;
74  typename std::vector<T>::const_iterator track;
75  for( track = tracks.begin(); track != tracks.end(); ++track ) {
77  mu = reco::Particle::LorentzVector(track->px(),track->py(),track->pz(),
78  sqrt(track->p()*track->p() + + 0.011163612));
79 
80  Double_t hitsTk(0), hitsMuon(0), ptError(0);
81  if ( const reco::Muon* myMu = dynamic_cast<const reco::Muon*>(&(*track)) ){
82  hitsTk = myMu->innerTrack()->hitPattern().numberOfValidTrackerHits();
83  hitsMuon = myMu->innerTrack()->hitPattern().numberOfValidMuonHits();
84  ptError = myMu->innerTrack()->ptError();
85  }
86  else if ( const pat::Muon* myMu = dynamic_cast<const pat::Muon*>(&(*track)) ) {
87  hitsTk = myMu->innerTrack()->hitPattern().numberOfValidTrackerHits();
88  hitsMuon = myMu->innerTrack()->hitPattern().numberOfValidMuonHits();
89  ptError = myMu->innerTrack()->ptError();
90  }
91  else if (const reco::Track* myMu = dynamic_cast<const reco::Track*>(&(*track))){
92  hitsTk = myMu->hitPattern().numberOfValidTrackerHits();
93  hitsMuon = myMu->hitPattern().numberOfValidMuonHits();
94  ptError = myMu->ptError();
95  }
96 
97  MuScleFitMuon muon(mu,track->charge(),ptError,hitsTk,hitsMuon,false);
98 
99  if (debug_>0) {
100  std::cout<<"[TestCorrection::fillMuonCollection] after MuScleFitMuon initialization"<<std::endl;
101  std::cout<<" muon = "<<muon<<std::endl;
102  }
103 
104  muons.push_back(muon);
105  }
106  return muons;
107  }
T sqrt(T t)
Definition: SSEVec.h:18
const int mu
Definition: Constants.h:22
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
Analysis-level muon class.
Definition: Muon.h:49
template<typename T >
std::vector<MuScleFitMuon> TestCorrection::fillMuonCollection ( const std::vector< T > &  tracks)
inlineprivate

Definition at line 73 of file TestCorrection.cc.

References correctMuon(), gather_cfg::cout, MuScleFitBase::debug_, RPCpg::mu, metsig::muon, electronCleaner_cfi::muons, mathSSE::sqrt(), and HiIsolationCommonParameters_cff::track.

Referenced by analyze().

74  {
75  std::vector<MuScleFitMuon> muons;
76  typename std::vector<T>::const_iterator track;
77  for( track = tracks.begin(); track != tracks.end(); ++track ) {
79  mu = reco::Particle::LorentzVector(track->px(),track->py(),track->pz(),
80  sqrt(track->p()*track->p() + + 0.011163612));
81 
82  Double_t hitsTk(0), hitsMuon(0), ptError(0);
83  if ( const reco::Muon* myMu = dynamic_cast<const reco::Muon*>(&(*track)) ){
84  hitsTk = myMu->innerTrack()->hitPattern().numberOfValidTrackerHits();
85  hitsMuon = myMu->innerTrack()->hitPattern().numberOfValidMuonHits();
86  ptError = myMu->innerTrack()->ptError();
87  }
88  else if ( const pat::Muon* myMu = dynamic_cast<const pat::Muon*>(&(*track)) ) {
89  hitsTk = myMu->innerTrack()->hitPattern().numberOfValidTrackerHits();
90  hitsMuon = myMu->innerTrack()->hitPattern().numberOfValidMuonHits();
91  ptError = myMu->innerTrack()->ptError();
92  }
93  else if (const reco::Track* myMu = dynamic_cast<const reco::Track*>(&(*track))){
94  hitsTk = myMu->hitPattern().numberOfValidTrackerHits();
95  hitsMuon = myMu->hitPattern().numberOfValidMuonHits();
96  ptError = myMu->ptError();
97  }
98 
99  MuScleFitMuon muon(mu,track->charge(),ptError,hitsTk,hitsMuon);
100 
101  if (debug_>0) {
102  std::cout<<"[TestCorrection::fillMuonCollection] after MuScleFitMuon initialization"<<std::endl;
103  std::cout<<" muon = "<<muon<<std::endl;
104  }
105 
106  muons.push_back(muon);
107  }
108  return muons;
109  }
T sqrt(T t)
Definition: SSEVec.h:18
const int mu
Definition: Constants.h:22
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
Analysis-level muon class.
Definition: Muon.h:49
virtual void TestCorrection::initialize ( const edm::EventSetup )
privatevirtual
void TestCorrection::initialize ( const edm::EventSetup )
privatevirtual

Definition at line 320 of file TestCorrection.cc.

References DEFINE_FWK_MODULE, and MuScleFitBase::readProbabilityDistributionsFromFile().

Referenced by analyze().

321 {
322  // Read the pdf from root file. They are used by massProb when finding the muon pair, needed
323  // for the mass histograms.
325 }
void readProbabilityDistributionsFromFile()
Read probability distributions from a local root file.

Member Data Documentation

std::auto_ptr< BackgroundFunction > TestCorrection::background_
private

Definition at line 128 of file TestCorrection.cc.

Referenced by fillMuonCollection(), and TestCorrection().

TH1F * TestCorrection::correctedPt_
private

Definition at line 121 of file TestCorrection.cc.

Referenced by analyze(), fillMuonCollection(), TestCorrection(), and ~TestCorrection().

TProfile * TestCorrection::correctedPtVsEta_
private

Definition at line 122 of file TestCorrection.cc.

Referenced by analyze(), fillMuonCollection(), TestCorrection(), and ~TestCorrection().

std::auto_ptr< MomentumScaleCorrector > TestCorrection::corrector_
private

Definition at line 126 of file TestCorrection.cc.

Referenced by correctMuon(), fillMuonCollection(), and TestCorrection().

int TestCorrection::eventCounter_
private

Definition at line 124 of file TestCorrection.cc.

Referenced by analyze(), fillMuonCollection(), TestCorrection(), and ~TestCorrection().

edm::EDGetTokenT<reco::MuonCollection> TestCorrection::glbMuonsToken_
private

Definition at line 130 of file TestCorrection.cc.

Referenced by analyze().

std::auto_ptr< ResolutionFunction > TestCorrection::resolution_
private

Definition at line 127 of file TestCorrection.cc.

Referenced by fillMuonCollection(), and TestCorrection().

edm::EDGetTokenT<reco::TrackCollection> TestCorrection::saMuonsToken_
private

Definition at line 131 of file TestCorrection.cc.

Referenced by analyze().

edm::EDGetTokenT<reco::TrackCollection> TestCorrection::tracksToken_
private

Definition at line 132 of file TestCorrection.cc.

Referenced by analyze().

TH1F * TestCorrection::uncorrectedPt_
private

Definition at line 119 of file TestCorrection.cc.

Referenced by analyze(), fillMuonCollection(), TestCorrection(), and ~TestCorrection().

TProfile * TestCorrection::uncorrectedPtVsEta_
private

Definition at line 120 of file TestCorrection.cc.

Referenced by analyze(), fillMuonCollection(), TestCorrection(), and ~TestCorrection().