CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastElectronSeedGenerator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EgammaElectronAlgos
4 // Class: FastElectronSeedGenerator.
5 //
13 //
14 // Original Author: Patrick Janot
15 //
16 //
18 
25 
28 
31 
34 
40 
45 
48 
49 #include <vector>
50 
51 //#define FAMOS_DEBUG
52 
54  const edm::ParameterSet &pset,
55  double pTMin,
56  const edm::InputTag& beamSpot)
57  :
58  dynamicphiroad_(pset.getParameter<bool>("dynamicPhiRoad")),
59  lowPtThreshold_(pset.getParameter<double>("LowPtThreshold")),
60  highPtThreshold_(pset.getParameter<double>("HighPtThreshold")),
61  sizeWindowENeg_(pset.getParameter<double>("SizeWindowENeg")),
62  phimin2_(pset.getParameter<double>("PhiMin2")),
63  phimax2_(pset.getParameter<double>("PhiMax2")),
64  deltaPhi1Low_(pset.getParameter<double>("DeltaPhi1Low")),
65  deltaPhi1High_(pset.getParameter<double>("DeltaPhi1High")),
66  deltaPhi2_(pset.getParameter<double>("DeltaPhi2")),
67  pTMin2(pTMin*pTMin),
68  myGSPixelMatcher(0),
69  fromTrackerSeeds_(pset.getParameter<bool>("fromTrackerSeeds")),
70  theUpdator(0), thePropagator(0),
71  // theMeasurementTracker(0),
72  // theNavigationSchool(0)
73  theSetup(0), theBeamSpot(beamSpot)
74 {
75 
76 #ifdef FAMOS_DEBUG
77  std::cout << "FromTrackerSeeds = " << fromTrackerSeeds_ << std::endl;
78 #endif
79 
80  // Instantiate the pixel hit matcher
81  searchInTIDTEC = pset.getParameter<bool>("searchInTIDTEC");
82  myGSPixelMatcher = new FastPixelHitMatcher(pset.getParameter<double>("ePhiMin1"),
83  pset.getParameter<double>("ePhiMax1"),
84  pset.getParameter<double>("pPhiMin1"),
85  pset.getParameter<double>("pPhiMax1"),
86  pset.getParameter<double>("PhiMin2"),
87  pset.getParameter<double>("PhiMax2"),
88  pset.getParameter<double>("z2MinB"),
89  pset.getParameter<double>("z2MaxB"),
90  pset.getParameter<double>("r2MinF"),
91  pset.getParameter<double>("r2MaxF"),
92  pset.getParameter<double>("rMinI"),
93  pset.getParameter<double>("rMaxI"),
94  pset.getParameter<bool>("searchInTIDTEC"));
95 
96 }
97 
99 
100  // delete theNavigationSchool;
101  delete myGSPixelMatcher;
102  // delete myMatchPos;
103  delete thePropagator;
104  delete theUpdator;
105 
106 }
107 
108 
110 
111  theSetup= &setup;
112 
114  setup.get<IdealMagneticFieldRecord>().get(pMF);
115  theMagField = &(*pMF);
116 
118  setup.get<TrackerDigiGeometryRecord>().get(geometry);
119  theTrackerGeometry = &(*geometry);
120 
122  setup.get<TrackerRecoGeometryRecord>().get( recoGeom );
123  theGeomSearchTracker = &(*recoGeom);
124 
126  setup.get<TrackerInteractionGeometryRecord>().get( interGeom );
127  theTrackerInteractionGeometry = &(*interGeom);
128 
130  setup.get<MagneticFieldMapRecord>().get(fieldMap);
131  theMagneticFieldMap = &(*fieldMap);
132 
134 
139 
140 }
141 
143  const reco::SuperClusterRefVector &sclRefs,
144  const SiTrackerGSMatchedRecHit2DCollection* theGSRecHits,
145  const edm::SimTrackContainer* theSimTracks,
148 
149  // Take the seed collection.
150  theInitialSeedColl=seeds;
151 
152  // Get the beam spot
153  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
154  e.getByLabel(theBeamSpot,recoBeamSpotHandle);
155 
156  // Get its position
157  BSPosition_ = recoBeamSpotHandle->position();
158  double sigmaZ=recoBeamSpotHandle->sigmaZ();
159  double sigmaZ0Error=recoBeamSpotHandle->sigmaZ0Error();
160  double sq=sqrt(sigmaZ*sigmaZ+sigmaZ0Error*sigmaZ0Error);
161  double zmin1 = BSPosition_.z()-3*sq;
162  double zmax1 = BSPosition_.z()+3*sq;
163 #ifdef FAMOS_DEBUG
164  std::cout << "Z Range for pixel matcher : " << zmin1 << " " << BSPosition_.z() << " " << zmax1 << std::endl;
165 #endif
166  myGSPixelMatcher->set1stLayerZRange(zmin1,zmax1);
167 
168  // A map of vector of pixel seeds, for each clusters
169  std::map<unsigned,std::vector<reco::ElectronSeed> > myPixelSeeds;
170 
171  // No seeding attempted if no hits !
172  if(theGSRecHits->size() == 0) return;
173 
174  if ( !fromTrackerSeeds_ ) {
175 
176  // The vector of simTrack Id's carrying GSRecHits
177  const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
178 
179  // Loop over the simTrack carrying GSRecHits
180  for ( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
181 
182  unsigned simTrackId = theSimTrackIds[tkId];
183  const SimTrack& theSimTrack = (*theSimTracks)[simTrackId];
184 
185  // Request a minimum pT for the sim track
186  if ( theSimTrack.momentum().perp2() < pTMin2 ) continue;
187 
188  // Request a minimum number of RecHits (total and in the pixel detector)
189  unsigned numberOfRecHits = 0;
190 
191  // The vector of rechits for seeding
192 
193  // 1) Cluster-pixel match seeding:
194  // Save a collection of Pixel +TEC +TID hits for seeding electrons
195  std::vector<unsigned> layerHit(6,static_cast<unsigned>(0));
196  // const SiTrackerGSMatchedRecHit2D *hit;
197  TrackerRecHit currentHit;
198  std::vector<TrackerRecHit> theHits;
199  TrajectorySeed theTrackerSeed;
200 
201  SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
202  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
203  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
206 
207  for ( iterRecHit = theRecHitRangeIteratorBegin;
208  iterRecHit != theRecHitRangeIteratorEnd;
209  ++iterRecHit) {
210  ++numberOfRecHits;
211 
212  currentHit = TrackerRecHit(&(*iterRecHit),theTrackerGeometry);
213  if ( ( currentHit.subDetId() <= 2 ) || // Pixel Hits
214  // Add TID/TEC (optional)
215  ( searchInTIDTEC &&
216  ( ( currentHit.subDetId() == 3 &&
217  currentHit.ringNumber() < 3 &&
218  currentHit.layerNumber() < 3 ) || // TID first two rings, first two layers
219  ( currentHit.subDetId() == 6 &&
220  currentHit.ringNumber() < 3 &&
221  currentHit.layerNumber() < 3 ) ) ) ) // TEC first two rings, first two layers
222  theHits.push_back(currentHit);
223  }
224 
225  // At least 3 hits
226  if ( numberOfRecHits < 3 ) continue;
227 
228  // At least 2 pixel hits
229  if ( theHits.size() < 2 ) continue;
230 
231  // Loop over clusters
232 
233  unsigned csize = sclRefs.size();
234  for (unsigned int i=0;i<csize;++i) {
235 
236  // Find the pixel seeds (actually only the best one is returned)
237  LogDebug ("run") << "new cluster, calling addAseedFromThisCluster";
238  addASeedToThisCluster(sclRefs[i],theHits,theTrackerSeed,myPixelSeeds[i]);
239 
240  }
241 
242  }
243  // 2) Check if the seed is in the a-priori seed collection
244  } else {
245 
246  // Loop over the tracker seed
247 #ifdef FAMOS_DEBUG
248  std::cout << "We have " << seeds->size() << " tracker seeds!" << std::endl;
249 #endif
250  for (unsigned int i=0;i<seeds->size();++i) {
251 
252  TrackerRecHit currentHit;
253  std::vector<TrackerRecHit> theHits;
254  const TrajectorySeed& theTrackerSeed = (*seeds)[i];
255  TrajectorySeed::range theSeedRange=theTrackerSeed.recHits();
256  TrajectorySeed::const_iterator theSeedRangeIteratorBegin = theSeedRange.first;
257  TrajectorySeed::const_iterator theSeedRangeIteratorEnd = theSeedRange.second;
258  TrajectorySeed::const_iterator theSeedItr = theSeedRangeIteratorBegin;
259 
260  for ( ; theSeedItr != theSeedRangeIteratorEnd; ++theSeedItr ) {
261  const SiTrackerGSMatchedRecHit2D * theSeedingRecHit =
262  (const SiTrackerGSMatchedRecHit2D*) (&(*theSeedItr));
263  currentHit = TrackerRecHit(theSeedingRecHit,theTrackerGeometry);
264  theHits.push_back(currentHit);
265  }
266 
267  // Loop over clusters
268  unsigned csize = sclRefs.size();
269  for (unsigned int i=0;i<csize;++i) {
270 
271  // Find the pixel seeds (actually only the best one is returned)
272 #ifdef FAMOS_DEBUG
273  std::cout << "new cluster, calling addAseedFromThisCluster" << std::endl;
274 #endif
275  addASeedToThisCluster(sclRefs[i],theHits,theTrackerSeed,myPixelSeeds[i]);
276 
277  }
278  // End loop over clusters
279  }
280  // End loop over seeds
281  }
282  // end else
283 
284  // Back to the expected collection
285 
286  std::map<unsigned,std::vector<reco::ElectronSeed> >::const_iterator px = myPixelSeeds.begin();
287  std::map<unsigned,std::vector<reco::ElectronSeed> >::const_iterator pxEnd = myPixelSeeds.end();
288  for ( ; px!=pxEnd; ++px ) {
289  unsigned nSeeds = (px->second).size();
290  for ( unsigned ipx = 0; ipx<nSeeds; ++ipx ) {
291  out.push_back((px->second)[ipx]);
292  reco::ElectronSeed is = px->second[ipx];
293  }
294  }
295 
296  LogDebug ("run") << ": For event "<<e.id();
297  LogDebug ("run") <<"Nr of superclusters: "<<sclRefs.size()
298  <<", no. of ElectronSeeds found = " << out.size();
299 #ifdef FAMOS_DEBUG
300  std::cout << ": For event "<<e.id() << std::endl;
301  std::cout <<"Nr of superclusters: "<<sclRefs.size()
302  <<", no. of ElectronSeeds found = " << out.size() << std::endl;
303 #endif
304 
305 }
306 
307 void
309  std::vector<TrackerRecHit>& theHits,
310  const TrajectorySeed& theTrackerSeed,
311  std::vector<reco::ElectronSeed>& result)
312 {
313 
314  float clusterEnergy = seedCluster->energy();
315  GlobalPoint clusterPos(seedCluster->position().x(),
316  seedCluster->position().y(),
317  seedCluster->position().z());
318  const GlobalPoint vertexPos(BSPosition_.x(),BSPosition_.y(),BSPosition_.z());
319 
320 #ifdef FAMOS_DEBUG
321  std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] "
322  << "new supercluster with energy: " << clusterEnergy << std::endl;
323  std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] "
324  << "and position: " << clusterPos << std::endl;
325  std::cout << "Vertex position : " << vertexPos << std::endl;
326 #endif
327 
328  //Here change the deltaPhi window of the first pixel layer in function of the seed pT
329  if (dynamicphiroad_) {
330  float clusterEnergyT = clusterEnergy*sin(seedCluster->position().theta()) ;
331 
332  float deltaPhi1 = 0.875/clusterEnergyT + 0.055;
333  if (clusterEnergyT < lowPtThreshold_) deltaPhi1= deltaPhi1Low_;
334  if (clusterEnergyT > highPtThreshold_) deltaPhi1= deltaPhi1High_;
335 
336  float ephimin1 = -deltaPhi1*sizeWindowENeg_ ;
337  float ephimax1 = deltaPhi1*(1.-sizeWindowENeg_);
338  float pphimin1 = -deltaPhi1*(1.-sizeWindowENeg_);
339  float pphimax1 = deltaPhi1*sizeWindowENeg_;
340 
341  float phimin2 = -deltaPhi2_/2.;
342  float phimax2 = deltaPhi2_/2.;
343 
344  myGSPixelMatcher->set1stLayer(ephimin1,ephimax1,pphimin1,pphimax1);
345  myGSPixelMatcher->set2ndLayer(phimin2,phimax2);
346 
347  }
348 
349 
350 
352 
353  // Find the best pixel pair compatible with the cluster
354  std::vector<std::pair<ConstRecHitPointer,ConstRecHitPointer> > compatPixelHits =
355  myGSPixelMatcher->compatibleHits(clusterPos, vertexPos, clusterEnergy, theHits);
356 
357  // The corresponding origin vertex
358  double vertexZ = myGSPixelMatcher->getVertex();
359  GlobalPoint theVertex(BSPosition_.x(),BSPosition_.y(),vertexZ);
360 
361  // Create the Electron pixel seed.
362  if (!compatPixelHits.empty() ) {
363 #ifdef FAMOS_DEBUG
364  std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] "
365  << " electron compatible hits found " << std::endl;
366 #endif
367  // Pixel-matching case: create the seed from scratch
368  if (!fromTrackerSeeds_) {
369 
370  std::vector<std::pair<ConstRecHitPointer,ConstRecHitPointer> >::iterator v;
371  for (v = compatPixelHits.begin(); v != compatPixelHits.end(); ++v ) {
372 
373  bool valid = prepareElTrackSeed(v->first,v->second, theVertex);
374  if (valid) {
377  result.push_back(s);
378  }
379 
380  }
381 
382  // Here we take instead the seed from a-priori seed collection
383  } else {
384 
385  reco::ElectronSeed s(theTrackerSeed);
387  result.push_back(s);
388  }
389 
390  }
391 
392 #ifdef FAMOS_DEBUG
393  else
394  std::cout << "[FastElectronSeedGenerator::seedsFromThisCluster] "
395  << " No electron compatible hits found " << std::endl;
396 #endif
397 
398 
399  // And return !
400  return ;
401 
402 }
403 
405  ConstRecHitPointer outerhit,
406  const GlobalPoint& vertexPos)
407 {
408 
409  // debug prints
410  LogDebug("") <<"[FastElectronSeedGenerator::prepareElTrackSeed] "
411  << "inner PixelHit x,y,z "<<innerhit->globalPosition();
412  LogDebug("") <<"[FastElectronSeedGenerator::prepareElTrackSeed] "
413  << "outer PixelHit x,y,z "<<outerhit->globalPosition();
414 
415  // make a spiral from the two hits and the vertex position
416  FastHelix helix(outerhit->globalPosition(),innerhit->globalPosition(),vertexPos,*theSetup);
417  if ( !helix.isValid()) return false;
418 
419  FreeTrajectoryState fts = helix.stateAtVertex();
420 
421  // Give infinite errors to start the fit (no pattern recognition here).
423  fts.setCurvilinearError(errorMatrix*100.);
424 
425  TrajectoryStateOnSurface propagatedState = thePropagator->propagate(fts,innerhit->det()->surface()) ;
426  if (!propagatedState.isValid()) return false;
427 
428  // The persitent trajectory state
429  pts_ = trajectoryStateTransform::persistentState(propagatedState, innerhit->geographicalId().rawId());
430 
431  // The corresponding rechits
432  recHits_.clear();
433  recHits_.push_back(innerhit->hit()->clone());
434  recHits_.push_back(outerhit->hit()->clone());
435 
436 
437  return true;
438 
439 }
440 
#define LogDebug(id)
void setCaloCluster(const CaloClusterRef &, unsigned char hitsMask=0, int subDet2=0, int subDet1=0, float hoe1=std::numeric_limits< float >::infinity(), float hoe2=std::numeric_limits< float >::infinity())
Definition: ElectronSeed.cc:71
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
void setupES(const edm::EventSetup &setup)
size_t size() const
return number of contained object
Definition: RangeMap.h:132
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
const edm::EventSetup * theSetup
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
const MagneticFieldMap * theMagneticFieldMap
range get(ID id, CMP comparator) const
Definition: RangeMap.h:81
const TrackerInteractionGeometry * theTrackerInteractionGeometry
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
PropagationDirection
TrajectorySeedCollection * theInitialSeedColl
unsigned int ringNumber() const
The Ring Number.
Definition: TrackerRecHit.h:84
void set1stLayer(float ephimin, float ephimax, float pphimin, float pphimax)
PropagatorWithMaterial * thePropagator
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
void push_back(D *&d)
Definition: OwnVector.h:273
void set1stLayerZRange(double zmin1, double zmax1)
std::vector< TrajectorySeed > TrajectorySeedCollection
FastElectronSeedGenerator(const edm::ParameterSet &pset, double pTMin, const edm::InputTag &beamSpot)
std::vector< std::pair< ConstRecHitPointer, ConstRecHitPointer > > compatibleHits(const GlobalPoint &xmeas, const GlobalPoint &vprim, float energy, std::vector< TrackerRecHit > &theHits)
recHitContainer::const_iterator const_iterator
T sqrt(T t)
Definition: SSEVec.h:46
void run(edm::Event &e, const reco::SuperClusterRefVector &sclRefs, const SiTrackerGSMatchedRecHit2DCollection *theGSRecHits, const edm::SimTrackContainer *theSimTracks, TrajectorySeedCollection *seeds, reco::ElectronSeedCollection &out)
tuple result
Definition: query.py:137
void clear()
Definition: OwnVector.h:370
std::pair< const_iterator, const_iterator > range
const GeometricSearchTracker * theGeomSearchTracker
std::vector< ElectronSeed > ElectronSeedCollection
collection of ElectronSeed objects
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
FastPixelHitMatcher * myGSPixelMatcher
tuple out
Definition: dbtoconf.py:99
virtual TrajectoryStateOnSurface propagate(const TrajectoryStateOnSurface &tsos, const Plane &plane) const
void setCurvilinearError(const CurvilinearTrajectoryError &err)
void addASeedToThisCluster(edm::Ref< reco::SuperClusterCollection > seedCluster, std::vector< TrackerRecHit > &theHits, const TrajectorySeed &theTrackerSeed, std::vector< reco::ElectronSeed > &result)
bool prepareElTrackSeed(ConstRecHitPointer outerhit, ConstRecHitPointer innerhit, const GlobalPoint &vertexPos)
const T & get() const
Definition: EventSetup.h:55
range recHits() const
unsigned int layerNumber() const
The Layer Number.
Definition: TrackerRecHit.h:81
ESHandle< TrackerGeometry > geometry
edm::EventID id() const
Definition: EventBase.h:56
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
unsigned int subDetId() const
The subdet Id.
Definition: TrackerRecHit.h:78
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
void setES(const MagneticFieldMap *aFieldMap, const TrackerGeometry *aTrackerGeometry, const GeometricSearchTracker *geomSearchTracker, const TrackerInteractionGeometry *interactionGeometry)
tuple cout
Definition: gather_cfg.py:121
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< ID > ids() const
indentifier vector
Definition: RangeMap.h:179
std::vector< SimTrack > SimTrackContainer
void set2ndLayer(float phimin, float phimax)
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
tuple size
Write out results.
mathSSE::Vec4< T > v
const TrackerGeometry * theTrackerGeometry