CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SubEventGenJetProducer.cc
Go to the documentation of this file.
1 
8 
9 using namespace std;
10 using namespace reco;
11 using namespace edm;
12 using namespace cms;
13 
14 namespace {
15  bool checkHydro(const reco::GenParticle * p){
16  const Candidate* m1 = p->mother();
17  while(m1){
18  int pdg = abs(m1->pdgId());
19  int st = m1->status();
20  LogDebug("SubEventMothers")<<"Pdg ID : "<<pdg<<endl;
21  if(st == 3 || pdg < 9 || pdg == 21){
22  LogDebug("SubEventMothers")<<"Sub-Collision Found! Pdg ID : "<<pdg<<endl;
23  return false;
24  }
25  const Candidate* m = m1->mother();
26  m1 = m;
27  if(!m1) LogDebug("SubEventMothers")<<"No Mother, particle is : "<<pdg<<" with status "<<st<<endl;
28  }
29  // return true;
30  return true; // Debugging - to be changed
31  }
32 }
33 
34 SubEventGenJetProducer::SubEventGenJetProducer(edm::ParameterSet const& conf):
35  VirtualJetProducer( conf )
36 {
37  // mapSrc_ = conf.getParameter<edm::InputTag>( "srcMap");
38  ignoreHydro_ = conf.getUntrackedParameter<bool>("ignoreHydro", true);
39  produces<reco::BasicJetCollection>();
40  // the subjet collections are set through the config file in the "jetCollInstanceName" field.
41 
42  input_cand_token_ = consumes<reco::CandidateView>(src_);
43 
44 }
45 
46 
48 {
49  std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
50  inEnd = inputs_.end(), i = inBegin;
51  for (; i != inEnd; ++i ) {
52  reco::CandidatePtr input = inputs_[i - inBegin];
53  if (edm::isNotFinite(input->pt())) continue;
54  if (input->et() <inputEtMin_) continue;
55  if (input->energy()<inputEMin_) continue;
56  if (isAnomalousTower(input)) continue;
57 
58  edm::Ptr<reco::Candidate> p = inputs_[i - inBegin];
59  const GenParticle * pref = dynamic_cast<const GenParticle *>(p.get());
60  int subevent = pref->collisionId();
61  LogDebug("SubEventContainers")<<"SubEvent is : "<<subevent<<endl;
62  LogDebug("SubEventContainers")<<"SubSize is : "<<subInputs_.size()<<endl;
63 
64  if(subevent >= (int)subInputs_.size()){
65  hydroTag_.resize(subevent+1, -1);
66  subInputs_.resize(subevent+1);
67  LogDebug("SubEventContainers")<<"SubSize is : "<<subInputs_.size()<<endl;
68  LogDebug("SubEventContainers")<<"HydroTagSize is : "<<hydroTag_.size()<<endl;
69  }
70 
71  LogDebug("SubEventContainers")<<"HydroTag is : "<<hydroTag_[subevent]<<endl;
72  if(hydroTag_[subevent] != 0) hydroTag_[subevent] = (int)checkHydro(pref);
73 
74  subInputs_[subevent].push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
75  input->energy()));
76 
77  subInputs_[subevent].back().set_user_index(i - inBegin);
78 
79  }
80 
81 }
82 
84  LogDebug("VirtualJetProducer") << "Entered produce\n";
85 
86  fjJets_.clear();
87  subInputs_.clear();
88  nSubParticles_.clear();
89  hydroTag_.clear();
90  inputs_.clear();
91 
92  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
94  iEvent.getByToken(input_cand_token_, inputsHandle);
95  for (size_t i = 0; i < inputsHandle->size(); ++i) {
96  inputs_.push_back(inputsHandle->ptrAt(i));
97  }
98  LogDebug("VirtualJetProducer") << "Got inputs\n";
99 
100  inputTowers();
101  // Convert candidates to fastjet::PseudoJets.
102  // Also correct to Primary Vertex. Will modify fjInputs_
103  // and use inputs_
104 
106 
107  std::auto_ptr<std::vector<GenJet> > jets(new std::vector<GenJet>() );
108  subJets_ = jets.get();
109 
110  LogDebug("VirtualJetProducer") << "Inputted towers\n";
111 
112  size_t nsub = subInputs_.size();
113 
114  for(size_t isub = 0; isub < nsub; ++isub){
115  if(ignoreHydro_ && hydroTag_[isub]) continue;
116  fjJets_.clear();
117  fjInputs_.clear();
118  fjInputs_ = subInputs_[isub];
119  runAlgorithm( iEvent, iSetup );
120  }
121 
122  //Finalize
123  LogDebug("SubEventJetProducer") << "Wrote jets\n";
124 
125  iEvent.put(jets);
126  return;
127 }
128 
130 {
131  // run algorithm
132  fjJets_.clear();
133 
134  fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
135  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
136 
137  using namespace reco;
138 
139  for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
140 
141  GenJet jet;
142  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
143 
144  std::vector<fastjet::PseudoJet> fjConstituents =
145  sorted_by_pt(fjClusterSeq_->constituents(fjJet));
146 
147  std::vector<CandidatePtr> constituents =
148  getConstituents(fjConstituents);
149 
150  double px=fjJet.px();
151  double py=fjJet.py();
152  double pz=fjJet.pz();
153  double E=fjJet.E();
154  double jetArea=0.0;
155  double pu=0.;
156 
157  writeSpecific( jet,
158  Particle::LorentzVector(px, py, pz, E),
159  vertex_,
160  constituents, iSetup);
161 
162  jet.setJetArea (jetArea);
163  jet.setPileup (pu);
164 
165  subJets_->push_back(jet);
166  }
167 }
168 
170 
171 
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
reco::Particle::Point vertex_
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
std::vector< fastjet::PseudoJet > fjJets_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
edm::EDGetTokenT< reco::CandidateView > input_cand_token_
virtual void setPileup(float fEnergy)
Set pileup energy contribution as calculated by algorithm.
Definition: Jet.h:108
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:143
virtual int status() const =0
status word
int collisionId() const
Definition: GenParticle.h:38
static std::string const input
Definition: EdmProvDump.cc:44
std::vector< reco::GenJet > * subJets_
virtual bool isAnomalousTower(reco::CandidatePtr input)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:103
std::vector< fastjet::PseudoJet > fjInputs_
int iEvent
Definition: GenABIO.cc:230
bool isNotFinite(T x)
Definition: isFinite.h:10
std::vector< edm::Ptr< reco::Candidate > > inputs_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
ClusterSequencePtr fjClusterSeq_
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:24
tuple conf
Definition: dbtoconf.py:185
virtual int pdgId() const =0
PDG identifier.
void produce(edm::Event &, const edm::EventSetup &)
std::vector< std::vector< fastjet::PseudoJet > > subInputs_
virtual const Candidate * mother(size_type=0) const
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
boost::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
math::PtEtaPhiELorentzVectorF LorentzVector
void runAlgorithm(edm::Event &, const edm::EventSetup &)
void writeSpecific(reco::CaloJet &jet, reco::Particle::LorentzVector const &p4, reco::Particle::Point const &point, std::vector< reco::CandidatePtr > const &constituents, edm::EventSetup const &c)
Definition: JetSpecific.cc:41