CMS 3D CMS Logo

ME0PreRecoGaussianModel.cc
Go to the documentation of this file.
5 
8 
9 #include "CLHEP/Random/RandFlat.h"
10 #include "CLHEP/Random/RandGaussQ.h"
11 #include "CLHEP/Random/RandPoissonQ.h"
12 
13 #include <cmath>
14 #include <utility>
15 #include <map>
16 
17 
18 const int bxwidth = 25; // [ns]
19 
21  ME0DigiPreRecoModel(config),
22  sigma_t(config.getParameter<double>("timeResolution")),
23  sigma_u(config.getParameter<double>("phiResolution")),
24  sigma_v(config.getParameter<double>("etaResolution")),
25  error_u(config.getParameter<double>("phiError")),
26  error_v(config.getParameter<double>("etaError")),
27  gaussianSmearing_(config.getParameter<bool>("gaussianSmearing")),
28  constPhiSmearing_(config.getParameter<bool>("constantPhiSpatialResolution")),
29  corr(config.getParameter<bool>("useCorrelation")),
30  etaproj(config.getParameter<bool>("useEtaProjectiveGEO")),
31  digitizeOnlyMuons_(config.getParameter<bool>("digitizeOnlyMuons")),
32  averageEfficiency_(config.getParameter<double>("averageEfficiency")),
33  // simulateIntrinsicNoise_(config.getParameter<bool>("simulateIntrinsicNoise")),
34  // averageNoiseRate_(config.getParameter<double>("averageNoiseRate")),
35  simulateElectronBkg_(config.getParameter<bool>("simulateElectronBkg")),
36  simulateNeutralBkg_(config.getParameter<bool>("simulateNeutralBkg")),
37  minBunch_(config.getParameter<int>("minBunch")),
38  maxBunch_(config.getParameter<int>("maxBunch")),
39  instLumi_(config.getParameter<double>("instLumi")),
40  rateFact_(config.getParameter<double>("rateFact"))
41 {
42  // polynomial parametrisation of neutral (n+g) and electron background
43  // This is the background for an Instantaneous Luminosity of L = 5E34 cm^-2 s^-1
44  neuBkg.push_back(899644.0); neuBkg.push_back(-30841.0); neuBkg.push_back(441.28);
45  neuBkg.push_back(-3.3405); neuBkg.push_back(0.0140588); neuBkg.push_back(-3.11473e-05); neuBkg.push_back(2.83736e-08);
46  eleBkg.push_back(4.68590e+05); eleBkg.push_back(-1.63834e+04); eleBkg.push_back(2.35700e+02);
47  eleBkg.push_back(-1.77706e+00); eleBkg.push_back(7.39960e-03); eleBkg.push_back(-1.61448e-05); eleBkg.push_back(1.44368e-08);
48 }
49 
51 {
52 }
53 
54 void ME0PreRecoGaussianModel::simulateSignal(const ME0EtaPartition* roll, const edm::PSimHitContainer& simHits, CLHEP::HepRandomEngine* engine)
55 {
56 for (const auto & hit: simHits)
57 {
58  // Digitize only Muons?
59  if (std::abs(hit.particleType()) != 13 && digitizeOnlyMuons_) continue;
60  // Digitize only in [minBunch,maxBunch] window
61  // window is: [(2n+1)*bxw/2, (2n+3)*bxw/2], n = [minBunch, maxBunch]
62  if(hit.timeOfFlight() < (2*minBunch_+1)*bxwidth*1.0/2 || hit.timeOfFlight() > (2*maxBunch_+3)*bxwidth*1.0/2) continue;
63  // is GEM efficient?
64  if (CLHEP::RandFlat::shoot(engine) > averageEfficiency_) continue;
65  // create digi
66  auto entry = hit.entryPoint();
67  double x=0.0, y=0.0;
68 
69  double sigma_u_new = sigma_u;
70  if(constPhiSmearing_) sigma_u_new = correctSigmaU(roll, entry.y());
71 
72  if(gaussianSmearing_) { // Gaussian Smearing
73  x=CLHEP::RandGaussQ::shoot(engine, entry.x(), sigma_u_new);
74  y=CLHEP::RandGaussQ::shoot(engine, entry.y(), sigma_v);
75  }
76  else { // Uniform Smearing ... use the sigmas as boundaries
77  x=entry.x()+(CLHEP::RandFlat::shoot(engine)-0.5)*sigma_u_new;
78  y=entry.y()+(CLHEP::RandFlat::shoot(engine)-0.5)*sigma_v;
79  }
80  double ex=sigma_u_new;
81  double ey=sigma_v;
82  double corr=0.;
83  double tof=CLHEP::RandGaussQ::shoot(engine, hit.timeOfFlight(), sigma_t);
84  int pdgid = hit.particleType();
85  if (ex == 0) ex = error_u;//errors cannot be zero
86  if (ey == 0) ey = error_v;
87 
88  int evtId = hit.eventId().event();
89  int bx = hit.eventId().bunchCrossing();
90  int procType = hit.processType();
91  int res = 1;
92  if(!(evtId == 0 && bx == 0 && procType == 0)) res = 2;
93 
94  ME0DigiPreReco digi(x,y,ex,ey,corr,tof,pdgid,res);
95  digi_.insert(digi);
96 
97  edm::LogVerbatim("ME0PreRecoGaussianModel") << "[ME0PreRecoDigi :: simulateSignal] :: simhit in "<<roll->id()<<" at loc x = "<<std::setw(8)<<entry.x()<<" [cm]"
98  << " loc y = "<<std::setw(8)<<entry.y()<<" [cm] time = "<<std::setw(8)<<hit.timeOfFlight()<<" [ns] pdgid = "<<std::showpos<<std::setw(4)<<pdgid;
99  edm::LogVerbatim("ME0PreRecoGaussianModel") << "[ME0PreRecoDigi :: simulateSignal] :: digi in "<<roll->id()<<" at loc x = "<<std::setw(8)<<x<<" [cm] loc y = "<<std::setw(8)<<y<<" [cm]"
100  <<" time = "<<std::setw(8)<<tof<<" [ns]";
101  edm::LogVerbatim("ME0PreRecoGaussianModel") << "[ME0PreRecoDigi :: simulateSignal] :: digi in "<<roll->id()<<" with DX = "<<std::setw(8)<<(entry.x()-x)<<" [cm]"
102  <<" DY = "<<std::setw(8)<<(entry.y()-y)<<" [cm] DT = "<<std::setw(8)<<(hit.timeOfFlight()-tof)<<" [ns]";
103 
104 }
105 }
106 
107 void ME0PreRecoGaussianModel::simulateNoise(const ME0EtaPartition* roll, CLHEP::HepRandomEngine* engine)
108 {
109  double trArea(0.0);
110  const ME0DetId me0Id(roll->id());
111 
112  // Extract detailed information from the Strip Topology:
113  // base_bottom, base_top, height, strips, pads
114  // note that (0,0) is in the middle of the roll ==> all param are at all half length
115  if (me0Id.region() == 0) { throw cms::Exception("Geometry") << "Asking TrapezoidalStripTopology from a ME0 will fail"; } // not sure we really need this
116  const TrapezoidalStripTopology* top_(dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology())));
117 
118  auto& parameters(roll->specs()->parameters());
119  double bottomLength(parameters[0]); bottomLength = 2*bottomLength; // bottom is largest length, so furtest away from beamline
120  double topLength(parameters[1]); topLength = 2*topLength; // top is shortest length, so closest to beamline
121  double height(parameters[2]); height = 2*height;
122  double myTanPhi = (topLength - bottomLength) / (height * 2);
123  double rollRadius = top_->radius();
124  trArea = height * (topLength + bottomLength) / 2.0;
125 
126  // Divide the detector area in different strips
127  // take smearing in y-coord as height for each strip
128  double initialHeight = sigma_v;
129  if(sigma_v < 1.0) initialHeight = 1.0;
130  double heightIt = initialHeight;
131  int heightbins = height/heightIt; // round down
132 
133  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: roll with id = "<<roll->id();
134  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: extracting parameters from the TrapezoidalStripTopology";
135  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: bottom = "<<bottomLength<<" [cm] top = "<<topLength<<" [cm] height = "<<height<<" [cm]"
136  <<" area = "<<trArea<<" [cm^2] Rmid = "<<rollRadius<<" [cm] => Rmin = "<<rollRadius-height*1.0/2.0<<" [cm] Rmax = "<<rollRadius+height*1.0/2.0<<" [cm]";
137  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: heightbins = "<<heightbins;
138 
139  for(int hx=0; hx<heightbins; ++hx) {
140  double bottomIt = bottomLength + hx *2*tan(10./180*3.14)*heightIt;
141  double topIt = bottomLength + (hx+1)*2*tan(10./180*3.14)*heightIt;
142  if(hx==heightbins-1) {
143  topIt = topLength; // last bin ... make strip a bit larger to cover entire roll
144  heightIt = height-hx*heightIt;
145  }
146  double areaIt = heightIt*(bottomIt+topIt)*1.0/2;
147 
148  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: height = "<<std::setw(12)<<heightIt<<" [cm] bottom = "<<std::setw(12)<<bottomIt<<" [cm]"
149  << " top = "<<std::setw(12)<<topIt<<" [cm] area = "<<std::setw(12)<<areaIt<<" [cm^2] || sin(10) = "<<sin(10./180*3.14);
150 
151  double myRandY = CLHEP::RandFlat::shoot(engine);
152  double y0_rand = (hx+myRandY)*heightIt; // Y coord, measured from the bottom of the roll
153  if(hx==heightbins-1) y0_rand = hx*initialHeight + myRandY*heightIt;
154  double yy_rand = (y0_rand-height*1.0/2); // Y coord, measured from the middle of the roll, which is the Y coord in Local Coords
155  double yy_glob = rollRadius + yy_rand; // R coord in Global Coords
156  // max length in x for given y coordinate (cfr trapezoidal eta partition)
157  double xMax = topLength/2.0 - (height/2.0 - yy_rand) * myTanPhi;
158 
159  double sigma_u_new = sigma_u;
160  if(constPhiSmearing_) sigma_u_new = correctSigmaU(roll, yy_rand);
161 
162 
163  // 1) Intrinsic Noise ... Not implemented right now
164  // ------------------------------------------------
165  // if (simulateIntrinsicNoise_)
166  // {
167  // }
168 
169  // 2) Background Noise
170  // ----------------------------
171 
172  // 2a) electron background
173  // -----------------------
174  if (simulateElectronBkg_) {
175  // Extract / Calculate the Average Electron Rate
176  // for the given global Y coord from Parametrization
177  double averageElectronRatePerRoll = 0.0;
178  double yy_helper = 1.0;
179  for(int j=0; j<7; ++j) { averageElectronRatePerRoll += eleBkg[j]*yy_helper; yy_helper *= yy_glob; }
180 
181  // Scale up/down for desired instantaneous lumi (reference is 5E34, double from config is in units of 1E34)
182  averageElectronRatePerRoll *= instLumi_*rateFact_*1.0/5;
183 
184  // Rate [Hz/cm^2] * Nbx * 25*10^-9 [s] * Area [cm] = # hits in this roll in this bx
185  const double averageElecRate(averageElectronRatePerRoll * (maxBunch_-minBunch_+1)*(bxwidth*1.0e-9) * areaIt);
186 
187  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: elebkg]["<<roll->id().rawId()<<"]" /* "] :: BX = "<<std::showpos<<bx*/
188  << " evaluation of Background Hit Rate at this coord :: "<<std::setw(12)<<averageElectronRatePerRoll<<" [Hz/cm^2]"
189  << " x 9 x 25*10^-9 [s] x Area (of strip = "<<std::setw(12)<<areaIt<<" [cm^2]) ==> "<<std::setw(12)<<averageElecRate<<" [hits]";
190 
191  bool ele_eff = (CLHEP::RandFlat::shoot(engine)<averageElecRate)?1:0;
192 
193  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: elebkg]["<<roll->id().rawId()<<"] :: myRandY = "<<std::setw(12)<<myRandY<<" => local y = "<<std::setw(12)<<yy_rand<<" [cm]"
194  <<" => global y (global R) = "<<std::setw(12)<<yy_glob<<" [cm] || Probability = "<<std::setw(12)<<averageElecRate
195  <<" => efficient? "<<ele_eff<<std::endl;
196 
197  if(ele_eff) {
198  //calculate xx_rand at a given yy_rand
199  double myRandX = CLHEP::RandFlat::shoot(engine);
200  double xx_rand = 2 * xMax * (myRandX - 0.5);
201  double ex = sigma_u_new;
202  double ey = sigma_v;
203  double corr = 0.;
204  // extract random BX
205  double myrandBX = CLHEP::RandFlat::shoot(engine);
206  int bx = int((maxBunch_-minBunch_+1)*myrandBX)+minBunch_;
207  // extract random time in this BX
208  double myrandT = CLHEP::RandFlat::shoot(engine);
209  double minBXtime = (bx-0.5)*bxwidth; // double maxBXtime = (bx+0.5)*bxwidth;
210  double time = myrandT*bxwidth+minBXtime;
211  double myrandP = CLHEP::RandFlat::shoot(engine);
212  int pdgid = 0;
213  if (myrandP <= 0.5) pdgid = -11; // electron
214  else pdgid = 11; // positron
215  if (ex == 0) ex = error_u;//errors cannot be zero
216  if (ey == 0) ey = error_v;
217  ME0DigiPreReco digi(xx_rand, yy_rand, ex, ey, corr, time, pdgid, 0);
218  digi_.insert(digi);
219  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: elebkg]["<<roll->id().rawId()<<"] =====> electron hit in "<<roll->id()<<" pdgid = "<<pdgid<<" bx = "<<bx
220  <<" ==> digitized"
221  <<" at loc x = "<<xx_rand<<" loc y = "<<yy_rand<<" time = "<<time<<" [ns]";
222  }
223  } // end if electron bkg
224 
225  // 2b) neutral (n+g) background
226  // ----------------------------
227  if (simulateNeutralBkg_) {
228  // Extract / Calculate the Average Neutral Rate
229  // for the given global Y coord from Parametrization
230  double averageNeutralRatePerRoll = 0.0;
231  double yy_helper = 1.0;
232  for(int j=0; j<7; ++j) { averageNeutralRatePerRoll += neuBkg[j]*yy_helper; yy_helper *= yy_glob; }
233 
234  // Scale up/down for desired instantaneous lumi (reference is 5E34, double from config is in units of 1E34)
235  averageNeutralRatePerRoll *= instLumi_*rateFact_*1.0/5;
236 
237  // Rate [Hz/cm^2] * Nbx * 25*10^-9 [s] * Area [cm] = # hits in this roll
238  const double averageNeutrRate(averageNeutralRatePerRoll * (maxBunch_-minBunch_+1)*(bxwidth*1.0e-9) * areaIt);
239 
240  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: neubkg]["<<roll->id().rawId()<<"]" /* "] :: BX = "<<std::showpos<<bx*/
241  << " evaluation of Background Hit Rate at this coord :: "<<std::setw(12)<<averageNeutralRatePerRoll<<" [Hz/cm^2]"
242  << " x 9 x 25*10^-9 [s] x Area (of strip = "<<std::setw(12)<<areaIt<<" [cm^2]) ==> "<<std::setw(12)<<averageNeutrRate<<" [hits]";
243 
244  bool neu_eff = (CLHEP::RandFlat::shoot(engine)<averageNeutrRate)?1:0;
245 
246  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: neubkg]["<<roll->id().rawId()<<"] :: myRandY = "<<std::setw(12)<<myRandY<<" => local y = "<<std::setw(12)<<yy_rand<<" [cm]"
247  <<" => global y (global R) = "<<std::setw(12)<<yy_glob<<" [cm] || Probability = "<<std::setw(12)<<averageNeutrRate
248  <<" => efficient? "<<neu_eff<<std::endl;
249 
250  if(neu_eff) {
251  //calculate xx_rand at a given yy_rand
252  double myRandX = CLHEP::RandFlat::shoot(engine);
253  double xx_rand = 2 * xMax * (myRandX - 0.5);
254  double ex = sigma_u_new;
255  double ey = sigma_v;
256  double corr = 0.;
257  // extract random BX
258  double myrandBX= CLHEP::RandFlat::shoot(engine);
259  int bx = int((maxBunch_-minBunch_+1)*myrandBX)+minBunch_;
260  // extract random time in this BX
261  double myrandT = CLHEP::RandFlat::shoot(engine);
262  double minBXtime = (bx-0.5)*bxwidth;
263  double time = myrandT*bxwidth+minBXtime;
264  int pdgid = 0;
265  double myrandP = CLHEP::RandFlat::shoot(engine);
266  if (myrandP <= 0.08) pdgid = 2112; // neutrons: GEM sensitivity for neutrons: 0.08%
267  else pdgid = 22; // photons: GEM sensitivity for photons: 1.04% ==> neutron fraction = (0.08 / 1.04) = 0.077 = 0.08
268  if (ex == 0) ex = error_u;//errors cannot be zero
269  if (ey == 0) ey = error_v;
270  ME0DigiPreReco digi(xx_rand, yy_rand, ex, ey, corr, time, pdgid, 0);
271  digi_.insert(digi);
272  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: neubkg]["<<roll->id().rawId()<<"] ======> neutral hit in "<<roll->id()<<" pdgid = "<<pdgid<<" bx = "<<bx
273  <<" ==> digitized"
274  <<" at loc x = "<<xx_rand<<" loc y = "<<yy_rand<<" time = "<<time<<" [ns]";
275  }
276 
277  } // end if neutral bkg
278  } // end loop over strips (= pseudo rolls)
279 }
280 
282  const TrapezoidalStripTopology* top_(dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology())));
283  auto& parameters(roll->specs()->parameters());
284  double height(parameters[2]); // height = height from Center of Roll
285  double rollRadius = top_->radius(); // rollRadius = Radius at Center of Roll
286  double Rmax = rollRadius+height; // MaxRadius = Radius at top of Roll
287  double Rx = rollRadius+y; // y in [-height,+height]
288  double sigma_u_new = Rx/Rmax*sigma_u;
289  return sigma_u_new;
290 }
const ME0EtaPartitionSpecs * specs() const
std::vector< double > neuBkg
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const ME0Specs & parameters() const
Definition: config.py:1
Definition: Electron.h:4
double correctSigmaU(const ME0EtaPartition *, double)
const int bxwidth
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
ME0DetId id() const
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
JetCorrectorParameters corr
Definition: classes.h:5
void simulateSignal(const ME0EtaPartition *, const edm::PSimHitContainer &, CLHEP::HepRandomEngine *) override
std::set< ME0DigiPreReco > digi_
void simulateNoise(const ME0EtaPartition *, CLHEP::HepRandomEngine *) override
const Topology & topology() const
ME0PreRecoGaussianModel(const edm::ParameterSet &)
const double Rmax[kNumberCalorimeter]
std::vector< PSimHit > PSimHitContainer
std::vector< double > eleBkg