CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
BjetAnalysis Class Reference
Inheritance diagram for BjetAnalysis:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
 BjetAnalysis (const edm::ParameterSet &cfg)
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

double chi2Cut_
 
double etaMax_
 
double etaMin_
 
TH1D * h_GlbMuChi2_
 
TH1D * h_GlbMuDxy_
 
TH1D * h_GlbMuNofHitsGlbMu_
 
TH1D * h_TrkMuNofHitsGlbMu_
 
double massMax_
 
double massMin_
 
std::vector< unsigned int > matched_
 
int nHitCut_
 
double ptMin_
 
edm::InputTag src_
 
double trkIso_
 
std::vector< unsigned int > unMatched_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- 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::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 6 of file BjetAnalyzer.cc.

Constructor & Destructor Documentation

BjetAnalysis::BjetAnalysis ( const edm::ParameterSet cfg)

Definition at line 37 of file BjetAnalyzer.cc.

References h_GlbMuChi2_, h_GlbMuDxy_, h_GlbMuNofHitsGlbMu_, h_TrkMuNofHitsGlbMu_, TFileDirectory::make(), and TFileDirectory::mkdir().

37  :
38  src_(cfg.getParameter<InputTag>("src")),
39  ptMin_(cfg.getUntrackedParameter<double>("ptMin")),
40  massMin_(cfg.getUntrackedParameter<double>("massMin")),
41  massMax_(cfg.getUntrackedParameter<double>("massMax")),
42  etaMin_(cfg.getUntrackedParameter<double>("etaMin")),
43  etaMax_(cfg.getUntrackedParameter<double>("etaMax")),
44  trkIso_(cfg.getUntrackedParameter<double>("trkIso")),
45  chi2Cut_(cfg.getUntrackedParameter<double>("chi2Cut")),
46  nHitCut_(cfg.getUntrackedParameter<int>("nHitCut"))
47 {
49  TFileDirectory trackEffDir = fs->mkdir("QualityOfGlbMu");
50  h_GlbMuNofHitsGlbMu_= trackEffDir.make<TH1D>("# of Hits of GlobalMuon", "# of Hits of GlobalMuon", 100, 0, 100);
51  h_TrkMuNofHitsGlbMu_= trackEffDir.make<TH1D>("# of Hits of TrackerMuon", "# of Hits of TrackerMuon", 100, 0, 100);
52  h_GlbMuChi2_= trackEffDir.make<TH1D>("chi2 of GlobalMuon", "chi2 of GlobalMuon", 100,0,10);
53  h_GlbMuDxy_= trackEffDir.make<TH1D>("Dxy of GlobalMuon", "Dxy of GlobalMuon", 1000,-5.,5.);
54 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
double massMax_
Definition: BjetAnalyzer.cc:14
double chi2Cut_
Definition: BjetAnalyzer.cc:14
TH1D * h_GlbMuDxy_
Definition: BjetAnalyzer.cc:16
double massMin_
Definition: BjetAnalyzer.cc:14
double trkIso_
Definition: BjetAnalyzer.cc:14
double etaMax_
Definition: BjetAnalyzer.cc:14
TH1D * h_GlbMuChi2_
Definition: BjetAnalyzer.cc:16
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
T * make() const
make new ROOT object
TH1D * h_TrkMuNofHitsGlbMu_
Definition: BjetAnalyzer.cc:16
edm::InputTag src_
Definition: BjetAnalyzer.cc:12
double etaMin_
Definition: BjetAnalyzer.cc:14
TH1D * h_GlbMuNofHitsGlbMu_
Definition: BjetAnalyzer.cc:16

Member Function Documentation

void BjetAnalysis::analyze ( const edm::Event evt,
const edm::EventSetup  
)
virtual

Implements edm::EDAnalyzer.

Definition at line 61 of file BjetAnalyzer.cc.

References gather_cfg::cout, edm::Event::getByLabel(), i, edm::Handle< T >::product(), and edm::AssociationVector< KeyRefProd, CVal, KeyRef, SizeType, KeyReferenceHelper >::size().

61  {
62 
63 
64  // Get b tag information
66  evt.getByLabel("trackCountingHighEffBJetTags", bTagHandle);
67  const reco::JetTagCollection & bTags = *(bTagHandle.product());
68 
69  // Loop over jets and study b tag info.
70  for (unsigned int i = 0; i != bTags.size(); ++i) {
71  cout<<" Jet "<< i
72  <<" has b tag discriminator (trackCountingHighEffBJetTags)= "<<bTags[i].second
73  << " and jet Pt = "<<bTags[i].first->pt()<<endl;
74  }
75 
76  // Get b tag information
78  evt.getByLabel("jetProbabilityBJetTags", bTagHandle2);
79  const reco::JetTagCollection & bTags2 = *(bTagHandle2.product());
80 
81  // Loop over jets and study b tag info.
82  for (unsigned int i = 0; i != bTags2.size(); ++i) {
83  cout<<" Jet "<< i
84  <<" has b tag discriminator (jetProbabilityBJetTags) = "<<bTags2[i].second
85  << " and jet Pt = "<<bTags2[i].first->pt()<<endl;
86  }
87 
88 
89  // Get b tag information
91  evt.getByLabel("jetBProbabilityBJetTags", bTagHandle3);
92  const reco::JetTagCollection & bTags3 = *(bTagHandle3.product());
93 
94  // Loop over jets and study b tag info.
95  for (unsigned int i = 0; i != bTags3.size(); ++i) {
96  cout<<" Jet "<< i
97  <<" has b tag discriminator (jetBProbabilityBJetTags) = "<<bTags3[i].second
98  << " and jet Pt = "<<bTags3[i].first->pt()<<endl;
99  }
100 
101 
102 
103 
104 
105 
106 
107 
108 
109 
110 }
int i
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
T const * product() const
Definition: Handle.h:74
tuple cout
Definition: gather_cfg.py:121
size_type size() const

Member Data Documentation

double BjetAnalysis::chi2Cut_
private

Definition at line 14 of file BjetAnalyzer.cc.

double BjetAnalysis::etaMax_
private

Definition at line 14 of file BjetAnalyzer.cc.

double BjetAnalysis::etaMin_
private

Definition at line 14 of file BjetAnalyzer.cc.

TH1D * BjetAnalysis::h_GlbMuChi2_
private

Definition at line 16 of file BjetAnalyzer.cc.

Referenced by BjetAnalysis().

TH1D * BjetAnalysis::h_GlbMuDxy_
private

Definition at line 16 of file BjetAnalyzer.cc.

Referenced by BjetAnalysis().

TH1D* BjetAnalysis::h_GlbMuNofHitsGlbMu_
private

Definition at line 16 of file BjetAnalyzer.cc.

Referenced by BjetAnalysis().

TH1D * BjetAnalysis::h_TrkMuNofHitsGlbMu_
private

Definition at line 16 of file BjetAnalyzer.cc.

Referenced by BjetAnalysis().

double BjetAnalysis::massMax_
private

Definition at line 14 of file BjetAnalyzer.cc.

double BjetAnalysis::massMin_
private

Definition at line 14 of file BjetAnalyzer.cc.

std::vector<unsigned int> BjetAnalysis::matched_
private

Definition at line 13 of file BjetAnalyzer.cc.

int BjetAnalysis::nHitCut_
private

Definition at line 15 of file BjetAnalyzer.cc.

double BjetAnalysis::ptMin_
private

Definition at line 14 of file BjetAnalyzer.cc.

edm::InputTag BjetAnalysis::src_
private

Definition at line 12 of file BjetAnalyzer.cc.

double BjetAnalysis::trkIso_
private

Definition at line 14 of file BjetAnalyzer.cc.

std::vector<unsigned int> BjetAnalysis::unMatched_
private

Definition at line 13 of file BjetAnalyzer.cc.