CMS 3D CMS Logo

GEMSimpleModel.cc
Go to the documentation of this file.
6 #include "CLHEP/Random/RandFlat.h"
7 #include "CLHEP/Random/RandPoissonQ.h"
8 #include "CLHEP/Random/RandGaussQ.h"
9 #include <cmath>
10 #include <utility>
11 #include <map>
12 
14 GEMDigiModel(config)
15 , averageEfficiency_(config.getParameter<double> ("averageEfficiency"))
16 , averageShapingTime_(config.getParameter<double> ("averageShapingTime"))
17 , timeResolution_(config.getParameter<double> ("timeResolution"))
18 , timeJitter_(config.getParameter<double> ("timeJitter"))
19 , averageNoiseRate_(config.getParameter<double> ("averageNoiseRate"))
20 , signalPropagationSpeed_(config.getParameter<double> ("signalPropagationSpeed"))
21 , bxwidth_(config.getParameter<int> ("bxwidth"))
22 , minBunch_(config.getParameter<int> ("minBunch"))
23 , maxBunch_(config.getParameter<int> ("maxBunch"))
24 , digitizeOnlyMuons_(config.getParameter<bool> ("digitizeOnlyMuons"))
25 , doBkgNoise_(config.getParameter<bool> ("doBkgNoise"))
26 , doNoiseCLS_(config.getParameter<bool> ("doNoiseCLS"))
27 , fixedRollRadius_(config.getParameter<bool> ("fixedRollRadius"))
28 , simulateElectronBkg_(config.getParameter<bool> ("simulateElectronBkg"))
29 , instLumi_(config.getParameter<double>("instLumi"))
30 , rateFact_(config.getParameter<double>("rateFact"))
31 , referenceInstLumi_(config.getParameter<double>("referenceInstLumi"))
32 , resolutionX_(config.getParameter<double>("resolutionX"))
33 , GE11ElecBkgParam0_(config.getParameter<double>("GE11ElecBkgParam0"))
34 , GE11ElecBkgParam1_(config.getParameter<double>("GE11ElecBkgParam1"))
35 , GE11ElecBkgParam2_(config.getParameter<double>("GE11ElecBkgParam2"))
36 , GE21ElecBkgParam0_(config.getParameter<double>("GE21ElecBkgParam0"))
37 , GE21ElecBkgParam1_(config.getParameter<double>("GE21ElecBkgParam1"))
38 , GE21ElecBkgParam2_(config.getParameter<double>("GE21ElecBkgParam2"))
39 , GE11ModNeuBkgParam0_(config.getParameter<double>("GE11ModNeuBkgParam0"))
40 , GE11ModNeuBkgParam1_(config.getParameter<double>("GE11ModNeuBkgParam1"))
41 , GE11ModNeuBkgParam2_(config.getParameter<double>("GE11ModNeuBkgParam2"))
42 , GE21ModNeuBkgParam0_(config.getParameter<double>("GE11ModNeuBkgParam0"))
43 , GE21ModNeuBkgParam1_(config.getParameter<double>("GE11ModNeuBkgParam1"))
44 , GE21ModNeuBkgParam2_(config.getParameter<double>("GE11ModNeuBkgParam2"))
45 {
46 }
47 
49 {
50 }
51 
53 {
54  return;
55 }
56 
57 void GEMSimpleModel::simulateSignal(const GEMEtaPartition* roll, const edm::PSimHitContainer& simHits, CLHEP::HepRandomEngine* engine)
58 {
60  detectorHitMap_.clear();
64  bool digiMuon = false;
65  bool digiElec = false;
66  for (const auto& hit : simHits)
67  {
68  if (std::abs(hit.particleType()) != 13 && digitizeOnlyMuons_)
69  continue;
70  double elecEff = 0.;
71  double partMom = hit.pabs();
72  double checkMuonEff = CLHEP::RandFlat::shoot(engine, 0., 1.);
73  double checkElecEff = CLHEP::RandFlat::shoot(engine, 0., 1.);
74  if (std::abs(hit.particleType()) == 13 && checkMuonEff < averageEfficiency_)
75  digiMuon = true;
76  if (std::abs(hit.particleType()) != 13) //consider all non muon particles with gem efficiency to electrons
77  {
78  if (partMom <= 1.95e-03)
79  elecEff = 1.7e-05 * std::exp(2.1 * partMom * 1000.);
80  if (partMom > 1.95e-03 && partMom < 10.e-03)
81  elecEff = 1.34 * log(7.96e-01 * partMom * 1000. - 5.75e-01)
82  / (1.34 + log(7.96e-01 * partMom * 1000. - 5.75e-01));
83  if (partMom > 10.e-03)
84  elecEff = 1.;
85  if (checkElecEff < elecEff)
86  digiElec = true;
87  }
88  if (!(digiMuon || digiElec))
89  continue;
90  const int bx(getSimHitBx(&hit, engine));
91  const std::vector<std::pair<int, int> >& cluster(simulateClustering(roll, &hit, bx, engine));
92  for (const auto & digi : cluster)
93  {
94  detectorHitMap_.emplace(digi,&hit);
95  strips_.emplace(digi);
96  }
97  }
98 }
99 
100 int GEMSimpleModel::getSimHitBx(const PSimHit* simhit, CLHEP::HepRandomEngine* engine)
101 {
102  int bx = -999;
103  const LocalPoint simHitPos(simhit->localPosition());
104  const float tof(simhit->timeOfFlight());
105  // random Gaussian time correction due to electronics jitter
106  float randomJitterTime = CLHEP::RandGaussQ::shoot(engine, 0., timeJitter_);
107  const GEMDetId id(simhit->detUnitId());
108  const GEMEtaPartition* roll(geometry_->etaPartition(id));
109  if (!roll)
110  {
111  throw cms::Exception("Geometry")<< "GEMSimpleModel::getSimHitBx() - GEM simhit id does not match any GEM roll id: " << id << "\n";
112  return 999;
113  }
114  if (roll->id().region() == 0)
115  {
116  throw cms::Exception("Geometry") << "GEMSimpleModel::getSimHitBx() - this GEM id is from barrel, which cannot happen: " << roll->id() << "\n";
117  }
118  const double cspeed = 299792458; // signal propagation speed in vacuum in [m/s]
119  const int nstrips = roll->nstrips();
120  float middleStrip = nstrips/2.;
121  const LocalPoint& middleOfRoll = roll->centreOfStrip(middleStrip);
122  const GlobalPoint& globMiddleRol = roll->toGlobal(middleOfRoll);
123  double muRadius = sqrt(globMiddleRol.x()*globMiddleRol.x() + globMiddleRol.y()*globMiddleRol.y() +globMiddleRol.z()*globMiddleRol.z());
124  double timeCalibrationOffset_ = (muRadius *1e+9)/(cspeed*1e+2); //[ns]
125 
126  const TrapezoidalStripTopology* top(dynamic_cast<const TrapezoidalStripTopology*> (&(roll->topology())));
127  const float halfStripLength(0.5 * top->stripLength());
128  const float distanceFromEdge(halfStripLength - simHitPos.y());
129 
130  // signal propagation speed in material in [cm/ns]
131  double signalPropagationSpeedTrue = signalPropagationSpeed_ * cspeed * 1e-7; // 1e+2 * 1e-9;
132 
133  // average time for the signal to propagate from the SimHit to the top of a strip
134  const float averagePropagationTime(distanceFromEdge / signalPropagationSpeedTrue);
135  // random Gaussian time correction due to the finite timing resolution of the detector
136  float randomResolutionTime = CLHEP::RandGaussQ::shoot(engine, 0., timeResolution_);
137  const float simhitTime(tof + averageShapingTime_ + randomResolutionTime + averagePropagationTime + randomJitterTime);
138  float referenceTime = 0.;
139  referenceTime = timeCalibrationOffset_ + halfStripLength / signalPropagationSpeedTrue + averageShapingTime_;
140  const float timeDifference(simhitTime - referenceTime);
141  // assign the bunch crossing
142  bx = static_cast<int> (std::round((timeDifference) / bxwidth_));
143 
144  // check time
145  const bool debug(false);
146  if (debug)
147  {
148  std::cout << "checktime " << "bx = " << bx << "\tdeltaT = " << timeDifference << "\tsimT = " << simhitTime
149  << "\trefT = " << referenceTime << "\ttof = " << tof << "\tavePropT = " << averagePropagationTime
150  << "\taveRefPropT = " << halfStripLength / signalPropagationSpeedTrue << std::endl;
151  }
152  return bx;
153 }
154 
155 void GEMSimpleModel::simulateNoise(const GEMEtaPartition* roll, CLHEP::HepRandomEngine* engine)
156 {
157  if (!doBkgNoise_)
158  return;
159  const GEMDetId& gemId(roll->id());
160  const int nstrips(roll->nstrips());
161  double trArea(0.0);
162  double trStripArea(0.0);
163  if (gemId.region() == 0)
164  {
165  throw cms::Exception("Geometry")
166  << "GEMSynchronizer::simulateNoise() - this GEM id is from barrel, which cannot happen.";
167  }
168  const TrapezoidalStripTopology* top_(dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology())));
169  const float striplength(top_->stripLength());
170  trStripArea = (roll->pitch()) * striplength;
171  trArea = trStripArea * nstrips;
172  const int nBxing(maxBunch_ - minBunch_ + 1);
173  const float rollRadius(fixedRollRadius_ ? top_->radius() :
174  top_->radius() + CLHEP::RandFlat::shoot(engine, -1.*top_->stripLength()/2., top_->stripLength()/2.));
175 
176  //calculate noise from model
177  double averageNeutralNoiseRatePerRoll = 0.;
178  double averageNoiseElectronRatePerRoll = 0.;
179  double averageNoiseRatePerRoll = 0.;
180  if (gemId.station() == 1)
181  {
182  //simulate neutral background for GE1/1
183  averageNeutralNoiseRatePerRoll = (GE11ModNeuBkgParam0_
184  + GE11ModNeuBkgParam1_ * rollRadius
185  + GE11ModNeuBkgParam2_ * rollRadius * rollRadius); //simulate electron background for GE1/1
187  averageNoiseElectronRatePerRoll = (GE11ElecBkgParam0_
188  + GE11ElecBkgParam1_ * rollRadius
189  + GE11ElecBkgParam2_ * rollRadius * rollRadius);
190 
191  // Scale up/down for desired instantaneous lumi (reference is 5E34, double from config is in units of 1E34)
192  averageNoiseRatePerRoll = averageNeutralNoiseRatePerRoll + averageNoiseElectronRatePerRoll;
193  averageNoiseRatePerRoll *= instLumi_*rateFact_*1.0/referenceInstLumi_;
194  }
195  if (gemId.station() == 2)
196  {
197  //simulate neutral background for GE2/1
198  averageNeutralNoiseRatePerRoll = (GE21ModNeuBkgParam0_
199  + GE21ModNeuBkgParam1_ * rollRadius
200  + GE21ModNeuBkgParam2_ * rollRadius * rollRadius);
201  //simulate electron background for GE2/1
203  averageNoiseElectronRatePerRoll = (GE21ElecBkgParam0_
204  + GE21ElecBkgParam1_ * rollRadius
205  + GE21ElecBkgParam2_ * rollRadius * rollRadius);
206 
207  // Scale up/down for desired instantaneous lumi (reference is 5E34, double from config is in units of 1E34)
208  averageNoiseRatePerRoll = averageNeutralNoiseRatePerRoll + averageNoiseElectronRatePerRoll;
209  averageNoiseRatePerRoll *= instLumi_*rateFact_*1.0/referenceInstLumi_;
210  }
211 
212  //simulate intrinsic noise
214  {
215  const double aveIntrinsicNoisePerStrip(averageNoiseRate_ * nBxing * bxwidth_ * trStripArea * 1.0e-9);
216  for(int j = 0; j < nstrips; ++j)
217  {
218  CLHEP::RandPoissonQ randPoissonQ(*engine, aveIntrinsicNoisePerStrip);
219  const int n_intrHits(randPoissonQ.fire());
220 
221  for (int k = 0; k < n_intrHits; k++ )
222  {
223  const int time_hit(static_cast<int>(CLHEP::RandFlat::shoot(engine, nBxing)) + minBunch_);
224  std::pair<int, int> digi(k+1,time_hit);
225  strips_.emplace(digi);
226  }
227  }
228  }//end simulate intrinsic noise
229 
230  //simulate bkg contribution
231  const double averageNoise(averageNoiseRatePerRoll * nBxing * bxwidth_ * trArea * 1.0e-9);
232  CLHEP::RandPoissonQ randPoissonQ(*engine, averageNoise);
233  const int n_hits(randPoissonQ.fire());
234  for (int i = 0; i < n_hits; ++i)
235  {
236  const int centralStrip(static_cast<int> (CLHEP::RandFlat::shoot(engine, 1, nstrips)));
237  const int time_hit(static_cast<int>(CLHEP::RandFlat::shoot(engine, nBxing)) + minBunch_);
238  if (doNoiseCLS_)
239  {
240  std::vector < std::pair<int, int> > cluster_;
241  cluster_.clear();
242  cluster_.emplace_back(centralStrip, time_hit);
243  int clusterSize((CLHEP::RandFlat::shoot(engine)) <= 0.53 ? 1 : 2);
244  if (clusterSize == 2)
245  {
246  if(CLHEP::RandFlat::shoot(engine) < 0.5)
247  {
248  if (CLHEP::RandFlat::shoot(engine) < averageEfficiency_ && (centralStrip - 1 > 0))
249  cluster_.emplace_back(centralStrip - 1, time_hit);
250  }
251  else
252  {
253  if (CLHEP::RandFlat::shoot(engine) < averageEfficiency_ && (centralStrip + 1 <= nstrips))
254  cluster_.emplace_back(centralStrip + 1, time_hit);
255  }
256  }
257  for (const auto& digi : cluster_)
258  {
259  strips_.emplace(digi);
260  }
261  } //end doNoiseCLS_
262  else
263  {
264  strips_.emplace(centralStrip, time_hit);
265  }
266  }
267  return;
268 }
269 
270 
271 std::vector<std::pair<int, int> > GEMSimpleModel::simulateClustering(
272  const GEMEtaPartition* roll,
273  const PSimHit* simHit,
274  const int bx,
275  CLHEP::HepRandomEngine* engine) {
276 
277  const LocalPoint & hit_entry(simHit->entryPoint());
278  const LocalPoint & hit_exit(simHit->exitPoint());
279 
280  LocalPoint start_point, end_point;
281  if(hit_entry.x() < hit_exit.x()) {
282  start_point = hit_entry;
283  end_point = hit_exit;
284  } else {
285  start_point = hit_exit;
286  end_point = hit_entry;
287  }
288 
289  // Add Gaussian noise to the points towards outside.
290  float smeared_start_x = start_point.x() - std::abs(CLHEP::RandGaussQ::shoot(engine, 0, resolutionX_));
291  float smeared_end_x = end_point.x() + std::abs(CLHEP::RandGaussQ::shoot(engine, 0, resolutionX_));
292 
293  LocalPoint smeared_start_point(smeared_start_x, start_point.y(), start_point.z());
294  LocalPoint smeared_end_point(smeared_end_x, end_point.y(), end_point.z());
295 
296  int cluster_start = roll->strip(smeared_start_point);
297  int cluster_end = roll->strip(smeared_end_point);
298 
299  std::vector< std::pair<int, int> > cluster;
300  for (int strip = cluster_start; strip <= cluster_end; strip++) {
301  cluster.emplace_back(strip, bx);
302  }
303 
304  return cluster;
305 }
double signalPropagationSpeed_
double GE21ElecBkgParam1_
double GE11ModNeuBkgParam2_
edm::DetSet< GEMDigiSimLink > GEMDigiSimLinks
Definition: GEMDigiModel.h:37
double GE11ElecBkgParam2_
double averageNoiseRate_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
T y() const
Definition: PV3DBase.h:63
int getSimHitBx(const PSimHit *, CLHEP::HepRandomEngine *)
double GE11ElecBkgParam0_
void setup() override
Definition: config.py:1
GEMDigiSimLinks theGemDigiSimLinks_
Definition: GEMDigiModel.h:80
const Topology & topology() const override
const GEMEtaPartition * etaPartition(GEMDetId id) const
Return a GEMEtaPartition given its id.
Definition: GEMGeometry.cc:99
GEMDetId id() const
double GE21ElecBkgParam2_
int nstrips() const
number of readout strips in partition
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
float timeOfFlight() const
Definition: PSimHit.h:69
bool simulateIntrinsicNoise_
~GEMSimpleModel() override
Local3DPoint localPosition() const
Definition: PSimHit.h:44
T sqrt(T t)
Definition: SSEVec.h:18
double GE21ModNeuBkgParam1_
double averageEfficiency_
double GE21ModNeuBkgParam2_
double GE11ModNeuBkgParam1_
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::pair< int, int > > simulateClustering(const GEMEtaPartition *, const PSimHit *, const int, CLHEP::HepRandomEngine *) override
double GE21ModNeuBkgParam0_
double GE11ModNeuBkgParam0_
GEMSimpleModel(const edm::ParameterSet &)
float pitch() const
double GE21ElecBkgParam0_
int k[5][pyjets_maxn]
const double referenceInstLumi_
bool simulateElectronBkg_
#define debug
Definition: HDRShower.cc:19
float strip(const LocalPoint &lp) const
double GE11ElecBkgParam1_
void simulateNoise(const GEMEtaPartition *, CLHEP::HepRandomEngine *) override
const GEMGeometry * geometry_
Definition: GEMDigiModel.h:63
double averageShapingTime_
std::set< std::pair< int, int > > strips_
Definition: GEMDigiModel.h:65
void clear()
Definition: DetSet.h:69
DetectorHitMap detectorHitMap_
Definition: GEMDigiModel.h:78
double timeResolution_
std::vector< PSimHit > PSimHitContainer
StripDigiSimLinks stripDigiSimLinks_
Definition: GEMDigiModel.h:79
edm::DetSet< StripDigiSimLink > StripDigiSimLinks
Definition: GEMDigiModel.h:36
T x() const
Definition: PV3DBase.h:62
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
void simulateSignal(const GEMEtaPartition *, const edm::PSimHitContainer &, CLHEP::HepRandomEngine *) override
unsigned int detUnitId() const
Definition: PSimHit.h:93