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