CMS 3D CMS Logo

IsolatedPixelTrackCandidateL1TProducer.cc
Go to the documentation of this file.
1 /* \class IsolatedPixelTrackCandidateL1TProducer
2  *
3  *
4  */
5 
6 #include <vector>
7 #include <memory>
8 #include <algorithm>
9 #include <cmath>
10 
17 // L1Extra
20 // l1t
29 //vertices
32 
33 // Framework
44 
45 // Math
46 #include "Math/GenVector/VectorUtil.h"
47 #include "Math/GenVector/PxPyPzE4D.h"
49 
50 //magF
54 
55 //for ECAL geometry
61 
62 //#define EDM_ML_DEBUG
63 
65 public:
68 
69  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
70 
71  void beginRun(const edm::Run&, const edm::EventSetup&) override;
72  void produce(edm::Event& evt, const edm::EventSetup& es) override;
73 
74  double getDistInCM(double eta1, double phi1, double eta2, double phi2);
75  std::pair<double, double> GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ);
76 
77 private:
78  struct seedAtEC {
79  seedAtEC(unsigned int i, bool f, double et, double fi) : index(i), ok(f), eta(et), phi(fi) {}
80  unsigned int index;
81  bool ok;
82  double eta, phi;
83  };
84 
88  const std::vector<edm::EDGetTokenT<reco::TrackCollection> > toks_pix_;
91 
93  const double prelimCone_;
95  const double vtxCutSeed_;
96  const double vtxCutIsol_;
97  const double tauAssocCone_;
98  const double tauUnbiasCone_;
99  const double minPTrackValue_;
101  const double ebEtaBoundary_;
102 
103  // these are read from the EventSetup, cannot be const
104  double rEB_;
105  double zEE_;
106  double bfVal_;
107 };
108 
110  : tok_hlt_(consumes<trigger::TriggerFilterObjectWithRefs>(config.getParameter<edm::InputTag>("L1GTSeedLabel"))),
111  tok_l1_(consumes<l1t::TauBxCollection>(config.getParameter<edm::InputTag>("L1eTauJetsSource"))),
112  tok_vert_(consumes<reco::VertexCollection>(config.getParameter<edm::InputTag>("VertexLabel"))),
113  toks_pix_(
114  edm::vector_transform(config.getParameter<std::vector<edm::InputTag> >("PixelTracksSources"),
115  [this](edm::InputTag const& tag) { return consumes<reco::TrackCollection>(tag); })),
116  tok_bFieldH_(esConsumes<MagneticField, IdealMagneticFieldRecord, edm::Transition::BeginRun>()),
117  tok_geom_(esConsumes<CaloGeometry, CaloGeometryRecord, edm::Transition::BeginRun>()),
118  bfield_(config.getParameter<std::string>("MagFieldRecordName")),
119  prelimCone_(config.getParameter<double>("ExtrapolationConeSize")),
120  pixelIsolationConeSizeAtEC_(config.getParameter<double>("PixelIsolationConeSizeAtEC")),
121  vtxCutSeed_(config.getParameter<double>("MaxVtxDXYSeed")),
122  vtxCutIsol_(config.getParameter<double>("MaxVtxDXYIsol")),
123  tauAssocCone_(config.getParameter<double>("tauAssociationCone")),
124  tauUnbiasCone_(config.getParameter<double>("tauUnbiasCone")),
125  minPTrackValue_(config.getParameter<double>("minPTrack")),
126  maxPForIsolationValue_(config.getParameter<double>("maxPTrackForIsolation")),
127  ebEtaBoundary_(config.getParameter<double>("EBEtaBoundary")),
128  rEB_(-1),
129  zEE_(-1),
130  bfVal_(0) {
131  // Register the product
132  produces<reco::IsolatedPixelTrackCandidateCollection>();
133 }
134 
136 
139  std::vector<edm::InputTag> tracksrc = {edm::InputTag("hltPixelTracks")};
140  desc.add<edm::InputTag>("L1eTauJetsSource", edm::InputTag("hltGtStage2Digis", "Tau"));
141  desc.add<double>("tauAssociationCone", 0.0);
142  desc.add<double>("tauUnbiasCone", 1.2);
143  desc.add<std::vector<edm::InputTag> >("PixelTracksSources", tracksrc);
144  desc.add<double>("ExtrapolationConeSize", 1.0);
145  desc.add<double>("PixelIsolationConeSizeAtEC", 40);
146  desc.add<edm::InputTag>("L1GTSeedLabel", edm::InputTag("hltL1sV0SingleJet60"));
147  desc.add<double>("MaxVtxDXYSeed", 101.0);
148  desc.add<double>("MaxVtxDXYIsol", 101.0);
149  desc.add<edm::InputTag>("VertexLabel", edm::InputTag("hltTrimmedPixelVertices"));
150  desc.add<std::string>("MagFieldRecordName", "VolumeBasedMagneticField");
151  desc.add<double>("minPTrack", 5.0);
152  desc.add<double>("maxPTrackForIsolation", 3.0);
153  desc.add<double>("EBEtaBoundary", 1.479);
154  descriptions.add("isolPixelTrackProdL1T", desc);
155 }
156 
158  const CaloGeometry* pG = &theEventSetup.getData(tok_geom_);
159 
160  const double rad(dynamic_cast<const EcalBarrelGeometry*>(pG->getSubdetectorGeometry(DetId::Ecal, EcalBarrel))
161  ->avgRadiusXYFrontFaceCenter());
162  const double zz(dynamic_cast<const EcalEndcapGeometry*>(pG->getSubdetectorGeometry(DetId::Ecal, EcalEndcap))
163  ->avgAbsZFrontFaceCenter());
164 
165  rEB_ = rad;
166  zEE_ = zz;
167 
168  const MagneticField* vbfField = &theEventSetup.getData(tok_bFieldH_);
169  const VolumeBasedMagneticField* vbfCPtr = dynamic_cast<const VolumeBasedMagneticField*>(&(*vbfField));
170  GlobalVector BField = vbfCPtr->inTesla(GlobalPoint(0, 0, 0));
171  bfVal_ = BField.mag();
172 #ifdef EDM_ML_DEBUG
173  edm::LogVerbatim("IsoTrack") << "rEB " << rEB_ << " zEE " << zEE_ << " B " << bfVal_;
174 #endif
175 }
176 
178  auto trackCollection = std::make_unique<reco::IsolatedPixelTrackCandidateCollection>();
179 
180  //create vector of refs from input collections
181  std::vector<reco::TrackRef> pixelTrackRefs;
182 
183  for (unsigned int iPix = 0; iPix < toks_pix_.size(); iPix++) {
184  const edm::Handle<reco::TrackCollection>& iPixCol = theEvent.getHandle(toks_pix_[iPix]);
185  for (reco::TrackCollection::const_iterator pit = iPixCol->begin(); pit != iPixCol->end(); pit++) {
186  pixelTrackRefs.push_back(reco::TrackRef(iPixCol, pit - iPixCol->begin()));
187  }
188  }
189 
190  const edm::Handle<l1t::TauBxCollection>& l1eTauJets = theEvent.getHandle(tok_l1_);
191 
193 
194  double ptTriggered = -10;
195  double etaTriggered = -100;
196  double phiTriggered = -100;
197 
199 
200  std::vector<edm::Ref<l1t::TauBxCollection> > l1tauobjref;
201  std::vector<edm::Ref<l1t::JetBxCollection> > l1jetobjref;
202 
203  // l1trigobj->getObjects(trigger::TriggerTau, l1tauobjref);
204  l1trigobj->getObjects(trigger::TriggerL1Tau, l1tauobjref);
205  // l1trigobj->getObjects(trigger::TriggerJet, l1jetobjref);
206  l1trigobj->getObjects(trigger::TriggerL1Jet, l1jetobjref);
207 
208  for (const auto& p : l1tauobjref) {
209  if (p->pt() > ptTriggered) {
210  ptTriggered = p->pt();
211  phiTriggered = p->phi();
212  etaTriggered = p->eta();
213  }
214  }
215  for (const auto& p : l1jetobjref) {
216  if (p->pt() > ptTriggered) {
217  ptTriggered = p->pt();
218  phiTriggered = p->phi();
219  etaTriggered = p->eta();
220  }
221  }
222 #ifdef EDM_ML_DEBUG
223  edm::LogVerbatim("IsoTrack") << "Sizes " << l1tauobjref.size() << ":" << l1jetobjref.size() << " Trig " << ptTriggered
224  << ":" << etaTriggered << ":" << phiTriggered;
225 #endif
226  double drMaxL1Track_ = tauAssocCone_;
227 #ifdef EDM_ML_DEBUG
228  int ntr = 0;
229 #endif
230  std::vector<seedAtEC> VecSeedsatEC;
231  //loop to select isolated tracks
232  for (unsigned iS = 0; iS < pixelTrackRefs.size(); iS++) {
233  bool vtxMatch = false;
234  //associate to vertex (in Z)
235  reco::VertexCollection::const_iterator vitSel;
236  double minDZ = 1000;
237  bool found(false);
238  for (reco::VertexCollection::const_iterator vit = pVert->begin(); vit != pVert->end(); vit++) {
239  if (std::abs(pixelTrackRefs[iS]->dz(vit->position())) < minDZ) {
240  minDZ = std::abs(pixelTrackRefs[iS]->dz(vit->position()));
241  vitSel = vit;
242  found = true;
243  }
244  }
245  //cut on dYX:
246  if (found) {
247  if (std::abs(pixelTrackRefs[iS]->dxy(vitSel->position())) < vtxCutSeed_)
248  vtxMatch = true;
249  } else {
250  vtxMatch = true;
251  }
252 #ifdef EDM_ML_DEBUG
253  edm::LogVerbatim("IsoTrack") << "minZD " << minDZ << " Found " << found << ":" << vtxMatch;
254 #endif
255  //select tracks not matched to triggered L1 jet
256  double R = reco::deltaR(etaTriggered, phiTriggered, pixelTrackRefs[iS]->eta(), pixelTrackRefs[iS]->phi());
257 #ifdef EDM_ML_DEBUG
258  edm::LogVerbatim("IsoTrack") << "Distance to L1 " << R << ":" << tauUnbiasCone_ << " Result "
259  << (R < tauUnbiasCone_);
260 #endif
261  if (R < tauUnbiasCone_)
262  continue;
263 
264  //check taujet matching
265  bool tmatch = false;
267  for (l1t::TauBxCollection::const_iterator tj = l1eTauJets->begin(); tj != l1eTauJets->end(); tj++) {
268  if (reco::deltaR(pixelTrackRefs[iS]->momentum().eta(),
269  pixelTrackRefs[iS]->momentum().phi(),
270  tj->momentum().eta(),
271  tj->momentum().phi()) > drMaxL1Track_)
272  continue;
273  selj = tj;
274  tmatch = true;
275  } //loop over L1 tau
276 #ifdef EDM_ML_DEBUG
277  edm::LogVerbatim("IsoTrack") << "tMatch " << tmatch;
278 #endif
279  //propagate seed track to ECAL surface:
280  std::pair<double, double> seedCooAtEC;
281  // in case vertex is found:
282  if (found)
283  seedCooAtEC = GetEtaPhiAtEcal(pixelTrackRefs[iS]->eta(),
284  pixelTrackRefs[iS]->phi(),
285  pixelTrackRefs[iS]->pt(),
286  pixelTrackRefs[iS]->charge(),
287  vitSel->z());
288  //in case vertex is not found:
289  else
290  seedCooAtEC = GetEtaPhiAtEcal(pixelTrackRefs[iS]->eta(),
291  pixelTrackRefs[iS]->phi(),
292  pixelTrackRefs[iS]->pt(),
293  pixelTrackRefs[iS]->charge(),
294  0);
295  seedAtEC seed(iS, (tmatch || vtxMatch), seedCooAtEC.first, seedCooAtEC.second);
296  VecSeedsatEC.push_back(seed);
297 #ifdef EDM_ML_DEBUG
298  edm::LogVerbatim("IsoTrack") << "Seed " << seedCooAtEC.first << seedCooAtEC.second;
299 #endif
300  }
301  for (unsigned int i = 0; i < VecSeedsatEC.size(); i++) {
302  unsigned int iSeed = VecSeedsatEC[i].index;
303  if (!VecSeedsatEC[i].ok)
304  continue;
305  if (pixelTrackRefs[iSeed]->p() < minPTrackValue_)
306  continue;
308  for (l1t::TauBxCollection::const_iterator tj = l1eTauJets->begin(); tj != l1eTauJets->end(); tj++) {
309  if (reco::deltaR(pixelTrackRefs[iSeed]->momentum().eta(),
310  pixelTrackRefs[iSeed]->momentum().phi(),
311  tj->momentum().eta(),
312  tj->momentum().phi()) > drMaxL1Track_)
313  continue;
314  selj = tj;
315  } //loop over L1 tau
316  double maxP = 0;
317  double sumP = 0;
318  for (unsigned int j = 0; j < VecSeedsatEC.size(); j++) {
319  if (i == j)
320  continue;
321  unsigned int iSurr = VecSeedsatEC[j].index;
322  //define preliminary cone around seed track impact point from which tracks will be extrapolated:
323  if (reco::deltaR(pixelTrackRefs[iSeed]->eta(),
324  pixelTrackRefs[iSeed]->phi(),
325  pixelTrackRefs[iSurr]->eta(),
326  pixelTrackRefs[iSurr]->phi()) > prelimCone_)
327  continue;
328  double minDZ2(1000);
329  bool found(false);
330  reco::VertexCollection::const_iterator vitSel2;
331  for (reco::VertexCollection::const_iterator vit = pVert->begin(); vit != pVert->end(); vit++) {
332  if (std::abs(pixelTrackRefs[iSurr]->dz(vit->position())) < minDZ2) {
333  minDZ2 = std::abs(pixelTrackRefs[iSurr]->dz(vit->position()));
334  vitSel2 = vit;
335  found = true;
336  }
337  }
338  //cut ot dXY:
339  if (found && std::abs(pixelTrackRefs[iSurr]->dxy(vitSel2->position())) > vtxCutIsol_)
340  continue;
341  //calculate distance at ECAL surface and update isolation:
342  if (getDistInCM(VecSeedsatEC[i].eta, VecSeedsatEC[i].phi, VecSeedsatEC[j].eta, VecSeedsatEC[j].phi) <
344  sumP += pixelTrackRefs[iSurr]->p();
345  if (pixelTrackRefs[iSurr]->p() > maxP)
346  maxP = pixelTrackRefs[iSurr]->p();
347  }
348  }
351  pixelTrackRefs[iSeed], l1t::TauRef(l1eTauJets, selj - l1eTauJets->begin()), maxP, sumP);
352  newCandidate.setEtaPhiEcal(VecSeedsatEC[i].eta, VecSeedsatEC[i].phi);
353  trackCollection->push_back(newCandidate);
354 #ifdef EDM_ML_DEBUG
355  ntr++;
356 #endif
357  }
358  }
359 #ifdef EDM_ML_DEBUG
360  edm::LogVerbatim("IsoTrack") << "Number of Isolated Track " << ntr;
361 #endif
362  // put the product in the event
363  theEvent.put(std::move(trackCollection));
364 }
365 
366 double IsolatedPixelTrackCandidateL1TProducer::getDistInCM(double eta1, double phi1, double eta2, double phi2) {
367  double Rec;
368  double theta1 = 2 * atan(exp(-eta1));
369  double theta2 = 2 * atan(exp(-eta2));
370  if (std::abs(eta1) < 1.479)
371  Rec = rEB_; //radius of ECAL barrel
372  else if (std::abs(eta1) > 1.479 && std::abs(eta1) < 7.0)
373  Rec = tan(theta1) * zEE_; //distance from IP to ECAL endcap
374  else
375  return 1000;
376 
377  //|vect| times tg of acos(scalar product)
378  double angle =
379  acos((sin(theta1) * sin(theta2) * (sin(phi1) * sin(phi2) + cos(phi1) * cos(phi2)) + cos(theta1) * cos(theta2)));
380  if (angle < M_PI_2)
381  return std::abs((Rec / sin(theta1)) * tan(angle));
382  else
383  return 1000;
384 }
385 
387  double etaIP, double phiIP, double pT, int charge, double vtxZ) {
388  double deltaPhi = 0;
389  double etaEC = 100;
390  double phiEC = 100;
391 
392  double Rcurv = 9999999;
393  if (bfVal_ != 0)
394  Rcurv = pT * 33.3 * 100 / (bfVal_ * 10); //r(m)=pT(GeV)*33.3/B(kG)
395 
396  double ecDist = zEE_; //distance to ECAL andcap from IP (cm), 317 - ecal (not preshower), preshower -300
397  double ecRad = rEB_; //radius of ECAL barrel (cm)
398  double theta = 2 * atan(exp(-etaIP));
399  double zNew = 0;
400  if (theta > M_PI_2)
401  theta = M_PI - theta;
402  if (std::abs(etaIP) < ebEtaBoundary_) {
403  if ((0.5 * ecRad / Rcurv) > 1) {
404  etaEC = 10000;
405  deltaPhi = 0;
406  } else {
407  deltaPhi = -charge * asin(0.5 * ecRad / Rcurv);
408  double alpha1 = 2 * asin(0.5 * ecRad / Rcurv);
409  double z = ecRad / tan(theta);
410  if (etaIP > 0)
411  zNew = z * (Rcurv * alpha1) / ecRad + vtxZ; //new z-coordinate of track
412  else
413  zNew = -z * (Rcurv * alpha1) / ecRad + vtxZ; //new z-coordinate of track
414  double zAbs = std::abs(zNew);
415  if (zAbs < ecDist) {
416  etaEC = -log(tan(0.5 * atan(ecRad / zAbs)));
417  deltaPhi = -charge * asin(0.5 * ecRad / Rcurv);
418  }
419  if (zAbs > ecDist) {
420  zAbs = (std::abs(etaIP) / etaIP) * ecDist;
421  double Zflight = std::abs(zAbs - vtxZ);
422  double alpha = (Zflight * ecRad) / (z * Rcurv);
423  double Rec = 2 * Rcurv * sin(alpha / 2);
424  deltaPhi = -charge * alpha / 2;
425  etaEC = -log(tan(0.5 * atan(Rec / ecDist)));
426  }
427  }
428  } else {
429  zNew = (std::abs(etaIP) / etaIP) * ecDist;
430  double Zflight = std::abs(zNew - vtxZ);
431  double Rvirt = std::abs(Zflight * tan(theta));
432  double Rec = 2 * Rcurv * sin(Rvirt / (2 * Rcurv));
433  deltaPhi = -(charge) * (Rvirt / (2 * Rcurv));
434  etaEC = -log(tan(0.5 * atan(Rec / ecDist)));
435  }
436 
437  if (zNew < 0)
438  etaEC = -etaEC;
439  phiEC = phiIP + deltaPhi;
440 
441  if (phiEC < -M_PI)
442  phiEC += M_2_PI;
443  if (phiEC > M_PI)
444  phiEC -= M_2_PI;
445 
446  std::pair<double, double> retVal(etaEC, phiEC);
447  return retVal;
448 }
451 
Log< level::Info, true > LogVerbatim
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
void produce(edm::Event &evt, const edm::EventSetup &es) override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
GlobalVector inTesla(const GlobalPoint &g) const override
Field value ad specified global point, in Tesla.
void beginRun(const edm::Run &, const edm::EventSetup &) override
double getDistInCM(double eta1, double phi1, double eta2, double phi2)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
#define M_PI_2
delete x;
Definition: CaloConfig.h:22
Definition: config.py:1
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
BXVector< Tau > TauBxCollection
Definition: Tau.h:10
const edm::EDGetTokenT< reco::VertexCollection > tok_vert_
const_iterator begin(int bx) const
std::vector< T >::const_iterator const_iterator
Definition: BXVector.h:18
const edm::EDGetTokenT< l1t::TauBxCollection > tok_l1_
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
void setEtaPhiEcal(double eta, double phi)
eta, phi at ECAL surface
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::pair< double, double > GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ)
trackCollection
Definition: JetHT_cfg.py:51
#define M_PI
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const_iterator end(int bx) const
fixed size matrix
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:552
HLT enums.
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > tok_bFieldH_
Geom::Theta< T > theta() const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tok_hlt_