CMS 3D CMS Logo

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)
28  LogDebug("SubEventMothers") << "No Mother, particle is : " << pdg << " with status " << st << endl;
29  }
30  // return true;
31  return true; // Debugging - to be changed
32  }
33 } // namespace
34 
35 SubEventGenJetProducer::SubEventGenJetProducer(edm::ParameterSet const& conf) : VirtualJetProducer(conf) {
36  // mapSrc_ = conf.getParameter<edm::InputTag>( "srcMap");
37  ignoreHydro_ = conf.getUntrackedParameter<bool>("ignoreHydro", true);
38  produces<reco::BasicJetCollection>();
39  // the subjet collections are set through the config file in the "jetCollInstanceName" field.
40 
41  input_cand_token_ = consumes<reco::CandidateView>(src_);
42 }
43 
45  std::vector<edm::Ptr<reco::Candidate>>::const_iterator inBegin = inputs_.begin(), inEnd = inputs_.end(), i = inBegin;
46  for (; i != inEnd; ++i) {
47  reco::CandidatePtr input = inputs_[i - inBegin];
48  if (edm::isNotFinite(input->pt()))
49  continue;
50  if (input->et() < inputEtMin_)
51  continue;
52  if (input->energy() < inputEMin_)
53  continue;
54  if (isAnomalousTower(input))
55  continue;
56 
57  edm::Ptr<reco::Candidate> p = inputs_[i - inBegin];
58  const GenParticle* pref = dynamic_cast<const GenParticle*>(p.get());
59  int subevent = pref->collisionId();
60  LogDebug("SubEventContainers") << "SubEvent is : " << subevent << endl;
61  LogDebug("SubEventContainers") << "SubSize is : " << subInputs_.size() << endl;
62 
63  if (subevent >= (int)subInputs_.size()) {
64  hydroTag_.resize(subevent + 1, -1);
65  subInputs_.resize(subevent + 1);
66  LogDebug("SubEventContainers") << "SubSize is : " << subInputs_.size() << endl;
67  LogDebug("SubEventContainers") << "HydroTagSize is : " << hydroTag_.size() << endl;
68  }
69 
70  LogDebug("SubEventContainers") << "HydroTag is : " << hydroTag_[subevent] << endl;
71  if (hydroTag_[subevent] != 0)
72  hydroTag_[subevent] = (int)checkHydro(pref);
73 
74  subInputs_[subevent].push_back(fastjet::PseudoJet(input->px(), input->py(), input->pz(), input->energy()));
75 
76  subInputs_[subevent].back().set_user_index(i - inBegin);
77  }
78 }
79 
81  LogDebug("VirtualJetProducer") << "Entered produce\n";
82 
83  fjJets_.clear();
84  subInputs_.clear();
85  nSubParticles_.clear();
86  hydroTag_.clear();
87  inputs_.clear();
88 
89  // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
91  iEvent.getByToken(input_cand_token_, inputsHandle);
92  for (size_t i = 0; i < inputsHandle->size(); ++i) {
93  inputs_.push_back(inputsHandle->ptrAt(i));
94  }
95  LogDebug("VirtualJetProducer") << "Got inputs\n";
96 
97  inputTowers();
98  // Convert candidates to fastjet::PseudoJets.
99  // Also correct to Primary Vertex. Will modify fjInputs_
100  // and use inputs_
101 
103 
104  auto jets = std::make_unique<std::vector<GenJet>>();
105  subJets_ = jets.get();
106 
107  LogDebug("VirtualJetProducer") << "Inputted towers\n";
108 
109  size_t nsub = subInputs_.size();
110 
111  for (size_t isub = 0; isub < nsub; ++isub) {
112  if (ignoreHydro_ && hydroTag_[isub])
113  continue;
114  fjJets_.clear();
115  fjInputs_.clear();
116  fjInputs_ = subInputs_[isub];
117  runAlgorithm(iEvent, iSetup);
118  }
119 
120  //Finalize
121  LogDebug("SubEventJetProducer") << "Wrote jets\n";
122 
123  iEvent.put(std::move(jets));
124  return;
125 }
126 
128  // run algorithm
129  fjJets_.clear();
130 
131  fjClusterSeq_ = ClusterSequencePtr(new fastjet::ClusterSequence(fjInputs_, *fjJetDefinition_));
132  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
133 
134  using namespace reco;
135 
136  for (unsigned int ijet = 0; ijet < fjJets_.size(); ++ijet) {
137  GenJet jet;
138  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
139 
140  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJet));
141 
142  std::vector<CandidatePtr> constituents = getConstituents(fjConstituents);
143 
144  double px = fjJet.px();
145  double py = fjJet.py();
146  double pz = fjJet.pz();
147  double E = fjJet.E();
148  double jetArea = 0.0;
149  double pu = 0.;
150 
151  writeSpecific(jet, Particle::LorentzVector(px, py, pz, E), vertex_, constituents, iSetup);
152 
153  jet.setJetArea(jetArea);
154  jet.setPileup(pu);
155 
156  subJets_->push_back(jet);
157  }
158 }
159 
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
virtual double pz() const =0
z coordinate of momentum vector
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
reco::Particle::Point vertex_
void runAlgorithm(edm::Event &, const edm::EventSetup &) override
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:525
edm::EDGetTokenT< reco::CandidateView > input_cand_token_
Ptr< value_type > ptrAt(size_type i) const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
virtual void setPileup(float fEnergy)
Set pileup energy contribution as calculated by algorithm.
Definition: Jet.h:106
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:139
std::shared_ptr< fastjet::ClusterSequence > ClusterSequencePtr
size_type size() const
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
int collisionId() const
Definition: GenParticle.h:36
void produce(edm::Event &, const edm::EventSetup &) override
static std::string const input
Definition: EdmProvDump.cc:48
std::vector< reco::GenJet > * subJets_
virtual int status() const =0
status word
virtual bool isAnomalousTower(reco::CandidatePtr input)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:101
std::vector< fastjet::PseudoJet > fjInputs_
virtual double et() const =0
transverse energy
virtual double energy() const =0
energy
virtual double py() const =0
y coordinate of momentum vector
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< edm::Ptr< reco::Candidate > > inputs_
ClusterSequencePtr fjClusterSeq_
virtual int pdgId() const =0
PDG identifier.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23
Namespace of DDCMS conversion namespace.
virtual double pt() const =0
transverse momentum
std::vector< std::vector< fastjet::PseudoJet > > subInputs_
fixed size matrix
HLT enums.
virtual double px() const =0
x coordinate of momentum vector
def move(src, dest)
Definition: eostools.py:511
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...
math::PtEtaPhiELorentzVectorF LorentzVector
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:34