CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
IsolatedPixelTrackCandidateProducer Class Reference

#include <IsolatedPixelTrackCandidateProducer.h>

Inheritance diagram for IsolatedPixelTrackCandidateProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

virtual void beginRun (const edm::Run &, const edm::EventSetup &) override
 
double getDistInCM (double eta1, double phi1, double eta2, double phi2)
 
std::pair< double, double > GetEtaPhiAtEcal (const edm::EventSetup &iSetup, double etaIP, double phiIP, double pT, int charge, double vtxZ)
 
 IsolatedPixelTrackCandidateProducer (const edm::ParameterSet &ps)
 
virtual void produce (edm::Event &evt, const edm::EventSetup &es) override
 
 ~IsolatedPixelTrackCandidateProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

std::string bfield_
 
double ebEtaBoundary_
 
double maxPForIsolationValue_
 
double minPTrackValue_
 
edm::ParameterSet parameters
 
double pixelIsolationConeSizeAtEC_
 
std::vector< edm::InputTagpixelTracksSources_
 
double prelimCone_
 
double rEB_
 
double tauAssocCone_
 
double tauUnbiasCone_
 
edm::EDGetTokenT
< trigger::TriggerFilterObjectWithRefs
tok_hlt_
 
edm::EDGetTokenT
< l1extra::L1JetParticleCollection
tok_l1_
 
edm::EDGetTokenT
< reco::VertexCollection
tok_vert_
 
std::vector< edm::EDGetTokenT
< reco::TrackCollection > > 
toks_pix_
 
double vtxCutIsol_
 
double vtxCutSeed_
 
double zEE_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 29 of file IsolatedPixelTrackCandidateProducer.h.

Constructor & Destructor Documentation

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

Definition at line 52 of file IsolatedPixelTrackCandidateProducer.cc.

References bfield_, ebEtaBoundary_, edm::ParameterSet::getParameter(), i, maxPForIsolationValue_, minPTrackValue_, pixelIsolationConeSizeAtEC_, pixelTracksSources_, prelimCone_, rEB_, AlCaHLTBitMon_QueryRunRegistry::string, tauAssocCone_, tauUnbiasCone_, tok_hlt_, tok_l1_, tok_vert_, toks_pix_, vtxCutIsol_, vtxCutSeed_, and zEE_.

52  {
53 
54  tok_l1_ = consumes<l1extra::L1JetParticleCollection>(config.getParameter<edm::InputTag>("L1eTauJetsSource"));
55  tauAssocCone_ = config.getParameter<double>("tauAssociationCone");
56  tauUnbiasCone_ = config.getParameter<double>("tauUnbiasCone");
57  pixelTracksSources_ = config.getParameter<std::vector<edm::InputTag> >("PixelTracksSources");
58 
59  const unsigned nLabels = pixelTracksSources_.size();
60  for ( unsigned i=0; i != nLabels; i++ )
61  toks_pix_.push_back(consumes<reco::TrackCollection>(pixelTracksSources_[i]));
62 
63  prelimCone_ = config.getParameter<double>("ExtrapolationConeSize");
64  pixelIsolationConeSizeAtEC_ = config.getParameter<double>("PixelIsolationConeSizeAtEC");
65  tok_hlt_ = consumes<trigger::TriggerFilterObjectWithRefs>(config.getParameter<edm::InputTag>("L1GTSeedLabel"));
66  vtxCutSeed_ = config.getParameter<double>("MaxVtxDXYSeed");
67  vtxCutIsol_ = config.getParameter<double>("MaxVtxDXYIsol");
68  tok_vert_ = consumes<reco::VertexCollection>(config.getParameter<edm::InputTag>("VertexLabel"));
69  bfield_ = config.getParameter<std::string>("MagFieldRecordName");
70  minPTrackValue_ = config.getParameter<double>("minPTrack");
71  maxPForIsolationValue_ = config.getParameter<double>("maxPTrackForIsolation");
72  ebEtaBoundary_ = config.getParameter<double>("EBEtaBoundary");
73  rEB_ = zEE_ = -1;
74 
75  // Register the product
76  produces< reco::IsolatedPixelTrackCandidateCollection >();
77 }
int i
Definition: DBlmapReader.cc:9
std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
edm::EDGetTokenT< reco::VertexCollection > tok_vert_
edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_l1_
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tok_hlt_
IsolatedPixelTrackCandidateProducer::~IsolatedPixelTrackCandidateProducer ( )

Definition at line 79 of file IsolatedPixelTrackCandidateProducer.cc.

79  {
80 
81 }

Member Function Documentation

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

Reimplemented from edm::EDProducer.

Definition at line 83 of file IsolatedPixelTrackCandidateProducer.cc.

References DetId::Ecal, EcalBarrel, EcalEndcap, edm::EventSetup::get(), rEB_, and zEE_.

84 {
85 
87  theEventSetup.get<CaloGeometryRecord>().get(pG);
88 
89  const double rad (dynamic_cast<const EcalBarrelGeometry*>( pG->getSubdetectorGeometry(DetId::Ecal, EcalBarrel ))->avgRadiusXYFrontFaceCenter() ) ;
90 
91  const double zz (dynamic_cast<const EcalEndcapGeometry*>( pG->getSubdetectorGeometry(DetId::Ecal, EcalEndcap ))->avgAbsZFrontFaceCenter() ) ;
92 
93  rEB_=rad;
94  zEE_=zz;
95 
96 }
const T & get() const
Definition: EventSetup.h:55
double IsolatedPixelTrackCandidateProducer::getDistInCM ( double  eta1,
double  phi1,
double  eta2,
double  phi2 
)

Definition at line 266 of file IsolatedPixelTrackCandidateProducer.cc.

References angle(), funct::cos(), create_public_lumi_plots::exp, funct::sin(), and funct::tan().

Referenced by produce().

267 {
268  double Rec;
269  double theta1=2*atan(exp(-eta1));
270  double theta2=2*atan(exp(-eta2));
271  if (fabs(eta1)<1.479) Rec=129; //radius of ECAL barrel
272  else if (fabs(eta1)>1.479&&fabs(eta1)<7.0) Rec=tan(theta1)*317; //distance from IP to ECAL endcap
273  else return 1000;
274 
275  //|vect| times tg of acos(scalar product)
276  double angle=acos((sin(theta1)*sin(theta2)*(sin(phi1)*sin(phi2)+cos(phi1)*cos(phi2))+cos(theta1)*cos(theta2)));
277  if (angle<acos(-1)/2) return fabs((Rec/sin(theta1))*tan(angle));
278  else return 1000;
279 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
std::pair< double, double > IsolatedPixelTrackCandidateProducer::GetEtaPhiAtEcal ( const edm::EventSetup iSetup,
double  etaIP,
double  phiIP,
double  pT,
int  charge,
double  vtxZ 
)

Definition at line 283 of file IsolatedPixelTrackCandidateProducer.cc.

References alpha, SQLiteEnsembleGenerator_cfg::BField, DeDxDiscriminatorTools::charge(), SiPixelRawToDigiRegional_cfi::deltaPhi, ebEtaBoundary_, create_public_lumi_plots::exp, edm::EventSetup::get(), VolumeBasedMagneticField::inTesla(), fff_deleter::log, PV3DBase< T, PVType, FrameType >::mag(), rEB_, funct::sin(), funct::tan(), theta(), detailsBasic3DVector::z, and zEE_.

Referenced by produce().

284 {
286  iSetup.get<IdealMagneticFieldRecord>().get(vbfField);
287  const VolumeBasedMagneticField* vbfCPtr = dynamic_cast<const VolumeBasedMagneticField*>(&(*vbfField));
288  GlobalVector BField=vbfCPtr->inTesla(GlobalPoint(0,0,0));
289  //test
290  //int curvSgn=int(BField.z()/fabs(BField.z()));
291 
292  double bfVal=BField.mag();
293 
294  double deltaPhi=0;
295  double etaEC = 100;
296  double phiEC = 100;
297 
298  double Rcurv = 9999999;
299  if (bfVal!=0) Rcurv=pT*33.3*100/(bfVal*10); //r(m)=pT(GeV)*33.3/B(kG)
300 
301  double ecDist = zEE_; //distance to ECAL andcap from IP (cm), 317 - ecal (not preshower), preshower -300
302  double ecRad = rEB_; //radius of ECAL barrel (cm)
303  double theta=2*atan(exp(-etaIP));
304  double zNew=0;
305  if (theta>0.5*acos(-1)) theta=acos(-1)-theta;
306  if (fabs(etaIP)<ebEtaBoundary_)
307  {
308  if ((0.5*ecRad/Rcurv)>1)
309  {
310  etaEC=10000;
311  deltaPhi=0;
312  }
313  else
314  {
315  deltaPhi =-charge*asin(0.5*ecRad/Rcurv);
316  double alpha1 = 2*asin(0.5*ecRad/Rcurv);
317  double z = ecRad/tan(theta);
318  if (etaIP>0) zNew = z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
319  else zNew =-z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
320  double zAbs=fabs(zNew);
321  if (zAbs<ecDist)
322  {
323  etaEC = -log(tan(0.5*atan(ecRad/zAbs)));
324  deltaPhi = -charge*asin(0.5*ecRad/Rcurv);
325  }
326  if (zAbs>ecDist)
327  {
328  zAbs = (fabs(etaIP)/etaIP)*ecDist;
329  double Zflight = fabs(zAbs-vtxZ);
330  double alpha = (Zflight*ecRad)/(z*Rcurv);
331  double Rec = 2*Rcurv*sin(alpha/2);
332  deltaPhi =-charge*alpha/2;
333  etaEC =-log(tan(0.5*atan(Rec/ecDist)));
334  }
335  }
336  }
337  else
338  {
339  zNew = (fabs(etaIP)/etaIP)*ecDist;
340  double Zflight = fabs(zNew-vtxZ);
341  double Rvirt = fabs(Zflight*tan(theta));
342  double Rec = 2*Rcurv*sin(Rvirt/(2*Rcurv));
343  deltaPhi =-(charge)*(Rvirt/(2*Rcurv));
344  etaEC =-log(tan(0.5*atan(Rec/ecDist)));
345  }
346 
347  if (zNew<0) etaEC=-etaEC;
348  phiEC = phiIP+deltaPhi;
349 
350  if (phiEC<-acos(-1)) phiEC = 2*acos(-1)+phiEC;
351  if (phiEC>acos(-1)) phiEC =-2*acos(-1)+phiEC;
352 
353  std::pair<double,double> retVal(etaEC,phiEC);
354  return retVal;
355 }
float alpha
Definition: AMPTWrapper.h:95
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
double charge(const std::vector< uint8_t > &Ampls)
float float float z
T mag() const
Definition: PV3DBase.h:67
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
GlobalVector inTesla(const GlobalPoint &g) const
Field value ad specified global point, in Tesla.
const T & get() const
Definition: EventSetup.h:55
void IsolatedPixelTrackCandidateProducer::produce ( edm::Event evt,
const edm::EventSetup es 
)
overridevirtual

Implements edm::EDProducer.

Definition at line 98 of file IsolatedPixelTrackCandidateProducer.cc.

References DeDxDiscriminatorTools::charge(), deltaR(), eta(), edm::Event::getByToken(), getDistInCM(), GetEtaPhiAtEcal(), maxPForIsolationValue_, minPTrackValue_, AlCaHLTBitMon_ParallelJobs::p, phi, pixelIsolationConeSizeAtEC_, prelimCone_, RecoTauCleanerPlugins::pt, edm::Event::put(), dttmaxenums::R, tauAssocCone_, tauUnbiasCone_, tok_hlt_, tok_l1_, tok_vert_, toks_pix_, trigger::TriggerL1CenJet, trigger::TriggerL1ForJet, trigger::TriggerL1TauJet, vtxCutIsol_, and vtxCutSeed_.

Referenced by JSONExport.JsonExport::export(), HTMLExport.HTMLExport::export(), and HTMLExport.HTMLExportStatic::export().

98  {
99 
101 
102  //create vector of refs from input collections
103  std::vector<reco::TrackRef> pixelTrackRefs;
104 
105  for (unsigned int iPix=0; iPix<toks_pix_.size(); iPix++)
106  {
108  theEvent.getByToken(toks_pix_[iPix],iPixCol);
109  for (reco::TrackCollection::const_iterator pit=iPixCol->begin(); pit!=iPixCol->end(); pit++)
110  {
111  pixelTrackRefs.push_back(reco::TrackRef(iPixCol,pit-iPixCol->begin()));
112  }
113  }
114 
116  theEvent.getByToken(tok_l1_,l1eTauJets);
117 
119  theEvent.getByToken(tok_vert_,pVert);
120 
121  double ptTriggered = -10;
122  double etaTriggered = -100;
123  double phiTriggered = -100;
124 
126  theEvent.getByToken(tok_hlt_, l1trigobj);
127 
128  std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1tauobjref;
129  std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1jetobjref;
130  std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1forjetobjref;
131 
132  l1trigobj->getObjects(trigger::TriggerL1TauJet, l1tauobjref);
133  l1trigobj->getObjects(trigger::TriggerL1CenJet, l1jetobjref);
134  l1trigobj->getObjects(trigger::TriggerL1ForJet, l1forjetobjref);
135 
136  for (unsigned int p=0; p<l1tauobjref.size(); p++)
137  {
138  if (l1tauobjref[p]->pt()>ptTriggered)
139  {
140  ptTriggered = l1tauobjref[p]->pt();
141  phiTriggered = l1tauobjref[p]->phi();
142  etaTriggered = l1tauobjref[p]->eta();
143  }
144  }
145  for (unsigned int p=0; p<l1jetobjref.size(); p++)
146  {
147  if (l1jetobjref[p]->pt()>ptTriggered)
148  {
149  ptTriggered = l1jetobjref[p]->pt();
150  phiTriggered = l1jetobjref[p]->phi();
151  etaTriggered = l1jetobjref[p]->eta();
152  }
153  }
154  for (unsigned int p=0; p<l1forjetobjref.size(); p++)
155  {
156  if (l1forjetobjref[p]->pt()>ptTriggered)
157  {
158  ptTriggered=l1forjetobjref[p]->pt();
159  phiTriggered=l1forjetobjref[p]->phi();
160  etaTriggered=l1forjetobjref[p]->eta();
161  }
162  }
163 
164  double minPTrack_ = minPTrackValue_;
165  double drMaxL1Track_ = tauAssocCone_;
166 
167  int ntr = 0;
168 
169  //loop to select isolated tracks
170  for (unsigned iSeed=0; iSeed<pixelTrackRefs.size(); iSeed++)
171  {
172  if(pixelTrackRefs[iSeed]->p()<minPTrack_) continue;
173 
174  bool good = false;
175  bool vtxMatch = false;
176 
177  //associate to vertex (in Z)
178  reco::VertexCollection::const_iterator vitSel;
179  double minDZ = 100;
180  for (reco::VertexCollection::const_iterator vit=pVert->begin(); vit!=pVert->end(); vit++)
181  {
182  if (fabs(pixelTrackRefs[iSeed]->dz(vit->position()))<minDZ)
183  {
184  minDZ = fabs(pixelTrackRefs[iSeed]->dz(vit->position()));
185  vitSel = vit;
186  }
187  }
188  //cut on dYX:
189  if (minDZ!=100&&fabs(pixelTrackRefs[iSeed]->dxy(vitSel->position()))<vtxCutSeed_) vtxMatch=true;
190  if (minDZ==100) vtxMatch=true;
191 
192  //select tracks not matched to triggered L1 jet
193  double R=deltaR(etaTriggered, phiTriggered, pixelTrackRefs[iSeed]->eta(), pixelTrackRefs[iSeed]->phi());
194  if (R<tauUnbiasCone_) continue;
195 
196  //check taujet matching
197  bool tmatch=false;
198  l1extra::L1JetParticleCollection::const_iterator selj;
199  for (l1extra::L1JetParticleCollection::const_iterator tj=l1eTauJets->begin(); tj!=l1eTauJets->end(); tj++)
200  {
201  if(ROOT::Math::VectorUtil::DeltaR(pixelTrackRefs[iSeed]->momentum(),tj->momentum()) > drMaxL1Track_) continue;
202  selj = tj;
203  tmatch = true;
204  } //loop over L1 tau
205 
206  //propagate seed track to ECAL surface:
207  std::pair<double,double> seedCooAtEC;
208  // in case vertex is found:
209  if (minDZ!=100) seedCooAtEC=GetEtaPhiAtEcal(theEventSetup, pixelTrackRefs[iSeed]->eta(), pixelTrackRefs[iSeed]->phi(), pixelTrackRefs[iSeed]->pt(), pixelTrackRefs[iSeed]->charge(), vitSel->z());
210  //in case vertex is not found:
211  else seedCooAtEC=GetEtaPhiAtEcal(theEventSetup, pixelTrackRefs[iSeed]->eta(), pixelTrackRefs[iSeed]->phi(), pixelTrackRefs[iSeed]->pt(), pixelTrackRefs[iSeed]->charge(), 0);
212 
213  //calculate isolation
214  double maxP = 0;
215  double sumP = 0;
216  for (unsigned iSurr=0; iSurr<pixelTrackRefs.size(); iSurr++)
217  {
218  if(iSeed==iSurr) continue;
219  //define preliminary cone around seed track impact point from which tracks will be extrapolated:
220  if (deltaR(seedCooAtEC.first, seedCooAtEC.second, pixelTrackRefs[iSurr]->eta(), pixelTrackRefs[iSurr]->phi())>prelimCone_) continue;
221  //associate to vertex (in Z):
222  double minDZ2=100;
223  reco::VertexCollection::const_iterator vitSel2;
224  for (reco::VertexCollection::const_iterator vit=pVert->begin(); vit!=pVert->end(); vit++)
225  {
226  if (fabs(pixelTrackRefs[iSurr]->dz(vit->position()))<minDZ2)
227  {
228  minDZ2 = fabs(pixelTrackRefs[iSurr]->dz(vit->position()));
229  vitSel2 = vit;
230  }
231  }
232  //cut ot dXY:
233  if (minDZ2!=100&&fabs(pixelTrackRefs[iSurr]->dxy(vitSel2->position()))>vtxCutIsol_) continue;
234  //propagate to ECAL surface:
235  std::pair<double,double> cooAtEC;
236  // in case vertex is found:
237  if (minDZ2!=100) cooAtEC=GetEtaPhiAtEcal(theEventSetup, pixelTrackRefs[iSurr]->eta(), pixelTrackRefs[iSurr]->phi(), pixelTrackRefs[iSurr]->pt(), pixelTrackRefs[iSurr]->charge(), vitSel2->z());
238  // in case vertex is not found:
239  else cooAtEC=GetEtaPhiAtEcal(theEventSetup, pixelTrackRefs[iSurr]->eta(), pixelTrackRefs[iSurr]->phi(), pixelTrackRefs[iSurr]->pt(), pixelTrackRefs[iSurr]->charge(), 0);
240 
241  //calculate distance at ECAL surface and update isolation:
242  if (getDistInCM(seedCooAtEC.first, seedCooAtEC.second, cooAtEC.first, cooAtEC.second)<pixelIsolationConeSizeAtEC_)
243  {
244  sumP+=pixelTrackRefs[iSurr]->p();
245  if(pixelTrackRefs[iSurr]->p()>maxP) maxP=pixelTrackRefs[iSurr]->p();
246  }
247  }
248 
249  if (tmatch||vtxMatch) good=true;
250 
251  if (good&&maxP<maxPForIsolationValue_)
252  {
253  reco::IsolatedPixelTrackCandidate newCandidate(pixelTrackRefs[iSeed], l1extra::L1JetParticleRef(l1eTauJets,selj-l1eTauJets->begin()), maxP, sumP);
254  trackCollection->push_back(newCandidate);
255  ntr++;
256  }
257  }//loop over pixel tracks
258 
259  // put the product in the event
260  std::auto_ptr< reco::IsolatedPixelTrackCandidateCollection > outCollection(trackCollection);
261  theEvent.put(outCollection);
262 
263 }
std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
edm::EDGetTokenT< reco::VertexCollection > tok_vert_
T eta() const
double charge(const std::vector< uint8_t > &Ampls)
std::pair< double, double > GetEtaPhiAtEcal(const edm::EventSetup &iSetup, double etaIP, double phiIP, double pT, int charge, double vtxZ)
edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_l1_
std::vector< IsolatedPixelTrackCandidate > IsolatedPixelTrackCandidateCollection
collectin of IsolatedPixelTrackCandidate objects
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double getDistInCM(double eta1, double phi1, double eta2, double phi2)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tok_hlt_
Definition: DDAxes.h:10

Member Data Documentation

std::string IsolatedPixelTrackCandidateProducer::bfield_
private
double IsolatedPixelTrackCandidateProducer::ebEtaBoundary_
private
double IsolatedPixelTrackCandidateProducer::maxPForIsolationValue_
private
double IsolatedPixelTrackCandidateProducer::minPTrackValue_
private
edm::ParameterSet IsolatedPixelTrackCandidateProducer::parameters
private
double IsolatedPixelTrackCandidateProducer::pixelIsolationConeSizeAtEC_
private
std::vector<edm::InputTag> IsolatedPixelTrackCandidateProducer::pixelTracksSources_
private
double IsolatedPixelTrackCandidateProducer::prelimCone_
private
double IsolatedPixelTrackCandidateProducer::rEB_
private
double IsolatedPixelTrackCandidateProducer::tauAssocCone_
private
double IsolatedPixelTrackCandidateProducer::tauUnbiasCone_
private
edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> IsolatedPixelTrackCandidateProducer::tok_hlt_
private
edm::EDGetTokenT<l1extra::L1JetParticleCollection> IsolatedPixelTrackCandidateProducer::tok_l1_
private
edm::EDGetTokenT<reco::VertexCollection> IsolatedPixelTrackCandidateProducer::tok_vert_
private
std::vector<edm::EDGetTokenT<reco::TrackCollection> > IsolatedPixelTrackCandidateProducer::toks_pix_
private
double IsolatedPixelTrackCandidateProducer::vtxCutIsol_
private
double IsolatedPixelTrackCandidateProducer::vtxCutSeed_
private
double IsolatedPixelTrackCandidateProducer::zEE_
private