CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
IsolatedParticlesGeneratedJets.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: IsolatedParticlesGeneratedJets
4 // Class: IsolatedParticlesGeneratedJets
5 //
13 //
14 // Original Author: Seema Sharma
15 // Created: Thu Mar 4 18:52:02 CST 2010
16 //
17 //
18 
21 
23 
24  debug = iConfig.getUntrackedParameter<bool> ("Debug", false);
25  tok_jets_ = consumes<reco::GenJetCollection>(iConfig.getParameter<edm::InputTag>("JetSource"));
26  tok_parts_ = consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("ParticleSource"));
27 }
28 
29 
31 
32 }
33 
35 
36  //using namespace edm;
38 
39  //=== genJet information
41  iEvent.getByToken(tok_jets_, genJets);
42 
43  //=== genJet information
45  iEvent.getByToken(tok_parts_, genParticles);
46 
47  JetMatchingTools jetMatching (iEvent, consumesCollector());
48  std::vector <std::vector <const reco::GenParticle*> > genJetConstituents (genJets->size());
49 
50  int njets = 0;
51  for (unsigned iGenJet = 0; iGenJet < genJets->size(); ++iGenJet) {
52  const reco::GenJet& genJet = (*genJets) [iGenJet];
53 
54  double genJetE = genJet.energy();
55  double genJetPt = genJet.pt();
56  double genJetEta = genJet.eta();
57  double genJetPhi = genJet.phi();
58 
59  if( genJetPt> 30.0 && std::abs(genJetEta)<3.0 ) {
60 
61  njets++;
62 
63  std::vector <const reco::GenParticle*> genJetConstituents = jetMatching.getGenParticles ((*genJets) [iGenJet]);
64 
65  std::vector<double> v_trkP, v_trkPt, v_trkEta, v_trkPhi, v_trkPdg, v_trkCharge;
66 
67  if(debug) std::cout<<"Jet(pt,Eta,Phi) "<<genJetPt<<" "<<genJetEta<<" "<<genJetPhi <<std::endl;
68  for(unsigned int ic=0; ic<genJetConstituents.size(); ic++) {
69 
70  if(debug) {
71  std::cout << "p,pt,eta,phi "<<genJetConstituents[ic]->p()<<" "<<genJetConstituents[ic]->pt()
72  <<" "<<genJetConstituents[ic]->eta()<<" "<<genJetConstituents[ic]->phi()
73  <<std::endl;
74  }
75 
76  v_trkP.push_back(genJetConstituents[ic]->p());
77  v_trkPt.push_back(genJetConstituents[ic]->pt());
78  v_trkEta.push_back(genJetConstituents[ic]->eta());
79  v_trkPhi.push_back(genJetConstituents[ic]->phi());
80  v_trkPdg.push_back(genJetConstituents[ic]->pdgId());
81  v_trkCharge.push_back(genJetConstituents[ic]->charge());
82 
83  } //loop over genjet constituents
84 
85  t_gjetE ->push_back(genJetE );
86  t_gjetPt ->push_back(genJetPt );
87  t_gjetEta ->push_back(genJetEta);
88  t_gjetPhi ->push_back(genJetPhi);
89 
90  t_jetTrkP ->push_back(v_trkP );
91  t_jetTrkPt ->push_back(v_trkPt );
92  t_jetTrkEta ->push_back(v_trkEta);
93  t_jetTrkPhi ->push_back(v_trkPhi);
94  t_jetTrkPdg ->push_back(v_trkPdg);
95  t_jetTrkCharge ->push_back(v_trkCharge);
96 
97  } // if jetPt>30
98 
99  } //loop over genjets
100 
101  t_gjetN->push_back(njets);
102 
103  unsigned int indx = 0;
104  for(reco::GenParticleCollection::const_iterator ig = genParticles->begin(); ig!= genParticles->end(); ++ig,++indx) {
105 
106  if (debug)
107  std::cout << "Track " << indx << " Status " << ig->status() << " charge "
108  << ig->charge() << " pdgId " << ig->pdgId() << " mass "
109  << ig->mass() << " P " << ig->momentum() << " E "
110  << ig->energy() << " Origin " << ig->vertex() << std::endl;
111  }
112 
113 
114  tree->Fill();
115 }
116 
118 
119  BookHistograms();
120 
121 }
122 
124  t_gjetN ->clear();
125  t_gjetE ->clear();
126  t_gjetPt ->clear();
127  t_gjetEta ->clear();
128  t_gjetPhi ->clear();
129 
130  t_jetTrkP ->clear();
131  t_jetTrkPt ->clear();
132  t_jetTrkEta ->clear();
133  t_jetTrkPhi ->clear();
134  t_jetTrkPdg ->clear();
135  t_jetTrkCharge ->clear();
136 }
137 
139 
140  tree = fs->make<TTree>("tree", "tree");
141 
142  t_gjetN = new std::vector<int> ();
143  t_gjetE = new std::vector<double>();
144  t_gjetPt = new std::vector<double>();
145  t_gjetEta = new std::vector<double>();
146  t_gjetPhi = new std::vector<double>();
147 
148  t_jetTrkP = new std::vector<std::vector<double> >();
149  t_jetTrkPt = new std::vector<std::vector<double> >();
150  t_jetTrkEta = new std::vector<std::vector<double> >();
151  t_jetTrkPhi = new std::vector<std::vector<double> >();
152  t_jetTrkPdg = new std::vector<std::vector<double> >();
153  t_jetTrkCharge = new std::vector<std::vector<double> >();
154 
155  tree->Branch("t_gjetN", "vector<int>", &t_gjetN);
156  tree->Branch("t_gjetE", "vector<double>", &t_gjetE);
157  tree->Branch("t_gjetPt", "vector<double>", &t_gjetPt);
158  tree->Branch("t_gjetEta", "vector<double>", &t_gjetEta);
159  tree->Branch("t_gjetPhi", "vector<double>", &t_gjetPhi);
160 
161  tree->Branch("t_jetTrkP", "vector<vector<double> >", &t_jetTrkP);
162  tree->Branch("t_jetTrkPt", "vector<vector<double> >", &t_jetTrkPt);
163  tree->Branch("t_jetTrkEta", "vector<vector<double> >", &t_jetTrkEta);
164  tree->Branch("t_jetTrkPhi", "vector<vector<double> >", &t_jetTrkPhi);
165  tree->Branch("t_jetTrkPdg", "vector<vector<double> >", &t_jetTrkPdg);
166  tree->Branch("t_jetTrkCharge", "vector<vector<double> >", &t_jetTrkCharge);
167 
168 }
169 
171 }
172 
173 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< std::vector< double > > * t_jetTrkEta
virtual float pt() const
transverse momentum
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual float phi() const
momentum azimuthal angle
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< reco::GenParticleCollection > tok_parts_
T eta() const
std::vector< std::vector< double > > * t_jetTrkPdg
double charge(const std::vector< uint8_t > &Ampls)
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:230
std::vector< std::vector< double > > * t_jetTrkPhi
virtual float eta() const
momentum pseudorapidity
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:24
std::vector< std::vector< double > > * t_jetTrkPt
std::vector< std::vector< double > > * t_jetTrkCharge
std::vector< const reco::GenParticle * > getGenParticles(const reco::CaloJet &fJet, bool fVerbose=true)
GenParticles for CaloJet.
edm::EDGetTokenT< reco::GenJetCollection > tok_jets_
tuple cout
Definition: gather_cfg.py:121
std::vector< std::vector< double > > * t_jetTrkP
virtual void analyze(const edm::Event &, const edm::EventSetup &)
IsolatedParticlesGeneratedJets(const edm::ParameterSet &)
Definition: DDAxes.h:10