CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 // system include files
21 #include <memory>
22 #include <vector>
23 
24 // user include files
27 
30 
37 #include "fastjet/contrib/SoftKiller.hh"
38 
39 
40 
41 //
42 // class declaration
43 //
44 
46 public:
48  typedef std::vector<LorentzVector> LorentzVectorCollection;
49  typedef std::vector< reco::PFCandidate > PFOutputCollection;
50 
51  explicit SoftKillerProducer(const edm::ParameterSet&);
53 
54 private:
55  virtual void produce(edm::Event&, const edm::EventSetup&) override;
56 
58 
59 
60  double Rho_EtaMax_;
61  double rParam_;
62 };
63 
64 //
65 // constants, enums and typedefs
66 //
67 
68 
69 //
70 // static data member definitions
71 //
72 
73 //
74 // constructors and destructor
75 //
77  Rho_EtaMax_( iConfig.getParameter<double>("Rho_EtaMax") ),
78  rParam_ ( iConfig.getParameter<double>("rParam") )
79 {
80 
82  = consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("PFCandidates"));
83 
84  produces<edm::ValueMap<LorentzVector> > ("SoftKillerP4s");
85  produces< PFOutputCollection >();
86 
87 }
88 
89 
91 {
92 
93 }
94 
95 
96 //
97 // member functions
98 //
99 
100 // ------------ method called to produce the data ------------
101 void
103 {
104 
105  std::auto_ptr< PFOutputCollection > pOutput( new PFOutputCollection );
106 
107  // get PF Candidates
109  iEvent.getByToken( tokenPFCandidates_, pfCandidates);
110 
111  std::vector<fastjet::PseudoJet> fjInputs;
112  for ( auto i = pfCandidates->begin(),
113  ibegin = pfCandidates->begin(),
114  iend = pfCandidates->end(); i != iend; ++i ) {
115  fjInputs.push_back( fastjet::PseudoJet( i->px(), i->py(), i->pz(), i->energy() ) );
116  fjInputs.back().set_user_index( i - ibegin );
117  }
118 
119  // soft killer:
120  fastjet::contrib::SoftKiller soft_killer(Rho_EtaMax_, rParam_);
121 
122  double pt_threshold = 0.;
123  std::vector<fastjet::PseudoJet> soft_killed_event;
124  soft_killer.apply(fjInputs, soft_killed_event, pt_threshold);
125 
126  std::auto_ptr<edm::ValueMap<LorentzVector> > p4SKOut(new edm::ValueMap<LorentzVector>());
128 
129  static const reco::PFCandidate dummySinceTranslateIsNotStatic;
130 
131  // To satisfy the value map, the size of the "killed" collection needs to be the
132  // same size as the input collection, so if the constituent is killed, just set E = 0
133  for ( auto j = fjInputs.begin(),
134  jend = fjInputs.end(); j != jend; ++j ) {
135  const reco::Candidate& cand = pfCandidates->at(j->user_index());
136  auto id = dummySinceTranslateIsNotStatic.translatePdgIdToType(cand.pdgId());
137  const reco::PFCandidate *pPF = dynamic_cast<const reco::PFCandidate*>(&cand);
138  reco::PFCandidate pCand( pPF ? *pPF : reco::PFCandidate(cand.charge(), cand.p4(), id) );
139  auto val = j->user_index();
140  auto skmatch = find_if( soft_killed_event.begin(), soft_killed_event.end(), [&val](fastjet::PseudoJet const & i){return i.user_index() == val;} );
141  LorentzVector pVec;
142  if ( skmatch != soft_killed_event.end() ) {
143  pVec.SetPxPyPzE(skmatch->px(),skmatch->py(),skmatch->pz(),skmatch->E());
144  } else {
145  pVec.SetPxPyPzE( 0., 0., 0., 0.);
146  }
147  pCand.setP4(pVec);
148  skP4s.push_back( pVec );
149  pOutput->push_back(pCand);
150  }
151 
152  //Compute the modified p4s
153  edm::ValueMap<LorentzVector>::Filler p4SKFiller(*p4SKOut);
154  p4SKFiller.insert(pfCandidates,skP4s.begin(), skP4s.end() );
155  p4SKFiller.fill();
156 
157  iEvent.put(p4SKOut,"SoftKillerP4s");
158  iEvent.put( pOutput );
159 
160 }
161 
162 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
virtual void setP4(const LorentzVector &p4)
set 4-momentum
math::XYZTLorentzVector LorentzVector
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
virtual void produce(edm::Event &, const edm::EventSetup &) override
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
std::vector< LorentzVector > LorentzVectorCollection
int j
Definition: DBlmapReader.cc:9
virtual int charge() const =0
electric charge
SoftKillerProducer(const edm::ParameterSet &)
virtual int pdgId() const =0
PDG identifier.
std::vector< reco::PFCandidate > PFOutputCollection
ParticleType translatePdgIdToType(int pdgid) const
Definition: PFCandidate.cc:224
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
edm::EDGetTokenT< reco::CandidateView > tokenPFCandidates_
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector