CMS 3D CMS Logo

HLTTauRefProducer.cc
Go to the documentation of this file.
2 #include "TLorentzVector.h"
3 // TAU includes
8 // ELECTRON includes
15 // MUON includes
18 #include "TLorentzVector.h"
20 //CaloTower includes
24 #include "Math/GenVector/VectorUtil.h"
25 
26 using namespace edm;
27 using namespace reco;
28 using namespace std;
29 
31 {
32  //One Parameter Set per Collection
33  {
34  auto const& pfTau = iConfig.getUntrackedParameter<edm::ParameterSet>("PFTaus");
35  PFTaus_ = consumes<reco::PFTauCollection>(pfTau.getUntrackedParameter<InputTag>("PFTauProducer"));
36  auto discs = pfTau.getUntrackedParameter<vector<InputTag>>("PFTauDiscriminators");
37  for (edm::InputTag& tag: discs) {
38  PFTauDis_.push_back(consumes<reco::PFTauDiscriminator>(tag));
39  }
40  doPFTaus_ = pfTau.getUntrackedParameter<bool>("doPFTaus",false);
41  ptMinPFTau_= pfTau.getUntrackedParameter<double>("ptMin",15.);
42  etaMinPFTau_= pfTau.getUntrackedParameter<double>("etaMin",-2.5);
43  etaMaxPFTau_= pfTau.getUntrackedParameter<double>("etaMax",2.5);
44  phiMinPFTau_= pfTau.getUntrackedParameter<double>("phiMin",-3.15);
45  phiMaxPFTau_= pfTau.getUntrackedParameter<double>("phiMax",3.15);
46  }
47 
48  {
49  auto const& electrons = iConfig.getUntrackedParameter<edm::ParameterSet>("Electrons");
50  Electrons_ = consumes<reco::GsfElectronCollection>(electrons.getUntrackedParameter<InputTag>("ElectronCollection"));
51  doElectrons_ = electrons.getUntrackedParameter<bool>("doElectrons",false);
52  e_ctfTrackCollectionSrc_ = electrons.getUntrackedParameter<InputTag>("TrackCollection");
53  e_ctfTrackCollection_ = consumes<reco::TrackCollection>(e_ctfTrackCollectionSrc_);
54  ptMinElectron_= electrons.getUntrackedParameter<double>("ptMin",15.);
55  e_doTrackIso_ = electrons.getUntrackedParameter<bool>("doTrackIso",false);
56  e_trackMinPt_= electrons.getUntrackedParameter<double>("ptMinTrack",1.5);
57  e_lipCut_= electrons.getUntrackedParameter<double>("lipMinTrack",1.5);
58  e_minIsoDR_= electrons.getUntrackedParameter<double>("InnerConeDR",0.02);
59  e_maxIsoDR_= electrons.getUntrackedParameter<double>("OuterConeDR",0.6);
60  e_isoMaxSumPt_= electrons.getUntrackedParameter<double>("MaxIsoVar",0.02);
61  }
62 
63  {
64  auto const& muons = iConfig.getUntrackedParameter<edm::ParameterSet>("Muons");
65  Muons_ = consumes<reco::MuonCollection>(muons.getUntrackedParameter<InputTag>("MuonCollection"));
66  doMuons_ = muons.getUntrackedParameter<bool>("doMuons",false);
67  ptMinMuon_= muons.getUntrackedParameter<double>("ptMin",15.);
68  }
69 
70  {
71  auto const& jets = iConfig.getUntrackedParameter<edm::ParameterSet>("Jets");
72  Jets_ = consumes<reco::CaloJetCollection>(jets.getUntrackedParameter<InputTag>("JetCollection"));
73  doJets_ = jets.getUntrackedParameter<bool>("doJets");
74  ptMinJet_= jets.getUntrackedParameter<double>("etMin");
75  }
76 
77  {
78  auto const& towers = iConfig.getUntrackedParameter<edm::ParameterSet>("Towers");
79  Towers_ = consumes<CaloTowerCollection>(towers.getUntrackedParameter<InputTag>("TowerCollection"));
80  doTowers_ = towers.getUntrackedParameter<bool>("doTowers");
81  ptMinTower_= towers.getUntrackedParameter<double>("etMin");
82  towerIsol_= towers.getUntrackedParameter<double>("towerIsolation");
83  }
84 
85  {
86  auto const& photons = iConfig.getUntrackedParameter<edm::ParameterSet>("Photons");
87  Photons_ = consumes<reco::PhotonCollection>(photons.getUntrackedParameter<InputTag>("PhotonCollection"));
88  doPhotons_ = photons.getUntrackedParameter<bool>("doPhotons");
89  ptMinPhoton_= photons.getUntrackedParameter<double>("etMin");
90  photonEcalIso_= photons.getUntrackedParameter<double>("ECALIso");
91  }
92 
93  {
94  auto const& met = iConfig.getUntrackedParameter<edm::ParameterSet>("MET");
95  MET_ = consumes<reco::CaloMETCollection>(met.getUntrackedParameter<InputTag>("METCollection"));
96  doMET_ = met.getUntrackedParameter<bool>("doMET",false);
97  ptMinMET_= met.getUntrackedParameter<double>("ptMin",15.);
98  }
99 
100  etaMin_ = iConfig.getUntrackedParameter<double>("EtaMin",-2.5);
101  etaMax_ = iConfig.getUntrackedParameter<double>("EtaMax",2.5);
102  phiMin_ = iConfig.getUntrackedParameter<double>("PhiMin",-3.15);
103  phiMax_ = iConfig.getUntrackedParameter<double>("PhiMax",3.15);
104 
105  //recoCollections
106  produces<LorentzVectorCollection>("PFTaus");
107  produces<LorentzVectorCollection>("Electrons");
108  produces<LorentzVectorCollection>("Muons");
109  produces<LorentzVectorCollection>("Jets");
110  produces<LorentzVectorCollection>("Photons");
111  produces<LorentzVectorCollection>("Towers");
112  produces<LorentzVectorCollection>("MET");
113 
114 }
115 
117 {
118  if(doPFTaus_)
119  doPFTaus(iEvent);
120  if(doElectrons_)
121  doElectrons(iEvent);
122  if(doMuons_)
123  doMuons(iEvent);
124  if(doJets_)
125  doJets(iEvent);
126  if(doPhotons_)
127  doPhotons(iEvent);
128  if(doTowers_)
129  doTowers(iEvent);
130  if(doMET_)
131  doMET(iEvent);
132 }
133 
134 void
136 {
137  auto product_PFTaus = make_unique<LorentzVectorCollection>();
138 
140  if (iEvent.getByToken(PFTaus_,pftaus)) {
141  for (unsigned int i=0; i<pftaus->size(); ++i) {
142  auto const& pftau = (*pftaus)[i];
143  if (pftau.pt()>ptMinPFTau_&&
144  pftau.eta()>etaMinPFTau_&&pftau.eta()<etaMaxPFTau_&&
145  pftau.phi()>phiMinPFTau_&&pftau.phi()<phiMaxPFTau_) {
146  reco::PFTauRef thePFTau {pftaus,i};
147  bool passAll {true};
148  for (edm::EDGetTokenT<reco::PFTauDiscriminator> const& token: PFTauDis_) {
150  if (iEvent.getByToken(token, pftaudis)) {
151  if ((*pftaudis)[thePFTau] < 0.5) {
152  passAll = false;
153  break;
154  }
155  }
156  else{
157  passAll = false;
158  break;
159  }
160  }
161  if (passAll) {
162  product_PFTaus->emplace_back(pftau.px(), pftau.py(), pftau.pz(), pftau.energy());
163  }
164  }
165  }
166  }
167  iEvent.put(move(product_PFTaus),"PFTaus");
168 }
169 
170 
171 void
173 {
174  auto product_Electrons = make_unique<LorentzVectorCollection>();
175 
177  if (!iEvent.getByToken(e_ctfTrackCollection_, pCtfTracks)) {
178  edm::LogInfo("") << "Error! Can't get " << e_ctfTrackCollectionSrc_.label() << " by label. ";
179  iEvent.put(move(product_Electrons),"Electrons");
180  return;
181  }
182 
184  if (iEvent.getByToken(Electrons_,electrons)) {
185  for (size_t i=0 ; i<electrons->size(); ++i) {
186  edm::Ref<reco::GsfElectronCollection> electronRef (electrons,i);
187  auto const& electron = (*electrons)[i];
188  if (electron.pt()>ptMinElectron_&&fabs(electron.eta())<etaMax_) {
189  if (e_doTrackIso_) {
190  double sum_of_pt_ele {};
191  for (auto const& tr : *pCtfTracks) {
192  double const lip {electron.gsfTrack()->dz() - tr.dz()};
193  if (tr.pt() > e_trackMinPt_ && fabs(lip) < e_lipCut_) {
194  double dphi {fabs(tr.phi()-electron.trackMomentumAtVtx().phi())};
195  if (dphi>acos(-1.)) {
196  dphi=2*acos(-1.)-dphi;
197  }
198  double const deta {fabs(tr.eta()-electron.trackMomentumAtVtx().eta())};
199  double const dr_ctf_ele {sqrt(deta*deta+dphi*dphi)};
200  if((dr_ctf_ele>e_minIsoDR_) && (dr_ctf_ele<e_maxIsoDR_)){
201  double const cft_pt_2 {tr.pt()*tr.pt()};
202  sum_of_pt_ele += cft_pt_2;
203  }
204  }
205  }
206  double const isolation_value_ele {sum_of_pt_ele/(electron.trackMomentumAtVtx().Rho()*electron.trackMomentumAtVtx().Rho())};
207  if (isolation_value_ele<e_isoMaxSumPt_) {
208  product_Electrons->emplace_back(electron.px(),electron.py(),electron.pz(),electron.energy());
209  }
210  }
211  else {
212  product_Electrons->emplace_back(electron.px(),electron.py(),electron.pz(),electron.energy());
213  }
214  }
215  }
216  }
217  iEvent.put(move(product_Electrons),"Electrons");
218 }
219 
220 void
222 {
223  auto product_Muons = make_unique<LorentzVectorCollection>();
224 
226  if (iEvent.getByToken(Muons_,muons)) {
227  for (auto const& muon : *muons) {
228  if (muon.pt()>ptMinMuon_ && muon.eta()>etaMin_&&muon.eta()<etaMax_&&muon.phi()>phiMin_&&muon.phi()<phiMax_) {
229  product_Muons->emplace_back(muon.px(),muon.py(),muon.pz(),muon.energy());
230  }
231  }
232  }
233  iEvent.put(move(product_Muons),"Muons");
234 }
235 
236 
237 void
239 {
240  auto product_Jets = make_unique<LorentzVectorCollection>();
241 
243  if (iEvent.getByToken(Jets_,jets)) {
244  for (auto const& jet : *jets) {
245  if (jet.et()>ptMinJet_ && jet.eta()>etaMin_&&jet.eta()<etaMax_&&jet.phi()>phiMin_&&jet.phi()<phiMax_) {
246  product_Jets->emplace_back(jet.px(),jet.py(),jet.pz(),jet.energy());
247  }
248  }
249  }
250  iEvent.put(move(product_Jets),"Jets");
251 }
252 
253 void
255 {
256  auto product_Towers = make_unique<LorentzVectorCollection>();
257 
259  if (iEvent.getByToken(Towers_,towers)) {
260  for (auto const& tower1 : *towers) {
261  if (tower1.pt()>ptMinTower_ && tower1.eta()>etaMin_&&tower1.eta()<etaMax_&&tower1.phi()>phiMin_&&tower1.phi()<phiMax_) {
262  //calculate isolation
263  double isolET {};
264  for (auto const& tower2 : *towers) {
265  if (ROOT::Math::VectorUtil::DeltaR(tower1.p4(),tower2.p4())<0.5) {
266  isolET+=tower2.pt();
267  }
268  isolET-=tower1.pt();
269  }
270  if (isolET<towerIsol_) {
271  product_Towers->emplace_back(tower1.px(),tower1.py(),tower1.pz(),tower1.energy());
272  }
273  }
274  }
275  }
276  iEvent.put(move(product_Towers),"Towers");
277 }
278 
279 
280 void
282 {
283  auto product_Gammas = make_unique<LorentzVectorCollection>();
284 
286  if (iEvent.getByToken(Photons_,photons)) {
287  for (auto const& photon : *photons) {
288  if (photon.ecalRecHitSumEtConeDR04()<photonEcalIso_ &&
289  photon.et()>ptMinPhoton_ && photon.eta()>etaMin_&&photon.eta()<etaMax_&&photon.phi()>phiMin_&&photon.phi()<phiMax_) {
290  product_Gammas->emplace_back(photon.px(),photon.py(),photon.pz(),photon.energy());
291  }
292  }
293  }
294  iEvent.put(move(product_Gammas),"Photons");
295 }
296 
297 void
299 {
300  auto product_MET = make_unique<LorentzVectorCollection>();
301 
303  if(iEvent.getByToken(MET_,met) && !met->empty()){
304  auto const& metMom = met->front().p4();
305  product_MET->emplace_back(metMom.Px(), metMom.Py(), 0, metMom.Pt());
306  }
307  iEvent.put(move(product_MET),"MET");
308 }
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
void produce(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
void doElectrons(edm::Event &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void doMET(edm::Event &) const
HLTTauRefProducer(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:18
vector< PseudoJet > jets
void doTowers(edm::Event &) const
met
===> hadronic RAZOR
void doPFTaus(edm::Event &) const
fixed size matrix
HLT enums.
void doMuons(edm::Event &) const
void doJets(edm::Event &) const
void doPhotons(edm::Event &) const
def move(src, dest)
Definition: eostools.py:511