CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TauJetSelectorForHLTTrackSeeding.cc
Go to the documentation of this file.
2 
4 
6  : ptMinCaloJet_(iConfig.getParameter<double>("ptMinCaloJet")),
7  etaMinCaloJet_(iConfig.getParameter<double>("etaMinCaloJet")),
8  etaMaxCaloJet_(iConfig.getParameter<double>("etaMaxCaloJet")),
9  tauConeSize_(iConfig.getParameter<double>("tauConeSize")),
10  isolationConeSize_(iConfig.getParameter<double>("isolationConeSize")),
11  fractionMinCaloInTauCone_(iConfig.getParameter<double>("fractionMinCaloInTauCone")),
12  fractionMaxChargedPUInCaloCone_(iConfig.getParameter<double>("fractionMaxChargedPUInCaloCone")),
13  ptTrkMaxInCaloCone_(iConfig.getParameter<double>("ptTrkMaxInCaloCone")),
14  nTrkMaxInCaloCone_(iConfig.getParameter<int>("nTrkMaxInCaloCone")) {
15  //now do what ever initialization is needed
16  inputTrackJetToken_ = consumes<reco::TrackJetCollection>(iConfig.getParameter<edm::InputTag>("inputTrackJetTag"));
17  inputCaloJetToken_ = consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("inputCaloJetTag"));
18  inputTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTrackTag"));
19 
20  produces<reco::TrackJetCollection>();
21 }
22 
24  // do anything here that needs to be done at desctruction time
25  // (e.g. close files, deallocate resources etc.)
26 }
27 
28 //
29 // member functions
30 //
31 
32 // ------------ method called on each new Event ------------
35  const edm::EventSetup& iSetup) const {
36  std::unique_ptr<reco::TrackJetCollection> augmentedTrackJets(new reco::TrackJetCollection);
37 
39  iEvent.getByToken(inputTrackJetToken_, trackjets);
40 
41  for (reco::TrackJetCollection::const_iterator trackjet = trackjets->begin(); trackjet != trackjets->end();
42  trackjet++) {
43  augmentedTrackJets->push_back(*trackjet);
44  }
45 
47  iEvent.getByToken(inputTrackToken_, tracks);
48 
50  iEvent.getByToken(inputCaloJetToken_, calojets);
51 
52  const double tauConeSize2 = tauConeSize_ * tauConeSize_;
53  const double isolationConeSize2 = isolationConeSize_ * isolationConeSize_;
54 
55  for (reco::CaloJetCollection::const_iterator calojet = calojets->begin(); calojet != calojets->end(); calojet++) {
56  if (calojet->pt() < ptMinCaloJet_)
57  continue;
58  double etaJet = calojet->eta();
59  double phiJet = calojet->phi();
60  if (etaJet < etaMinCaloJet_)
61  continue;
62  if (etaJet > etaMaxCaloJet_)
63  continue;
64 
65  std::vector<CaloTowerPtr> const& theTowers = calojet->getCaloConstituents();
66  double ptIn = 0.;
67  double ptOut = 0.;
68  for (unsigned int itwr = 0; itwr < theTowers.size(); ++itwr) {
69  double etaTwr = theTowers[itwr]->eta() - etaJet;
70  double phiTwr = deltaPhi(theTowers[itwr]->phi(), phiJet);
71  double deltaR2 = etaTwr * etaTwr + phiTwr * phiTwr;
72  //std::cout << "Tower eta/phi/et : " << etaTwr << " " << phiTwr << " " << theTowers[itwr]->pt() << std::endl;
73  if (deltaR2 < tauConeSize2) {
74  ptIn += theTowers[itwr]->pt();
75  } else if (deltaR2 < isolationConeSize2) {
76  ptOut += theTowers[itwr]->pt();
77  }
78  }
79  double ptTot = ptIn + ptOut;
80  double fracIn = ptIn / ptTot;
81 
82  // We are looking for isolated tracks
83  if (fracIn < fractionMinCaloInTauCone_)
84  continue;
85 
86  int ntrk = 0;
87  double ptTrk = 0.;
88 
89  for (reco::TrackJetCollection::const_iterator trackjet = trackjets->begin(); trackjet != trackjets->end();
90  trackjet++) {
91  for (unsigned itr = 0; itr < trackjet->numberOfTracks(); ++itr) {
92  edm::Ptr<reco::Track> track = trackjet->track(itr);
93  double trackEta = track->eta() - etaJet;
94  double trackPhi = deltaPhi(track->phi(), phiJet);
95  double deltaR2 = trackEta * trackEta + trackPhi * trackPhi;
96  if (deltaR2 < isolationConeSize2) {
97  ntrk++;
98  ptTrk += track->pt();
99  }
100  }
101  }
102  // We are looking for calojets without signal tracks already in
103  if (ntrk > nTrkMaxInCaloCone_)
104  continue;
105  if (ptTrk > ptTrkMaxInCaloCone_)
106  continue;
107 
108  int ntrk2 = 0;
109  double ptTrk2 = 0.;
110 
111  for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); track++) {
112  double trackEta = track->eta() - etaJet;
113  double trackPhi = deltaPhi(track->phi(), phiJet);
114  double deltaR2 = trackEta * trackEta + trackPhi * trackPhi;
115  if (deltaR2 < isolationConeSize2) {
116  ntrk2++;
117  ptTrk2 += track->pt();
118  }
119  }
120  // We are looking for signal jets, not PU jets
121  double fractionChargedPU = ptTrk2 / calojet->pt();
122  if (fractionChargedPU > fractionMaxChargedPUInCaloCone_)
123  continue;
124  /*
125  std::cout << "Calo Jet " << calojet->pt() << " " << calojet->eta()
126  << " " << ptIn << " " << ptOut << " " << fracIn
127  << " " << ptTrk << " " << ntrk
128  << " " << fractionChargedPU
129  << std::endl;
130  */
131  math::XYZTLorentzVector p4(calojet->p4());
132  math::XYZPoint vertex(calojet->vertex());
133  augmentedTrackJets->push_back(reco::TrackJet(p4, vertex));
134  }
135 
136  iEvent.put(std::move(augmentedTrackJets));
137 }
138 
139 // ------------ method called once each job just before starting event loop ------------
141 
142 // ------------ method called once each job just after ending the event loop ------------
144 
147  desc.add<edm::InputTag>("inputTrackJetTag", edm::InputTag("hltAntiKT5TrackJetsIter0"))
148  ->setComment("Oryginal TrackJet collection");
149  desc.add<edm::InputTag>("inputCaloJetTag", edm::InputTag("hltAntiKT5CaloJetsPFEt5"))
150  ->setComment("CaloJet collection");
151  desc.add<edm::InputTag>("inputTrackTag", edm::InputTag("hltPFlowTrackSelectionHighPurity"))
152  ->setComment("Track collection for track isolation");
153  desc.add<double>("ptMinCaloJet", 5.0)->setComment("Min. pT of CaloJet");
154  desc.add<double>("etaMinCaloJet", -2.7)->setComment("Min. eta of CaloJet");
155  desc.add<double>("etaMaxCaloJet", +2.7)->setComment("Max. eta of CaloJet");
156  desc.add<double>("tauConeSize", 0.2)->setComment("Radius of tau signal cone");
157  desc.add<double>("isolationConeSize", 0.5)->setComment("Radius of isolation cone");
158  desc.add<double>("fractionMinCaloInTauCone", 0.8)->setComment("Min. fraction of calo energy in tau signal cone");
159  desc.add<double>("fractionMaxChargedPUInCaloCone", 0.2)
160  ->setComment("Max. fraction of PU charged energy Sum(Pt_trks)/Pt_jet in signal cone");
161  desc.add<double>("ptTrkMaxInCaloCone", 1.0)->setComment("Max. sum of pT of tracks in isolation cone");
162  desc.add<int>("nTrkMaxInCaloCone", 0)->setComment("Max. no. of tracks in isolation cone");
163  descriptions.setComment(
164  "This module produces a collection of TrackJets to seed iterative tracking.\nIt is done by enriching the input "
165  "TrackJet collection with CaloJets passing isolation critera - tau candidates");
166  descriptions.add("tauJetSelectorForHLTTrackSeeding", desc);
167 }
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
edm::EDGetTokenT< reco::TrackJetCollection > inputTrackJetToken_
edm::EDGetTokenT< reco::TrackCollection > inputTrackToken_
auto const & tracks
cannot be loose
TauJetSelectorForHLTTrackSeeding(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Jets made out of tracks.
Definition: TrackJet.h:24
void setComment(std::string const &value)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< TrackJet > TrackJetCollection
collection of TrackJet objects
edm::EDGetTokenT< reco::CaloJetCollection > inputCaloJetToken_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override