CMS 3D CMS Logo

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