CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Attributes
IsolatedPixelTrackCandidateProducer Class Reference
Inheritance diagram for IsolatedPixelTrackCandidateProducer:
edm::stream::EDProducer<>

Classes

struct  seedAtEC
 

Public Member Functions

void beginRun (const edm::Run &, const edm::EventSetup &) override
 
double getDistInCM (double eta1, double phi1, double eta2, double phi2)
 
std::pair< double, double > GetEtaPhiAtEcal (double etaIP, double phiIP, double pT, int charge, double vtxZ)
 
 IsolatedPixelTrackCandidateProducer (const edm::ParameterSet &ps)
 
void produce (edm::Event &evt, const edm::EventSetup &es) override
 
 ~IsolatedPixelTrackCandidateProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

const std::string bfield_
 
double bfVal_
 
const double ebEtaBoundary_
 
const double maxPForIsolationValue_
 
const double minPTrackValue_
 
const double pixelIsolationConeSizeAtEC_
 
const double prelimCone_
 
double rEB_
 
const double tauAssocCone_
 
const double tauUnbiasCone_
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordtok_bFieldH_
 
const edm::ESGetToken< CaloGeometry, CaloGeometryRecordtok_geom_
 
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefstok_hlt_
 
const edm::EDGetTokenT< l1extra::L1JetParticleCollectiontok_l1_
 
const edm::EDGetTokenT< reco::VertexCollectiontok_vert_
 
const std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
 
const double vtxCutIsol_
 
const double vtxCutSeed_
 
double zEE_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 63 of file IsolatedPixelTrackCandidateProducer.cc.

Constructor & Destructor Documentation

◆ IsolatedPixelTrackCandidateProducer()

IsolatedPixelTrackCandidateProducer::IsolatedPixelTrackCandidateProducer ( const edm::ParameterSet ps)

Definition at line 108 of file IsolatedPixelTrackCandidateProducer.cc.

References makeGlobalPositionRcd_cfg::tag.

109  : tok_hlt_(consumes<trigger::TriggerFilterObjectWithRefs>(config.getParameter<edm::InputTag>("L1GTSeedLabel"))),
110  tok_l1_(consumes<l1extra::L1JetParticleCollection>(config.getParameter<edm::InputTag>("L1eTauJetsSource"))),
111  tok_vert_(consumes<reco::VertexCollection>(config.getParameter<edm::InputTag>("VertexLabel"))),
112  toks_pix_(
113  edm::vector_transform(config.getParameter<std::vector<edm::InputTag> >("PixelTracksSources"),
114  [this](edm::InputTag const& tag) { return consumes<reco::TrackCollection>(tag); })),
115  tok_bFieldH_(esConsumes<MagneticField, IdealMagneticFieldRecord, edm::Transition::BeginRun>()),
116  tok_geom_(esConsumes<CaloGeometry, CaloGeometryRecord, edm::Transition::BeginRun>()),
117  bfield_(config.getParameter<std::string>("MagFieldRecordName")),
118  prelimCone_(config.getParameter<double>("ExtrapolationConeSize")),
119  pixelIsolationConeSizeAtEC_(config.getParameter<double>("PixelIsolationConeSizeAtEC")),
120  vtxCutSeed_(config.getParameter<double>("MaxVtxDXYSeed")),
121  vtxCutIsol_(config.getParameter<double>("MaxVtxDXYIsol")),
122  tauAssocCone_(config.getParameter<double>("tauAssociationCone")),
123  tauUnbiasCone_(config.getParameter<double>("tauUnbiasCone")),
124  minPTrackValue_(config.getParameter<double>("minPTrack")),
125  maxPForIsolationValue_(config.getParameter<double>("maxPTrackForIsolation")),
126  ebEtaBoundary_(config.getParameter<double>("EBEtaBoundary")),
127  rEB_(-1),
128  zEE_(-1),
129  bfVal_(0) {
130  // Register the product
131  produces<reco::IsolatedPixelTrackCandidateCollection>();
132 }
const std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tok_hlt_
const edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_l1_
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
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > tok_bFieldH_
const edm::EDGetTokenT< reco::VertexCollection > tok_vert_

◆ ~IsolatedPixelTrackCandidateProducer()

IsolatedPixelTrackCandidateProducer::~IsolatedPixelTrackCandidateProducer ( )
override

Definition at line 134 of file IsolatedPixelTrackCandidateProducer.cc.

134 {}

Member Function Documentation

◆ beginRun()

void IsolatedPixelTrackCandidateProducer::beginRun ( const edm::Run run,
const edm::EventSetup theEventSetup 
)
override

Definition at line 156 of file IsolatedPixelTrackCandidateProducer.cc.

References ProducerSetup_cfi::BField, bfVal_, DetId::Ecal, EcalBarrel, EcalEndcap, edm::EventSetup::getData(), CaloGeometry::getSubdetectorGeometry(), VolumeBasedMagneticField::inTesla(), rEB_, tok_bFieldH_, tok_geom_, zEE_, and geometryCSVtoXML::zz.

156  {
157  const CaloGeometry* pG = &theEventSetup.getData(tok_geom_);
158 
159  const double rad(dynamic_cast<const EcalBarrelGeometry*>(pG->getSubdetectorGeometry(DetId::Ecal, EcalBarrel))
160  ->avgRadiusXYFrontFaceCenter());
161  const double zz(dynamic_cast<const EcalEndcapGeometry*>(pG->getSubdetectorGeometry(DetId::Ecal, EcalEndcap))
162  ->avgAbsZFrontFaceCenter());
163 
164  rEB_ = rad;
165  zEE_ = zz;
166 
167  const MagneticField* vbfField = &theEventSetup.getData(tok_bFieldH_);
168  const VolumeBasedMagneticField* vbfCPtr = dynamic_cast<const VolumeBasedMagneticField*>(&(*vbfField));
169  GlobalVector BField = vbfCPtr->inTesla(GlobalPoint(0, 0, 0));
170  bfVal_ = BField.mag();
171 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tok_geom_
GlobalVector inTesla(const GlobalPoint &g) const override
Field value ad specified global point, in Tesla.
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > tok_bFieldH_
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34

◆ fillDescriptions()

void IsolatedPixelTrackCandidateProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 136 of file IsolatedPixelTrackCandidateProducer.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, AlCaHLTBitMon_QueryRunRegistry::string, and SiPixelMonitorTrack_cfi::tracksrc.

136  {
138  std::vector<edm::InputTag> tracksrc = {edm::InputTag("hltPixelTracks")};
139  desc.add<edm::InputTag>("L1eTauJetsSource", edm::InputTag("hltCaloStage2Digis", "Tau"));
140  desc.add<double>("tauAssociationCone", 0.0);
141  desc.add<double>("tauUnbiasCone", 1.2);
142  desc.add<std::vector<edm::InputTag> >("PixelTracksSources", tracksrc);
143  desc.add<double>("ExtrapolationConeSize", 1.0);
144  desc.add<double>("PixelIsolationConeSizeAtEC", 40);
145  desc.add<edm::InputTag>("L1GTSeedLabel", edm::InputTag("hltL1sIsoTrack"));
146  desc.add<double>("MaxVtxDXYSeed", 101.0);
147  desc.add<double>("MaxVtxDXYIsol", 101.0);
148  desc.add<edm::InputTag>("VertexLabel", edm::InputTag("hltTrimmedPixelVertices"));
149  desc.add<std::string>("MagFieldRecordName", "VolumeBasedMagneticField");
150  desc.add<double>("minPTrack", 5.0);
151  desc.add<double>("maxPTrackForIsolation", 3.0);
152  desc.add<double>("EBEtaBoundary", 1.479);
153  descriptions.add("isolPixelTrackProd", desc);
154 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ getDistInCM()

double IsolatedPixelTrackCandidateProducer::getDistInCM ( double  eta1,
double  phi1,
double  eta2,
double  phi2 
)

Definition at line 322 of file IsolatedPixelTrackCandidateProducer.cc.

References funct::abs(), angle(), funct::cos(), HLT_2023v12_cff::eta1, HLT_2023v12_cff::eta2, JetChargeProducer_cfi::exp, M_PI_2, rEB_, funct::sin(), funct::tan(), and zEE_.

Referenced by produce().

322  {
323  double Rec;
324  double theta1 = 2 * atan(exp(-eta1));
325  double theta2 = 2 * atan(exp(-eta2));
326  if (std::abs(eta1) < 1.479)
327  Rec = rEB_; //radius of ECAL barrel
328  else if (std::abs(eta1) > 1.479 && std::abs(eta1) < 7.0)
329  Rec = tan(theta1) * zEE_; //distance from IP to ECAL endcap
330  else
331  return 1000;
332 
333  //|vect| times tg of acos(scalar product)
334  double angle =
335  acos((sin(theta1) * sin(theta2) * (sin(phi1) * sin(phi2) + cos(phi1) * cos(phi2)) + cos(theta1) * cos(theta2)));
336  if (angle < M_PI_2)
337  return std::abs((Rec / sin(theta1)) * tan(angle));
338  else
339  return 1000;
340 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define M_PI_2
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11

◆ GetEtaPhiAtEcal()

std::pair< double, double > IsolatedPixelTrackCandidateProducer::GetEtaPhiAtEcal ( double  etaIP,
double  phiIP,
double  pT,
int  charge,
double  vtxZ 
)

Definition at line 342 of file IsolatedPixelTrackCandidateProducer.cc.

References funct::abs(), simBeamSpotPI::alpha, bfVal_, ALCARECOTkAlJpsiMuMu_cff::charge, SiPixelRawToDigiRegional_cfi::deltaPhi, ebEtaBoundary_, JetChargeProducer_cfi::exp, dqm-mbProfile::log, M_PI, M_PI_2, pv::pT, rEB_, funct::sin(), funct::tan(), theta(), z, and zEE_.

Referenced by produce().

343  {
344  double deltaPhi = 0;
345  double etaEC = 100;
346  double phiEC = 100;
347 
348  double Rcurv = 9999999;
349  if (bfVal_ != 0)
350  Rcurv = pT * 33.3 * 100 / (bfVal_ * 10); //r(m)=pT(GeV)*33.3/B(kG)
351 
352  double ecDist = zEE_; //distance to ECAL andcap from IP (cm), 317 - ecal (not preshower), preshower -300
353  double ecRad = rEB_; //radius of ECAL barrel (cm)
354  double theta = 2 * atan(exp(-etaIP));
355  double zNew = 0;
356  if (theta > M_PI_2)
357  theta = M_PI - theta;
358  if (std::abs(etaIP) < ebEtaBoundary_) {
359  if ((0.5 * ecRad / Rcurv) > 1) {
360  etaEC = 10000;
361  deltaPhi = 0;
362  } else {
363  deltaPhi = -charge * asin(0.5 * ecRad / Rcurv);
364  double alpha1 = 2 * asin(0.5 * ecRad / Rcurv);
365  double z = ecRad / tan(theta);
366  if (etaIP > 0)
367  zNew = z * (Rcurv * alpha1) / ecRad + vtxZ; //new z-coordinate of track
368  else
369  zNew = -z * (Rcurv * alpha1) / ecRad + vtxZ; //new z-coordinate of track
370  double zAbs = std::abs(zNew);
371  if (zAbs < ecDist) {
372  etaEC = -log(tan(0.5 * atan(ecRad / zAbs)));
373  deltaPhi = -charge * asin(0.5 * ecRad / Rcurv);
374  }
375  if (zAbs > ecDist) {
376  zAbs = (std::abs(etaIP) / etaIP) * ecDist;
377  double Zflight = std::abs(zAbs - vtxZ);
378  double alpha = (Zflight * ecRad) / (z * Rcurv);
379  double Rec = 2 * Rcurv * sin(alpha / 2);
380  deltaPhi = -charge * alpha / 2;
381  etaEC = -log(tan(0.5 * atan(Rec / ecDist)));
382  }
383  }
384  } else {
385  zNew = (std::abs(etaIP) / etaIP) * ecDist;
386  double Zflight = std::abs(zNew - vtxZ);
387  double Rvirt = std::abs(Zflight * tan(theta));
388  double Rec = 2 * Rcurv * sin(Rvirt / (2 * Rcurv));
389  deltaPhi = -(charge) * (Rvirt / (2 * Rcurv));
390  etaEC = -log(tan(0.5 * atan(Rec / ecDist)));
391  }
392 
393  if (zNew < 0)
394  etaEC = -etaEC;
395  phiEC = phiIP + deltaPhi;
396 
397  if (phiEC < -M_PI)
398  phiEC += M_2_PI;
399  if (phiEC > M_PI)
400  phiEC -= M_2_PI;
401 
402  std::pair<double, double> retVal(etaEC, phiEC);
403  return retVal;
404 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
#define M_PI_2
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define M_PI
Geom::Theta< T > theta() const

◆ produce()

void IsolatedPixelTrackCandidateProducer::produce ( edm::Event evt,
const edm::EventSetup es 
)
override

Definition at line 173 of file IsolatedPixelTrackCandidateProducer.cc.

References funct::abs(), ALCARECOTkAlJpsiMuMu_cff::charge, reco::deltaR(), PVValHelper::dxy, PVValHelper::dz, PVValHelper::eta, newFWLiteAna::found, getDistInCM(), GetEtaPhiAtEcal(), edm::Event::getHandle(), mps_fire::i, dqmiolumiharvest::j, RecoMuonValidator_cfi::maxP, maxPForIsolationValue_, minPTrackValue_, eostools::move(), convertSQLiteXML::ok, AlCaHLTBitMon_ParallelJobs::p, phi, pixelIsolationConeSizeAtEC_, prelimCone_, DiDispStaMuonMonitor_cfi::pt, edm::Event::put(), fileCollector::seed, reco::IsolatedPixelTrackCandidate::setEtaPhiEcal(), tauAssocCone_, tok_l1_, tok_vert_, toks_pix_, JetHT_cfg::trackCollection, vtxCutIsol_, and vtxCutSeed_.

173  {
174  auto trackCollection = std::make_unique<reco::IsolatedPixelTrackCandidateCollection>();
175 
176  //create vector of refs from input collections
177  std::vector<reco::TrackRef> pixelTrackRefs;
178 #ifdef EDM_ML_DEBUG
179  edm::LogVerbatim("HcalIsoTrack") << "IsolatedPixelTrakCandidate: with" << toks_pix_.size()
180  << " candidates to start with\n";
181 #endif
182  for (unsigned int iPix = 0; iPix < toks_pix_.size(); iPix++) {
183  const edm::Handle<reco::TrackCollection>& iPixCol = theEvent.getHandle(toks_pix_[iPix]);
184  for (reco::TrackCollection::const_iterator pit = iPixCol->begin(); pit != iPixCol->end(); pit++) {
185  pixelTrackRefs.push_back(reco::TrackRef(iPixCol, pit - iPixCol->begin()));
186  }
187  }
188 
189  const edm::Handle<l1extra::L1JetParticleCollection>& l1eTauJets = theEvent.getHandle(tok_l1_);
190 
191  const edm::Handle<reco::VertexCollection>& pVert = theEvent.getHandle(tok_vert_);
192 
193  double drMaxL1Track_ = tauAssocCone_;
194 #ifdef EDM_ML_DEBUG
195  int ntr = 0;
196 #endif
197  std::vector<seedAtEC> VecSeedsatEC;
198  //loop to select isolated tracks
199  for (unsigned iS = 0; iS < pixelTrackRefs.size(); iS++) {
200  bool vtxMatch = false;
201  //associate to vertex (in Z)
202  reco::VertexCollection::const_iterator vitSel;
203  double minDZ = 1000;
204  bool found(false);
205  for (reco::VertexCollection::const_iterator vit = pVert->begin(); vit != pVert->end(); vit++) {
206  if (std::abs(pixelTrackRefs[iS]->dz(vit->position())) < minDZ) {
207  minDZ = std::abs(pixelTrackRefs[iS]->dz(vit->position()));
208  vitSel = vit;
209  found = true;
210  }
211  }
212  //cut on dYX:
213  if (found) {
214  if (std::abs(pixelTrackRefs[iS]->dxy(vitSel->position())) < vtxCutSeed_)
215  vtxMatch = true;
216  } else {
217  vtxMatch = true;
218  }
219 
220  //check taujet matching
221  bool tmatch = false;
222  l1extra::L1JetParticleCollection::const_iterator selj;
223  for (l1extra::L1JetParticleCollection::const_iterator tj = l1eTauJets->begin(); tj != l1eTauJets->end(); tj++) {
224  if (reco::deltaR(pixelTrackRefs[iS]->momentum().eta(),
225  pixelTrackRefs[iS]->momentum().phi(),
226  tj->momentum().eta(),
227  tj->momentum().phi()) > drMaxL1Track_)
228  continue;
229  selj = tj;
230  tmatch = true;
231  } //loop over L1 tau
232 
233  //propagate seed track to ECAL surface:
234  std::pair<double, double> seedCooAtEC;
235  // in case vertex is found:
236  if (found)
237  seedCooAtEC = GetEtaPhiAtEcal(pixelTrackRefs[iS]->eta(),
238  pixelTrackRefs[iS]->phi(),
239  pixelTrackRefs[iS]->pt(),
240  pixelTrackRefs[iS]->charge(),
241  vitSel->z());
242  //in case vertex is not found:
243  else
244  seedCooAtEC = GetEtaPhiAtEcal(pixelTrackRefs[iS]->eta(),
245  pixelTrackRefs[iS]->phi(),
246  pixelTrackRefs[iS]->pt(),
247  pixelTrackRefs[iS]->charge(),
248  0);
249  seedAtEC seed(iS, (tmatch || vtxMatch), seedCooAtEC.first, seedCooAtEC.second);
250  VecSeedsatEC.push_back(seed);
251  }
252 #ifdef EDM_ML_DEBUG
253  edm::LogVerbatim("HcalIsoTrack") << "IsolatedPixelTrakCandidate: " << VecSeedsatEC.size()
254  << " seeds after propagation\n";
255 #endif
256 
257  for (unsigned int i = 0; i < VecSeedsatEC.size(); i++) {
258  unsigned int iSeed = VecSeedsatEC[i].index;
259  if (!VecSeedsatEC[i].ok)
260  continue;
261  if (pixelTrackRefs[iSeed]->p() < minPTrackValue_)
262  continue;
263  l1extra::L1JetParticleCollection::const_iterator selj;
264  for (l1extra::L1JetParticleCollection::const_iterator tj = l1eTauJets->begin(); tj != l1eTauJets->end(); tj++) {
265  if (reco::deltaR(pixelTrackRefs[iSeed]->momentum().eta(),
266  pixelTrackRefs[iSeed]->momentum().phi(),
267  tj->momentum().eta(),
268  tj->momentum().phi()) > drMaxL1Track_)
269  continue;
270  selj = tj;
271  } //loop over L1 tau
272  double maxP = 0;
273  double sumP = 0;
274  for (unsigned int j = 0; j < VecSeedsatEC.size(); j++) {
275  if (i == j)
276  continue;
277  unsigned int iSurr = VecSeedsatEC[j].index;
278  //define preliminary cone around seed track impact point from which tracks will be extrapolated:
279  if (reco::deltaR(pixelTrackRefs[iSeed]->eta(),
280  pixelTrackRefs[iSeed]->phi(),
281  pixelTrackRefs[iSurr]->eta(),
282  pixelTrackRefs[iSurr]->phi()) > prelimCone_)
283  continue;
284  double minDZ2(1000);
285  bool found(false);
286  reco::VertexCollection::const_iterator vitSel2;
287  for (reco::VertexCollection::const_iterator vit = pVert->begin(); vit != pVert->end(); vit++) {
288  if (std::abs(pixelTrackRefs[iSurr]->dz(vit->position())) < minDZ2) {
289  minDZ2 = std::abs(pixelTrackRefs[iSurr]->dz(vit->position()));
290  vitSel2 = vit;
291  found = true;
292  }
293  }
294  //cut ot dXY:
295  if (found && std::abs(pixelTrackRefs[iSurr]->dxy(vitSel2->position())) > vtxCutIsol_)
296  continue;
297  //calculate distance at ECAL surface and update isolation:
298  if (getDistInCM(VecSeedsatEC[i].eta, VecSeedsatEC[i].phi, VecSeedsatEC[j].eta, VecSeedsatEC[j].phi) <
300  sumP += pixelTrackRefs[iSurr]->p();
301  if (pixelTrackRefs[iSurr]->p() > maxP)
302  maxP = pixelTrackRefs[iSurr]->p();
303  }
304  }
307  pixelTrackRefs[iSeed], l1extra::L1JetParticleRef(l1eTauJets, selj - l1eTauJets->begin()), maxP, sumP);
308  newCandidate.setEtaPhiEcal(VecSeedsatEC[i].eta, VecSeedsatEC[i].phi);
309  trackCollection->push_back(newCandidate);
310 #ifdef EDM_ML_DEBUG
311  ntr++;
312 #endif
313  }
314  }
315  // put the product in the event
316  theEvent.put(std::move(trackCollection));
317 #ifdef EDM_ML_DEBUG
318  edm::LogVerbatim("HcalIsoTrack") << "IsolatedPixelTrackCandidate: Final # of candiates " << ntr << "\n";
319 #endif
320 }
std::pair< double, double > GetEtaPhiAtEcal(double etaIP, double phiIP, double pT, int charge, double vtxZ)
Log< level::Info, true > LogVerbatim
const std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
const edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_l1_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
trackCollection
Definition: JetHT_cfg.py:51
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
const edm::EDGetTokenT< reco::VertexCollection > tok_vert_
double getDistInCM(double eta1, double phi1, double eta2, double phi2)
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ bfield_

const std::string IsolatedPixelTrackCandidateProducer::bfield_
private

Definition at line 91 of file IsolatedPixelTrackCandidateProducer.cc.

◆ bfVal_

double IsolatedPixelTrackCandidateProducer::bfVal_
private

Definition at line 105 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by beginRun(), and GetEtaPhiAtEcal().

◆ ebEtaBoundary_

const double IsolatedPixelTrackCandidateProducer::ebEtaBoundary_
private

Definition at line 100 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by GetEtaPhiAtEcal().

◆ maxPForIsolationValue_

const double IsolatedPixelTrackCandidateProducer::maxPForIsolationValue_
private

Definition at line 99 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ minPTrackValue_

const double IsolatedPixelTrackCandidateProducer::minPTrackValue_
private

Definition at line 98 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ pixelIsolationConeSizeAtEC_

const double IsolatedPixelTrackCandidateProducer::pixelIsolationConeSizeAtEC_
private

Definition at line 93 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ prelimCone_

const double IsolatedPixelTrackCandidateProducer::prelimCone_
private

Definition at line 92 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ rEB_

double IsolatedPixelTrackCandidateProducer::rEB_
private

Definition at line 103 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by beginRun(), getDistInCM(), and GetEtaPhiAtEcal().

◆ tauAssocCone_

const double IsolatedPixelTrackCandidateProducer::tauAssocCone_
private

Definition at line 96 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ tauUnbiasCone_

const double IsolatedPixelTrackCandidateProducer::tauUnbiasCone_
private

Definition at line 97 of file IsolatedPixelTrackCandidateProducer.cc.

◆ tok_bFieldH_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> IsolatedPixelTrackCandidateProducer::tok_bFieldH_
private

Definition at line 88 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by beginRun().

◆ tok_geom_

const edm::ESGetToken<CaloGeometry, CaloGeometryRecord> IsolatedPixelTrackCandidateProducer::tok_geom_
private

Definition at line 89 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by beginRun().

◆ tok_hlt_

const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> IsolatedPixelTrackCandidateProducer::tok_hlt_
private

Definition at line 84 of file IsolatedPixelTrackCandidateProducer.cc.

◆ tok_l1_

const edm::EDGetTokenT<l1extra::L1JetParticleCollection> IsolatedPixelTrackCandidateProducer::tok_l1_
private

Definition at line 85 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ tok_vert_

const edm::EDGetTokenT<reco::VertexCollection> IsolatedPixelTrackCandidateProducer::tok_vert_
private

Definition at line 86 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ toks_pix_

const std::vector<edm::EDGetTokenT<reco::TrackCollection> > IsolatedPixelTrackCandidateProducer::toks_pix_
private

Definition at line 87 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ vtxCutIsol_

const double IsolatedPixelTrackCandidateProducer::vtxCutIsol_
private

Definition at line 95 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ vtxCutSeed_

const double IsolatedPixelTrackCandidateProducer::vtxCutSeed_
private

Definition at line 94 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by produce().

◆ zEE_

double IsolatedPixelTrackCandidateProducer::zEE_
private

Definition at line 104 of file IsolatedPixelTrackCandidateProducer.cc.

Referenced by beginRun(), getDistInCM(), and GetEtaPhiAtEcal().