CMS 3D CMS Logo

TrackerSimHitProducer.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <memory>
3 
4 // framework
9 
10 // tracking
16 
17 // fastsim
24 
25 // data formats
33 
34 // other
37 
39 // Author: L. Vanelderen, S. Kurz
40 // Date: 29 May 2017
42 
43 typedef std::pair<const GeomDet*, TrajectoryStateOnSurface> DetWithState;
44 
45 namespace fastsim {
47 
53  public:
56 
58  ~TrackerSimHitProducer() override { ; }
59 
61 
67  void interact(Particle& particle,
68  const SimplifiedGeometry& layer,
69  std::vector<std::unique_ptr<Particle>>& secondaries,
70  const RandomEngineAndDistribution& random) override;
71 
73  void registerProducts(edm::ProducesCollector) const override;
74 
76  void storeProducts(edm::Event& iEvent) override;
77 
79 
89  std::pair<double, std::unique_ptr<PSimHit>> createHitOnDetector(const TrajectoryStateOnSurface& particle,
90  int pdgId,
91  double layerThickness,
92  double eLoss,
93  int simTrackId,
94  const GeomDet& detector,
95  GlobalPoint& refPos);
96 
97  private:
98  const double
100  std::unique_ptr<edm::PSimHitContainer> simHitContainer_;
101  double minMomentum_;
103  };
104 } // namespace fastsim
105 
107  : fastsim::InteractionModel(name), onSurfaceTolerance_(0.01), simHitContainer_(new edm::PSimHitContainer) {
108  // Set the minimal momentum
109  minMomentum_ = cfg.getParameter<double>("minMomentumCut");
110  // - if not set, particles from outside the beampipe with a negative speed in R direction are propagated but no SimHits
111  // - particle with positive (negative) z and negative (positive) speed in z direction: no SimHits
112  // -> this is not neccesary since a track reconstruction is not possible in this case anyways
113  doHitsFromInboundParticles_ = cfg.getParameter<bool>("doHitsFromInboundParticles");
114 }
115 
117  producesCollector.produces<edm::PSimHitContainer>("TrackerHits");
118 }
119 
121  iEvent.put(std::move(simHitContainer_), "TrackerHits");
122  simHitContainer_.reset(new edm::PSimHitContainer);
123 }
124 
126  const SimplifiedGeometry& layer,
127  std::vector<std::unique_ptr<Particle>>& secondaries,
128  const RandomEngineAndDistribution& random) {
129  // the energy deposit in the layer
130  double energyDeposit = particle.getEnergyDeposit();
131  particle.setEnergyDeposit(0);
132 
133  //
134  // 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)
135  //
136  if (!doHitsFromInboundParticles_) {
137  if (particle.isLooper()) {
138  return;
139  }
140  if (particle.momentum().X() * particle.position().X() + particle.momentum().Y() * particle.position().Y() < 0) {
141  particle.setLooper();
142  return;
143  }
144  if (layer.isForward() && ((particle.position().Z() > 0 && particle.momentum().Z() < 0) ||
145  (particle.position().Z() < 0 && particle.momentum().Z() > 0))) {
146  particle.setLooper();
147  return;
148  }
149  }
150 
151  //
152  // check that layer has tracker modules
153  //
154  if (!layer.getDetLayer()) {
155  return;
156  }
157 
158  //
159  // no material
160  //
161  if (layer.getThickness(particle.position(), particle.momentum()) < 1E-10) {
162  return;
163  }
164 
165  //
166  // only charged particles
167  //
168  if (particle.charge() == 0) {
169  return;
170  }
171 
172  //
173  // save hit only if momentum higher than threshold
174  //
175  if (particle.momentum().Perp2() < minMomentum_ * minMomentum_) {
176  return;
177  }
178 
179  // create the trajectory of the particle
181  GlobalPoint position(particle.position().X(), particle.position().Y(), particle.position().Z());
182  GlobalVector momentum(particle.momentum().Px(), particle.momentum().Py(), particle.momentum().Pz());
183  auto plane = layer.getDetLayer()->surface().tangentPlane(position);
184  TrajectoryStateOnSurface trajectory(
185  GlobalTrajectoryParameters(position, momentum, TrackCharge(particle.charge()), &magneticField), *plane);
186 
187  // find detectors compatible with the particle's trajectory
190  std::vector<DetWithState> compatibleDetectors = layer.getDetLayer()->compatibleDets(trajectory, propagator, est);
191 
193  // You have to sort the simHits in the order they occur!
195 
196  // The old algorithm (sorting by distance to IP) doesn't seem to make sense to me (what if particle moves inwards?)
197 
198  // Detector layers have to be sorted by proximity to particle.position
199  // Propagate particle backwards a bit to make sure it's outside any components
200  std::map<double, std::unique_ptr<PSimHit>> distAndHits;
201  // Position relative to which the hits should be sorted
202  GlobalPoint positionOutside(particle.position().x() - particle.momentum().x() / particle.momentum().P() * 10.,
203  particle.position().y() - particle.momentum().y() / particle.momentum().P() * 10.,
204  particle.position().z() - particle.momentum().z() / particle.momentum().P() * 10.);
205 
206  // FastSim: cheat tracking -> assign hits to closest charged daughter if particle decays
207  int pdgId = particle.pdgId();
208  int simTrackId =
209  particle.getMotherSimTrackIndex() >= 0 ? particle.getMotherSimTrackIndex() : particle.simTrackIndex();
210 
211  // loop over the compatible detectors
212  for (const auto& detectorWithState : compatibleDetectors) {
213  const GeomDet& detector = *detectorWithState.first;
214  const TrajectoryStateOnSurface& particleState = detectorWithState.second;
215 
216  // if the detector has no components
217  if (detector.isLeaf()) {
218  std::pair<double, std::unique_ptr<PSimHit>> hitPair = createHitOnDetector(particleState,
219  pdgId,
220  layer.getThickness(particle.position()),
221  energyDeposit,
222  simTrackId,
223  detector,
224  positionOutside);
225  if (hitPair.second) {
226  distAndHits.insert(distAndHits.end(), std::move(hitPair));
227  }
228  } else {
229  // if the detector has components
230  for (const auto component : detector.components()) {
231  std::pair<double, std::unique_ptr<PSimHit>> hitPair =
232  createHitOnDetector(particleState,
233  pdgId,
234  layer.getThickness(particle.position()),
235  energyDeposit,
236  simTrackId,
237  *component,
238  positionOutside);
239  if (hitPair.second) {
240  distAndHits.insert(distAndHits.end(), std::move(hitPair));
241  }
242  }
243  }
244  }
245 
246  // Fill simHitContainer
247  for (std::map<double, std::unique_ptr<PSimHit>>::const_iterator it = distAndHits.begin(); it != distAndHits.end();
248  it++) {
249  simHitContainer_->push_back(*(it->second));
250  }
251 }
252 
253 // 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)
254 std::pair<double, std::unique_ptr<PSimHit>> fastsim::TrackerSimHitProducer::createHitOnDetector(
255  const TrajectoryStateOnSurface& particle,
256  int pdgId,
257  double layerThickness,
258  double eLoss,
259  int simTrackId,
260  const GeomDet& detector,
261  GlobalPoint& refPos) {
262  // determine position and momentum of particle in the coordinate system of the detector
263  LocalPoint localPosition;
264  LocalVector localMomentum;
265 
266  // if the particle is close enough, no further propagation is needed
267  if (fabs(detector.toLocal(particle.globalPosition()).z()) < onSurfaceTolerance_) {
268  localPosition = particle.localPosition();
269  localMomentum = particle.localMomentum();
270  }
271  // else, propagate
272  else {
273  // find crossing of particle with detector layer
275  particle.globalMomentum().basicVector(),
276  particle.transverseCurvature(),
277  anyDirection);
278  std::pair<bool, double> path = crossing.pathLength(detector.surface());
279  // case propagation succeeds
280  if (path.first) {
281  localPosition = detector.toLocal(GlobalPoint(crossing.position(path.second)));
282  localMomentum = detector.toLocal(GlobalVector(crossing.direction(path.second)));
283  localMomentum = localMomentum.unit() * particle.localMomentum().mag();
284  }
285  // case propagation fails
286  else {
287  return std::pair<double, std::unique_ptr<PSimHit>>(0, std::unique_ptr<PSimHit>(nullptr));
288  }
289  }
290 
291  // find entry and exit point of particle in detector
292  const Plane& detectorPlane = detector.surface();
293  double halfThick = 0.5 * detectorPlane.bounds().thickness();
294  double pZ = localMomentum.z();
295  LocalPoint entry = localPosition + (-halfThick / pZ) * localMomentum;
296  LocalPoint exit = localPosition + (halfThick / pZ) * localMomentum;
297  double tof = particle.globalPosition().mag() / fastsim::Constants::speedOfLight; // in nanoseconds
298 
299  // make sure the simhit is physically on the module
300  double boundX = detectorPlane.bounds().width() / 2.;
301  double boundY = detectorPlane.bounds().length() / 2.;
302  // Special treatment for TID and TEC trapeziodal modules
303  unsigned subdet = DetId(detector.geographicalId()).subdetId();
304  if (subdet == 4 || subdet == 6)
305  boundX *= 1. - localPosition.y() / detectorPlane.position().perp();
306  if (fabs(localPosition.x()) > boundX || fabs(localPosition.y()) > boundY) {
307  return std::pair<double, std::unique_ptr<PSimHit>>(0, std::unique_ptr<PSimHit>(nullptr));
308  }
309 
310  // Create the hit: the energy loss rescaled to the module thickness
311  // Total thickness is in radiation lengths, 1 radlen = 9.36 cm
312  // Sensitive module thickness is about 30 microns larger than
313  // the module thickness itself
314  eLoss *= (2. * halfThick - 0.003) / (9.36 * layerThickness);
315 
316  // Position of the hit in global coordinates
317  GlobalPoint hitPos(detector.surface().toGlobal(localPosition));
318 
319  return std::pair<double, std::unique_ptr<PSimHit>>(
320  (hitPos - refPos).mag(),
321  std::unique_ptr<PSimHit>(new PSimHit(entry,
322  exit,
323  localMomentum.mag(),
324  tof,
325  eLoss,
326  pdgId,
327  detector.geographicalId().rawId(),
328  simTrackId,
329  localMomentum.theta(),
330  localMomentum.phi())));
331 }
332 
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:99
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
HLT_2018_cff.magneticField
magneticField
Definition: HLT_2018_cff.py:348
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:116
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:254
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:106
GeometricSearchDet::compatibleDets
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
Definition: GeometricSearchDet.cc:35
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
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
fastsim::SimplifiedGeometry::getMagneticFieldZ
virtual const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const =0
Return magnetic field (field only has Z component!) on the layer.
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
fastsim::SimplifiedGeometry::getThickness
virtual const double getThickness(const math::XYZTLorentzVector &position) const =0
Return thickness of the layer at a given position.
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:100
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
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
DetWithState
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
Definition: TrackerSimHitProducer.cc:43
fastsim::TrackerSimHitProducer::~TrackerSimHitProducer
~TrackerSimHitProducer() override
Default destructor.
Definition: TrackerSimHitProducer.cc:58
Surface::tangentPlane
virtual ConstReferenceCountingPointer< TangentPlane > tangentPlane(const GlobalPoint &) const =0
Bounds::thickness
virtual float thickness() const =0
TrajectoryStateOnSurface::localPosition
LocalPoint localPosition() const
Definition: TrajectoryStateOnSurface.h:74
HelixArbitraryPlaneCrossing.h
edm::ParameterSet
Definition: ParameterSet.h:36
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:120
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:125
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
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
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:102
SimplifiedGeometry.h
PV3DBase::basicVector
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
looper.cfg
cfg
Definition: looper.py:297
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
fastsim::SimplifiedGeometry::getDetLayer
const DetLayer * getDetLayer() const
Return pointer to the assigned active layer (if any).
Definition: SimplifiedGeometry.h:127
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:52
LocalPoint.h
Plane
Definition: Plane.h:16
GeometricSearchDet::surface
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
DetLayer.h
GlobalVector.h
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
Particle
Definition: Particle.py:1
fastsim::SimplifiedGeometry::isForward
virtual bool isForward() const =0
Returns false/true depending if the object is a (non-abstract) barrel/forward layer.
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:53
ParameterSet.h
PSimHit
Definition: PSimHit.h:15
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:101
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