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