CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  corr(config.getParameter<bool>("useCorrelation")),
26  etaproj(config.getParameter<bool>("useEtaProjectiveGEO")),
27  digitizeOnlyMuons_(config.getParameter<bool>("digitizeOnlyMuons")),
28  gaussianSmearing_(config.getParameter<bool>("gaussianSmearing")),
29  averageEfficiency_(config.getParameter<double>("averageEfficiency")),
30  // simulateIntrinsicNoise_(config.getParameter<bool>("simulateIntrinsicNoise")),
31  // averageNoiseRate_(config.getParameter<double>("averageNoiseRate")),
32  simulateElectronBkg_(config.getParameter<bool>("simulateElectronBkg")),
33  simulateNeutralBkg_(config.getParameter<bool>("simulateNeutralBkg")),
34  minBunch_(config.getParameter<int>("minBunch")),
35  maxBunch_(config.getParameter<int>("maxBunch"))
36 {
37  // polynomial parametrisation of neutral (n+g) and electron background
38  neuBkg.push_back(899644.0); neuBkg.push_back(-30841.0); neuBkg.push_back(441.28);
39  neuBkg.push_back(-3.3405); neuBkg.push_back(0.0140588); neuBkg.push_back(-3.11473e-05); neuBkg.push_back(2.83736e-08);
40  eleBkg.push_back(4.68590e+05); eleBkg.push_back(-1.63834e+04); eleBkg.push_back(2.35700e+02);
41  eleBkg.push_back(-1.77706e+00); eleBkg.push_back(7.39960e-03); eleBkg.push_back(-1.61448e-05); eleBkg.push_back(1.44368e-08);
42 }
43 
45 {
46 }
47 
48 void ME0PreRecoGaussianModel::simulateSignal(const ME0EtaPartition* roll, const edm::PSimHitContainer& simHits, CLHEP::HepRandomEngine* engine)
49 {
50 for (const auto & hit: simHits)
51 {
52  // Digitize only Muons?
53  if (std::abs(hit.particleType()) != 13 && digitizeOnlyMuons_) continue;
54  // Digitize only in [minBunch,maxBunch] window
55  // window is: [(2n-1)*bxw/2, (2n+1)*bxw/2], n = [minBunch, maxBunch]
56  if(hit.timeOfFlight() < (2*minBunch_-1)*bxwidth*1.0/2 || hit.timeOfFlight() > (2*maxBunch_+1)*bxwidth*1.0/2) continue;
57  // is GEM efficient?
58  if (CLHEP::RandFlat::shoot(engine) > averageEfficiency_) continue;
59  // create digi
60  auto entry = hit.entryPoint();
61  double x=0.0, y=0.0;
62  if(gaussianSmearing_) { // Gaussian Smearing
63  x=CLHEP::RandGaussQ::shoot(engine, entry.x(), sigma_u);
64  y=CLHEP::RandGaussQ::shoot(engine, entry.y(), sigma_v);
65  }
66  else { // Uniform Smearing ... use the sigmas as boundaries
67  x=entry.x()+(CLHEP::RandFlat::shoot(engine)-0.5)*sigma_u;
68  y=entry.y()+(CLHEP::RandFlat::shoot(engine)-0.5)*sigma_v;
69  }
70  double ex=sigma_u;
71  double ey=sigma_v;
72  double corr=0.;
73  double tof=CLHEP::RandGaussQ::shoot(engine, hit.timeOfFlight(), sigma_t);
74  int pdgid = hit.particleType();
75  ME0DigiPreReco digi(x,y,ex,ey,corr,tof,pdgid);
76  digi_.insert(digi);
77 
78  edm::LogVerbatim("ME0PreRecoGaussianModel") << "[ME0PreRecoDigi :: simulateSignal] :: simhit in "<<roll->id()<<" at loc x = "<<std::setw(8)<<entry.x()<<" [cm]"
79  << " loc y = "<<std::setw(8)<<entry.y()<<" [cm] time = "<<std::setw(8)<<hit.timeOfFlight()<<" [ns] pdgid = "<<std::showpos<<std::setw(4)<<pdgid;
80  edm::LogVerbatim("ME0PreRecoGaussianModel") << "[ME0PreRecoDigi :: simulateSignal] :: digi in "<<roll->id()<<" at loc x = "<<std::setw(8)<<x<<" [cm] loc y = "<<std::setw(8)<<y<<" [cm]"
81  <<" time = "<<std::setw(8)<<tof<<" [ns]";
82  edm::LogVerbatim("ME0PreRecoGaussianModel") << "[ME0PreRecoDigi :: simulateSignal] :: digi in "<<roll->id()<<" with DX = "<<std::setw(8)<<(entry.x()-x)<<" [cm]"
83  <<" DY = "<<std::setw(8)<<(entry.y()-y)<<" [cm] DT = "<<std::setw(8)<<(hit.timeOfFlight()-tof)<<" [ns]";
84 
85 }
86 }
87 
88 void ME0PreRecoGaussianModel::simulateNoise(const ME0EtaPartition* roll, CLHEP::HepRandomEngine* engine)
89 {
90  double trArea(0.0);
91  const ME0DetId me0Id(roll->id());
92 
93  // Extract detailed information from the Strip Topology:
94  // base_bottom, base_top, height, strips, pads
95  // note that (0,0) is in the middle of the roll ==> all param are at all half length
96  if (me0Id.region() == 0) { throw cms::Exception("Geometry") << "Asking TrapezoidalStripTopology from a ME0 will fail"; } // not sure we really need this
97  const TrapezoidalStripTopology* top_(dynamic_cast<const TrapezoidalStripTopology*>(&(roll->topology())));
98 
99  auto& parameters(roll->specs()->parameters());
100  double bottomLength(parameters[0]); bottomLength = 2*bottomLength; // bottom is largest length, so furtest away from beamline
101  double topLength(parameters[1]); topLength = 2*topLength; // top is shortest length, so closest to beamline
102  double height(parameters[2]); height = 2*height;
103  double myTanPhi = (topLength - bottomLength) / (height * 2);
104  double rollRadius = top_->radius();
105  trArea = height * (topLength + bottomLength) / 2.0;
106 
107  // Divide the detector area in different strips
108  // take smearing in y-coord as height for each strip
109  double heightIt = sigma_v;
110  int heightbins = height/heightIt; // round down
111 
112  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: roll with id = "<<roll->id();
113  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: extracting parameters from the TrapezoidalStripTopology";
114  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: bottom = "<<bottomLength<<" [cm] top = "<<topLength<<" [cm] height = "<<height<<" [cm]"
115  <<" area = "<<trArea<<" [cm^2] Rmid = "<<rollRadius<<" [cm] => Rmin = "<<rollRadius-height*1.0/2.0<<" [cm] Rmax = "<<rollRadius+height*1.0/2.0<<" [cm]";
116  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: heightbins = "<<heightbins;
117 
118  for(int hx=0; hx<heightbins; ++hx) {
119  double bottomIt = bottomLength + hx *2*tan(10./180*3.14)*heightIt;
120  double topIt = bottomLength + (hx+1)*2*tan(10./180*3.14)*heightIt;
121  if(hx==heightbins-1) {
122  topIt = topLength; // last bin ... make strip a bit larger to cover entire roll
123  heightIt = height-hx*heightIt;
124  }
125  double areaIt = heightIt*(bottomIt+topIt)*1.0/2;
126 
127  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: sNoise]["<<roll->id().rawId()<<"] :: height = "<<std::setw(12)<<heightIt<<" [cm] bottom = "<<std::setw(12)<<bottomIt<<" [cm]"
128  << " top = "<<std::setw(12)<<topIt<<" [cm] area = "<<std::setw(12)<<areaIt<<" [cm^2] || sin(10) = "<<sin(10./180*3.14);
129 
130  double myRandY = CLHEP::RandFlat::shoot(engine);
131  double y0_rand = (hx+myRandY)*heightIt; // Y coord, measured from the bottom of the roll
132  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
133  double yy_glob = rollRadius + yy_rand; // R coord in Global Coords
134  // max length in x for given y coordinate (cfr trapezoidal eta partition)
135  double xMax = topLength/2.0 - (height/2.0 - yy_rand) * myTanPhi;
136 
137  // 1) Intrinsic Noise ... Not implemented right now
138  // ------------------------------------------------
139  // if (simulateIntrinsicNoise_)
140  // {
141  // }
142 
143  // 2) Background Noise
144  // ----------------------------
145 
146  // 2a) electron background
147  // -----------------------
148  if (simulateElectronBkg_) {
149  // Extract / Calculate the Average Electron Rate
150  // for the given global Y coord from Parametrization
151  double averageElectronRatePerRoll = 0.0;
152  double yy_helper = 1.0;
153  for(int j=0; j<7; ++j) { averageElectronRatePerRoll += eleBkg[j]*yy_helper; yy_helper *= yy_glob; }
154 
155  // Rate [Hz/cm^2] * Nbx * 25*10^-9 [s] * Area [cm] = # hits in this roll in this bx
156  const double averageElecRate(averageElectronRatePerRoll * (maxBunch_-minBunch_+1)*(bxwidth*1.0e-9) * areaIt);
157 
158  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: elebkg]["<<roll->id().rawId()<<"]" /* "] :: BX = "<<std::showpos<<bx*/
159  << " evaluation of Background Hit Rate at this coord :: "<<std::setw(12)<<averageElectronRatePerRoll<<" [Hz/cm^2]"
160  << " x 9 x 25*10^-9 [s] x Area (of strip = "<<std::setw(12)<<areaIt<<" [cm^2]) ==> "<<std::setw(12)<<averageElecRate<<" [hits]";
161 
162  bool ele_eff = (CLHEP::RandFlat::shoot(engine)<averageElecRate)?1:0;
163 
164  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: elebkg]["<<roll->id().rawId()<<"] :: myRandY = "<<std::setw(12)<<myRandY<<" => local y = "<<std::setw(12)<<yy_rand<<" [cm]"
165  <<" => global y (global R) = "<<std::setw(12)<<yy_glob<<" [cm] || Probability = "<<std::setw(12)<<averageElecRate
166  <<" => efficient? "<<ele_eff<<std::endl;
167 
168  if(ele_eff) {
169  //calculate xx_rand at a given yy_rand
170  double myRandX = CLHEP::RandFlat::shoot(engine);
171  double xx_rand = 2 * xMax * (myRandX - 0.5);
172  double ex = sigma_u;
173  double ey = sigma_v;
174  double corr = 0.;
175  // extract random BX
176  double myrandBX = CLHEP::RandFlat::shoot(engine);
177  int bx = int((maxBunch_-minBunch_+1)*myrandBX)+minBunch_;
178  // extract random time in this BX
179  double myrandT = CLHEP::RandFlat::shoot(engine);
180  double minBXtime = (bx-0.5)*bxwidth; // double maxBXtime = (bx+0.5)*bxwidth;
181  double time = myrandT*bxwidth+minBXtime;
182  double myrandP = CLHEP::RandFlat::shoot(engine);
183  int pdgid = 0;
184  if (myrandP <= 0.5) pdgid = -11; // electron
185  else pdgid = 11; // positron
186  ME0DigiPreReco digi(xx_rand, yy_rand, ex, ey, corr, time, pdgid);
187  digi_.insert(digi);
188  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: elebkg]["<<roll->id().rawId()<<"] =====> electron hit in "<<roll->id()<<" pdgid = "<<pdgid<<" bx = "<<bx
189  <<" ==> digitized"
190  <<" at loc x = "<<xx_rand<<" loc y = "<<yy_rand<<" time = "<<time<<" [ns]";
191  }
192  } // end if electron bkg
193 
194  // 2b) neutral (n+g) background
195  // ----------------------------
196  if (simulateNeutralBkg_) {
197  // Extract / Calculate the Average Neutral Rate
198  // for the given global Y coord from Parametrization
199  double averageNeutralRatePerRoll = 0.0;
200  double yy_helper = 1.0;
201  for(int j=0; j<7; ++j) { averageNeutralRatePerRoll += neuBkg[j]*yy_helper; yy_helper *= yy_glob; }
202 
203  // Rate [Hz/cm^2] * Nbx * 25*10^-9 [s] * Area [cm] = # hits in this roll
204  const double averageNeutrRate(averageNeutralRatePerRoll * (maxBunch_-minBunch_+1)*(bxwidth*1.0e-9) * areaIt);
205 
206  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: neubkg]["<<roll->id().rawId()<<"]" /* "] :: BX = "<<std::showpos<<bx*/
207  << " evaluation of Background Hit Rate at this coord :: "<<std::setw(12)<<averageNeutralRatePerRoll<<" [Hz/cm^2]"
208  << " x 9 x 25*10^-9 [s] x Area (of strip = "<<std::setw(12)<<areaIt<<" [cm^2]) ==> "<<std::setw(12)<<averageNeutrRate<<" [hits]";
209 
210  bool neu_eff = (CLHEP::RandFlat::shoot(engine)<averageNeutrRate)?1:0;
211 
212  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: neubkg]["<<roll->id().rawId()<<"] :: myRandY = "<<std::setw(12)<<myRandY<<" => local y = "<<std::setw(12)<<yy_rand<<" [cm]"
213  <<" => global y (global R) = "<<std::setw(12)<<yy_glob<<" [cm] || Probability = "<<std::setw(12)<<averageNeutrRate
214  <<" => efficient? "<<neu_eff<<std::endl;
215 
216  if(neu_eff) {
217  //calculate xx_rand at a given yy_rand
218  double myRandX = CLHEP::RandFlat::shoot(engine);
219  double xx_rand = 2 * xMax * (myRandX - 0.5);
220  double ex = sigma_u;
221  double ey = sigma_v;
222  double corr = 0.;
223  // extract random BX
224  double myrandBX= CLHEP::RandFlat::shoot(engine);
225  int bx = int((maxBunch_-minBunch_+1)*myrandBX)+minBunch_;
226  // extract random time in this BX
227  double myrandT = CLHEP::RandFlat::shoot(engine);
228  double minBXtime = (bx-0.5)*bxwidth;
229  double time = myrandT*bxwidth+minBXtime;
230  int pdgid = 0;
231  double myrandP = CLHEP::RandFlat::shoot(engine);
232  if (myrandP <= 0.08) pdgid = 2112; // neutrons: GEM sensitivity for neutrons: 0.08%
233  else pdgid = 22; // photons: GEM sensitivity for photons: 1.04% ==> neutron fraction = (0.08 / 1.04) = 0.077 = 0.08
234  ME0DigiPreReco digi(xx_rand, yy_rand, ex, ey, corr, time, pdgid);
235  digi_.insert(digi);
236  edm::LogVerbatim("ME0PreRecoGaussianModelNoise") << "[ME0PreRecoDigi :: neubkg]["<<roll->id().rawId()<<"] ======> neutral hit in "<<roll->id()<<" pdgid = "<<pdgid<<" bx = "<<bx
237  <<" ==> digitized"
238  <<" at loc x = "<<xx_rand<<" loc y = "<<yy_rand<<" time = "<<time<<" [ns]";
239  }
240 
241  } // end if neutral bkg
242  } // end loop over strips (= pseudo rolls)
243 }
244 
const ME0EtaPartitionSpecs * specs() const
std::vector< double > neuBkg
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const ME0Specs & parameters() const
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_
tuple simHits
Definition: trackerHits.py:16
void simulateNoise(const ME0EtaPartition *, CLHEP::HepRandomEngine *) override
const Topology & topology() const
ME0PreRecoGaussianModel(const edm::ParameterSet &)
std::vector< PSimHit > PSimHitContainer
std::vector< double > eleBkg