CMS 3D CMS Logo

SoftKillerProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CommonTools/PileupAlgos
4 // Class: SoftKillerProducer
5 //
13 //
14 // Original Author: Salvatore Rappoccio
15 // Created: Wed, 22 Oct 2014 15:14:20 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 #include <vector>
22 
23 // user include files
26 
29 
36 #include "fastjet/contrib/SoftKiller.hh"
37 
38 //
39 // class declaration
40 //
41 
43 public:
45  typedef std::vector<LorentzVector> LorentzVectorCollection;
46  typedef std::vector<reco::PFCandidate> PFOutputCollection;
47 
48  explicit SoftKillerProducer(const edm::ParameterSet&);
49  ~SoftKillerProducer() override;
50 
51 private:
52  void produce(edm::Event&, const edm::EventSetup&) override;
53 
55 
56  double Rho_EtaMax_;
57  double rParam_;
58 };
59 
60 //
61 // constants, enums and typedefs
62 //
63 
64 //
65 // static data member definitions
66 //
67 
68 //
69 // constructors and destructor
70 //
72  : Rho_EtaMax_(iConfig.getParameter<double>("Rho_EtaMax")), rParam_(iConfig.getParameter<double>("rParam")) {
73  tokenPFCandidates_ = consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("PFCandidates"));
74 
75  produces<edm::ValueMap<LorentzVector> >("SoftKillerP4s");
76  produces<PFOutputCollection>();
77 }
78 
80 
81 //
82 // member functions
83 //
84 
85 // ------------ method called to produce the data ------------
87  std::unique_ptr<PFOutputCollection> pOutput(new PFOutputCollection);
88 
89  // get PF Candidates
92 
93  std::vector<fastjet::PseudoJet> fjInputs;
94  for (auto i = pfCandidates->begin(), ibegin = pfCandidates->begin(), iend = pfCandidates->end(); i != iend; ++i) {
95  fjInputs.push_back(fastjet::PseudoJet(i->px(), i->py(), i->pz(), i->energy()));
96  fjInputs.back().set_user_index(i - ibegin);
97  }
98 
99  // soft killer:
100  fastjet::contrib::SoftKiller soft_killer(Rho_EtaMax_, rParam_);
101 
102  double pt_threshold = 0.;
103  std::vector<fastjet::PseudoJet> soft_killed_event;
104  soft_killer.apply(fjInputs, soft_killed_event, pt_threshold);
105 
106  std::unique_ptr<edm::ValueMap<LorentzVector> > p4SKOut(new edm::ValueMap<LorentzVector>());
108 
109  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
110 
111  // To satisfy the value map, the size of the "killed" collection needs to be the
112  // same size as the input collection, so if the constituent is killed, just set E = 0
113  for (auto j = fjInputs.begin(), jend = fjInputs.end(); j != jend; ++j) {
114  const reco::Candidate& cand = pfCandidates->at(j->user_index());
115  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(cand.pdgId());
116  const reco::PFCandidate* pPF = dynamic_cast<const reco::PFCandidate*>(&cand);
117  reco::PFCandidate pCand(pPF ? *pPF : reco::PFCandidate(cand.charge(), cand.p4(), id));
118  auto val = j->user_index();
119  auto skmatch = find_if(soft_killed_event.begin(), soft_killed_event.end(), [&val](fastjet::PseudoJet const& i) {
120  return i.user_index() == val;
121  });
122  LorentzVector pVec;
123  if (skmatch != soft_killed_event.end()) {
124  pVec.SetPxPyPzE(skmatch->px(), skmatch->py(), skmatch->pz(), skmatch->E());
125  } else {
126  pVec.SetPxPyPzE(0., 0., 0., 0.);
127  }
128  pCand.setP4(pVec);
129  skP4s.push_back(pVec);
130  pOutput->push_back(pCand);
131  }
132 
133  //Compute the modified p4s
134  edm::ValueMap<LorentzVector>::Filler p4SKFiller(*p4SKOut);
135  p4SKFiller.insert(pfCandidates, skP4s.begin(), skP4s.end());
136  p4SKFiller.fill();
137 
138  iEvent.put(std::move(p4SKOut), "SoftKillerP4s");
139  iEvent.put(std::move(pOutput));
140 }
141 
142 //define this as a plug-in
SoftKillerProducer::LorentzVectorCollection
std::vector< LorentzVector > LorentzVectorCollection
Definition: SoftKillerProducer.cc:45
zmumugammaAnalyzer_cfi.pfCandidates
pfCandidates
Definition: zmumugammaAnalyzer_cfi.py:11
edm::helper::Filler::insert
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
mps_fire.i
i
Definition: mps_fire.py:355
SoftKillerProducer::PFOutputCollection
std::vector< reco::PFCandidate > PFOutputCollection
Definition: SoftKillerProducer.cc:46
PFCandidate.h
edm::EDGetTokenT
Definition: EDGetToken.h:33
SoftKillerProducer::SoftKillerProducer
SoftKillerProducer(const edm::ParameterSet &)
Definition: SoftKillerProducer.cc:71
edm::helper::Filler::fill
void fill()
Definition: ValueMap.h:65
SoftKillerProducer::rParam_
double rParam_
Definition: SoftKillerProducer.cc:57
EDProducer.h
edm::Handle
Definition: AssociativeIterator.h:50
MakerMacros.h
SoftKillerProducer::~SoftKillerProducer
~SoftKillerProducer() override
Definition: SoftKillerProducer.cc:79
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SoftKillerProducer::LorentzVector
math::XYZTLorentzVector LorentzVector
Definition: SoftKillerProducer.cc:44
edm::ParameterSet
Definition: ParameterSet.h:36
Event.h
cand
Definition: decayParser.h:34
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:57
SoftKillerProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: SoftKillerProducer.cc:86
SoftKillerProducer
Definition: SoftKillerProducer.cc:42
SoftKillerProducer::tokenPFCandidates_
edm::EDGetTokenT< reco::CandidateView > tokenPFCandidates_
Definition: SoftKillerProducer.cc:54
reco::Candidate
Definition: Candidate.h:27
SoftKillerProducer::Rho_EtaMax_
double Rho_EtaMax_
Definition: SoftKillerProducer.cc:56
ValueMap.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
reco::LeafCandidate::setP4
void setP4(const LorentzVector &p4) final
set 4-momentum
Definition: LeafCandidate.h:158
heppy_batch.val
val
Definition: heppy_batch.py:351
eostools.move
def move(src, dest)
Definition: eostools.py:511
Frameworkfwd.h
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:31
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
edm::ValueMap
Definition: ValueMap.h:107
reco::PFCandidate::translatePdgIdToType
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:209
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
Candidate.h
edm::helper::Filler
Definition: ValueMap.h:22
View.h
ParameterSet.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
PFCandidateFwd.h