CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
ElectronSeedGenerator Class Reference

#include <ElectronSeedGenerator.h>

Classes

struct  Tokens
 

Public Types

typedef TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
 
typedef edm::OwnVector< TrackingRecHitPRecHitContainer
 
typedef TransientTrackingRecHit::RecHitContainer RecHitContainer
 
typedef TransientTrackingRecHit::RecHitPointer RecHitPointer
 

Public Member Functions

 ElectronSeedGenerator (const edm::ParameterSet &, const Tokens &, edm::ConsumesCollector &&)
 
void run (edm::Event &, const reco::SuperClusterRefVector &, const std::vector< const TrajectorySeedCollection *> &seedsV, reco::ElectronSeedCollection &)
 
void setupES (const edm::EventSetup &setup)
 

Private Member Functions

void seedsFromThisCluster (edm::Ref< reco::SuperClusterCollection > seedCluster, reco::BeamSpot const &beamSpot, std::vector< reco::Vertex > const *vertices, reco::ElectronSeedCollection &out)
 

Private Attributes

const edm::EDGetTokenT< reco::BeamSpotbeamSpotTag_
 
const float deltaPhi1High_
 
const float deltaPhi1Low_
 
const float deltaPhi2B_
 
const float deltaPhi2F_
 
const float deltaZ1WithVertex_
 
const double dPhi1Coef1_
 
const double dPhi1Coef2_
 
const bool dynamicPhiRoad_
 
const float highPtThresh_
 
const std::vector< const TrajectorySeedCollection * > * initialSeedCollectionVector_ = nullptr
 
const float lowPtThresh_
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagFieldToken_
 
edm::ESWatcher< IdealMagneticFieldRecordmagneticFieldWatcher_
 
PixelHitMatcher matcher_
 
const float nSigmasDeltaZ1_
 
const float phiMax2B_
 
const float phiMax2F_
 
const float phiMin2B_
 
const float phiMin2F_
 
const float sizeWindowENeg_
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtrackerGeometryToken_
 
edm::ESWatcher< TrackerDigiGeometryRecordtrackerGeometryWatcher_
 
const bool useRecoVertex_
 
const edm::EDGetTokenT< std::vector< reco::Vertex > > verticesTag_
 

Detailed Description

Class to generate the trajectory seed from two hits in the pixel detector which have been found compatible with an ECAL cluster.

Author
U.Berthon, C.Charlot, LLR Palaiseau
Version
1st Version May 30, 2006

Description: Top algorithm producing ElectronSeeds, ported from ORCA

Implementation: future redesign...

Definition at line 46 of file ElectronSeedGenerator.h.

Member Typedef Documentation

◆ ConstRecHitPointer

Definition at line 54 of file ElectronSeedGenerator.h.

◆ PRecHitContainer

Definition at line 53 of file ElectronSeedGenerator.h.

◆ RecHitContainer

Definition at line 56 of file ElectronSeedGenerator.h.

◆ RecHitPointer

Definition at line 55 of file ElectronSeedGenerator.h.

Constructor & Destructor Documentation

◆ ElectronSeedGenerator()

ElectronSeedGenerator::ElectronSeedGenerator ( const edm::ParameterSet pset,
const Tokens ts,
edm::ConsumesCollector &&  cc 
)

Definition at line 128 of file ElectronSeedGenerator.cc.

References gpuPixelDoublets::cc.

131  : dynamicPhiRoad_(pset.getParameter<bool>("dynamicPhiRoad")),
132  verticesTag_(ts.token_vtx),
133  beamSpotTag_(ts.token_bs),
134  magFieldToken_{cc.esConsumes()},
135  trackerGeometryToken_{cc.esConsumes()},
136  lowPtThresh_(pset.getParameter<double>("LowPtThreshold")),
137  highPtThresh_(pset.getParameter<double>("HighPtThreshold")),
138  nSigmasDeltaZ1_(pset.getParameter<double>("nSigmasDeltaZ1")),
139  deltaZ1WithVertex_(pset.getParameter<double>("deltaZ1WithVertex")),
140  sizeWindowENeg_(pset.getParameter<double>("SizeWindowENeg")),
141  deltaPhi1Low_(pset.getParameter<double>("DeltaPhi1Low")),
142  deltaPhi1High_(pset.getParameter<double>("DeltaPhi1High")),
143  // so that deltaPhi1 = dPhi1Coef1_ + dPhi1Coef2_/clusterEnergyT
146  // use of reco vertex
147  useRecoVertex_(pset.getParameter<bool>("useRecoVertex")),
148  // new B/F configurables
149  deltaPhi2B_(pset.getParameter<double>("DeltaPhi2B")),
150  deltaPhi2F_(pset.getParameter<double>("DeltaPhi2F")),
151  phiMin2B_(-pset.getParameter<double>("PhiMax2B")),
152  phiMin2F_(-pset.getParameter<double>("PhiMax2F")),
153  phiMax2B_(pset.getParameter<double>("PhiMax2B")),
154  phiMax2F_(pset.getParameter<double>("PhiMax2F")),
155  matcher_(pset.getParameter<double>("ePhiMin1"),
156  pset.getParameter<double>("ePhiMax1"),
157  phiMin2B_,
158  phiMax2B_,
159  phiMin2F_,
160  phiMax2F_,
161  pset.getParameter<double>("z2MaxB"),
162  pset.getParameter<double>("r2MaxF"),
163  pset.getParameter<double>("rMaxI"),
164  useRecoVertex_) {}
const edm::EDGetTokenT< std::vector< reco::Vertex > > verticesTag_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeometryToken_
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotTag_

Member Function Documentation

◆ run()

void ElectronSeedGenerator::run ( edm::Event e,
const reco::SuperClusterRefVector sclRefs,
const std::vector< const TrajectorySeedCollection *> &  seedsV,
reco::ElectronSeedCollection out 
)

Definition at line 174 of file ElectronSeedGenerator.cc.

References pwdgSkimBPark_cfi::beamSpot, beamSpotTag_, MillePedeFileConverter_cfg::e, mps_fire::i, initialSeedCollectionVector_, LogDebug, MillePedeFileConverter_cfg::out, seedsFromThisCluster(), edm::RefVector< C, T, F >::size(), useRecoVertex_, AlignmentTracksFromVertexSelector_cfi::vertices, and verticesTag_.

177  {
179 
180  // get the beamspot from the Event:
181  auto const &beamSpot = e.get(beamSpotTag_);
182 
183  // if required get the vertices
184  std::vector<reco::Vertex> const *vertices = nullptr;
185  if (useRecoVertex_)
186  vertices = &e.get(verticesTag_);
187 
188  for (unsigned int i = 0; i < sclRefs.size(); ++i) {
189  // Find the seeds
190  LogDebug("ElectronSeedGenerator") << "new cluster, calling seedsFromThisCluster";
192  }
193 
194  LogDebug("ElectronSeedGenerator") << ": For event " << e.id();
195  LogDebug("ElectronSeedGenerator") << "Nr of superclusters after filter: " << sclRefs.size()
196  << ", no. of ElectronSeeds found = " << out.size();
197 }
void seedsFromThisCluster(edm::Ref< reco::SuperClusterCollection > seedCluster, reco::BeamSpot const &beamSpot, std::vector< reco::Vertex > const *vertices, reco::ElectronSeedCollection &out)
const edm::EDGetTokenT< std::vector< reco::Vertex > > verticesTag_
const std::vector< const TrajectorySeedCollection * > * initialSeedCollectionVector_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotTag_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
#define LogDebug(id)

◆ seedsFromThisCluster()

void ElectronSeedGenerator::seedsFromThisCluster ( edm::Ref< reco::SuperClusterCollection seedCluster,
reco::BeamSpot const &  beamSpot,
std::vector< reco::Vertex > const *  vertices,
reco::ElectronSeedCollection out 
)
private

Definition at line 199 of file ElectronSeedGenerator.cc.

References pwdgSkimBPark_cfi::beamSpot, deltaPhi1High_, deltaPhi1Low_, deltaPhi2B_, deltaPhi2F_, deltaZ1WithVertex_, dPhi1Coef1_, dPhi1Coef2_, dynamicPhiRoad_, ele_convert(), PVValHelper::eta, highPtThresh_, initialSeedCollectionVector_, lowPtThresh_, matcher_, nSigmasDeltaZ1_, MillePedeFileConverter_cfg::out, PixelHitMatcher::set1stLayer(), PixelHitMatcher::set1stLayerZRange(), PixelHitMatcher::set2ndLayer(), beamSpotPI::sigmaZ, sizeWindowENeg_, mathSSE::sqrt(), useRecoVertex_, bphysicsOniaDQM_cfi::vertex, and AlignmentTracksFromVertexSelector_cfi::vertices.

Referenced by run().

202  {
203  float clusterEnergy = seedCluster->energy();
204  GlobalPoint clusterPos(seedCluster->position().x(), seedCluster->position().y(), seedCluster->position().z());
205  reco::ElectronSeed::CaloClusterRef caloCluster(seedCluster);
206 
207  if (dynamicPhiRoad_) {
208  float clusterEnergyT = clusterEnergy / cosh(EleRelPoint(clusterPos, beamSpot.position()).eta());
209 
210  float deltaPhi1;
211  if (clusterEnergyT < lowPtThresh_) {
212  deltaPhi1 = deltaPhi1Low_;
213  } else if (clusterEnergyT > highPtThresh_) {
214  deltaPhi1 = deltaPhi1High_;
215  } else {
216  deltaPhi1 = dPhi1Coef1_ + dPhi1Coef2_ / clusterEnergyT;
217  }
218 
219  matcher_.set1stLayer(-deltaPhi1 * sizeWindowENeg_, deltaPhi1 * (1. - sizeWindowENeg_));
221  }
222 
223  if (!useRecoVertex_) // here use the beam spot position
224  {
225  double sigmaZ = beamSpot.sigmaZ();
226  double sigmaZ0Error = beamSpot.sigmaZ0Error();
227  double sq = sqrt(sigmaZ * sigmaZ + sigmaZ0Error * sigmaZ0Error);
228  double myZmin1 = beamSpot.position().z() - nSigmasDeltaZ1_ * sq;
229  double myZmax1 = beamSpot.position().z() + nSigmasDeltaZ1_ * sq;
230 
231  GlobalPoint vertexPos;
232  ele_convert(beamSpot.position(), vertexPos);
233 
234  matcher_.set1stLayerZRange(myZmin1, myZmax1);
235 
236  // try electron
237  auto elePixelSeeds = matcher_(*initialSeedCollectionVector_, clusterPos, vertexPos, clusterEnergy, -1.);
238  seedsFromTrajectorySeeds(elePixelSeeds, caloCluster, out, false);
239  // try positron
240  auto posPixelSeeds = matcher_(*initialSeedCollectionVector_, clusterPos, vertexPos, clusterEnergy, 1.);
241  seedsFromTrajectorySeeds(posPixelSeeds, caloCluster, out, true);
242 
243  } else if (vertices) // here we use the reco vertices
244  {
245  for (auto const &vertex : *vertices) {
246  GlobalPoint vertexPos(vertex.position().x(), vertex.position().y(), vertex.position().z());
247  double myZmin1, myZmax1;
248  if (vertexPos.z() == beamSpot.position().z()) { // in case vetex not found
249  double sigmaZ = beamSpot.sigmaZ();
250  double sigmaZ0Error = beamSpot.sigmaZ0Error();
251  double sq = sqrt(sigmaZ * sigmaZ + sigmaZ0Error * sigmaZ0Error);
252  myZmin1 = beamSpot.position().z() - nSigmasDeltaZ1_ * sq;
253  myZmax1 = beamSpot.position().z() + nSigmasDeltaZ1_ * sq;
254  } else { // a vertex has been recoed
255  myZmin1 = vertex.position().z() - deltaZ1WithVertex_;
256  myZmax1 = vertex.position().z() + deltaZ1WithVertex_;
257  }
258 
259  matcher_.set1stLayerZRange(myZmin1, myZmax1);
260 
261  // try electron
262  auto elePixelSeeds = matcher_(*initialSeedCollectionVector_, clusterPos, vertexPos, clusterEnergy, -1.);
263  seedsFromTrajectorySeeds(elePixelSeeds, caloCluster, out, false);
264  // try positron
265  auto posPixelSeeds = matcher_(*initialSeedCollectionVector_, clusterPos, vertexPos, clusterEnergy, 1.);
266  seedsFromTrajectorySeeds(posPixelSeeds, caloCluster, out, true);
267  }
268  }
269 }
T z() const
Definition: PV3DBase.h:61
void set2ndLayer(float dummyphi2minB, float dummyphi2maxB, float dummyphi2minF, float dummyphi2maxF)
const std::vector< const TrajectorySeedCollection * > * initialSeedCollectionVector_
void set1stLayerZRange(float zmin1, float zmax1)
T sqrt(T t)
Definition: SSEVec.h:19
void set1stLayer(float dummyphi1min, float dummyphi1max)
void ele_convert(const Type1 &obj1, Type2 &obj2)

◆ setupES()

void ElectronSeedGenerator::setupES ( const edm::EventSetup setup)

Definition at line 166 of file ElectronSeedGenerator.cc.

References edm::ESWatcher< T >::check(), magFieldToken_, magneticFieldWatcher_, matcher_, PixelHitMatcher::setES(), singleTopDQM_cfi::setup, trackerGeometryToken_, and trackerGeometryWatcher_.

166  {
167  auto newMagField = magneticFieldWatcher_.check(setup);
168  auto newTrackerGeom = trackerGeometryWatcher_.check(setup);
169  if (newMagField || newTrackerGeom) {
171  }
172 }
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerGeometryToken_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
edm::ESWatcher< TrackerDigiGeometryRecord > trackerGeometryWatcher_
edm::ESWatcher< IdealMagneticFieldRecord > magneticFieldWatcher_
void setES(MagneticField const &, TrackerGeometry const &trackerGeometry)

Member Data Documentation

◆ beamSpotTag_

const edm::EDGetTokenT<reco::BeamSpot> ElectronSeedGenerator::beamSpotTag_
private

Definition at line 75 of file ElectronSeedGenerator.h.

Referenced by run().

◆ deltaPhi1High_

const float ElectronSeedGenerator::deltaPhi1High_
private

Definition at line 88 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ deltaPhi1Low_

const float ElectronSeedGenerator::deltaPhi1Low_
private

Definition at line 87 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ deltaPhi2B_

const float ElectronSeedGenerator::deltaPhi2B_
private

Definition at line 98 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ deltaPhi2F_

const float ElectronSeedGenerator::deltaPhi2F_
private

Definition at line 99 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ deltaZ1WithVertex_

const float ElectronSeedGenerator::deltaZ1WithVertex_
private

Definition at line 84 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ dPhi1Coef1_

const double ElectronSeedGenerator::dPhi1Coef1_
private

Definition at line 92 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ dPhi1Coef2_

const double ElectronSeedGenerator::dPhi1Coef2_
private

Definition at line 91 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ dynamicPhiRoad_

const bool ElectronSeedGenerator::dynamicPhiRoad_
private

Definition at line 72 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ highPtThresh_

const float ElectronSeedGenerator::highPtThresh_
private

Definition at line 82 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ initialSeedCollectionVector_

const std::vector<const TrajectorySeedCollection*>* ElectronSeedGenerator::initialSeedCollectionVector_ = nullptr
private

Definition at line 94 of file ElectronSeedGenerator.h.

Referenced by run(), and seedsFromThisCluster().

◆ lowPtThresh_

const float ElectronSeedGenerator::lowPtThresh_
private

Definition at line 81 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ magFieldToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> ElectronSeedGenerator::magFieldToken_
private

Definition at line 76 of file ElectronSeedGenerator.h.

Referenced by setupES().

◆ magneticFieldWatcher_

edm::ESWatcher<IdealMagneticFieldRecord> ElectronSeedGenerator::magneticFieldWatcher_
private

Definition at line 78 of file ElectronSeedGenerator.h.

Referenced by setupES().

◆ matcher_

PixelHitMatcher ElectronSeedGenerator::matcher_
private

Definition at line 106 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster(), and setupES().

◆ nSigmasDeltaZ1_

const float ElectronSeedGenerator::nSigmasDeltaZ1_
private

Definition at line 83 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ phiMax2B_

const float ElectronSeedGenerator::phiMax2B_
private

Definition at line 103 of file ElectronSeedGenerator.h.

◆ phiMax2F_

const float ElectronSeedGenerator::phiMax2F_
private

Definition at line 104 of file ElectronSeedGenerator.h.

◆ phiMin2B_

const float ElectronSeedGenerator::phiMin2B_
private

Definition at line 101 of file ElectronSeedGenerator.h.

◆ phiMin2F_

const float ElectronSeedGenerator::phiMin2F_
private

Definition at line 102 of file ElectronSeedGenerator.h.

◆ sizeWindowENeg_

const float ElectronSeedGenerator::sizeWindowENeg_
private

Definition at line 85 of file ElectronSeedGenerator.h.

Referenced by seedsFromThisCluster().

◆ trackerGeometryToken_

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> ElectronSeedGenerator::trackerGeometryToken_
private

Definition at line 77 of file ElectronSeedGenerator.h.

Referenced by setupES().

◆ trackerGeometryWatcher_

edm::ESWatcher<TrackerDigiGeometryRecord> ElectronSeedGenerator::trackerGeometryWatcher_
private

Definition at line 79 of file ElectronSeedGenerator.h.

Referenced by setupES().

◆ useRecoVertex_

const bool ElectronSeedGenerator::useRecoVertex_
private

Definition at line 96 of file ElectronSeedGenerator.h.

Referenced by run(), and seedsFromThisCluster().

◆ verticesTag_

const edm::EDGetTokenT<std::vector<reco::Vertex> > ElectronSeedGenerator::verticesTag_
private

Definition at line 74 of file ElectronSeedGenerator.h.

Referenced by run().