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
TrajectorySeedProducer Class Reference

#include <TrajectorySeedProducer.h>

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

Public Member Functions

virtual void beginRun (edm::Run const &run, const edm::EventSetup &es)
 
bool compatibleWithBeamAxis (const GlobalPoint &gpos1, const GlobalPoint &gpos2, double error, bool forward) const
 
const SeedingNode
< TrackingLayer > * 
insertHit (const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
 method inserts hit into the tree structure at an empty position. More...
 
bool isHitOnLayer (const TrajectorySeedHitCandidate &trackerRecHit, const TrackingLayer &layer) const
 
virtual std::vector< unsigned int > iterateHits (unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
 method tries to insert all hits into the tree structure. More...
 
bool pass2HitsCuts (const TrajectorySeedHitCandidate &hit1, const TrajectorySeedHitCandidate &hit2) const
 
bool passHitTuplesCuts (const SeedingNode< TrackingLayer > &seedingNode, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, const std::vector< int > &hitIndicesInTree, const TrajectorySeedHitCandidate &currentTrackerHit) const
 method checks if a TrajectorySeedHitCandidate fulfills the quality requirements. More...
 
virtual bool passSimTrackQualityCuts (const SimTrack &theSimTrack, const SimVertex &theSimVertex) const
 method checks if a SimTrack fulfills the quality requirements. More...
 
virtual void produce (edm::Event &e, const edm::EventSetup &es)
 
 TrajectorySeedProducer (const edm::ParameterSet &conf)
 
virtual ~TrajectorySeedProducer ()
 
- 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
 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

SeedingTree< TrackingLayer_seedingTree
 
unsigned int absMinRecHits
 
math::XYZPoint beamspotPosition
 
edm::EDGetTokenT< reco::BeamSpotbeamSpotToken
 
edm::InputTag hitProducer
 
const MagneticFieldmagneticField
 
const MagneticFieldMapmagneticFieldMap
 
double maxD0
 
double maxZ0
 
unsigned int minRecHits
 
unsigned int numberOfHits
 
double originHalfLength
 
double originpTMin
 
double originRadius
 
std::string outputSeedCollectionName
 
double pTMin
 
edm::EDGetTokenT
< SiTrackerGSMatchedRecHit2DCollection
recHitToken
 
edm::EDGetTokenT
< reco::VertexCollection
recoVertexToken
 
bool seedCleaning
 
std::vector< std::vector
< TrackingLayer > > 
seedingLayers
 
edm::EDGetTokenT
< edm::SimTrackContainer
simTrackToken
 
edm::EDGetTokenT
< edm::SimVertexContainer
simVertexToken
 
bool skipPVCompatibility
 
std::vector< edm::EDGetTokenT
< std::vector< int > > > 
skipSimTrackIdTokens
 
edm::InputTag theBeamSpot
 
PropagatorWithMaterialthePropagator
 
const TrackerGeometrytrackerGeometry
 
const TrackerTopologytrackerTopology
 
const reco::VertexCollectionvertices
 
double zVertexConstraint
 

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
 
- 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::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 35 of file TrajectorySeedProducer.h.

Constructor & Destructor Documentation

TrajectorySeedProducer::TrajectorySeedProducer ( const edm::ParameterSet conf)

Definition at line 44 of file TrajectorySeedProducer.cc.

References _seedingTree, absMinRecHits, beamSpotToken, edm::EDConsumerBase::consumes(), TrackingLayer::createFromString(), edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hitProducer, HLT_25ns14e33_v1_cff::InputTag, SeedingTree< DATA >::insert(), relval_steps::k, geometryCSVtoXML::line, maxD0, maxZ0, minRecHits, eostools::move(), numberOfHits, originHalfLength, originpTMin, originRadius, outputSeedCollectionName, HLT_25ns14e33_v1_cff::primaryVertex, pTMin, recHitToken, recoVertexToken, seedCleaning, seedingLayers, simTrackToken, simVertexToken, skipPVCompatibility, skipSimTrackIdTokens, AlCaHLTBitMon_QueryRunRegistry::string, theBeamSpot, and zVertexConstraint.

44  :
45  thePropagator(nullptr),
46  vertices(nullptr) //TODO:: what else should be initialized properly?
47 {
49  if (conf.exists("outputSeedCollectionName"))
50  {
51  outputSeedCollectionName=conf.getParameter<std::string>("outputSeedCollectionName");
52  }
53  // The input tag for the beam spot
54  theBeamSpot = conf.getParameter<edm::InputTag>("beamSpot");
55 
56  // The name of the TrajectorySeed Collections
57  produces<TrajectorySeedCollection>(outputSeedCollectionName);
58 
59  // The smallest true pT for a track candidate
60  pTMin = conf.getParameter<double>("pTMin");
61 
62 
63  // The smallest number of Rec Hits for a track candidate
64  minRecHits = conf.getParameter<unsigned int>("minRecHits");
65 
66  //TODO: REMOVE
67  // Set the overall number hits to be checked
69 
70  // The smallest true impact parameters (d0 and z0) for a track candidate
71  maxD0 = conf.getParameter<double>("maxD0");
72 
73  maxZ0 = conf.getParameter<double>("maxZ0");
74 
75  // The name of the hit producer
76  hitProducer = conf.getParameter<edm::InputTag>("HitProducer");
77 
78  // The cuts for seed cleaning
79  seedCleaning = conf.getParameter<bool>("seedCleaning");
80 
81  // Number of hits needed for a seed
82  numberOfHits = conf.getParameter<unsigned int>("numberOfHits");
83 
84  // read Layers
85  std::vector<std::string> layerStringList = conf.getParameter<std::vector<std::string>>("layerList");
86 
87  for(auto it=layerStringList.cbegin(); it < layerStringList.cend(); ++it)
88  {
89  std::vector<TrackingLayer> trackingLayerList;
90  std::string line = *it;
92  while (pos != std::string::npos)
93  {
94  pos=line.find("+");
95  std::string layer = line.substr(0, pos);
97 
98  trackingLayerList.push_back(layerSpec);
99  line=line.substr(pos+1,std::string::npos);
100  }
101  _seedingTree.insert(trackingLayerList);
102  seedingLayers.push_back(std::move(trackingLayerList));
103  }
104 
105  originRadius = conf.getParameter<double>("originRadius");
106 
107  //collections to read in
108  std::vector<edm::InputTag> defaultSkip;
109  std::vector<edm::InputTag> skipSimTrackIdTags = conf.getUntrackedParameter<std::vector<edm::InputTag> >("skipSimTrackIdTags",defaultSkip);
110  for ( unsigned int k=0; k<skipSimTrackIdTags.size(); ++k ) {
111  skipSimTrackIdTokens.push_back(consumes<std::vector<int> >(skipSimTrackIdTags[k]));
112  }
113 
114  originHalfLength = conf.getParameter<double>("originHalfLength");
115 
116  originpTMin = conf.getParameter<double>("originpTMin");
117 
118  edm::InputTag primaryVertex = conf.getParameter<edm::InputTag>("primaryVertex");
119 
120  zVertexConstraint = conf.getParameter<double>("zVertexConstraint");
121 
122 
123 
124 
125  skipPVCompatibility=false;
126  if (conf.exists("skipPVCompatibility"))
127  {
128  skipPVCompatibility = conf.getParameter<bool>("skipPVCompatibility");
129  }
130 
131 
132  // consumes
133  beamSpotToken = consumes<reco::BeamSpot>(theBeamSpot);
134  edm::InputTag("famosSimHits");
135  simTrackToken = consumes<edm::SimTrackContainer>(edm::InputTag("famosSimHits"));
136  simVertexToken = consumes<edm::SimVertexContainer>(edm::InputTag("famosSimHits"));
137  recHitToken = consumes<SiTrackerGSMatchedRecHit2DCollection>(hitProducer);
138 
139  recoVertexToken=consumes<reco::VertexCollection>(primaryVertex);
140 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const reco::VertexCollection * vertices
bool exists(std::string const &parameterName) const
checks if a parameter exists
uint16_t size_type
edm::EDGetTokenT< SiTrackerGSMatchedRecHit2DCollection > recHitToken
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
PropagatorWithMaterial * thePropagator
bool insert(const std::vector< DATA > &dataList)
Definition: SeedingTree.h:179
def move
Definition: eostools.py:508
static TrackingLayer createFromString(std::string layerSpecification)
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken
std::vector< edm::EDGetTokenT< std::vector< int > > > skipSimTrackIdTokens
edm::EDGetTokenT< edm::SimTrackContainer > simTrackToken
std::vector< std::vector< TrackingLayer > > seedingLayers
edm::EDGetTokenT< reco::VertexCollection > recoVertexToken
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken
SeedingTree< TrackingLayer > _seedingTree
TrajectorySeedProducer::~TrajectorySeedProducer ( )
virtual

Definition at line 143 of file TrajectorySeedProducer.cc.

References thePropagator.

143  {
144 
145  if(thePropagator) delete thePropagator;
146 }
PropagatorWithMaterial * thePropagator

Member Function Documentation

void TrajectorySeedProducer::beginRun ( edm::Run const &  run,
const edm::EventSetup es 
)
virtual

Reimplemented from edm::stream::EDProducerBase.

Definition at line 150 of file TrajectorySeedProducer.cc.

References alongMomentum, edm::EventSetup::get(), magneticField, magneticFieldMap, thePropagator, trackerGeometry, and trackerTopology.

151 {
152  edm::ESHandle<MagneticField> magneticFieldHandle;
153  edm::ESHandle<TrackerGeometry> trackerGeometryHandle;
154  edm::ESHandle<MagneticFieldMap> magneticFieldMapHandle;
155  edm::ESHandle<TrackerTopology> trackerTopologyHandle;
156 
157  es.get<IdealMagneticFieldRecord>().get(magneticFieldHandle);
158  es.get<TrackerDigiGeometryRecord>().get(trackerGeometryHandle);
159  es.get<MagneticFieldMapRecord>().get(magneticFieldMapHandle);
160  es.get<IdealGeometryRecord>().get(trackerTopologyHandle);
161 
162  magneticField = &(*magneticFieldHandle);
163  trackerGeometry = &(*trackerGeometryHandle);
164  magneticFieldMap = &(*magneticFieldMapHandle);
165  trackerTopology = &(*trackerTopologyHandle);
166 
168 }
const MagneticField * magneticField
PropagatorWithMaterial * thePropagator
const TrackerTopology * trackerTopology
const T & get() const
Definition: EventSetup.h:55
const TrackerGeometry * trackerGeometry
const MagneticFieldMap * magneticFieldMap
bool TrajectorySeedProducer::compatibleWithBeamAxis ( const GlobalPoint gpos1,
const GlobalPoint gpos2,
double  error,
bool  forward 
) const

Check that the seed is compatible with a track coming from within a cylinder of radius originRadius, with a decent pT.

Check that the seed is compatible with a track coming from within a cylinder of radius originRadius, with a decent pT, and propagate to the distance of closest approach, for the appropriate charge

Definition at line 524 of file TrajectorySeedProducer.cc.

References beamspotPosition, relativeConstraints::error, magneticFieldMap, bookConverter::max, min(), originHalfLength, originpTMin, originRadius, BaseParticlePropagator::propagateToBeamCylinder(), seedCleaning, mathSSE::sqrt(), vertices, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), RawParticle::Z(), and zVertexConstraint.

Referenced by pass2HitsCuts().

530 {
531 
532  const double x0 = beamspotPosition.X();
533  const double y0 = beamspotPosition.Y();
534  const double z0 = beamspotPosition.Z();
535 
536  if ( !seedCleaning )
537  {
538  return true;
539  }
540 
541  // The hits 1 and 2 positions, in HepLorentzVector's
542  XYZTLorentzVector thePos1(gpos1.x(),gpos1.y(),gpos1.z(),0.);
543  XYZTLorentzVector thePos2(gpos2.x(),gpos2.y(),gpos2.z(),0.);
544 
545  // Create new particles that pass through the second hit with pT = ptMin
546  // and charge = +/-1
547 
548  // The momentum direction is by default joining the two hits
549  XYZTLorentzVector theMom2 = (thePos2-thePos1);
550 
551  // The corresponding RawParticle, with an (irrelevant) electric charge
552  // (The charge is determined in the next step)
553  ParticlePropagator myPart(theMom2,thePos2,1.,magneticFieldMap);
554 
558  bool intersect = myPart.propagateToBeamCylinder(thePos1,originRadius*1.);
559  if ( !intersect )
560  {
561  return false;
562  }
563 
564  // Check if the constraints are satisfied
565  // 1. pT at cylinder with radius originRadius
566  if ( myPart.Pt() < originpTMin )
567  {
568  return false;
569  }
570 
571  // 2. Z compatible with beam spot size
572  if ( fabs(myPart.Z()-z0) > originHalfLength )
573  {
574  return false;
575  }
576 
577 
578  // 3. Z compatible with one of the primary vertices (always the case if no primary vertex)
579  if (!vertices)
580  {
581  return true;
582  }
583  unsigned int nVertices = vertices->size();
584  if ( !nVertices || zVertexConstraint < 0. )
585  {
586  return true;
587  }
588  // Radii of the two hits with respect to the beam spot position
589  double R1 = std::sqrt ( (thePos1.X()-x0)*(thePos1.X()-x0) + (thePos1.Y()-y0)*(thePos1.Y()-y0) );
590  double R2 = std::sqrt ( (thePos2.X()-x0)*(thePos2.X()-x0) + (thePos2.Y()-y0)*(thePos2.Y()-y0) );
591  // Loop on primary vertices
592 
593  //TODO: Check if pTMin is correctly used (old code stored pTMin^2 in pTMin)
594 
595  for ( unsigned iv=0; iv<nVertices; ++iv )
596  {
597  // Z position of the primary vertex
598  double zV = (*vertices)[iv].z();
599  // Constraints on the inner hit
600  double checkRZ1 = forward ?
601  (thePos1.Z()-zV+zVertexConstraint) / (thePos2.Z()-zV+zVertexConstraint) * R2 :
602  -zVertexConstraint + R1/R2*(thePos2.Z()-zV+zVertexConstraint);
603  double checkRZ2 = forward ?
604  (thePos1.Z()-zV-zVertexConstraint)/(thePos2.Z()-zV-zVertexConstraint) * R2 :
605  +zVertexConstraint + R1/R2*(thePos2.Z()-zV-zVertexConstraint);
606  double checkRZmin = std::min(checkRZ1,checkRZ2)-3.*error;
607  double checkRZmax = std::max(checkRZ1,checkRZ2)+3.*error;
608  // Check if the innerhit is within bounds
609  bool compat = forward ?
610  checkRZmin < R1 && R1 < checkRZmax :
611  checkRZmin < thePos1.Z()-zV && thePos1.Z()-zV < checkRZmax;
612  // If it is, just return ok
613  if ( compat )
614  {
615  return compat;
616  }
617  }
618  // Otherwise, return not ok
619  return false;
620 
621 }
const reco::VertexCollection * vertices
T y() const
Definition: PV3DBase.h:63
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
T min(T a, T b)
Definition: MathUtil.h:58
T x() const
Definition: PV3DBase.h:62
const MagneticFieldMap * magneticFieldMap
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
const SeedingNode< TrackingLayer > * TrajectorySeedProducer::insertHit ( const std::vector< TrajectorySeedHitCandidate > &  trackerRecHits,
std::vector< int > &  hitIndicesInTree,
const SeedingNode< TrackingLayer > *  node,
unsigned int  trackerHit 
) const

method inserts hit into the tree structure at an empty position.

Parameters
trackerRecHitslist of all TrackerRecHits.
hitIndicesInTreehit indices which translates the tree node to the hits in trackerRecHits. Empty positions are identified with '-1'.
nodewhere to look for an empty position. Important for recursive tree traversing (Breadth-first). Starts with the root.
trackerHithit which is tested.
Returns
pointer if this hit is inserted at a leaf which means that a seed has been found. Returns 'nullptr' otherwise.

Definition at line 232 of file TrajectorySeedProducer.cc.

References SeedingNode< DATA >::getChild(), SeedingNode< DATA >::getChildrenSize(), SeedingNode< DATA >::getData(), SeedingNode< DATA >::getIndex(), SeedingNode< DATA >::getParent(), isHitOnLayer(), python.Node::node, passHitTuplesCuts(), and fileCollector::seed.

Referenced by iterateHits().

237 {
238  if (!node->getParent() || hitIndicesInTree[node->getParent()->getIndex()]>=0)
239  {
240  if (hitIndicesInTree[node->getIndex()]<0)
241  {
242  const TrajectorySeedHitCandidate& currentTrackerHit = trackerRecHits[trackerHit];
243  if (!isHitOnLayer(currentTrackerHit,node->getData()))
244  {
245  return nullptr;
246  }
247  if (!passHitTuplesCuts(*node,trackerRecHits,hitIndicesInTree,currentTrackerHit))
248  {
249  return nullptr;
250  }
251  hitIndicesInTree[node->getIndex()]=trackerHit;
252  if (node->getChildrenSize()==0)
253  {
254  return node;
255  }
256  return nullptr;
257  }
258  else
259  {
260  for (unsigned int ichild = 0; ichild<node->getChildrenSize(); ++ichild)
261  {
262  const SeedingNode<TrackingLayer>* seed = insertHit(trackerRecHits,hitIndicesInTree,node->getChild(ichild),trackerHit);
263  if (seed)
264  {
265  return seed;
266  }
267  }
268  }
269  }
270  return nullptr;
271 }
bool passHitTuplesCuts(const SeedingNode< TrackingLayer > &seedingNode, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, const std::vector< int > &hitIndicesInTree, const TrajectorySeedHitCandidate &currentTrackerHit) const
method checks if a TrajectorySeedHitCandidate fulfills the quality requirements.
const DATA & getData() const
Definition: SeedingTree.h:136
unsigned int getChildrenSize() const
Definition: SeedingTree.h:121
tuple node
Definition: Node.py:50
const SeedingNode * getParent() const
Definition: SeedingTree.h:112
bool isHitOnLayer(const TrajectorySeedHitCandidate &trackerRecHit, const TrackingLayer &layer) const
const SeedingNode< TrackingLayer > * insertHit(const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
method inserts hit into the tree structure at an empty position.
unsigned int getIndex() const
Definition: SeedingTree.h:107
const SeedingNode< DATA > * getChild(unsigned int ichild) const
Definition: SeedingTree.h:126
bool TrajectorySeedProducer::isHitOnLayer ( const TrajectorySeedHitCandidate trackerRecHit,
const TrackingLayer layer 
) const
inline

Definition at line 170 of file TrajectorySeedProducer.h.

References TrajectorySeedHitCandidate::getTrackingLayer().

Referenced by insertHit().

171  {
172  return layer==trackerRecHit.getTrackingLayer();
173  }
const TrackingLayer & getTrackingLayer() const
std::vector< unsigned int > TrajectorySeedProducer::iterateHits ( unsigned int  start,
const std::vector< TrajectorySeedHitCandidate > &  trackerRecHits,
std::vector< int >  hitIndicesInTree,
bool  processSkippedHits 
) const
virtual

method tries to insert all hits into the tree structure.

Parameters
startindex where to begin insertion. Important for recursion.
trackerRecHitslist of all TrackerRecHits.
hitIndicesInTreehit indices which translates the tree node to the hits in trackerRecHits.
currentTrackerHithit which is tested.
Returns
list of hit indices which form a found seed. Returns empty list if no seed was found.

Definition at line 274 of file TrajectorySeedProducer.cc.

References _seedingTree, SeedingNode< DATA >::getDepth(), SeedingNode< DATA >::getIndex(), SeedingNode< DATA >::getParent(), SeedingTree< DATA >::getRoot(), insertHit(), and SeedingTree< DATA >::numberOfRoots().

Referenced by produce().

280 {
281  for (unsigned int irecHit = start; irecHit<trackerRecHits.size(); ++irecHit)
282  {
283  unsigned int currentHitIndex=irecHit;
284 
285  for (unsigned int inext=currentHitIndex+1; inext< trackerRecHits.size(); ++inext)
286  {
287  //if multiple hits are on the same layer -> follow all possibilities by recusion
288  if (trackerRecHits[currentHitIndex].getTrackingLayer()==trackerRecHits[inext].getTrackingLayer())
289  {
290 
291  if (processSkippedHits)
292  {
293  //recusively call the method again with hit 'inext' but skip all following on the same layer though 'processSkippedHits=false'
294  std::vector<unsigned int> seedHits = iterateHits(
295  inext,
296  trackerRecHits,
297  hitIndicesInTree,
298  false
299  );
300  if (seedHits.size()>0)
301  {
302  return seedHits;
303  }
304 
305 
306  }
307  irecHit+=1;
308  }
309  else
310  {
311  break;
312  }
313  }
314 
315  processSkippedHits=true;
316 
317  const SeedingNode<TrackingLayer>* seedNode = nullptr;
318  for (unsigned int iroot=0; seedNode==nullptr && iroot<_seedingTree.numberOfRoots(); ++iroot)
319  {
320  seedNode=insertHit(trackerRecHits,hitIndicesInTree,_seedingTree.getRoot(iroot), currentHitIndex);
321  }
322  if (seedNode)
323  {
324  std::vector<unsigned int> seedIndices(seedNode->getDepth()+1);
325  while (seedNode)
326  {
327  seedIndices[seedNode->getDepth()]=hitIndicesInTree[seedNode->getIndex()];
328  seedNode=seedNode->getParent();
329  }
330  return seedIndices;
331  }
332 
333  }
334 
335  return std::vector<unsigned int>();
336 
337 }
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
unsigned int getDepth() const
Definition: SeedingTree.h:84
const SeedingNode * getParent() const
Definition: SeedingTree.h:112
unsigned int numberOfRoots() const
Definition: SeedingTree.h:222
virtual std::vector< unsigned int > iterateHits(unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
method tries to insert all hits into the tree structure.
const SeedingNode< DATA > * getRoot(unsigned int i) const
Definition: SeedingTree.h:232
const SeedingNode< TrackingLayer > * insertHit(const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > &hitIndicesInTree, const SeedingNode< TrackingLayer > *node, unsigned int trackerHit) const
method inserts hit into the tree structure at an empty position.
unsigned int getIndex() const
Definition: SeedingTree.h:107
SeedingTree< TrackingLayer > _seedingTree
bool TrajectorySeedProducer::pass2HitsCuts ( const TrajectorySeedHitCandidate hit1,
const TrajectorySeedHitCandidate hit2 
) const

Definition at line 214 of file TrajectorySeedProducer.cc.

References compatibleWithBeamAxis(), relativeConstraints::error, TrajectorySeedHitCandidate::globalPosition(), TrajectorySeedHitCandidate::isForward(), TrajectorySeedHitCandidate::largerError(), skipPVCompatibility, and mathSSE::sqrt().

Referenced by passHitTuplesCuts().

215 {
216  bool compatible=false;
218  {
219  compatible = true;
220  }
221  else
222  {
223  GlobalPoint gpos1 = hit1.globalPosition();
224  GlobalPoint gpos2 = hit2.globalPosition();
225  bool forward = hit1.isForward();
226  double error = std::sqrt(hit1.largerError()+hit2.largerError());
227  compatible = compatibleWithBeamAxis(gpos1,gpos2,error,forward);
228  }
229  return compatible;
230 }
bool isForward() const
Is it a forward hit ?
GlobalPoint globalPosition() const
The global position.
T sqrt(T t)
Definition: SSEVec.h:48
bool compatibleWithBeamAxis(const GlobalPoint &gpos1, const GlobalPoint &gpos2, double error, bool forward) const
bool TrajectorySeedProducer::passHitTuplesCuts ( const SeedingNode< TrackingLayer > &  seedingNode,
const std::vector< TrajectorySeedHitCandidate > &  trackerRecHits,
const std::vector< int > &  hitIndicesInTree,
const TrajectorySeedHitCandidate currentTrackerHit 
) const
inline

method checks if a TrajectorySeedHitCandidate fulfills the quality requirements.

Parameters
seedingNodetree node at which the hit will be inserted.
trackerRecHitslist of all TrackerRecHits.
hitIndicesInTreehit indices which translates the tree node to the hits in trackerRecHits.
currentTrackerHithit which is tested.
Returns
true if a hit fulfills the requirements.

Definition at line 109 of file TrajectorySeedProducer.h.

References SeedingNode< DATA >::getDepth(), SeedingNode< DATA >::getIndex(), SeedingNode< DATA >::getParent(), and pass2HitsCuts().

Referenced by insertHit().

115  {
116  switch (seedingNode.getDepth())
117  {
118  case 0:
119  {
120  return true;
121  /* example for 1 hits
122  const TrajectorySeedHitCandidate& hit1 = currentTrackerHit;
123  return pass1HitsCuts(hit1,trackingAlgorithmId);
124  */
125  }
126 
127  case 1:
128  {
129  const SeedingNode<TrackingLayer>* parentNode = &seedingNode;
130  parentNode = parentNode->getParent();
131  const TrajectorySeedHitCandidate& hit1 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
132  const TrajectorySeedHitCandidate& hit2 = currentTrackerHit;
133 
134  return pass2HitsCuts(hit1,hit2);
135  }
136  case 2:
137  {
138  return true;
139  /* example for 3 hits
140  const SeedingNode<LayerSpec>* parentNode = &seedingNode;
141  parentNode = parentNode->getParent();
142  const TrajectorySeedHitCandidate& hit2 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
143  parentNode = parentNode->getParent();
144  const TrajectorySeedHitCandidate& hit1 = trackerRecHits[hitIndicesInTree[parentNode->getIndex()]];
145  const TrajectorySeedHitCandidate& hit3 = currentTrackerHit;
146  return pass3HitsCuts(hit1,hit2,hit3,trackingAlgorithmId);
147  */
148  }
149  }
150  return true;
151  }
unsigned int getDepth() const
Definition: SeedingTree.h:84
const SeedingNode * getParent() const
Definition: SeedingTree.h:112
bool pass2HitsCuts(const TrajectorySeedHitCandidate &hit1, const TrajectorySeedHitCandidate &hit2) const
unsigned int getIndex() const
Definition: SeedingTree.h:107
bool TrajectorySeedProducer::passSimTrackQualityCuts ( const SimTrack theSimTrack,
const SimVertex theSimVertex 
) const
virtual

method checks if a SimTrack fulfills the quality requirements.

Parameters
theSimTrackthe SimTrack to be tested.
theSimVertexthe associated SimVertex of the SimTrack.
Returns
true if a track fulfills the requirements.

Definition at line 172 of file TrajectorySeedProducer.cc.

References beamspotPosition, CoreSimTrack::charge(), maxD0, maxZ0, CoreSimTrack::momentum(), CoreSimVertex::position(), pTMin, RawParticle::setCharge(), BaseParticlePropagator::xyImpactParameter(), and BaseParticlePropagator::zImpactParameter().

Referenced by produce().

173 {
174 
175  //require min pT of the simtrack
176  if ( theSimTrack.momentum().Perp2() < pTMin*pTMin)
177  {
178  return false;
179  }
180 
181  //require impact parameter of the simtrack
183  RawParticle(
185  theSimTrack.momentum().px(),
186  theSimTrack.momentum().py(),
187  theSimTrack.momentum().pz(),
188  theSimTrack.momentum().e()
189  ),
191  theSimVertex.position().x(),
192  theSimVertex.position().y(),
193  theSimVertex.position().z(),
194  theSimVertex.position().t())
195  ),
196  0.,0.,4.
197  );
198  theParticle.setCharge(theSimTrack.charge());
199  const double x0 = beamspotPosition.X();
200  const double y0 = beamspotPosition.Y();
201  const double z0 = beamspotPosition.Z();
202  if ( theParticle.xyImpactParameter(x0,y0) > maxD0 )
203  {
204  return false;
205  }
206  if ( fabs( theParticle.zImpactParameter(x0,y0) - z0 ) > maxZ0)
207  {
208  return false;
209  }
210  return true;
211 }
void setCharge(float q)
set the MEASURED charge
Definition: RawParticle.cc:139
double zImpactParameter(double x0=0, double y0=0.) const
Longitudinal impact parameter.
float charge() const
charge
Definition: CoreSimTrack.cc:18
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:26
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
double xyImpactParameter(double x0=0., double y0=0.) const
Transverse impact parameter.
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
void TrajectorySeedProducer::produce ( edm::Event e,
const edm::EventSetup es 
)
virtual

should check if corresponds to m

Implements edm::stream::EDProducerBase.

Definition at line 340 of file TrajectorySeedProducer.cc.

References _seedingTree, alongMomentum, edm::OwnVector< T, P >::back(), beamspotPosition, beamSpotToken, CoreSimTrack::charge(), RecoTauCleanerPlugins::charge, TrackingRecHit::clone(), TrackingRecHit::geographicalId(), edm::Event::getByToken(), SeedingTree< DATA >::getSingleSet(), TrajectorySeedHitCandidate::getTrackingLayer(), TrajectoryStateOnSurface::globalMomentum(), TrackingRecHit::hit(), i, TrackerGeometry::idToDet(), TrajectorySeedHitCandidate::isOnTheSameLayer(), TrajectoryStateOnSurface::isValid(), iterateHits(), j, relval_steps::k, TrajectoryStateOnSurface::localError(), TrajectoryStateOnSurface::localParameters(), visualization-live-secondInstance_cfg::m, magneticField, LocalTrajectoryError::matrix(), minRecHits, CoreSimTrack::momentum(), eostools::move(), SeedingTree< DATA >::numberOfNodes(), convertSQLitetoXML_cfg::output, outputSeedCollectionName, passSimTrackQualityCuts(), PV3DBase< T, PVType, FrameType >::perp(), position, Propagator::propagate(), edm::OwnVector< T, P >::push_back(), edm::Event::put(), DetId::rawId(), HLT_25ns14e33_v3_cff::recHits, recHitToken, recoVertexToken, simTrackToken, simVertexToken, skipSimTrackIdTokens, GeomDetEnumerators::subDetId, GeomDet::surface(), TrajectoryStateOnSurface::surfaceSide(), thePropagator, trackerGeometry, trackerTopology, vertices, SimTrack::vertIndex(), x, detailsBasic3DVector::y, and detailsBasic3DVector::z.

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

341 {
342  // unsigned nTrackCandidates = 0;
343  PTrajectoryStateOnDet initialState;
344 
345  // First, the tracks to be removed
346  std::set<unsigned int> skipSimTrackIds;
347  for ( unsigned int i=0; i<skipSimTrackIdTokens.size(); ++i ) {
348  edm::Handle<std::vector<int> > skipSimTrackIds_temp;
349  e.getByToken(skipSimTrackIdTokens[i],skipSimTrackIds_temp);
350  for ( unsigned int j=0; j<skipSimTrackIds_temp->size(); ++j ) {
351  unsigned int mySimTrackId = (*skipSimTrackIds_temp)[j];
352  skipSimTrackIds.insert((unsigned int)mySimTrackId);
353  }
354  }
355 
356  // Beam spot
357  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
358  e.getByToken(beamSpotToken,recoBeamSpotHandle);
359  beamspotPosition = recoBeamSpotHandle->position();
360 
361  // SimTracks and SimVertices
363  e.getByToken(simTrackToken,theSimTracks);
364 
366  e.getByToken(simVertexToken,theSimVtx);
367 
368  // edm::Handle<SiTrackerGSRecHit2DCollection> theGSRecHits;
370  e.getByToken(recHitToken, theGSRecHits);
371 
372  // Primary vertices
374  if (e.getByToken(recoVertexToken,theRecVtx))
375  {
376 
377  //this can be nullptr if the PV compatiblity should not be tested against
378  vertices = &(*theRecVtx);
379  }
380 
381 
382  // Output - gets moved, no delete needed
383  std::auto_ptr<TrajectorySeedCollection> output{new TrajectorySeedCollection()};
384 
385  //if no hits -> directly write empty collection
386  if(theGSRecHits->size() == 0)
387  {
389  return;
390  }
391  for (SiTrackerGSMatchedRecHit2DCollection::id_iterator itSimTrackId=theGSRecHits->id_begin(); itSimTrackId!=theGSRecHits->id_end(); ++itSimTrackId )
392  {
393 
394  const unsigned int currentSimTrackId = *itSimTrackId;
395 
396  if(skipSimTrackIds.find(currentSimTrackId)!=skipSimTrackIds.end()) continue;
397 
398  const SimTrack& theSimTrack = (*theSimTracks)[currentSimTrackId];
399 
400  int vertexIndex = theSimTrack.vertIndex();
401  if (vertexIndex<0)
402  {
403  //tracks are required to be associated to a vertex
404  continue;
405  }
406  const SimVertex& theSimVertex = (*theSimVtx)[vertexIndex];
407 
408  if (!this->passSimTrackQualityCuts(theSimTrack,theSimVertex))
409  {
410  continue;
411 
412  }
413  SiTrackerGSMatchedRecHit2DCollection::range recHitRange = theGSRecHits->get(currentSimTrackId);
414 
415  TrajectorySeedHitCandidate previousTrackerHit;
416  TrajectorySeedHitCandidate currentTrackerHit;
417  unsigned int numberOfNonEqualHits=0;
418 
419  std::vector<TrajectorySeedHitCandidate> trackerRecHits;
420  for (SiTrackerGSMatchedRecHit2DCollection::const_iterator itRecHit = recHitRange.first; itRecHit!=recHitRange.second; ++itRecHit)
421  {
422  const SiTrackerGSMatchedRecHit2D& vec = *itRecHit;
423  previousTrackerHit=currentTrackerHit;
424 
426 
427  if (!currentTrackerHit.isOnTheSameLayer(previousTrackerHit))
428  {
429  ++numberOfNonEqualHits;
430  }
431 
432 
433  if (_seedingTree.getSingleSet().find(currentTrackerHit.getTrackingLayer())!=_seedingTree.getSingleSet().end())
434  {
435  //add only the hits which are actually on the requested layers
436  trackerRecHits.push_back(std::move(currentTrackerHit));
437  }
438 
439  }
440  if ( numberOfNonEqualHits < minRecHits) continue;
441 
442  std::vector<int> hitIndicesInTree(_seedingTree.numberOfNodes(),-1);
443  //A SeedingNode is associated by its index to this list. The list stores the indices of the hits in 'trackerRecHits'
444  /* example
445  SeedingNode | hit index | hit
446  -------------------------------------------------------------------------------
447  index= 0: [BPix1] | hitIndicesInTree[0] (=1) | trackerRecHits[1]
448  index= 1: -- [BPix2] | hitIndicesInTree[1] (=3) | trackerRecHits[3]
449  index= 2: -- -- [BPix3] | hitIndicesInTree[2] (=4) | trackerRecHits[4]
450  index= 3: -- -- [FPix1_pos] | hitIndicesInTree[3] (=6) | trackerRecHits[6]
451  index= 4: -- -- [FPix1_neg] | hitIndicesInTree[4] (=7) | trackerRecHits[7]
452 
453  The implementation has been chosen such that the tree only needs to be build once upon construction.
454  */
455 
456 
457  std::vector<unsigned int> seedHitNumbers = iterateHits(0,trackerRecHits,hitIndicesInTree,true);
458 
459  if (seedHitNumbers.size()>0)
460  {
461 
462 
464  for ( unsigned ihit=0; ihit<seedHitNumbers.size(); ++ihit )
465  {
466  TrackingRecHit* aTrackingRecHit = trackerRecHits[seedHitNumbers[ihit]].hit()->clone();
467  recHits.push_back(aTrackingRecHit);
468  }
469 
470 
471  GlobalPoint position((*theSimVtx)[vertexIndex].position().x(),
472  (*theSimVtx)[vertexIndex].position().y(),
473  (*theSimVtx)[vertexIndex].position().z());
474 
475  GlobalVector momentum(theSimTrack.momentum().x(),theSimTrack.momentum().y(),theSimTrack.momentum().z());
476  float charge = theSimTrack.charge();
477  GlobalTrajectoryParameters initialParams(position,momentum,(int)charge,magneticField);
479  //this line help the fit succeed in the case of pixelless tracks (4th and 5th iteration)
480  //for the future: probably the best thing is to use the mini-kalmanFilter
481  if(trackerRecHits[seedHitNumbers[0]].subDetId() !=1 ||trackerRecHits[seedHitNumbers[0]].subDetId() !=2)
482  {
483  errorMatrix = errorMatrix * 0.0000001;
484  }
485  CurvilinearTrajectoryError initialError(errorMatrix);
486  FreeTrajectoryState initialFTS(initialParams, initialError);
487 
488 
489  const GeomDet* initialLayer = trackerGeometry->idToDet( recHits.back().geographicalId() );
490  const TrajectoryStateOnSurface initialTSOS = thePropagator->propagate(initialFTS,initialLayer->surface()) ;
491 
492 
493  if (!initialTSOS.isValid())
494  {
495  break; //continues with the next seeding algorithm
496  }
497 
498  const AlgebraicSymMatrix55& m = initialTSOS.localError().matrix();
499  int dim = 5;
500  float localErrors[15];
501  int k = 0;
502  for (int i=0; i<dim; ++i)
503  {
504  for (int j=0; j<=i; ++j)
505  {
506  localErrors[k++] = m(i,j);
507  }
508  }
509  int surfaceSide = static_cast<int>(initialTSOS.surfaceSide());
510  initialState = PTrajectoryStateOnDet( initialTSOS.localParameters(),initialTSOS.globalMomentum().perp(),localErrors, recHits.back().geographicalId().rawId(), surfaceSide);
511  output->push_back(TrajectorySeed(initialState, recHits, PropagationDirection::alongMomentum));
512 
513 
514 
515  }
516  } //end loop over simtracks
517 
518 
520 }
virtual FreeTrajectoryState propagate(const FreeTrajectoryState &ftsStart, const GlobalPoint &pDest) const final
Definition: Propagator.h:119
int i
Definition: DBlmapReader.cc:9
reference back()
Definition: OwnVector.h:327
const MagneticField * magneticField
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
T perp() const
Definition: PV3DBase.h:72
const LocalTrajectoryParameters & localParameters() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
const SingleSet & getSingleSet() const
Definition: SeedingTree.h:204
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
const reco::VertexCollection * vertices
identifier iterator
Definition: RangeMap.h:136
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
float charge() const
charge
Definition: CoreSimTrack.cc:18
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
edm::EDGetTokenT< SiTrackerGSMatchedRecHit2DCollection > recHitToken
float float float z
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
PropagatorWithMaterial * thePropagator
unsigned int numberOfNodes() const
Definition: SeedingTree.h:227
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
void push_back(D *&d)
Definition: OwnVector.h:280
const TrackerTopology * trackerTopology
std::vector< TrajectorySeed > TrajectorySeedCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
unsigned int subDetId[18]
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
def move
Definition: eostools.py:508
virtual std::vector< unsigned int > iterateHits(unsigned int start, const std::vector< TrajectorySeedHitCandidate > &trackerRecHits, std::vector< int > hitIndicesInTree, bool processSkippedHits) const
method tries to insert all hits into the tree structure.
int j
Definition: DBlmapReader.cc:9
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
edm::EDGetTokenT< edm::SimVertexContainer > simVertexToken
virtual TrackingRecHit * clone() const =0
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
std::vector< edm::EDGetTokenT< std::vector< int > > > skipSimTrackIdTokens
edm::EDGetTokenT< edm::SimTrackContainer > simTrackToken
virtual bool passSimTrackQualityCuts(const SimTrack &theSimTrack, const SimVertex &theSimVertex) const
method checks if a SimTrack fulfills the quality requirements.
virtual TrackingRecHit const * hit() const
GlobalVector globalMomentum() const
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:22
static int position[264][3]
Definition: ReadPGInfo.cc:509
edm::EDGetTokenT< reco::VertexCollection > recoVertexToken
const TrackerGeometry * trackerGeometry
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken
DetId geographicalId() const
const TrackingLayer & getTrackingLayer() const
Definition: DDAxes.h:10
SeedingTree< TrackingLayer > _seedingTree
bool isOnTheSameLayer(const TrajectorySeedHitCandidate &other) const
Check if two hits are on the same layer of the same subdetector.
virtual const TrackerGeomDet * idToDet(DetId) const

Member Data Documentation

SeedingTree<TrackingLayer> TrajectorySeedProducer::_seedingTree
private

Definition at line 38 of file TrajectorySeedProducer.h.

Referenced by iterateHits(), produce(), and TrajectorySeedProducer().

unsigned int TrajectorySeedProducer::absMinRecHits
private

Definition at line 56 of file TrajectorySeedProducer.h.

Referenced by TrajectorySeedProducer().

math::XYZPoint TrajectorySeedProducer::beamspotPosition
private
edm::EDGetTokenT<reco::BeamSpot> TrajectorySeedProducer::beamSpotToken
private

Definition at line 77 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

edm::InputTag TrajectorySeedProducer::hitProducer
private

Definition at line 52 of file TrajectorySeedProducer.h.

Referenced by TrajectorySeedProducer().

const MagneticField* TrajectorySeedProducer::magneticField
private

Definition at line 41 of file TrajectorySeedProducer.h.

Referenced by beginRun(), and produce().

const MagneticFieldMap* TrajectorySeedProducer::magneticFieldMap
private

Definition at line 42 of file TrajectorySeedProducer.h.

Referenced by beginRun(), and compatibleWithBeamAxis().

double TrajectorySeedProducer::maxD0
private

Definition at line 49 of file TrajectorySeedProducer.h.

Referenced by passSimTrackQualityCuts(), and TrajectorySeedProducer().

double TrajectorySeedProducer::maxZ0
private

Definition at line 50 of file TrajectorySeedProducer.h.

Referenced by passSimTrackQualityCuts(), and TrajectorySeedProducer().

unsigned int TrajectorySeedProducer::minRecHits
private

Definition at line 51 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

unsigned int TrajectorySeedProducer::numberOfHits
private

Definition at line 57 of file TrajectorySeedProducer.h.

Referenced by TrajectorySeedProducer().

double TrajectorySeedProducer::originHalfLength
private

Definition at line 67 of file TrajectorySeedProducer.h.

Referenced by compatibleWithBeamAxis(), and TrajectorySeedProducer().

double TrajectorySeedProducer::originpTMin
private

Definition at line 68 of file TrajectorySeedProducer.h.

Referenced by compatibleWithBeamAxis(), and TrajectorySeedProducer().

double TrajectorySeedProducer::originRadius
private

Definition at line 66 of file TrajectorySeedProducer.h.

Referenced by compatibleWithBeamAxis(), and TrajectorySeedProducer().

std::string TrajectorySeedProducer::outputSeedCollectionName
private

Definition at line 59 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

double TrajectorySeedProducer::pTMin
private

Definition at line 48 of file TrajectorySeedProducer.h.

Referenced by passSimTrackQualityCuts(), and TrajectorySeedProducer().

edm::EDGetTokenT<SiTrackerGSMatchedRecHit2DCollection> TrajectorySeedProducer::recHitToken
private

Definition at line 80 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

edm::EDGetTokenT<reco::VertexCollection> TrajectorySeedProducer::recoVertexToken
private

Definition at line 81 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

bool TrajectorySeedProducer::seedCleaning
private

Definition at line 55 of file TrajectorySeedProducer.h.

Referenced by compatibleWithBeamAxis(), and TrajectorySeedProducer().

std::vector<std::vector<TrackingLayer> > TrajectorySeedProducer::seedingLayers
private

Definition at line 62 of file TrajectorySeedProducer.h.

Referenced by TrajectorySeedProducer().

edm::EDGetTokenT<edm::SimTrackContainer> TrajectorySeedProducer::simTrackToken
private

Definition at line 78 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

edm::EDGetTokenT<edm::SimVertexContainer> TrajectorySeedProducer::simVertexToken
private

Definition at line 79 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

bool TrajectorySeedProducer::skipPVCompatibility
private

Definition at line 72 of file TrajectorySeedProducer.h.

Referenced by pass2HitsCuts(), and TrajectorySeedProducer().

std::vector<edm::EDGetTokenT<std::vector<int> > > TrajectorySeedProducer::skipSimTrackIdTokens
private

Definition at line 82 of file TrajectorySeedProducer.h.

Referenced by produce(), and TrajectorySeedProducer().

edm::InputTag TrajectorySeedProducer::theBeamSpot
private

Definition at line 53 of file TrajectorySeedProducer.h.

Referenced by TrajectorySeedProducer().

PropagatorWithMaterial* TrajectorySeedProducer::thePropagator
private

Definition at line 46 of file TrajectorySeedProducer.h.

Referenced by beginRun(), produce(), and ~TrajectorySeedProducer().

const TrackerGeometry* TrajectorySeedProducer::trackerGeometry
private

Definition at line 43 of file TrajectorySeedProducer.h.

Referenced by beginRun(), and produce().

const TrackerTopology* TrajectorySeedProducer::trackerTopology
private

Definition at line 44 of file TrajectorySeedProducer.h.

Referenced by beginRun(), and produce().

const reco::VertexCollection* TrajectorySeedProducer::vertices
private

Definition at line 74 of file TrajectorySeedProducer.h.

Referenced by compatibleWithBeamAxis(), and produce().

double TrajectorySeedProducer::zVertexConstraint
private

Definition at line 70 of file TrajectorySeedProducer.h.

Referenced by compatibleWithBeamAxis(), and TrajectorySeedProducer().