CMS 3D CMS Logo

TrackerSimHitProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
3 #include <memory>
4 #include <vector>
5 
6 // framework
11 
12 // tracking
18 
19 // fastsim
26 
27 // data formats
35 
36 // other
39 
41 // Author: L. Vanelderen, S. Kurz
42 // Date: 29 May 2017
44 
45 typedef std::pair<const GeomDet*, TrajectoryStateOnSurface> DetWithState;
46 
47 namespace fastsim {
49 
55  public:
58 
60  ~TrackerSimHitProducer() override { ; }
61 
63 
69  void interact(Particle& particle,
71  std::vector<std::unique_ptr<Particle>>& secondaries,
72  const RandomEngineAndDistribution& random) override;
73 
75  void registerProducts(edm::ProducesCollector) const override;
76 
78  void storeProducts(edm::Event& iEvent) override;
79 
81 
91  std::pair<double, std::unique_ptr<PSimHit>> createHitOnDetector(const TrajectoryStateOnSurface& particle,
92  int pdgId,
93  double layerThickness,
94  double eLoss,
95  int simTrackId,
96  const GeomDet& detector,
97  GlobalPoint& refPos);
98 
99  private:
100  const double
102  std::unique_ptr<edm::PSimHitContainer> simHitContainer_;
103  double minMomentum_;
105  };
106 } // namespace fastsim
107 
109  : fastsim::InteractionModel(name), onSurfaceTolerance_(0.01), simHitContainer_(new edm::PSimHitContainer) {
110  // Set the minimal momentum
111  minMomentum_ = cfg.getParameter<double>("minMomentumCut");
112  // - if not set, particles from outside the beampipe with a negative speed in R direction are propagated but no SimHits
113  // - particle with positive (negative) z and negative (positive) speed in z direction: no SimHits
114  // -> this is not neccesary since a track reconstruction is not possible in this case anyways
115  doHitsFromInboundParticles_ = cfg.getParameter<bool>("doHitsFromInboundParticles");
116 }
117 
119  producesCollector.produces<edm::PSimHitContainer>("TrackerHits");
120 }
121 
123  iEvent.put(std::move(simHitContainer_), "TrackerHits");
124  simHitContainer_ = std::make_unique<edm::PSimHitContainer>();
125 }
126 
128  const SimplifiedGeometry& layer,
129  std::vector<std::unique_ptr<Particle>>& secondaries,
130  const RandomEngineAndDistribution& random) {
131  // the energy deposit in the layer
132  double energyDeposit = particle.getEnergyDeposit();
133  particle.setEnergyDeposit(0);
134 
135  //
136  // don't save hits from particle that did a loop or are inbound (coming from the outer part of the tracker, going towards the center)
137  //
138  if (!doHitsFromInboundParticles_) {
139  if (particle.isLooper()) {
140  return;
141  }
142  if (particle.momentum().X() * particle.position().X() + particle.momentum().Y() * particle.position().Y() < 0) {
143  particle.setLooper();
144  return;
145  }
146  if (layer.isForward() && ((particle.position().Z() > 0 && particle.momentum().Z() < 0) ||
147  (particle.position().Z() < 0 && particle.momentum().Z() > 0))) {
148  particle.setLooper();
149  return;
150  }
151  }
152 
153  //
154  // check that layer has tracker modules
155  //
156  if (!layer.getDetLayer()) {
157  return;
158  }
159 
160  //
161  // no material
162  //
163  if (layer.getThickness(particle.position(), particle.momentum()) < 1E-10) {
164  return;
165  }
166 
167  //
168  // only charged particles
169  //
170  if (particle.charge() == 0) {
171  return;
172  }
173 
174  //
175  // save hit only if momentum higher than threshold
176  //
177  if (particle.momentum().Perp2() < minMomentum_ * minMomentum_) {
178  return;
179  }
180 
181  // create the trajectory of the particle
182  UniformMagneticField magneticField(layer.getMagneticFieldZ(particle.position()));
183  GlobalPoint position(particle.position().X(), particle.position().Y(), particle.position().Z());
184  GlobalVector momentum(particle.momentum().Px(), particle.momentum().Py(), particle.momentum().Pz());
185  auto plane = layer.getDetLayer()->surface().tangentPlane(position);
186  TrajectoryStateOnSurface trajectory(
187  GlobalTrajectoryParameters(position, momentum, TrackCharge(particle.charge()), &magneticField), *plane);
188 
189  // find detectors compatible with the particle's trajectory
192  std::vector<DetWithState> compatibleDetectors = layer.getDetLayer()->compatibleDets(trajectory, propagator, est);
193 
195  // You have to sort the simHits in the order they occur!
197 
198  // The old algorithm (sorting by distance to IP) doesn't seem to make sense to me (what if particle moves inwards?)
199 
200  // Detector layers have to be sorted by proximity to particle.position
201  // Propagate particle backwards a bit to make sure it's outside any components
202  std::map<double, std::unique_ptr<PSimHit>> distAndHits;
203  // Position relative to which the hits should be sorted
204  GlobalPoint positionOutside(particle.position().x() - particle.momentum().x() / particle.momentum().P() * 10.,
205  particle.position().y() - particle.momentum().y() / particle.momentum().P() * 10.,
206  particle.position().z() - particle.momentum().z() / particle.momentum().P() * 10.);
207 
208  // FastSim: cheat tracking -> assign hits to closest charged daughter if particle decays
209  int pdgId = particle.pdgId();
210  int simTrackId =
211  particle.getMotherSimTrackIndex() >= 0 ? particle.getMotherSimTrackIndex() : particle.simTrackIndex();
212 
213  // loop over the compatible detectors
214  for (const auto& detectorWithState : compatibleDetectors) {
215  const GeomDet& detector = *detectorWithState.first;
216  const TrajectoryStateOnSurface& particleState = detectorWithState.second;
217 
218  // if the detector has no components
219  if (detector.isLeaf()) {
220  std::pair<double, std::unique_ptr<PSimHit>> hitPair = createHitOnDetector(particleState,
221  pdgId,
222  layer.getThickness(particle.position()),
223  energyDeposit,
224  simTrackId,
225  detector,
226  positionOutside);
227  if (hitPair.second) {
228  distAndHits.insert(distAndHits.end(), std::move(hitPair));
229  }
230  } else {
231  // if the detector has components
232  for (const auto component : detector.components()) {
233  std::pair<double, std::unique_ptr<PSimHit>> hitPair =
234  createHitOnDetector(particleState,
235  pdgId,
236  layer.getThickness(particle.position()),
237  energyDeposit,
238  simTrackId,
239  *component,
240  positionOutside);
241  if (hitPair.second) {
242  distAndHits.insert(distAndHits.end(), std::move(hitPair));
243  }
244  }
245  }
246  }
247 
248  // Fill simHitContainer
249  for (std::map<double, std::unique_ptr<PSimHit>>::const_iterator it = distAndHits.begin(); it != distAndHits.end();
250  it++) {
251  simHitContainer_->push_back(*(it->second));
252  }
253 }
254 
255 // Also returns distance to simHit since hits have to be ordered (in time) afterwards. Necessary to do explicit copy of TrajectoryStateOnSurface particle (not call by reference)
256 std::pair<double, std::unique_ptr<PSimHit>> fastsim::TrackerSimHitProducer::createHitOnDetector(
257  const TrajectoryStateOnSurface& particle,
258  int pdgId,
259  double layerThickness,
260  double eLoss,
261  int simTrackId,
262  const GeomDet& detector,
263  GlobalPoint& refPos) {
264  // determine position and momentum of particle in the coordinate system of the detector
265  LocalPoint localPosition;
266  LocalVector localMomentum;
267 
268  // if the particle is close enough, no further propagation is needed
269  if (fabs(detector.toLocal(particle.globalPosition()).z()) < onSurfaceTolerance_) {
270  localPosition = particle.localPosition();
271  localMomentum = particle.localMomentum();
272  }
273  // else, propagate
274  else {
275  // find crossing of particle with detector layer
277  particle.globalMomentum().basicVector(),
278  particle.transverseCurvature(),
279  anyDirection);
280  std::pair<bool, double> path = crossing.pathLength(detector.surface());
281  // case propagation succeeds
282  if (path.first) {
283  localPosition = detector.toLocal(GlobalPoint(crossing.position(path.second)));
284  localMomentum = detector.toLocal(GlobalVector(crossing.direction(path.second)));
285  localMomentum = localMomentum.unit() * particle.localMomentum().mag();
286  }
287  // case propagation fails
288  else {
289  return std::pair<double, std::unique_ptr<PSimHit>>(0, std::unique_ptr<PSimHit>(nullptr));
290  }
291  }
292 
293  // find entry and exit point of particle in detector
294  const Plane& detectorPlane = detector.surface();
295  double halfThick = 0.5 * detectorPlane.bounds().thickness();
296  double pZ = localMomentum.z();
297  LocalPoint entry = localPosition + (-halfThick / pZ) * localMomentum;
298  LocalPoint exit = localPosition + (halfThick / pZ) * localMomentum;
299  double tof = particle.globalPosition().mag() / fastsim::Constants::speedOfLight; // in nanoseconds
300 
301  // make sure the simhit is physically on the module
302  double boundX = detectorPlane.bounds().width() / 2.;
303  double boundY = detectorPlane.bounds().length() / 2.;
304  // Special treatment for TID and TEC trapeziodal modules
305  unsigned subdet = DetId(detector.geographicalId()).subdetId();
306  if (subdet == 4 || subdet == 6)
307  boundX *= 1. - localPosition.y() / detectorPlane.position().perp();
308  if (fabs(localPosition.x()) > boundX || fabs(localPosition.y()) > boundY) {
309  return std::pair<double, std::unique_ptr<PSimHit>>(0, std::unique_ptr<PSimHit>(nullptr));
310  }
311 
312  // Create the hit: the energy loss rescaled to the module thickness
313  // Total thickness is in radiation lengths, 1 radlen = 9.36 cm
314  // Sensitive module thickness is about 30 microns larger than
315  // the module thickness itself
316  eLoss *= (2. * halfThick - 0.003) / (9.36 * layerThickness);
317 
318  // Position of the hit in global coordinates
319  GlobalPoint hitPos(detector.surface().toGlobal(localPosition));
320 
321  return std::pair<double, std::unique_ptr<PSimHit>>((hitPos - refPos).mag(),
322  std::make_unique<PSimHit>(entry,
323  exit,
324  localMomentum.mag(),
325  tof,
326  eLoss,
327  pdgId,
328  detector.geographicalId().rawId(),
329  simTrackId,
330  localMomentum.theta(),
331  localMomentum.phi()));
332 }
333 
Vector3DBase
Definition: Vector3DBase.h:8
fastsim::TrackerSimHitProducer::onSurfaceTolerance_
const double onSurfaceTolerance_
Max distance between particle and active (sub-) module. Otherwise particle has to be propagated.
Definition: TrackerSimHitProducer.cc:101
TrajectoryStateOnSurface.h
fastsim::Constants::speedOfLight
static constexpr double speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:12
TrajectoryStateOnSurface::localMomentum
LocalVector localMomentum() const
Definition: TrajectoryStateOnSurface.h:75
fastsim::Particle::momentum
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:143
anyDirection
Definition: PropagationDirection.h:4
Bounds::width
virtual float width() const =0
TrackCharge
int TrackCharge
Definition: TrackCharge.h:4
GeomDet
Definition: GeomDet.h:27
InteractionModel.h
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
edm
HLT enums.
Definition: AlignableModifier.h:19
mps_splice.entry
entry
Definition: mps_splice.py:68
PV3DBase::theta
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
fastsim::TrackerSimHitProducer::registerProducts
void registerProducts(edm::ProducesCollector) const override
Register the SimHit collection.
Definition: TrackerSimHitProducer.cc:118
TrajectoryStateOnSurface::globalPosition
GlobalPoint globalPosition() const
Definition: TrajectoryStateOnSurface.h:65
PSimHitContainer.h
fastsim::TrackerSimHitProducer::createHitOnDetector
std::pair< double, std::unique_ptr< PSimHit > > createHitOnDetector(const TrajectoryStateOnSurface &particle, int pdgId, double layerThickness, double eLoss, int simTrackId, const GeomDet &detector, GlobalPoint &refPos)
Helper funtion to create the actual SimHit on a detector (sub-) module.
Definition: TrackerSimHitProducer.cc:256
fastsim::SimplifiedGeometry
Implementation of a generic detector layer (base class for forward/barrel layers).
Definition: SimplifiedGeometry.h:35
fastsim::Particle::getMotherSimTrackIndex
int getMotherSimTrackIndex() const
Get index of simTrack of mother particle (only makes sense if mother and daughter charged).
Definition: Particle.h:216
fastsim::TrackerSimHitProducer::TrackerSimHitProducer
TrackerSimHitProducer(const std::string &name, const edm::ParameterSet &cfg)
Constructor.
Definition: TrackerSimHitProducer.cc:108
TrajectoryStateOnSurface::transverseCurvature
double transverseCurvature() const
Definition: TrajectoryStateOnSurface.h:70
fastsim::Particle::simTrackIndex
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:153
GlobalVector
Global3DVector GlobalVector
Definition: GlobalVector.h:10
HLT_FULL_cff.magneticField
magneticField
Definition: HLT_FULL_cff.py:348
fastsim::Particle::setEnergyDeposit
void setEnergyDeposit(double energyDeposit)
Set the energy the particle deposited in the tracker layer that was last hit (ionization).
Definition: Particle.h:87
simKBmtfDigis_cfi.eLoss
eLoss
Definition: simKBmtfDigis_cfi.py:9
InteractionModelFactory.h
ProducesCollector.h
Particle.h
Bounds::length
virtual float length() const =0
fastsim::InteractionModel
Base class for any interaction model between a particle and a tracker layer.
Definition: InteractionModel.h:29
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
DetId
Definition: DetId.h:17
Plane.h
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
fastsim::Particle::getEnergyDeposit
double getEnergyDeposit() const
Return the energy the particle deposited in the tracker layer that was last hit (ionization).
Definition: Particle.h:188
InsideBoundsMeasurementEstimator.h
fastsim::Particle::setLooper
void setLooper()
This particle is about to do a loop in the tracker or the direction of the momentum is going inwards.
Definition: Particle.h:94
Vector3DBase::unit
Vector3DBase unit() const
Definition: Vector3DBase.h:54
PSimHit.h
TrackCandidateProducer_cfi.propagator
propagator
Definition: TrackCandidateProducer_cfi.py:17
fastsim::TrackerSimHitProducer::simHitContainer_
std::unique_ptr< edm::PSimHitContainer > simHitContainer_
The SimHit.
Definition: TrackerSimHitProducer.cc:102
Surface::bounds
const Bounds & bounds() const
Definition: Surface.h:87
HelixArbitraryPlaneCrossing
Definition: HelixArbitraryPlaneCrossing.h:10
GlobalTrajectoryParameters
Definition: GlobalTrajectoryParameters.h:15
GlobalPoint
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Point3DBase< float, GlobalTag >
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
fastsim::Particle::position
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:140
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
DetWithState
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
Definition: TrackerSimHitProducer.cc:45
fastsim::TrackerSimHitProducer::~TrackerSimHitProducer
~TrackerSimHitProducer() override
Default destructor.
Definition: TrackerSimHitProducer.cc:60
Bounds::thickness
virtual float thickness() const =0
TrajectoryStateOnSurface::localPosition
LocalPoint localPosition() const
Definition: TrajectoryStateOnSurface.h:74
HelixArbitraryPlaneCrossing.h
edm::ParameterSet
Definition: ParameterSet.h:47
Constants.h
Event.h
edmplugin::PluginFactory
Definition: PluginFactory.h:34
fastsim::Particle::pdgId
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:134
fastsim::TrackerSimHitProducer::storeProducts
void storeProducts(edm::Event &iEvent) override
Store the SimHit collection.
Definition: TrackerSimHitProducer.cc:122
fastsim::Particle::isLooper
double isLooper() const
Check if this particle is about to do a loop in the tracker or the direction of the momentum is going...
Definition: Particle.h:195
DetID.h
fastsim::TrackerSimHitProducer::interact
void interact(Particle &particle, const SimplifiedGeometry &layer, std::vector< std::unique_ptr< Particle >> &secondaries, const RandomEngineAndDistribution &random) override
Perform the interaction.
Definition: TrackerSimHitProducer.cc:127
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
iEvent
int iEvent
Definition: GenABIO.cc:224
fastsim::Particle::charge
double charge() const
Return charge of the particle.
Definition: Particle.h:137
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:117
edm::ProducesCollector::produces
ProductRegistryHelper::BranchAliasSetterT< ProductType > produces()
Definition: ProducesCollector.h:52
UniformMagneticField.h
AnalyticalPropagator.h
fastsim::TrackerSimHitProducer::doHitsFromInboundParticles_
bool doHitsFromInboundParticles_
If not set, incoming particles (negative speed relative to center of detector) don't create a SimHits...
Definition: TrackerSimHitProducer.cc:104
SimplifiedGeometry.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PV3DBase::basicVector
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
looper.cfg
cfg
Definition: looper.py:296
PV3DBase::mag
T mag() const
Definition: PV3DBase.h:64
InsideBoundsMeasurementEstimator
Definition: InsideBoundsMeasurementEstimator.h:6
GloballyPositioned::position
const PositionType & position() const
Definition: GloballyPositioned.h:36
AnalyticalPropagator
Definition: AnalyticalPropagator.h:22
TrajectoryStateOnSurface::globalMomentum
GlobalVector globalMomentum() const
Definition: TrajectoryStateOnSurface.h:66
mag
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Definition: Basic3DVectorLD.h:127
GeomDet.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
edm::ProducesCollector
Definition: ProducesCollector.h:43
fastsim::TrackerSimHitProducer
Produces SimHits in the tracker layers.
Definition: TrackerSimHitProducer.cc:54
LocalPoint.h
Plane
Definition: Plane.h:16
DetLayer.h
GlobalVector.h
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
Particle
Definition: Particle.py:1
hgcalTestNeighbor_cfi.detector
detector
Definition: hgcalTestNeighbor_cfi.py:6
castor_dqm_sourceclient_file_cfg.path
path
Definition: castor_dqm_sourceclient_file_cfg.py:37
fastsim
Definition: BarrelSimplifiedGeometry.h:15
genParticles_cff.map
map
Definition: genParticles_cff.py:11
edm::PSimHitContainer
std::vector< PSimHit > PSimHitContainer
Definition: PSimHitContainer.h:11
beamvalidation.exit
def exit(msg="")
Definition: beamvalidation.py:52
ParameterSet.h
edm::Event
Definition: Event.h:73
LocalVector.h
GlobalPoint.h
fastsim::TrackerSimHitProducer::minMomentum_
double minMomentum_
Set the minimal momentum of incoming particle.
Definition: TrackerSimHitProducer.cc:103
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
RandomEngineAndDistribution
Definition: RandomEngineAndDistribution.h:18
UniformMagneticField
Definition: UniformMagneticField.h:13