CMS 3D CMS Logo

SubEventGenJetProducer.cc
Go to the documentation of this file.
1 
4 
5 #include <memory>
6 
11 
12 using namespace std;
13 using namespace reco;
14 using namespace edm;
15 using namespace cms;
16 
17 namespace {
18  bool checkHydro(const reco::GenParticle* p) {
19  const Candidate* m1 = p->mother();
20  while (m1) {
21  int pdg = abs(m1->pdgId());
22  int st = m1->status();
23  LogDebug("SubEventMothers") << "Pdg ID : " << pdg << endl;
24  if (st == 3 || pdg < 9 || pdg == 21) {
25  LogDebug("SubEventMothers") << "Sub-Collision Found! Pdg ID : " << pdg << endl;
26  return false;
27  }
28  const Candidate* m = m1->mother();
29  m1 = m;
30  if (!m1)
31  LogDebug("SubEventMothers") << "No Mother, particle is : " << pdg << " with status " << st << endl;
32  }
33  // return true;
34  return true; // Debugging - to be changed
35  }
36 } // namespace
37 
38 SubEventGenJetProducer::SubEventGenJetProducer(edm::ParameterSet const& conf) : VirtualJetProducer(conf) {
39  // mapSrc_ = conf.getParameter<edm::InputTag>( "srcMap");
40  ignoreHydro_ = conf.getUntrackedParameter<bool>("ignoreHydro", true);
41  produces<reco::BasicJetCollection>();
42  // the subjet collections are set through the config file in the "jetCollInstanceName" field.
43 
44  input_cand_token_ = consumes<reco::CandidateView>(src_);
45 }
46 
48  std::vector<edm::Ptr<reco::Candidate>>::const_iterator inBegin = inputs_.begin(), inEnd = inputs_.end(), i = inBegin;
49  for (; i != inEnd; ++i) {
50  reco::CandidatePtr input = inputs_[i - inBegin];
51  if (edm::isNotFinite(input->pt()))
52  continue;
53  if (input->et() < inputEtMin_)
54  continue;
55  if (input->energy() < inputEMin_)
56  continue;
58  continue;
59 
60  edm::Ptr<reco::Candidate> p = inputs_[i - inBegin];
61  const GenParticle* pref = dynamic_cast<const GenParticle*>(p.get());
62  int subevent = pref->collisionId();
63  LogDebug("SubEventContainers") << "SubEvent is : " << subevent << endl;
64  LogDebug("SubEventContainers") << "SubSize is : " << subInputs_.size() << endl;
65 
66  if (subevent >= (int)subInputs_.size()) {
67  hydroTag_.resize(subevent + 1, -1);
68  subInputs_.resize(subevent + 1);
69  LogDebug("SubEventContainers") << "SubSize is : " << subInputs_.size() << endl;
70  LogDebug("SubEventContainers") << "HydroTagSize is : " << hydroTag_.size() << endl;
71  }
72 
73  LogDebug("SubEventContainers") << "HydroTag is : " << hydroTag_[subevent] << endl;
74  if (hydroTag_[subevent] != 0)
75  hydroTag_[subevent] = (int)checkHydro(pref);
76 
77  subInputs_[subevent].push_back(fastjet::PseudoJet(input->px(), input->py(), input->pz(), input->energy()));
78 
79  subInputs_[subevent].back().set_user_index(i - inBegin);
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  auto jets = std::make_unique<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])
116  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(std::move(jets));
127  return;
128 }
129 
131  // run algorithm
132  fjJets_.clear();
133 
134  fjClusterSeq_ = std::make_shared<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  GenJet jet;
141  const fastjet::PseudoJet& fjJet = fjJets_[ijet];
142 
143  std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJet));
144 
145  std::vector<CandidatePtr> constituents = getConstituents(fjConstituents);
146 
147  double px = fjJet.px();
148  double py = fjJet.py();
149  double pz = fjJet.pz();
150  double E = fjJet.E();
151  double jetArea = 0.0;
152  double pu = 0.;
153 
154  writeSpecific(jet, Particle::LorentzVector(px, py, pz, E), vertex_, constituents, iSetup);
155 
156  jet.setJetArea(jetArea);
157  jet.setPileup(pu);
158 
159  subJets_->push_back(jet);
160  }
161 }
162 
cms::SubEventGenJetProducer
Definition: SubEventGenJetProducer.h:18
reco::writeSpecific
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
mps_fire.i
i
Definition: mps_fire.py:355
input
static const std::string input
Definition: EdmProvDump.cc:48
reco::GenJet
Jets made from MC generator particles.
Definition: GenJet.h:23
reco::GenParticle::collisionId
int collisionId() const
Definition: GenParticle.h:36
cms::SubEventGenJetProducer::runAlgorithm
void runAlgorithm(edm::Event &, const edm::EventSetup &) override
Definition: SubEventGenJetProducer.cc:130
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
reco::GenParticle
Definition: GenParticle.h:21
cms::SubEventGenJetProducer::input_cand_token_
edm::EDGetTokenT< reco::CandidateView > input_cand_token_
Definition: SubEventGenJetProducer.h:38
multPhiCorr_741_25nsDY_cfi.py
py
Definition: multPhiCorr_741_25nsDY_cfi.py:12
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
cms::SubEventGenJetProducer::hydroTag_
std::vector< int > hydroTag_
Definition: SubEventGenJetProducer.h:28
reco::Candidate::status
virtual int status() const =0
status word
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
reco::Candidate::mother
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
edm::Handle
Definition: AssociativeIterator.h:50
GenParticle
Definition: GenParticle.py:1
VirtualJetProducer::vertex_
reco::Particle::Point vertex_
Definition: VirtualJetProducer.h:184
VirtualJetProducer::fjClusterSeq_
ClusterSequencePtr fjClusterSeq_
Definition: VirtualJetProducer.h:185
cms::SubEventGenJetProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: SubEventGenJetProducer.cc:83
GenParticle.h
VirtualJetProducer::inputEtMin_
double inputEtMin_
Definition: VirtualJetProducer.h:153
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
cms::SubEventGenJetProducer::inputTowers
void inputTowers() override
Definition: SubEventGenJetProducer.cc:47
cms::SubEventGenJetProducer::subJets_
std::vector< reco::GenJet > * subJets_
Definition: SubEventGenJetProducer.h:27
nsub
const int nsub
Definition: CMTRawAnalyzer.h:421
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
VirtualJetProducer
Definition: VirtualJetProducer.h:35
VirtualJetProducer::inputs_
std::vector< edm::Ptr< reco::Candidate > > inputs_
Definition: VirtualJetProducer.h:183
edm::View::size
size_type size() const
VirtualJetProducer::src_
edm::InputTag src_
Definition: VirtualJetProducer.h:148
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
SubEventGenJetProducer.h
VirtualJetProducer::isAnomalousTower
virtual bool isAnomalousTower(reco::CandidatePtr input)
Definition: VirtualJetProducer.cc:499
VirtualJetProducer::jetPtMin_
double jetPtMin_
Definition: VirtualJetProducer.h:155
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
cms::SubEventGenJetProducer::ignoreHydro_
bool ignoreHydro_
Definition: SubEventGenJetProducer.h:30
cms::SubEventGenJetProducer::nSubParticles_
std::vector< int > nSubParticles_
Definition: SubEventGenJetProducer.h:29
edm::EventSetup
Definition: EventSetup.h:57
reco::Candidate::pdgId
virtual int pdgId() const =0
PDG identifier.
edm::Ptr< Candidate >
reco::Candidate
Definition: Candidate.h:27
reco::JetExtendedAssociation::LorentzVector
math::PtEtaPhiELorentzVectorF LorentzVector
Definition: JetExtendedAssociation.h:25
multPhiCorr_741_25nsDY_cfi.px
px
Definition: multPhiCorr_741_25nsDY_cfi.py:10
cms::SubEventGenJetProducer::subInputs_
std::vector< std::vector< fastjet::PseudoJet > > subInputs_
Definition: SubEventGenJetProducer.h:26
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
VirtualJetProducer::getConstituents
virtual std::vector< reco::CandidatePtr > getConstituents(const std::vector< fastjet::PseudoJet > &fjConstituents)
Definition: VirtualJetProducer.cc:522
isFinite.h
metsig::jet
Definition: SignAlgoResolutions.h:47
VirtualJetProducer::fjJets_
std::vector< fastjet::PseudoJet > fjJets_
Definition: VirtualJetProducer.h:192
VirtualJetProducer::fjJetDefinition_
JetDefPtr fjJetDefinition_
Definition: VirtualJetProducer.h:186
JetSpecific.h
pdg
Definition: pdg_functions.h:28
Exception.h
muons2muons_cfi.pu
pu
Definition: muons2muons_cfi.py:31
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
VirtualJetProducer::fjInputs_
std::vector< fastjet::PseudoJet > fjInputs_
Definition: VirtualJetProducer.h:191
edm::Event
Definition: Event.h:73
edm::View::ptrAt
Ptr< value_type > ptrAt(size_type i) const
VirtualJetProducer::inputEMin_
double inputEMin_
Definition: VirtualJetProducer.h:154
cms
Namespace of DDCMS conversion namespace.
Definition: ProducerAnalyzer.cc:21