CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GEMSignalModel.cc
Go to the documentation of this file.
7 #include "CLHEP/Random/RandFlat.h"
8 #include "CLHEP/Random/RandPoissonQ.h"
9 #include "CLHEP/Random/RandGaussQ.h"
10 #include "CLHEP/Units/GlobalPhysicalConstants.h"
12 #include <cmath>
13 #include <utility>
14 #include <map>
15 
17  : GEMDigiModel(config),
18  averageEfficiency_(config.getParameter<double>("averageEfficiency")),
19  averageShapingTime_(config.getParameter<double>("averageShapingTime")),
20  timeResolution_(config.getParameter<double>("timeResolution")),
21  timeJitter_(config.getParameter<double>("timeJitter")),
22  signalPropagationSpeed_(config.getParameter<double>("signalPropagationSpeed")),
23  bx0filter_(config.getParameter<bool>("bx0filter")),
24  digitizeOnlyMuons_(config.getParameter<bool>("digitizeOnlyMuons")),
25  resolutionX_(config.getParameter<double>("resolutionX")),
26  muonPdgId(13),
27  cspeed(geant_units::operators::convertMmToCm(CLHEP::c_light)),
28  momConvFact(1000.),
29  elecMomCut1(1.96e-03),
30  elecMomCut2(10.e-03),
31  elecEffLowCoeff(1.7e-05),
32  elecEffLowParam0(2.1),
33  elecEffMidCoeff(1.34),
34  elecEffMidParam0(-5.75e-01),
35  elecEffMidParam1(7.96e-01) {}
36 
38 
41  CLHEP::HepRandomEngine* engine,
42  Strips& strips_,
43  DetectorHitMap& detectorHitMap_) {
44  bool digiMuon = false;
45  bool digiElec = false;
46  const GEMStripTopology* top(dynamic_cast<const GEMStripTopology*>(&(roll->topology())));
47  for (const auto& hit : simHits) {
48  if (std::abs(hit.particleType()) != muonPdgId && digitizeOnlyMuons_)
49  continue;
50  double elecEff = 0.;
51  double partMom = hit.pabs();
52  double checkMuonEff = CLHEP::RandFlat::shoot(engine, 0., 1.);
53  double checkElecEff = CLHEP::RandFlat::shoot(engine, 0., 1.);
54  if (std::abs(hit.particleType()) == muonPdgId && checkMuonEff < averageEfficiency_)
55  digiMuon = true;
56  if (std::abs(hit.particleType()) != muonPdgId) //consider all non muon particles with gem efficiency to electrons
57  {
58  if (partMom <= elecMomCut1)
59  elecEff = elecEffLowCoeff * std::exp(elecEffLowParam0 * partMom * momConvFact);
60  if (partMom > elecMomCut1 && partMom < elecMomCut2)
61  elecEff = elecEffMidCoeff * log(elecEffMidParam1 * partMom * momConvFact + elecEffMidParam0) /
62  (elecEffMidCoeff + log(elecEffMidParam1 * partMom * momConvFact + elecEffMidParam0));
63  if (partMom > elecMomCut2)
64  elecEff = 1.;
65  if (checkElecEff < elecEff)
66  digiElec = true;
67  }
68  if (!(digiMuon || digiElec))
69  continue;
70  const int bx(getSimHitBx(&hit, engine));
71  if (bx != 0 and bx0filter_)
72  continue;
73  const std::vector<std::pair<int, int> >& cluster(simulateClustering(top, &hit, bx, engine));
74  for (const auto& digi : cluster) {
75  detectorHitMap_.emplace(digi, &hit);
76  strips_.emplace(digi);
77  }
78  }
79 }
80 
81 int GEMSignalModel::getSimHitBx(const PSimHit* simhit, CLHEP::HepRandomEngine* engine) {
82  int bx = -999;
83  const LocalPoint simHitPos(simhit->localPosition());
84  const float tof(simhit->timeOfFlight());
85  // random Gaussian time correction due to electronics jitter
86  float randomJitterTime = CLHEP::RandGaussQ::shoot(engine, 0., timeJitter_);
87  const GEMDetId id(simhit->detUnitId());
88  const GEMEtaPartition* roll(geometry_->etaPartition(id));
89  if (!roll) {
90  throw cms::Exception("Geometry") << "GEMSignalModel::getSimHitBx() - GEM simhit id does not match any GEM roll id: "
91  << id << "\n";
92  return 999;
93  }
94  if (roll->id().region() == 0) {
95  throw cms::Exception("Geometry")
96  << "GEMSignalModel::getSimHitBx() - this GEM id is from barrel, which cannot happen: " << roll->id() << "\n";
97  }
98  const int nstrips = roll->nstrips();
99  float middleStrip = nstrips / 2.;
100  const LocalPoint& middleOfRoll = roll->centreOfStrip(middleStrip);
101  const GlobalPoint& globMiddleRol = roll->toGlobal(middleOfRoll);
102  double muRadius = sqrt(globMiddleRol.x() * globMiddleRol.x() + globMiddleRol.y() * globMiddleRol.y() +
103  globMiddleRol.z() * globMiddleRol.z());
104  double timeCalibrationOffset_ = muRadius / cspeed; //[ns]
105 
106  const GEMStripTopology* top(dynamic_cast<const GEMStripTopology*>(&(roll->topology())));
107  const float halfStripLength(0.5 * top->stripLength());
108  const float distanceFromEdge(halfStripLength - simHitPos.y());
109 
110  // signal propagation speed in material in [cm/ns]
111  double signalPropagationSpeedTrue = signalPropagationSpeed_ * cspeed;
112 
113  // average time for the signal to propagate from the SimHit to the top of a strip
114  const float averagePropagationTime(distanceFromEdge / signalPropagationSpeedTrue);
115  // random Gaussian time correction due to the finite timing resolution of the detector
116  float randomResolutionTime = CLHEP::RandGaussQ::shoot(engine, 0., timeResolution_);
117  const float simhitTime(tof + averageShapingTime_ + randomResolutionTime + averagePropagationTime + randomJitterTime);
118  float referenceTime = 0.;
119  referenceTime = timeCalibrationOffset_ + halfStripLength / signalPropagationSpeedTrue + averageShapingTime_;
120  const float timeDifference(simhitTime - referenceTime);
121  // assign the bunch crossing
122  bx = static_cast<int>(std::round((timeDifference) / 25.));
123 
124  // check time
125  LogDebug("GEMDigiProducer") << "checktime "
126  << "bx = " << bx << "\tdeltaT = " << timeDifference << "\tsimT = " << simhitTime
127  << "\trefT = " << referenceTime << "\ttof = " << tof
128  << "\tavePropT = " << averagePropagationTime
129  << "\taveRefPropT = " << halfStripLength / signalPropagationSpeedTrue << "\n";
130  return bx;
131 }
132 
133 std::vector<std::pair<int, int> > GEMSignalModel::simulateClustering(const GEMStripTopology* top,
134  const PSimHit* simHit,
135  const int bx,
136  CLHEP::HepRandomEngine* engine) {
137  const LocalPoint& hit_entry(simHit->entryPoint());
138  const LocalPoint& hit_exit(simHit->exitPoint());
139 
140  LocalPoint start_point, end_point;
141  if (hit_entry.x() < hit_exit.x()) {
142  start_point = hit_entry;
143  end_point = hit_exit;
144  } else {
145  start_point = hit_exit;
146  end_point = hit_entry;
147  }
148 
149  // Add Gaussian noise to the points towards outside.
150  float smeared_start_x = start_point.x() - std::abs(CLHEP::RandGaussQ::shoot(engine, 0, resolutionX_));
151  float smeared_end_x = end_point.x() + std::abs(CLHEP::RandGaussQ::shoot(engine, 0, resolutionX_));
152 
153  LocalPoint smeared_start_point(smeared_start_x, start_point.y(), start_point.z());
154  LocalPoint smeared_end_point(smeared_end_x, end_point.y(), end_point.z());
155 
156  int cluster_start = top->channel(smeared_start_point);
157  int cluster_end = top->channel(smeared_end_point);
158 
159  std::vector<std::pair<int, int> > cluster;
160  for (int strip = cluster_start; strip <= cluster_end; strip++) {
161  cluster.emplace_back(strip, bx);
162  }
163 
164  return cluster;
165 }
const double elecMomCut1
static std::vector< std::string > checklist log
double averageEfficiency_
const double elecEffMidParam1
uint16_t *__restrict__ id
const double elecEffMidCoeff
T y() const
Definition: PV3DBase.h:60
std::vector< std::pair< int, int > > simulateClustering(const GEMStripTopology *, const PSimHit *, const int, CLHEP::HepRandomEngine *)
int channel(const LocalPoint &) const override
const double cspeed
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
const GEMEtaPartition * etaPartition(GEMDetId id) const
Return a GEMEtaPartition given its id.
Definition: GEMGeometry.cc:77
double signalPropagationSpeed_
const double elecEffLowCoeff
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
float timeOfFlight() const
Definition: PSimHit.h:73
Local3DPoint localPosition() const
Definition: PSimHit.h:52
const double elecEffLowParam0
const double elecMomCut2
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void simulate(const GEMEtaPartition *, const edm::PSimHitContainer &, CLHEP::HepRandomEngine *, Strips &, DetectorHitMap &) override
const double elecEffMidParam0
double timeResolution_
int getSimHitBx(const PSimHit *, CLHEP::HepRandomEngine *)
const Topology & topology() const override
tuple simHits
Definition: trackerHits.py:16
std::multimap< std::pair< unsigned int, int >, const PSimHit *, std::less< std::pair< unsigned int, int > > > DetectorHitMap
Definition: GEMDigiModel.h:35
std::set< std::pair< int, int > > Strips
Definition: GEMDigiModel.h:31
const GEMGeometry * geometry_
Definition: GEMDigiModel.h:47
tuple config
parse the configuration file
constexpr NumType convertMmToCm(NumType millimeters)
Definition: GeantUnits.h:63
std::vector< PSimHit > PSimHitContainer
~GEMSignalModel() override
double averageShapingTime_
T x() const
Definition: PV3DBase.h:59
const int muonPdgId
GEMSignalModel(const edm::ParameterSet &)
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
const double momConvFact
unsigned int detUnitId() const
Definition: PSimHit.h:97
#define LogDebug(id)