CMS 3D CMS Logo

PATPuppiJetSpecificsProducer.cc
Go to the documentation of this file.
1 /*
2  * PATPuppiJetSpecificProducer
3  *
4  * Author: Andreas Hinzmann
5  *
6  * Compute weighted constituent multiplicites for PUPPI PAT jets
7  *
8  */
9 
16 
19 
23 
25 {
26 public:
27 
30  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
31  static void fillDescriptions(edm::ConfigurationDescriptions & descriptions);
32 
33 private:
34  // input collection
37 
38 };
39 
41 {
42  srcjets_ = cfg.getParameter<edm::InputTag>("src");
43  jets_token = consumes<edm::View<pat::Jet> >(srcjets_);
44 
45  produces<edm::ValueMap<float> >("puppiMultiplicity");
46  produces<edm::ValueMap<float> >("neutralPuppiMultiplicity");
47  produces<edm::ValueMap<float> >("neutralHadronPuppiMultiplicity");
48  produces<edm::ValueMap<float> >("photonPuppiMultiplicity");
49  produces<edm::ValueMap<float> >("HFHadronPuppiMultiplicity");
50  produces<edm::ValueMap<float> >("HFEMPuppiMultiplicity");
51 }
52 
54 {
56  evt.getByToken(jets_token, jets);
57 
58  std::vector<float> puppiMultiplicities;
59  std::vector<float> neutralPuppiMultiplicities;
60  std::vector<float> neutralHadronPuppiMultiplicities;
61  std::vector<float> photonPuppiMultiplicities;
62  std::vector<float> HFHadronPuppiMultiplicities;
63  std::vector<float> HFEMPuppiMultiplicities;
64 
65  for( auto const& c : *jets ) {
66  float puppiMultiplicity = 0;
67  float neutralPuppiMultiplicity = 0;
68  float neutralHadronPuppiMultiplicity = 0;
69  float photonPuppiMultiplicity = 0;
70  float HFHadronPuppiMultiplicity = 0;
71  float HFEMPuppiMultiplicity = 0;
72 
73  for (unsigned i = 0; i < c.numberOfDaughters(); i++) {
74  const pat::PackedCandidate &dau = static_cast<const pat::PackedCandidate &>(*c.daughter(i));
75  auto weight = dau.puppiWeight();
76  puppiMultiplicity += weight;
77  // This logic is taken from RecoJets/JetProducers/src/JetSpecific.cc
78  switch (std::abs(dau.pdgId())) {
79  case 130: //PFCandidate::h0 : // neutral hadron
80  neutralHadronPuppiMultiplicity += weight;
81  neutralPuppiMultiplicity += weight;
82  break;
83  case 22: //PFCandidate::gamma: // photon
84  photonPuppiMultiplicity += weight;
85  neutralPuppiMultiplicity += weight;
86  break;
87  case 1: // PFCandidate::h_HF : // hadron in HF
88  HFHadronPuppiMultiplicity += weight;
89  neutralPuppiMultiplicity += weight;
90  break;
91  case 2: //PFCandidate::egamma_HF : // electromagnetic in HF
92  HFEMPuppiMultiplicity += weight;
93  neutralPuppiMultiplicity += weight;
94  break;
95  }
96  }
97 
98  puppiMultiplicities.push_back(puppiMultiplicity);
99  neutralPuppiMultiplicities.push_back(neutralPuppiMultiplicity);
100  neutralHadronPuppiMultiplicities.push_back(neutralHadronPuppiMultiplicity);
101  photonPuppiMultiplicities.push_back(photonPuppiMultiplicity);
102  HFHadronPuppiMultiplicities.push_back(HFHadronPuppiMultiplicity);
103  HFEMPuppiMultiplicities.push_back(HFEMPuppiMultiplicity);
104  }
105 
106  std::unique_ptr<edm::ValueMap<float> > puppiMultiplicities_out(new edm::ValueMap<float>());
107  edm::ValueMap<float>::Filler puppiMultiplicities_filler(*puppiMultiplicities_out);
108  puppiMultiplicities_filler.insert(jets,puppiMultiplicities.begin(),puppiMultiplicities.end());
109  puppiMultiplicities_filler.fill();
110  evt.put(std::move(puppiMultiplicities_out),"puppiMultiplicity");
111 
112  std::unique_ptr<edm::ValueMap<float> > neutralPuppiMultiplicities_out(new edm::ValueMap<float>());
113  edm::ValueMap<float>::Filler neutralPuppiMultiplicities_filler(*neutralPuppiMultiplicities_out);
114  neutralPuppiMultiplicities_filler.insert(jets,neutralPuppiMultiplicities.begin(),neutralPuppiMultiplicities.end());
115  neutralPuppiMultiplicities_filler.fill();
116  evt.put(std::move(neutralPuppiMultiplicities_out),"neutralPuppiMultiplicity");
117 
118  std::unique_ptr<edm::ValueMap<float> > neutralHadronPuppiMultiplicities_out(new edm::ValueMap<float>());
119  edm::ValueMap<float>::Filler neutralHadronPuppiMultiplicities_filler(*neutralHadronPuppiMultiplicities_out);
120  neutralHadronPuppiMultiplicities_filler.insert(jets,neutralHadronPuppiMultiplicities.begin(),neutralHadronPuppiMultiplicities.end());
121  neutralHadronPuppiMultiplicities_filler.fill();
122  evt.put(std::move(neutralHadronPuppiMultiplicities_out),"neutralHadronPuppiMultiplicity");
123 
124  std::unique_ptr<edm::ValueMap<float> > photonPuppiMultiplicities_out(new edm::ValueMap<float>());
125  edm::ValueMap<float>::Filler photonPuppiMultiplicities_filler(*photonPuppiMultiplicities_out);
126  photonPuppiMultiplicities_filler.insert(jets,photonPuppiMultiplicities.begin(),photonPuppiMultiplicities.end());
127  photonPuppiMultiplicities_filler.fill();
128  evt.put(std::move(photonPuppiMultiplicities_out),"photonPuppiMultiplicity");
129 
130  std::unique_ptr<edm::ValueMap<float> > HFHadronPuppiMultiplicities_out(new edm::ValueMap<float>());
131  edm::ValueMap<float>::Filler HFHadronPuppiMultiplicities_filler(*HFHadronPuppiMultiplicities_out);
132  HFHadronPuppiMultiplicities_filler.insert(jets,HFHadronPuppiMultiplicities.begin(),HFHadronPuppiMultiplicities.end());
133  HFHadronPuppiMultiplicities_filler.fill();
134  evt.put(std::move(HFHadronPuppiMultiplicities_out),"HFHadronPuppiMultiplicity");
135 
136  std::unique_ptr<edm::ValueMap<float> > HFEMPuppiMultiplicities_out(new edm::ValueMap<float>());
137  edm::ValueMap<float>::Filler HFEMPuppiMultiplicities_filler(*HFEMPuppiMultiplicities_out);
138  HFEMPuppiMultiplicities_filler.insert(jets,HFEMPuppiMultiplicities.begin(),HFEMPuppiMultiplicities.end());
139  HFEMPuppiMultiplicities_filler.fill();
140  evt.put(std::move(HFEMPuppiMultiplicities_out),"HFEMPuppiMultiplicity");
141 }
142 
145  desc.add<edm::InputTag>("src", edm::InputTag("slimmedJets"));
146  descriptions.add("patPuppiJetSpecificProducer", desc);
147 }
148 
150 
float puppiWeight() const
Set both weights at once (with option for only full PUPPI)
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::View< pat::Jet > > jets_token
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
int pdgId() const override
PDG identifier.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
Definition: weight.py:1
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
PATPuppiJetSpecificProducer(const edm::ParameterSet &cfg)
vector< PseudoJet > jets
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
def move(src, dest)
Definition: eostools.py:510