CMS 3D CMS Logo

HPSPFTauProducer.cc
Go to the documentation of this file.
1 #include "HPSPFTauProducer.h"
3 #include <cmath> // std::fabs
4 
6  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
7  tauBuilder_(cfg),
8  useChargedPFCandSeeds_(cfg.getParameter<bool>("useChargedPFCandSeeds")),
9  minSeedChargedPFCandPt_(cfg.getParameter<double>("minSeedChargedPFCandPt")),
10  maxSeedChargedPFCandEta_(cfg.getParameter<double>("maxSeedChargedPFCandEta")),
11  maxSeedChargedPFCandDz_(cfg.getParameter<double>("maxSeedChargedPFCandDz")),
12  useJetSeeds_(cfg.getParameter<bool>("useJetSeeds")),
13  minSeedJetPt_(cfg.getParameter<double>("minSeedJetPt")),
14  maxSeedJetEta_(cfg.getParameter<double>("maxSeedJetEta")),
15  minPFTauPt_(cfg.getParameter<double>("minPFTauPt")),
16  maxPFTauEta_(cfg.getParameter<double>("maxPFTauEta")),
17  minLeadChargedPFCandPt_(cfg.getParameter<double>("minLeadChargedPFCandPt")),
18  maxLeadChargedPFCandEta_(cfg.getParameter<double>("maxLeadChargedPFCandEta")),
19  maxLeadChargedPFCandDz_(cfg.getParameter<double>("maxLeadChargedPFCandDz")),
20  maxChargedIso_(cfg.getParameter<double>("maxChargedIso")),
21  maxChargedRelIso_(cfg.getParameter<double>("maxChargedRelIso")),
22  deltaRCleaning_(cfg.getParameter<double>("deltaRCleaning")),
23  applyPreselection_(cfg.getParameter<bool>("applyPreselection")),
24  debug_(cfg.getUntrackedParameter<bool>("debug", false)) {
25  srcL1PFCands_ = cfg.getParameter<edm::InputTag>("srcL1PFCands");
26  tokenL1PFCands_ = consumes<l1t::PFCandidateCollection>(srcL1PFCands_);
27  srcL1Jets_ = cfg.getParameter<edm::InputTag>("srcL1Jets");
28  if (useJetSeeds_) {
29  tokenL1Jets_ = consumes<std::vector<reco::CaloJet>>(srcL1Jets_);
30  }
31  srcL1Vertices_ = cfg.getParameter<edm::InputTag>("srcL1Vertices");
32  if (!srcL1Vertices_.label().empty()) {
33  tokenL1Vertices_ = consumes<l1t::VertexWordCollection>(srcL1Vertices_);
34  }
36 
37  edm::ParameterSet cfg_signalQualityCuts = cfg.getParameter<edm::ParameterSet>("signalQualityCuts");
38  signalQualityCutsDzCutDisabled_ = readL1PFTauQualityCuts(cfg_signalQualityCuts, "disabled");
39  edm::ParameterSet cfg_isolationQualityCuts = cfg.getParameter<edm::ParameterSet>("isolationQualityCuts");
40  isolationQualityCutsDzCutDisabled_ = readL1PFTauQualityCuts(cfg_isolationQualityCuts, "disabled");
41 
42  produces<l1t::HPSPFTauCollection>();
43 }
44 
45 namespace {
46  bool isHigherPt_pfCandRef(const l1t::PFCandidateRef& l1PFCand1, const l1t::PFCandidateRef& l1PFCand2) {
47  return l1PFCand1->pt() > l1PFCand2->pt();
48  }
49 
50  bool isHigherPt_pfTau(const l1t::HPSPFTau& l1PFTau1, const l1t::HPSPFTau& l1PFTau2) {
51  return l1PFTau1.pt() > l1PFTau2.pt();
52  }
53 } // namespace
54 
56  std::unique_ptr<l1t::HPSPFTauCollection> l1PFTauCollectionCleaned(new l1t::HPSPFTauCollection());
57 
59  evt.getByToken(tokenL1PFCands_, l1PFCands);
60 
62  float primaryVertex_z = 0.;
63  if (!srcL1Vertices_.label().empty()) {
66  if (!vertices->empty()) {
68  primaryVertex_z = primaryVertex->z0();
69  }
70  }
71 
72  // build collection of selected PFCandidates
73  std::vector<l1t::PFCandidateRef> selectedL1PFCandsSignalQualityCuts;
74  std::vector<l1t::PFCandidateRef> selectedL1PFCandsSignalOrIsolationQualityCuts;
75  size_t numL1PFCands = l1PFCands->size();
76  for (size_t idxL1PFCand = 0; idxL1PFCand < numL1PFCands; ++idxL1PFCand) {
77  l1t::PFCandidateRef l1PFCand(l1PFCands, idxL1PFCand);
78  bool passesSignalQualityCuts = isSelected(signalQualityCutsDzCutDisabled_, *l1PFCand, primaryVertex_z);
79  bool passesIsolationQualityCuts = isSelected(isolationQualityCutsDzCutDisabled_, *l1PFCand, primaryVertex_z);
80  if (passesSignalQualityCuts) {
81  selectedL1PFCandsSignalQualityCuts.push_back(l1PFCand);
82  }
83  if (passesSignalQualityCuts || passesIsolationQualityCuts) {
84  selectedL1PFCandsSignalOrIsolationQualityCuts.push_back(l1PFCand);
85  }
86  }
87 
88  // sort PFCandidate collection by decreasing pT
89  std::sort(selectedL1PFCandsSignalQualityCuts.begin(), selectedL1PFCandsSignalQualityCuts.end(), isHigherPt_pfCandRef);
90  std::sort(selectedL1PFCandsSignalOrIsolationQualityCuts.begin(),
91  selectedL1PFCandsSignalOrIsolationQualityCuts.end(),
92  isHigherPt_pfCandRef);
93 
94  l1t::HPSPFTauCollection l1PFTauCollectionUncleaned;
95 
97  for (const auto& l1PFCand : selectedL1PFCandsSignalQualityCuts) {
98  if (l1PFCand->charge() != 0 && l1PFCand->pt() > minSeedChargedPFCandPt_ &&
99  std::fabs(l1PFCand->eta()) < maxSeedChargedPFCandEta_) {
100  bool isFromPrimaryVertex = false;
101  if (primaryVertex.get()) {
102  l1t::PFTrackRef l1PFTrack = l1PFCand->pfTrack();
103  double dz = std::fabs(l1PFTrack->vertex().z() - primaryVertex_z);
104  if (dz < maxSeedChargedPFCandDz_) {
105  isFromPrimaryVertex = true;
106  }
107  } else {
108  isFromPrimaryVertex = true;
109  }
110  if (isFromPrimaryVertex) {
111  tauBuilder_.reset();
112  tauBuilder_.setL1PFCandProductID(l1PFCands.id());
114  tauBuilder_.setL1PFTauSeed(l1PFCand);
115  tauBuilder_.addL1PFCandidates(selectedL1PFCandsSignalOrIsolationQualityCuts);
118  if (l1PFTau.pt() > isPFTauPt_)
119  l1PFTauCollectionUncleaned.push_back(l1PFTau);
120  }
121  }
122  }
123  }
124 
125  if (useJetSeeds_) {
127  evt.getByToken(tokenL1Jets_, l1Jets);
128 
129  size_t numL1Jets = l1Jets->size();
130  for (size_t idxL1Jet = 0; idxL1Jet < numL1Jets; ++idxL1Jet) {
131  reco::CaloJetRef l1Jet(l1Jets, idxL1Jet);
132  if (l1Jet->pt() > minSeedJetPt_ && std::fabs(l1Jet->eta()) < maxSeedJetEta_) {
133  tauBuilder_.reset();
134  tauBuilder_.setL1PFCandProductID(l1PFCands.id());
136  //tauBuilder_.setL1PFTauSeed(l1Jet);
137  tauBuilder_.setL1PFTauSeed(l1Jet, selectedL1PFCandsSignalQualityCuts);
138  tauBuilder_.addL1PFCandidates(selectedL1PFCandsSignalOrIsolationQualityCuts);
141  if (l1PFTau.pt() > isPFTauPt_)
142  l1PFTauCollectionUncleaned.push_back(l1PFTau);
143  }
144  }
145  }
146 
147  // sort PFTau candidate collection by decreasing pT
148  std::sort(l1PFTauCollectionUncleaned.begin(), l1PFTauCollectionUncleaned.end(), isHigherPt_pfTau);
149 
150  for (const auto& l1PFTau : l1PFTauCollectionUncleaned) {
151  if (applyPreselection_ &&
152  !(l1PFTau.pt() > minPFTauPt_ && std::fabs(l1PFTau.eta()) < maxPFTauEta_ &&
153  l1PFTau.leadChargedPFCand().isNonnull() && l1PFTau.leadChargedPFCand()->pt() > minLeadChargedPFCandPt_ &&
154  std::fabs(l1PFTau.leadChargedPFCand()->eta()) < maxLeadChargedPFCandEta_ &&
155  (srcL1Vertices_.label().empty() ||
156  (primaryVertex.isNonnull() && l1PFTau.leadChargedPFCand()->pfTrack().isNonnull() &&
157  std::fabs(l1PFTau.leadChargedPFCand()->pfTrack()->vertex().z() - primaryVertex->z0()) <
159  l1PFTau.sumChargedIso() < maxChargedIso_ && l1PFTau.sumChargedIso() < maxChargedRelIso_ * l1PFTau.pt()))
160  continue;
161 
162  bool isOverlap = false;
163  for (const auto& l1PFTau2 : *l1PFTauCollectionCleaned) {
164  double deltaEta = l1PFTau.eta() - l1PFTau2.eta();
165  double deltaPhi = l1PFTau.phi() - l1PFTau2.phi();
167  isOverlap = true;
168  }
169  }
170  if (!isOverlap) {
171  l1PFTauCollectionCleaned->push_back(l1PFTau);
172  }
173  }
174 
175  evt.put(std::move(l1PFTauCollectionCleaned));
176 }
177 
180 
182  // HPSPFTauProducerPF
184  desc.add<bool>("useJetSeeds", true);
185  desc.add<double>("minPFTauPt", 20.0);
186  desc.add<double>("minSignalConeSize", 0.05);
187  {
189  {
191  psd1.add<double>("minPt", 0.0);
192  psd0.add<edm::ParameterSetDescription>("neutralHadron", psd1);
193  }
194  {
196  psd1.add<double>("maxDz", 0.4);
197  psd1.add<double>("minPt", 0.0);
198  psd0.add<edm::ParameterSetDescription>("muon", psd1);
199  }
200  {
202  psd1.add<double>("maxDz", 0.4);
203  psd1.add<double>("minPt", 0.0);
204  psd0.add<edm::ParameterSetDescription>("electron", psd1);
205  }
206  {
208  psd1.add<double>("minPt", 0.0);
209  psd0.add<edm::ParameterSetDescription>("photon", psd1);
210  }
211  {
213  psd1.add<double>("maxDz", 0.4);
214  psd1.add<double>("minPt", 0.0);
215  psd0.add<edm::ParameterSetDescription>("chargedHadron", psd1);
216  }
217  desc.add<edm::ParameterSetDescription>("isolationQualityCuts", psd0);
218  }
219  desc.add<double>("stripSizePhi", 0.2);
220  desc.add<double>("minSeedJetPt", 30.0);
221  desc.add<double>("maxChargedRelIso", 1.0);
222  desc.add<double>("minSeedChargedPFCandPt", 5.0);
223  desc.add<edm::InputTag>("srcL1PFCands", edm::InputTag("l1tLayer1", "PF"));
224  desc.add<double>("stripSizeEta", 0.05);
225  desc.add<double>("maxLeadChargedPFCandEta", 2.4);
226  desc.add<double>("deltaRCleaning", 0.4);
227  desc.add<bool>("useStrips", true);
228  desc.add<double>("maxSeedChargedPFCandDz", 1000.0);
229  desc.add<double>("minLeadChargedPFCandPt", 1.0);
230  desc.add<double>("maxSeedChargedPFCandEta", 2.4);
231  desc.add<bool>("applyPreselection", false);
232  desc.add<double>("isolationConeSize", 0.4);
233  desc.add<edm::InputTag>("srcL1Vertices", edm::InputTag("l1tVertexFinderEmulator", "l1verticesEmulation"));
234  desc.add<double>("maxChargedIso", 1000.0);
235  {
237  {
239  psd1.add<double>("minPt", 0.0);
240  psd0.add<edm::ParameterSetDescription>("neutralHadron", psd1);
241  }
242  {
244  psd1.add<double>("maxDz", 0.4);
245  psd1.add<double>("minPt", 0.0);
246  psd0.add<edm::ParameterSetDescription>("muon", psd1);
247  }
248  {
250  psd1.add<double>("maxDz", 0.4);
251  psd1.add<double>("minPt", 0.0);
252  psd0.add<edm::ParameterSetDescription>("electron", psd1);
253  }
254  {
256  psd1.add<double>("minPt", 0.0);
257  psd0.add<edm::ParameterSetDescription>("photon", psd1);
258  }
259  {
261  psd1.add<double>("maxDz", 0.4);
262  psd1.add<double>("minPt", 0.0);
263  psd0.add<edm::ParameterSetDescription>("chargedHadron", psd1);
264  }
265  desc.add<edm::ParameterSetDescription>("signalQualityCuts", psd0);
266  }
267  desc.add<bool>("useChargedPFCandSeeds", true);
268  desc.add<double>("maxLeadChargedPFCandDz", 1000.0);
269  desc.add<double>("maxSeedJetEta", 2.4);
270  desc.add<std::string>("signalConeSize", "2.8/max(1., pt)");
271  desc.add<edm::InputTag>("srcL1Jets",
272  edm::InputTag("l1tPhase1JetCalibrator9x9trimmed", "Phase1L1TJetFromPfCandidates"));
273  desc.addUntracked<bool>("debug", false);
274  desc.add<double>("maxPFTauEta", 2.4);
275  desc.add<double>("maxSignalConeSize", 0.1);
276  descriptions.addWithDefaultLabel(desc);
277 }
278 
edm::InputTag srcL1Vertices_
const double isPFTauPt_
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
ProductID id() const
Definition: HandleBase.cc:29
double pt() const final
transverse momentum
double maxLeadChargedPFCandEta_
void setL1PFCandProductID(const edm::ProductID &l1PFCandProductID)
void addL1PFCandidates(const std::vector< l1t::PFCandidateRef > &l1PFCands)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:540
std::vector< HPSPFTau > HPSPFTauCollection
Definition: HPSPFTauFwd.h:7
std::string const & label() const
Definition: InputTag.h:36
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static const double deltaEta
Definition: CaloConstants.h:8
L1HPSPFTauBuilder tauBuilder_
edm::InputTag srcL1Jets_
double maxSeedChargedPFCandEta_
double minLeadChargedPFCandPt_
void produce(edm::Event &evt, const edm::EventSetup &es) override
void setL1PFTauSeed(const l1t::PFCandidateRef &l1PFCandSeed)
edm::Ref< VertexWordCollection > VertexWordRef
Definition: VertexWord.h:198
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double minSeedChargedPFCandPt_
std::vector< L1HPSPFTauQualityCut > readL1PFTauQualityCuts(const edm::ParameterSet &cfg, const std::string &dzCut, bool debug=false)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setVertex(const l1t::VertexWordRef &primaryVertex)
double maxSeedChargedPFCandDz_
bool isSelected(const std::vector< L1HPSPFTauQualityCut > &qualityCuts, const l1t::PFCandidate &pfCand, float_t primaryVertexZ)
std::vector< L1HPSPFTauQualityCut > isolationQualityCutsDzCutDisabled_
edm::EDGetTokenT< std::vector< reco::CaloJet > > tokenL1Jets_
l1t::HPSPFTau getL1PFTau() const
HPSPFTauProducer(const edm::ParameterSet &cfg)
edm::InputTag srcL1PFCands_
double maxLeadChargedPFCandDz_
std::vector< L1HPSPFTauQualityCut > signalQualityCutsDzCutDisabled_
edm::EDGetTokenT< l1t::PFCandidateCollection > tokenL1PFCands_
double phi() const final
momentum azimuthal angle
primaryVertex
hltOfflineBeamSpot for HLTMON
def move(src, dest)
Definition: eostools.py:511
edm::EDGetTokenT< l1t::VertexWordCollection > tokenL1Vertices_
double eta() const final
momentum pseudorapidity