CMS 3D CMS Logo

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

#include <IsolatedPixelTrackCandidateProducer.h>

Inheritance diagram for IsolatedPixelTrackCandidateProducer:
edm::stream::EDProducer<> edm::stream::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Classes

struct  seedAtEC
 

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 (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::stream::EDProducer<>
 EDProducer ()=default
 
- Public Member Functions inherited from edm::stream::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducerBase ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 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
std::vector< ConsumesInfoconsumesInfo () const
 
 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
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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::EDGetTokenT
< trigger::TriggerFilterObjectWithRefs
tok_hlt_
 
const edm::EDGetTokenT
< l1extra::L1JetParticleCollection
tok_l1_
 
const edm::EDGetTokenT
< reco::VertexCollection
tok_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<>
typedef CacheContexts< T...> CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T...> HasAbility
 
typedef
CacheTypes::LuminosityBlockCache 
LuminosityBlockCache
 
typedef
LuminosityBlockContextT
< LuminosityBlockCache,
RunCache, GlobalCache
LuminosityBlockContext
 
typedef
CacheTypes::LuminosityBlockSummaryCache 
LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache,
GlobalCache
RunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDProducerBase
typedef EDProducerAdaptorBase ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::stream::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- 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 41 of file IsolatedPixelTrackCandidateProducer.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file IsolatedPixelTrackCandidateProducer.cc.

References GlobalPosition_Frontier_DevDB_cff::tag.

41  :
42  tok_hlt_( consumes<trigger::TriggerFilterObjectWithRefs>(config.getParameter<edm::InputTag>("L1GTSeedLabel")) ),
43  tok_l1_( consumes<l1extra::L1JetParticleCollection>(config.getParameter<edm::InputTag>("L1eTauJetsSource")) ),
44  tok_vert_( consumes<reco::VertexCollection>(config.getParameter<edm::InputTag>("VertexLabel")) ),
46  config.getParameter<std::vector<edm::InputTag> >("PixelTracksSources"),
47  [this](edm::InputTag const & tag){return consumes<reco::TrackCollection>(tag);}) ),
const std::vector< edm::EDGetTokenT< reco::TrackCollection > > toks_pix_
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tok_hlt_
const edm::EDGetTokenT< l1extra::L1JetParticleCollection > tok_l1_
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::EDGetTokenT< reco::VertexCollection > tok_vert_
IsolatedPixelTrackCandidateProducer::~IsolatedPixelTrackCandidateProducer ( )

Definition at line 66 of file IsolatedPixelTrackCandidateProducer.cc.

66  {
67 
68 }

Member Function Documentation

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

Reimplemented from edm::stream::EDProducerBase.

Definition at line 70 of file IsolatedPixelTrackCandidateProducer.cc.

References SQLiteEnsembleGenerator_cfg::BField, bfVal_, DetId::Ecal, EcalBarrel, EcalEndcap, edm::EventSetup::get(), VolumeBasedMagneticField::inTesla(), PV3DBase< T, PVType, FrameType >::mag(), rEB_, and zEE_.

70  {
71 
73  theEventSetup.get<CaloGeometryRecord>().get(pG);
74 
75  const double rad (dynamic_cast<const EcalBarrelGeometry*>( pG->getSubdetectorGeometry(DetId::Ecal, EcalBarrel ))->avgRadiusXYFrontFaceCenter() ) ;
76  const double zz (dynamic_cast<const EcalEndcapGeometry*>( pG->getSubdetectorGeometry(DetId::Ecal, EcalEndcap ))->avgAbsZFrontFaceCenter() ) ;
77 
78  rEB_ = rad;
79  zEE_ = zz;
80 
82  theEventSetup.get<IdealMagneticFieldRecord>().get(vbfField);
83  const VolumeBasedMagneticField* vbfCPtr = dynamic_cast<const VolumeBasedMagneticField*>(&(*vbfField));
84  GlobalVector BField = vbfCPtr->inTesla(GlobalPoint(0,0,0));
85  bfVal_=BField.mag();
86 }
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T mag() const
Definition: PV3DBase.h:67
GlobalVector inTesla(const GlobalPoint &g) const
Field value ad specified global point, in Tesla.
const T & get() const
Definition: EventSetup.h:56
double IsolatedPixelTrackCandidateProducer::getDistInCM ( double  eta1,
double  phi1,
double  eta2,
double  phi2 
)

Definition at line 243 of file IsolatedPixelTrackCandidateProducer.cc.

References funct::abs(), angle(), funct::cos(), create_public_lumi_plots::exp, M_PI_2, rEB_, funct::sin(), funct::tan(), and zEE_.

Referenced by produce().

243  {
244  double Rec;
245  double theta1=2*atan(exp(-eta1));
246  double theta2=2*atan(exp(-eta2));
247  if (std::abs(eta1)<1.479) Rec=rEB_; //radius of ECAL barrel
248  else if (std::abs(eta1)>1.479&&std::abs(eta1)<7.0) Rec=tan(theta1)*zEE_; //distance from IP to ECAL endcap
249  else return 1000;
250 
251  //|vect| times tg of acos(scalar product)
252  double angle=acos((sin(theta1)*sin(theta2)*(sin(phi1)*sin(phi2)+cos(phi1)*cos(phi2))+cos(theta1)*cos(theta2)));
253  if (angle<M_PI_2) return std::abs((Rec/sin(theta1))*tan(angle));
254  else return 1000;
255 }
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
std::pair< double, double > IsolatedPixelTrackCandidateProducer::GetEtaPhiAtEcal ( double  etaIP,
double  phiIP,
double  pT,
int  charge,
double  vtxZ 
)

Definition at line 258 of file IsolatedPixelTrackCandidateProducer.cc.

References funct::abs(), alpha, bfVal_, RecoTauCleanerPlugins::charge, SiPixelRawToDigiRegional_cfi::deltaPhi, ebEtaBoundary_, create_public_lumi_plots::exp, cmsBatch::log, M_PI, M_PI_2, rEB_, funct::sin(), funct::tan(), theta(), z, and zEE_.

Referenced by produce().

258  {
259 
260  double deltaPhi=0;
261  double etaEC = 100;
262  double phiEC = 100;
263 
264  double Rcurv = 9999999;
265  if (bfVal_!=0) Rcurv=pT*33.3*100/(bfVal_*10); //r(m)=pT(GeV)*33.3/B(kG)
266 
267  double ecDist = zEE_; //distance to ECAL andcap from IP (cm), 317 - ecal (not preshower), preshower -300
268  double ecRad = rEB_; //radius of ECAL barrel (cm)
269  double theta = 2*atan(exp(-etaIP));
270  double zNew = 0;
271  if (theta>M_PI_2) theta=M_PI-theta;
272  if (std::abs(etaIP)<ebEtaBoundary_) {
273  if ((0.5*ecRad/Rcurv)>1) {
274  etaEC = 10000;
275  deltaPhi = 0;
276  } else {
277  deltaPhi =-charge*asin(0.5*ecRad/Rcurv);
278  double alpha1 = 2*asin(0.5*ecRad/Rcurv);
279  double z = ecRad/tan(theta);
280  if (etaIP>0) zNew = z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
281  else zNew =-z*(Rcurv*alpha1)/ecRad+vtxZ; //new z-coordinate of track
282  double zAbs=std::abs(zNew);
283  if (zAbs<ecDist) {
284  etaEC = -log(tan(0.5*atan(ecRad/zAbs)));
285  deltaPhi = -charge*asin(0.5*ecRad/Rcurv);
286  }
287  if (zAbs>ecDist) {
288  zAbs = (std::abs(etaIP)/etaIP)*ecDist;
289  double Zflight = std::abs(zAbs-vtxZ);
290  double alpha = (Zflight*ecRad)/(z*Rcurv);
291  double Rec = 2*Rcurv*sin(alpha/2);
292  deltaPhi =-charge*alpha/2;
293  etaEC =-log(tan(0.5*atan(Rec/ecDist)));
294  }
295  }
296  } else {
297  zNew = (std::abs(etaIP)/etaIP)*ecDist;
298  double Zflight = std::abs(zNew-vtxZ);
299  double Rvirt = std::abs(Zflight*tan(theta));
300  double Rec = 2*Rcurv*sin(Rvirt/(2*Rcurv));
301  deltaPhi =-(charge)*(Rvirt/(2*Rcurv));
302  etaEC =-log(tan(0.5*atan(Rec/ecDist)));
303  }
304 
305  if (zNew<0) etaEC=-etaEC;
306  phiEC = phiIP+deltaPhi;
307 
308  if (phiEC<-M_PI) phiEC += M_2_PI;
309  if (phiEC>M_PI) phiEC -= M_2_PI;
310 
311  std::pair<double,double> retVal(etaEC,phiEC);
312  return retVal;
313 }
float alpha
Definition: AMPTWrapper.h:95
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
#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
tuple log
Definition: cmsBatch.py:341
void IsolatedPixelTrackCandidateProducer::produce ( edm::Event evt,
const edm::EventSetup es 
)
overridevirtual

Implements edm::stream::EDProducerBase.

Definition at line 88 of file IsolatedPixelTrackCandidateProducer.cc.

References funct::abs(), RecoTauCleanerPlugins::charge, reco::deltaR(), eta, newFWLiteAna::found, edm::Event::getByToken(), getDistInCM(), GetEtaPhiAtEcal(), i, j, maxPForIsolationValue_, minPTrackValue_, convertSQLiteXML::ok, AlCaHLTBitMon_ParallelJobs::p, phi, pixelIsolationConeSizeAtEC_, prelimCone_, EnergyCorrector::pt, edm::Event::put(), dttmaxenums::R, fileCollector::seed, reco::IsolatedPixelTrackCandidate::setEtaPhiEcal(), tauAssocCone_, tauUnbiasCone_, tok_hlt_, tok_l1_, tok_vert_, toks_pix_, HLT_25ns14e33_v3_cff::trackCollection, trigger::TriggerL1CenJet, trigger::TriggerL1ForJet, trigger::TriggerL1TauJet, vtxCutIsol_, and vtxCutSeed_.

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

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

Member Data Documentation

const std::string IsolatedPixelTrackCandidateProducer::bfield_
private

Definition at line 68 of file IsolatedPixelTrackCandidateProducer.h.

double IsolatedPixelTrackCandidateProducer::bfVal_
private

Definition at line 82 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by beginRun(), and GetEtaPhiAtEcal().

const double IsolatedPixelTrackCandidateProducer::ebEtaBoundary_
private

Definition at line 77 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by GetEtaPhiAtEcal().

const double IsolatedPixelTrackCandidateProducer::maxPForIsolationValue_
private

Definition at line 76 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

const double IsolatedPixelTrackCandidateProducer::minPTrackValue_
private

Definition at line 75 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

const double IsolatedPixelTrackCandidateProducer::pixelIsolationConeSizeAtEC_
private

Definition at line 70 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

const double IsolatedPixelTrackCandidateProducer::prelimCone_
private

Definition at line 69 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

double IsolatedPixelTrackCandidateProducer::rEB_
private

Definition at line 80 of file IsolatedPixelTrackCandidateProducer.h.

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

const double IsolatedPixelTrackCandidateProducer::tauAssocCone_
private

Definition at line 73 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

const double IsolatedPixelTrackCandidateProducer::tauUnbiasCone_
private

Definition at line 74 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

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

Definition at line 63 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

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

Definition at line 64 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

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

Definition at line 65 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

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

Definition at line 66 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

const double IsolatedPixelTrackCandidateProducer::vtxCutIsol_
private

Definition at line 72 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

const double IsolatedPixelTrackCandidateProducer::vtxCutSeed_
private

Definition at line 71 of file IsolatedPixelTrackCandidateProducer.h.

Referenced by produce().

double IsolatedPixelTrackCandidateProducer::zEE_
private

Definition at line 81 of file IsolatedPixelTrackCandidateProducer.h.

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