CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
TopHypothesisFWLiteAnalyzer.cc File Reference
#include <memory>
#include <string>
#include <cstdlib>
#include <sstream>
#include <fstream>
#include <iostream>
#include "FWCore/FWLite/interface/FWLiteEnabler.h"
#include "AnalysisDataFormats/TopObjects/interface/TtSemiLeptonicEvent.h"
#include "TopQuarkAnalysis/Examples/bin/NiceStyle.cc"
#include "TopQuarkAnalysis/Examples/interface/RootSystem.h"
#include "TopQuarkAnalysis/Examples/interface/RootHistograms.h"
#include "TopQuarkAnalysis/Examples/interface/RootPostScript.h"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

int main ( int  argc,
char *  argv[] 
)

CALO JETS


PF JETS

Definition at line 16 of file TopHypothesisFWLiteAnalyzer.cc.

References assert(), ecal_dqm_sourceclient-live_cfg::cerr, gather_cfg::cout, FWLiteEnabler::enable(), genEvt_(), TtSemiLeptonicEvent::hadronicDecayTop(), TtSemiLeptonicEvent::hadronicDecayW(), TtEvent::isHypoAvailable(), TtEvent::isHypoValid(), TtSemiLeptonicEvent::leptonicDecayTop(), TtSemiLeptonicEvent::leptonicDecayW(), reco::Candidate::mass(), visualization-live-secondInstance_cfg::outDir, GetRecoTauVFromDQM_MC_cff::outFile, reco::Candidate::pt(), setNiceStyle(), and AlCaHLTBitMon_QueryRunRegistry::string.

17 {
18  if( argc<4 ){
19  // -------------------------------------------------
20  std::cerr << "ERROR: Wrong number of arguments!" << std::endl
21  << "Please specify: * file name" << std::endl
22  << " * process name" << std::endl
23  << " * HypoClassKey" << std::endl
24  << "Example: TopHypothesisFWLiteAnalyzer ttSemiLepEvtBuilder.root TEST kGeom" << std::endl;
25  // -------------------------------------------------
26  return -1;
27  }
28 
29  // get HypoClassKey
30  std::string hypoClassKey = argv[3];
31 
32  // load framework libraries
33  gSystem->Load( "libFWCoreFWLite" );
35 
36  // set nice style for histograms
37  setNiceStyle();
38 
39  // define some histograms
40  TH1F* hadWPt_ = new TH1F("hadWPt", "p_{t} (W_{had}) [GeV]", 100, 0., 500.);
41  TH1F* hadWMass_ = new TH1F("hadWMass", "M (W_{had}) [GeV]", 50, 0., 150.);
42  TH1F* hadTopPt_ = new TH1F("hadTopPt", "p_{t} (t_{had}) [GeV]", 100, 0., 500.);
43  TH1F* hadTopMass_= new TH1F("hadTopMass", "M (t_{had}) [GeV]", 50, 50., 250.);
44 
45  TH1F* lepWPt_ = new TH1F("lepWPt", "p_{t} (W_{lep}) [GeV]", 100, 0., 500.);
46  TH1F* lepWMass_ = new TH1F("lepWMass", "M (W_{lep}) [GeV]", 50, 0., 150.);
47  TH1F* lepTopPt_ = new TH1F("lepTopPt", "p_{t} (t_{lep}) [GeV]", 100, 0., 500.);
48  TH1F* lepTopMass_= new TH1F("lepTopMass", "M (t_{lep}) [GeV]", 50, 50., 250.);
49 
50  // -------------------------------------------------
51  std::cout << "open file: " << argv[1] << std::endl;
52  // -------------------------------------------------
53  TFile* inFile = TFile::Open(argv[1]);
54  TTree* events_= 0;
55  if( inFile ) inFile->GetObject("Events", events_);
56  if( events_==0 ){
57  // -------------------------------------------------
58  std::cerr << "ERROR: Unable to retrieve TTree Events!" << std::endl
59  << "Either wrong file name or the tree doesn't exist." << std::endl;
60  // -------------------------------------------------
61  return -1;
62  }
63 
64  // acess branch of ttSemiLepEvent
65  char decayName[50];
66  sprintf(decayName, "recoGenParticles_decaySubset__%s.obj", argv[2]);
67  TBranch* decay_ = events_->GetBranch( decayName ); // referred to from within TtGenEvent class
68  assert( decay_ != 0 );
69  char genEvtName[50];
70  sprintf(genEvtName, "TtGenEvent_genEvt__%s.obj", argv[2]);
71  TBranch* genEvt_ = events_->GetBranch( genEvtName ); // referred to from within TtSemiLeptonicEvent class
72  assert( genEvt_ != 0 );
73  char semiLepEvtName[50];
74  sprintf(semiLepEvtName, "TtSemiLeptonicEvent_ttSemiLepEvent__%s.obj", argv[2]);
75  TBranch* semiLepEvt_ = events_->GetBranch( semiLepEvtName );
76  assert( semiLepEvt_ != 0 );
77 
78  // loop over events and fill histograms
79  int nevt = events_->GetEntries();
80  TtSemiLeptonicEvent semiLepEvt;
81  // -------------------------------------------------
82  std::cout << "start looping " << nevt << " events..." << std::endl;
83  // -------------------------------------------------
84  for(int evt=0; evt<nevt; ++evt){
85  // set branch address
86  semiLepEvt_->SetAddress( &semiLepEvt );
87  // get event
88  decay_ ->GetEntry( evt );
89  genEvt_ ->GetEntry( evt );
90  semiLepEvt_->GetEntry( evt );
91  events_ ->GetEntry( evt, 0 );
92 
93  // -------------------------------------------------
94  if(evt>0 && !(evt%10)) std::cout << " processing event: " << evt << std::endl;
95  // -------------------------------------------------
96 
97  // fill histograms
98  if( !semiLepEvt.isHypoAvailable(hypoClassKey) ){
99  std::cerr << "NonValidHyp:: " << "Hypothesis not available for this event" << std::endl;
100  continue;
101  }
102  if( !semiLepEvt.isHypoValid(hypoClassKey) ){
103  std::cerr << "NonValidHyp::" << "Hypothesis not valid for this event" << std::endl;
104  continue;
105  }
106 
107  const reco::Candidate* hadTop = semiLepEvt.hadronicDecayTop(hypoClassKey);
108  const reco::Candidate* hadW = semiLepEvt.hadronicDecayW (hypoClassKey);
109  const reco::Candidate* lepTop = semiLepEvt.leptonicDecayTop(hypoClassKey);
110  const reco::Candidate* lepW = semiLepEvt.leptonicDecayW (hypoClassKey);
111 
112  if(hadTop && hadW && lepTop && lepW){
113  hadWPt_ ->Fill( hadW->pt() );
114  hadWMass_ ->Fill( hadW->mass() );
115  hadTopPt_ ->Fill( hadTop->pt() );
116  hadTopMass_->Fill( hadTop->mass());
117 
118  lepWPt_ ->Fill( lepW->pt() );
119  lepWMass_ ->Fill( lepW->mass() );
120  lepTopPt_ ->Fill( lepTop->pt() );
121  lepTopMass_->Fill( lepTop->mass());
122  }
123  }
124  // -------------------------------------------------
125  std::cout << "close file" << std::endl;
126  // -------------------------------------------------
127  inFile->Close();
128 
129  // save histograms to file
130  TFile outFile( "analyzeHypothesis.root", "recreate" );
131  // strip the leading "k" from the hypoClassKey to build directory name
132  TString outDir = "analyze" + std::string(hypoClassKey, 1, hypoClassKey.length());
133  outFile.mkdir(outDir);
134  outFile.cd(outDir);
135  hadWPt_ ->Write( );
136  hadWMass_ ->Write( );
137  hadTopPt_ ->Write( );
138  hadTopMass_->Write( );
139  lepWPt_ ->Write( );
140  lepWMass_ ->Write( );
141  lepTopPt_ ->Write( );
142  lepTopMass_->Write( );
143  outFile.Close();
144 
145  // free allocated space
146  delete hadWPt_;
147  delete hadWMass_;
148  delete hadTopPt_;
149  delete hadTopMass_;
150  delete lepWPt_;
151  delete lepWMass_;
152  delete lepTopPt_;
153  delete lepTopMass_;
154 
155  return 0;
156 }
const reco::Candidate * leptonicDecayW(const std::string &key, const unsigned &cmb=0) const
get leptonic W of the given hypothesis
void setNiceStyle()
Definition: NiceStyle.cc:3
bool isHypoAvailable(const std::string &key, const unsigned &cmb=0) const
Definition: TtEvent.h:60
const reco::Candidate * hadronicDecayW(const std::string &key, const unsigned &cmb=0) const
get hadronic W of the given hypothesis
virtual double pt() const =0
transverse momentum
virtual double mass() const =0
mass
assert(m_qm.get())
const reco::Candidate * hadronicDecayTop(const std::string &key, const unsigned &cmb=0) const
get hadronic top of the given hypothesis
Class derived from the TtEvent for the semileptonic decay channel.
static void enable()
enable automatic library loading
const reco::Candidate * leptonicDecayTop(const std::string &key, const unsigned &cmb=0) const
get leptonic top of the given hypothesis
tuple argc
Definition: dir2webdir.py:38
bool isHypoValid(const std::string &key, const unsigned &cmb=0) const
check if hypothesis &#39;cmb&#39; within the hypothesis class was valid; if not it lead to an empty Composite...
Definition: TtEvent.h:64
tuple cout
Definition: gather_cfg.py:121
genEvt_(cfg.getParameter< edm::InputTag >("genEvent"))