CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TestResolution.cc
Go to the documentation of this file.
1 #ifndef TESTRESOLUTION_CC
2 #define TESTRESOLUTION_CC
3 
4 #include "TestResolution.h"
5 
6 #include "TCanvas.h"
7 #include "TLegend.h"
8 
9 //
10 // constants, enums and typedefs
11 //
12 
13 //
14 // static data member definitions
15 //
16 
17 //
18 // constructors and destructor
19 //
21  theMuonLabel_( iConfig.getParameter<edm::InputTag>( "MuonLabel" ) ),
22  glbMuonsToken_(mayConsume<reco::MuonCollection>(theMuonLabel_)),
23  saMuonsToken_(mayConsume<reco::TrackCollection>(theMuonLabel_)),
24  tracksToken_(mayConsume<reco::TrackCollection>(theMuonLabel_)),
25  theMuonType_( iConfig.getParameter<int>( "MuonType" ) ),
26  theRootFileName_( iConfig.getUntrackedParameter<std::string>("OutputFileName") )
27 {
28  //now do what ever initialization is needed
29  outputFile_ = new TFile(theRootFileName_.c_str(), "RECREATE");
30  outputFile_->cd();
31  sigmaPt_ = new TProfile("sigmaPtOverPt", "sigmaPt/Pt vs muon Pt", 1000, 0, 100);
32  eventCounter_ = 0;
33  // Create the corrector and set the parameters
34  resolutionFunction_.reset(new ResolutionFunction( iConfig.getUntrackedParameter<std::string>("ResolutionsIdentifier") ) );
35  std::cout << "resolutionFunction_ = " << &*resolutionFunction_ << std::endl;
36 }
37 
38 
40 {
41  outputFile_->cd();
42  TCanvas canvas("sigmaPtOverPt","sigmaPt/Pt vs muon Pt", 1000, 800);
43  canvas.cd();
44  sigmaPt_->GetXaxis()->SetTitle("Pt(GeV)");
45 // TLegend * legend = new TLegend(0.7,0.71,0.98,1.);
46 // legend->SetTextSize(0.02);
47 // legend->SetFillColor(0); // Have a white background
48 // legend->AddEntry(uncorrectedPt_, "original pt");
49 // legend->AddEntry(correctedPt_, "corrected pt");
50  sigmaPt_->Draw();
51 // legend->Draw("same");
52 
53  canvas.Write();
54  sigmaPt_->Write();
55  outputFile_->Close();
56 
57  std::cout << "Total analyzed events = " << eventCounter_ << std::endl;
58 }
59 
60 
61 //
62 // member functions
63 //
64 
65 // ------------ method called to for each event ------------
67  using namespace edm;
68 
69  ++eventCounter_;
70  if ( eventCounter_%100 == 0 ) {
71  std::cout << "Event number " << eventCounter_ << std::endl;
72  }
73 
74  // Take the reco-muons, depending on the type selected in the cfg
75  // --------------------------------------------------------------
76 
77  std::vector<reco::LeafCandidate> muons;
78 
79  if (theMuonType_==1) { // GlobalMuons
81  iEvent.getByToken (glbMuonsToken_, glbMuons);
82  muons = fillMuonCollection(*glbMuons);
83  }
84  else if (theMuonType_==2) { // StandaloneMuons
86  iEvent.getByToken (saMuonsToken_, saMuons);
87  muons = fillMuonCollection(*saMuons);
88  }
89  else if (theMuonType_==3) { // Tracker tracks
91  iEvent.getByToken (tracksToken_, tracks);
92  muons = fillMuonCollection(*tracks);
93  }
94 
95  // Loop on the recMuons
96  std::vector<reco::LeafCandidate>::const_iterator recMuon = muons.begin();
97  for ( ; recMuon!=muons.end(); ++recMuon ) {
98 
99  // Fill the histogram with uncorrected pt values
100  sigmaPt_->Fill(resolutionFunction_->sigmaPt(*recMuon, 0), recMuon->pt());
101 
102  }
103 }
104 
105 //define this as a plug-in
107 
108 #endif // TESTRESOLUTION_CC
T getUntrackedParameter(std::string const &, T const &) const
TestResolution(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::TrackCollection > saMuonsToken_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::string theRootFileName_
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
def canvas
Definition: svgfig.py:481
edm::EDGetTokenT< reco::TrackCollection > tracksToken_
TFile * outputFile_
int iEvent
Definition: GenABIO.cc:230
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::auto_ptr< ResolutionFunction > resolutionFunction_
std::vector< reco::LeafCandidate > fillMuonCollection(const std::vector< T > &tracks)
tuple tracks
Definition: testEve_cfg.py:39
TProfile * sigmaPt_
tuple muons
Definition: patZpeak.py:38
edm::EDGetTokenT< reco::MuonCollection > glbMuonsToken_
tuple cout
Definition: gather_cfg.py:121