test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiTrivialInduceChargeOnStrips.cc
Go to the documentation of this file.
2 
10 
11 #include <algorithm>
12 #include <iostream>
13 
14 namespace {
15  struct Count {
16 #ifdef SISTRIP_COUNT
17  // note: this code is not thread safe, counts will be inaccurate if run with multiple threads
18  double ncall=0;
19  double ndep=0, ndep2=0, maxdep=0;
20  double nstr=0, nstr2=0;
21  double ncv=0, nval=0, nval2=0, maxv=0;
22  double dzero=0;
23  void dep(double d) { ncall++; ndep+=d; ndep2+=d*d;maxdep=std::max(d,maxdep);}
24  void str(double d) { nstr+=d; nstr2+=d*d;}
25  void val(double d) { ncv++; nval+=d; nval2+=d*d; maxv=std::max(d,maxv);}
26  void zero() { dzero++;}
27  ~Count() {
28  std::cout << "deposits " << ncall << " " << maxdep << " " << ndep/ncall << " " << std::sqrt(ndep2*ncall -ndep*ndep)/ncall << std::endl;
29  std::cout << "zeros " << dzero << std::endl;
30  std::cout << "strips " << nstr/ndep << " " << std::sqrt(nstr2*ndep -nstr*nstr)/ndep << std::endl;
31  std::cout << "vaules " << ncv << " " << maxv << " " << nval/ncv << " " << std::sqrt(nval2*ncv -nval*nval)/ncv << std::endl;
32  }
33 #else
34  void dep(double) {}
35  void str(double) {}
36  void val(double) {}
37  void zero() {}
38 #endif
39  };
40 
41  Count count;
42 }
43 
44 
45 namespace {
46  constexpr int Ntypes = 14;
47  // 0 1 2 3 4 5 6 7
48  const std::string type[Ntypes] = { "IB1", "IB2","OB1","OB2","W1a","W2a","W3a","W1b","W2b","W3b","W4","W5","W6","W7"};
49  enum { indexOfIB1=0, indexOfIB2=1, indexOfOB1=2, indexOfOB2=3, indexOfW1a=4, indexOfW1b=7};
50 
51  inline
52  std::vector<std::vector<float> >
53  fillSignalCoupling(const edm::ParameterSet& conf, int nTypes, const std::string* typeArray) {
54  std::vector<std::vector<float> > signalCoupling;
55  signalCoupling.reserve(nTypes);
56  std::string mode = conf.getParameter<bool>("APVpeakmode") ? "Peak" : "Dec";
57  for(int i=0; i<nTypes; ++i) {
58  auto dc = conf.getParameter<std::vector<double> >("CouplingConstant"+mode+typeArray[i]);
59  signalCoupling.emplace_back(dc.begin(),dc.end());
60  }
61  return signalCoupling;
62  }
63 
64  inline unsigned int indexOf(const std::string& t) { return std::find( type, type + Ntypes, t) - type;}
65 
66 
67  inline unsigned int typeOf(const StripGeomDetUnit& det, const TrackerTopology *tTopo) {
68  DetId id = det.geographicalId();
69  switch (det.specificType().subDetector()) {
70  case GeomDetEnumerators::TIB: {return (tTopo->tibLayer(id) < 3) ? indexOfIB1 : indexOfIB2;}
71  case GeomDetEnumerators::TOB: {return (tTopo->tobLayer(id) > 4) ? indexOfOB1 : indexOfOB2;}
72  case GeomDetEnumerators::TID: {return indexOfW1a -1 + tTopo->tidRing(id);} //fragile: relies on ordering of 'type'
73  case GeomDetEnumerators::TEC: {return indexOfW1b -1 + tTopo->tecRing(id);} //fragile: relies on ordering of 'type'
74  default: throw cms::Exception("Invalid subdetector") << id();
75  }
76  }
77 
78  inline double
79  chargeDeposited(size_t strip, size_t Nstrips, double amplitude, double chargeSpread, double chargePosition) {
80  double integralUpToStrip = (strip == 0) ? 0. : ( approx_erf( (strip-chargePosition)/chargeSpread/1.41421356237309515 ) );
81  double integralUpToNext = (strip+1 == Nstrips) ? 1. : ( approx_erf( (strip+1-chargePosition)/chargeSpread/1.41421356237309515 ) );
82 
83  return 0.5 * (integralUpToNext - integralUpToStrip) * amplitude;
84  }
85 
86 }
87 
88 
91  : signalCoupling(fillSignalCoupling(conf, Ntypes, type)), Nsigma(3.), geVperElectron(g) {
92 }
93 
94 void
97  const StripGeomDetUnit& det,
98  std::vector<float>& localAmplitudes,
99  size_t& recordMinAffectedStrip,
100  size_t& recordMaxAffectedStrip,
101  const TrackerTopology *tTopo) const {
102 
103 
104  induceVector(collection_points, det, localAmplitudes, recordMinAffectedStrip, recordMaxAffectedStrip, tTopo);
105 
106  /*
107  auto ominA=recordMinAffectedStrip, omaxA=recordMaxAffectedStrip;
108  std::vector<float> oampl(localAmplitudes);
109  induceOriginal(collection_points, det, oampl, ominA, omaxA, tTopo);
110 
111  // std::cout << "orig " << ominA << " " << omaxA << " ";
112  //for (auto a : oampl) std::cout << a << ",";
113  //std::cout << std::endl;
114 
115  auto minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
116  std::vector<float> ampl(localAmplitudes);
117  induceVector(collection_points, det, ampl, minA, maxA, tTopo);
118 
119  // std::cout << "vect " << minA << " " << maxA << " ";
120  //for (auto a : ampl) std::cout << a << ",";
121  //std::cout << std::endl;
122 
123  float diff=0;
124  for (size_t i=0; i!=ampl.size(); ++i) { diff = std::max(diff,ampl[i]>0 ? std::abs(ampl[i]-oampl[i])/ampl[i] : 0);}
125  if (diff> 1.e-4) {
126  std::cout << diff << std::endl;
127  std::cout << "orig " << ominA << " " << omaxA << " ";
128 // for (auto a : oampl) std::cout << a << ",";
129  std::cout << std::endl;
130  std::cout << "vect " << minA << " " << maxA << " ";
131 // for (auto a : ampl) std::cout << a << ",";
132  std::cout << std::endl;
133  }
134 
135  localAmplitudes.swap(ampl);
136  recordMinAffectedStrip=minA;
137  recordMaxAffectedStrip=maxA;
138  */
139 
140 }
141 
142 void
145  const StripGeomDetUnit& det,
146  std::vector<float>& localAmplitudes,
147  size_t& recordMinAffectedStrip,
148  size_t& recordMaxAffectedStrip,
149  const TrackerTopology *tTopo) const {
150 
151 
152  auto const & coupling = signalCoupling[typeOf(det,tTopo)];
153  const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
154  const int Nstrips = topology.nstrips();
155 
156  if (Nstrips == 0) return;
157 
158  const int NP = collection_points.size();
159  if(0==NP) return;
160 
161  constexpr int MaxN = 512;
162  // if NP too large split...
163 
164  for (int ip=0; ip<NP; ip+=MaxN) {
165  auto N = std::min(NP-ip,MaxN);
166 
167  count.dep(N);
168  float amplitude[N];
169  float chargePosition[N];
170  float chargeSpread[N];
171  int fromStrip[N];
172  int nStrip[N];
173 
174  // load not vectorize
175  //In strip coordinates:
176  for (int i=0; i!=N;++i) {
177  auto j = ip+i;
178  if (0==collection_points[j].amplitude()) count.zero();
179  chargePosition[i]=topology.strip(collection_points[j].position());
180  chargeSpread[i]= collection_points[j].sigma() / topology.localPitch(collection_points[j].position());
181  amplitude[i]=0.5f*collection_points[j].amplitude() / geVperElectron;
182  }
183 
184  // this vectorize
185  for (int i=0; i!=N;++i) {
186  fromStrip[i] = std::max( 0, int(std::floor( chargePosition[i] - Nsigma*chargeSpread[i])) );
187  nStrip[i] = std::min( Nstrips, int(std::ceil( chargePosition[i] + Nsigma*chargeSpread[i])) ) - fromStrip[i];
188  }
189  int tot=0;
190  for (int i=0; i!=N;++i) tot += nStrip[i];
191  tot+=N; // add last strip
192  count.val(tot);
193  float value[tot];
194 
195  // assign relative position (lower bound of strip) in value;
196  int kk=0;
197  for (int i=0; i!=N;++i) {
198  auto delta = 1.f/(std::sqrt(2.f)*chargeSpread[i]);
199  auto pos = delta*(float(fromStrip[i])-chargePosition[i]);
200 
201  // VI: before value[0] was not defined and value[tot] was filled
202  // to fix this the loop below was changed
203  for (int j=0;j<=nStrip[i]; ++j) {
204  value[kk] = pos+float(j)*delta;
205  ++kk;
206  }
207  }
208  assert(kk==tot);
209 
210  // main loop fully vectorized
211  for (int k=0;k!=tot; ++k)
212  value[k] = approx_erf(value[k]);
213 
214  // saturate 0 & NStrips strip to 0 and 1???
215  kk=0;
216  for (int i=0; i!=N;++i) {
217  if (0 == fromStrip[i]) value[kk]=0;
218  kk+=nStrip[i];
219  if (Nstrips == fromStrip[i]+nStrip[i]) value[kk]=1.f;
220  ++kk;
221  }
222  assert(kk==tot);
223 
224  // compute integral over strip (lower bound becomes the value)
225  for (int k=0;k!=tot-1; ++k)
226  value[k]-=value[k+1]; // this is negative!
227 
228 
229  float charge[Nstrips]; for (int i=0;i!=Nstrips; ++i) charge[i]=0;
230  kk=0;
231  for (int i=0; i!=N;++i){
232  for (int j=0;j!=nStrip[i]; ++j)
233  charge[fromStrip[i]+j]-= amplitude[i]*value[kk++];
234  ++kk; // skip last "strip"
235  }
236  assert(kk==tot);
237 
238 
240  int minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
241  int sc = coupling.size();
242  for (int i=0;i!=Nstrips; ++i) {
243  int strip = i;
244  if (0==charge[i]) continue;
245  auto affectedFromStrip = std::max( 0, strip - sc + 1);
246  auto affectedUntilStrip = std::min(Nstrips, strip + sc);
247  for (auto affectedStrip=affectedFromStrip; affectedStrip < affectedUntilStrip; ++affectedStrip)
248  localAmplitudes[affectedStrip] += charge[i] * coupling[std::abs(affectedStrip - strip)] ;
249 
250  if( affectedFromStrip < minA ) minA = affectedFromStrip;
251  if( affectedUntilStrip > maxA ) maxA = affectedUntilStrip;
252  }
253  recordMinAffectedStrip=minA;
254  recordMaxAffectedStrip=maxA;
255  } // end loop ip
256 
257 }
258 
259 void
262  const StripGeomDetUnit& det,
263  std::vector<float>& localAmplitudes,
264  size_t& recordMinAffectedStrip,
265  size_t& recordMaxAffectedStrip,
266  const TrackerTopology *tTopo) const {
267 
268 
269  auto const & coupling = signalCoupling[typeOf(det,tTopo)];
270  const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
271  size_t Nstrips = topology.nstrips();
272 
273  if (!collection_points.empty()) count.dep(collection_points.size());
274 
275  for (auto signalpoint = collection_points.begin(); signalpoint != collection_points.end(); signalpoint++ ) {
276 
277  //In strip coordinates:
278  double chargePosition = topology.strip(signalpoint->position());
279  double chargeSpread = signalpoint->sigma() / topology.localPitch(signalpoint->position());
280 
281  size_t fromStrip = size_t(std::max( 0, int(std::floor( chargePosition - Nsigma*chargeSpread))));
282  size_t untilStrip = size_t(std::min( Nstrips, size_t(std::ceil( chargePosition + Nsigma*chargeSpread) )));
283 
284  count.str(std::max(0,int(untilStrip)-int(fromStrip)));
285  for (size_t strip = fromStrip; strip < untilStrip; strip++) {
286 
287  double chargeDepositedOnStrip = chargeDeposited( strip, Nstrips, signalpoint->amplitude() / geVperElectron, chargeSpread, chargePosition);
288 
289  size_t affectedFromStrip = size_t(std::max( 0, int(strip - coupling.size() + 1)));
290  size_t affectedUntilStrip = size_t(std::min( Nstrips, strip + coupling.size()) );
291  for (size_t affectedStrip = affectedFromStrip; affectedStrip < affectedUntilStrip; affectedStrip++) {
292  localAmplitudes.at( affectedStrip ) += chargeDepositedOnStrip * coupling.at(abs( affectedStrip - strip )) ;
293  }
294 
295  if( affectedFromStrip < recordMinAffectedStrip ) recordMinAffectedStrip = affectedFromStrip;
296  if( affectedUntilStrip > recordMaxAffectedStrip ) recordMaxAffectedStrip = affectedUntilStrip;
297  }
298  }
299 
300 }
301 
dbl * delta
Definition: mlp_gen.cc:36
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
virtual int nstrips() const =0
int i
Definition: DBlmapReader.cc:9
unsigned int tibLayer(const DetId &id) const
unsigned int tidRing(const DetId &id) const
CaloTopology const * topology(0)
unsigned int tecRing(const DetId &id) const
ring id
assert(m_qm.get())
std::vector< SignalPoint > collection_type
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
int dzero(double a, double b, double &x0, double &rv, double eps, int mxf, F func)
Private version of the exponential integral.
Definition: VVIObj.cc:570
void induce(const SiChargeCollectionDrifter::collection_type &collection_points, const StripGeomDetUnit &det, std::vector< float > &localAmplitudes, size_t &recordMinAffectedStrip, size_t &recordMaxAffectedStrip, const TrackerTopology *tTopo) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
#define constexpr
virtual float strip(const LocalPoint &) const =0
float approx_erf(float x)
Definition: approx_erf.h:6
tuple d
Definition: ztail.py:151
virtual float localPitch(const LocalPoint &) const =0
virtual StripGeomDetType const & specificType() const
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
double f[11][100]
T min(T a, T b)
Definition: MathUtil.h:58
Definition: DetId.h:18
#define N
Definition: blowfish.cc:9
static int position[264][3]
Definition: ReadPGInfo.cc:509
tuple cout
Definition: gather_cfg.py:145
SiTrivialInduceChargeOnStrips(const edm::ParameterSet &conf, double g)
void induceVector(const SiChargeCollectionDrifter::collection_type &collection_points, const StripGeomDetUnit &det, std::vector< float > &localAmplitudes, size_t &recordMinAffectedStrip, size_t &recordMaxAffectedStrip, const TrackerTopology *tTopo) const
SubDetector subDetector() const
Definition: GeomDetType.h:22
unsigned int tobLayer(const DetId &id) const
void induceOriginal(const SiChargeCollectionDrifter::collection_type &collection_points, const StripGeomDetUnit &det, std::vector< float > &localAmplitudes, size_t &recordMinAffectedStrip, size_t &recordMaxAffectedStrip, const TrackerTopology *tTopo) const