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  double ncall=0;
17  double ndep=0, ndep2=0, maxdep=0;
18  double nstr=0, nstr2=0;
19  double ncv=0, nval=0, nval2=0, maxv=0;
20  double dzero=0;
21  void dep(double d) { ncall++; ndep+=d; ndep2+=d*d;maxdep=std::max(d,maxdep);}
22  void str(double d) { nstr+=d; nstr2+=d*d;}
23  void val(double d) { ncv++; nval+=d; nval2+=d*d; maxv=std::max(d,maxv);}
24  void zero() { dzero++;}
25  ~Count() {
26 #ifdef SISTRIP_COUNT
27  std::cout << "deposits " << ncall << " " << maxdep << " " << ndep/ncall << " " << std::sqrt(ndep2*ncall -ndep*ndep)/ncall << std::endl;
28  std::cout << "zeros " << dzero << std::endl;
29  std::cout << "strips " << nstr/ndep << " " << std::sqrt(nstr2*ndep -nstr*nstr)/ndep << std::endl;
30  std::cout << "vaules " << ncv << " " << maxv << " " << nval/ncv << " " << std::sqrt(nval2*ncv -nval*nval)/ncv << std::endl;
31 #endif
32  }
33  };
34 
35  Count count;
36 }
37 
38 
39 namespace {
40  constexpr int Ntypes = 14;
41  // 0 1 2 3 4 5 6 7
42  const std::string type[Ntypes] = { "IB1", "IB2","OB1","OB2","W1a","W2a","W3a","W1b","W2b","W3b","W4","W5","W6","W7"};
43  enum { indexOfIB1=0, indexOfIB2=1, indexOfOB1=2, indexOfOB2=3, indexOfW1a=4, indexOfW1b=7};
44 
45  inline
46  std::vector<std::vector<float> >
47  fillSignalCoupling(const edm::ParameterSet& conf, int nTypes, const std::string* typeArray) {
48  std::vector<std::vector<float> > signalCoupling;
49  signalCoupling.reserve(nTypes);
50  std::string mode = conf.getParameter<bool>("APVpeakmode") ? "Peak" : "Dec";
51  for(int i=0; i<nTypes; ++i) {
52  auto dc = conf.getParameter<std::vector<double> >("CouplingConstant"+mode+typeArray[i]);
53  signalCoupling.emplace_back(dc.begin(),dc.end());
54  }
55  return signalCoupling;
56  }
57 
58  inline unsigned int indexOf(const std::string& t) { return std::find( type, type + Ntypes, t) - type;}
59 
60 
61  inline unsigned int typeOf(const StripGeomDetUnit& det, const TrackerTopology *tTopo) {
62  DetId id = det.geographicalId();
63  switch (det.specificType().subDetector()) {
64  case GeomDetEnumerators::TIB: {return (tTopo->tibLayer(id) < 3) ? indexOfIB1 : indexOfIB2;}
65  case GeomDetEnumerators::TOB: {return (tTopo->tobLayer(id) > 4) ? indexOfOB1 : indexOfOB2;}
66  case GeomDetEnumerators::TID: {return indexOfW1a -1 + tTopo->tidRing(id);} //fragile: relies on ordering of 'type'
67  case GeomDetEnumerators::TEC: {return indexOfW1b -1 + tTopo->tecRing(id);} //fragile: relies on ordering of 'type'
68  default: throw cms::Exception("Invalid subdetector") << id();
69  }
70  }
71 
72  inline double
73  chargeDeposited(size_t strip, size_t Nstrips, double amplitude, double chargeSpread, double chargePosition) {
74  double integralUpToStrip = (strip == 0) ? 0. : ( approx_erf( (strip-chargePosition)/chargeSpread/1.41421356237309515 ) );
75  double integralUpToNext = (strip+1 == Nstrips) ? 1. : ( approx_erf( (strip+1-chargePosition)/chargeSpread/1.41421356237309515 ) );
76 
77  return 0.5 * (integralUpToNext - integralUpToStrip) * amplitude;
78  }
79 
80 }
81 
82 
85  : signalCoupling(fillSignalCoupling(conf, Ntypes, type)), Nsigma(3.), geVperElectron(g) {
86 }
87 
88 void
91  const StripGeomDetUnit& det,
92  std::vector<float>& localAmplitudes,
93  size_t& recordMinAffectedStrip,
94  size_t& recordMaxAffectedStrip,
95  const TrackerTopology *tTopo) const {
96 
97 
98  induceVector(collection_points, det, localAmplitudes, recordMinAffectedStrip, recordMaxAffectedStrip, tTopo);
99 
100  /*
101  auto ominA=recordMinAffectedStrip, omaxA=recordMaxAffectedStrip;
102  std::vector<float> oampl(localAmplitudes);
103  induceOriginal(collection_points, det, oampl, ominA, omaxA, tTopo);
104 
105  // std::cout << "orig " << ominA << " " << omaxA << " ";
106  //for (auto a : oampl) std::cout << a << ",";
107  //std::cout << std::endl;
108 
109  auto minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
110  std::vector<float> ampl(localAmplitudes);
111  induceVector(collection_points, det, ampl, minA, maxA, tTopo);
112 
113  // std::cout << "vect " << minA << " " << maxA << " ";
114  //for (auto a : ampl) std::cout << a << ",";
115  //std::cout << std::endl;
116 
117  float diff=0;
118  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);}
119  if (diff> 1.e-4) {
120  std::cout << diff << std::endl;
121  std::cout << "orig " << ominA << " " << omaxA << " ";
122 // for (auto a : oampl) std::cout << a << ",";
123  std::cout << std::endl;
124  std::cout << "vect " << minA << " " << maxA << " ";
125 // for (auto a : ampl) std::cout << a << ",";
126  std::cout << std::endl;
127  }
128 
129  localAmplitudes.swap(ampl);
130  recordMinAffectedStrip=minA;
131  recordMaxAffectedStrip=maxA;
132  */
133 
134 }
135 
136 void
139  const StripGeomDetUnit& det,
140  std::vector<float>& localAmplitudes,
141  size_t& recordMinAffectedStrip,
142  size_t& recordMaxAffectedStrip,
143  const TrackerTopology *tTopo) const {
144 
145 
146  auto const & coupling = signalCoupling[typeOf(det,tTopo)];
147  const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
148  const int Nstrips = topology.nstrips();
149 
150 
151  const int NP = collection_points.size();
152  if(0==NP) return;
153 
154  constexpr int MaxN = 512;
155  // if NP too large split...
156 
157  for (int ip=0; ip<NP; ip+=MaxN) {
158  auto N = std::min(NP-ip,MaxN);
159 
160  count.dep(N);
161  float amplitude[N];
162  float chargePosition[N];
163  float chargeSpread[N];
164  int fromStrip[N];
165  int nStrip[N];
166 
167  // load not vectorize
168  //In strip coordinates:
169  for (int i=0; i!=N;++i) {
170  auto j = ip+i;
171  if (0==collection_points[j].amplitude()) count.zero();
172  chargePosition[i]=topology.strip(collection_points[j].position());
173  chargeSpread[i]= collection_points[j].sigma() / topology.localPitch(collection_points[j].position());
174  amplitude[i]=0.5f*collection_points[j].amplitude() / geVperElectron;
175  }
176 
177  // this vectorize
178  for (int i=0; i!=N;++i) {
179  fromStrip[i] = std::max( 0, int(std::floor( chargePosition[i] - Nsigma*chargeSpread[i])) );
180  nStrip[i] = std::min( Nstrips, int(std::ceil( chargePosition[i] + Nsigma*chargeSpread[i])) ) - fromStrip[i];
181  }
182  int tot=0;
183  for (int i=0; i!=N;++i) tot += nStrip[i];
184  tot+=N; // add last strip
185  count.val(tot);
186  float value[tot];
187 
188  // assign relative position (lower bound of strip) in value;
189  int kk=0;
190  for (int i=0; i!=N;++i) {
191  auto delta = 1.f/(std::sqrt(2.f)*chargeSpread[i]);
192  auto pos = delta*(float(fromStrip[i])-chargePosition[i]);
193  for (int j=0;j<=nStrip[i]; ++j)
194  value[kk++] = pos+float(j)*delta;
195  }
196  assert(kk==tot);
197 
198  // main loop fully vectorized
199  for (int k=0;k!=tot; ++k)
200  value[k] = approx_erf(value[k]);
201 
202  // saturate 0 & NStrips strip to 0 and 1???
203  kk=0;
204  for (int i=0; i!=N;++i) {
205  if (0 == fromStrip[i]) value[kk]=0;
206  kk+=nStrip[i];
207  if (Nstrips == fromStrip[i]+nStrip[i]) value[kk]=1.f;
208  ++kk;
209  }
210  assert(kk==tot);
211 
212  // compute integral over strip (lower bound becomes the value)
213  for (int k=0;k!=tot-1; ++k)
214  value[k]-=value[k+1]; // this is negative!
215 
216 
217  float charge[Nstrips]; for (int i=0;i!=Nstrips; ++i) charge[i]=0;
218  kk=0;
219  for (int i=0; i!=N;++i){
220  for (int j=0;j!=nStrip[i]; ++j)
221  charge[fromStrip[i]+j]-= amplitude[i]*value[kk++];
222  ++kk; // skip last "strip"
223  }
224  assert(kk==tot);
225 
226 
228  int minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
229  int sc = coupling.size();
230  for (int i=0;i!=Nstrips; ++i) {
231  int strip = i;
232  if (0==charge[i]) continue;
233  auto affectedFromStrip = std::max( 0, strip - sc + 1);
234  auto affectedUntilStrip = std::min(Nstrips, strip + sc);
235  for (auto affectedStrip=affectedFromStrip; affectedStrip < affectedUntilStrip; ++affectedStrip)
236  localAmplitudes[affectedStrip] += charge[i] * coupling[std::abs(affectedStrip - strip)] ;
237 
238  if( affectedFromStrip < minA ) minA = affectedFromStrip;
239  if( affectedUntilStrip > maxA ) maxA = affectedUntilStrip;
240  }
241  recordMinAffectedStrip=minA;
242  recordMaxAffectedStrip=maxA;
243  } // end loop ip
244 
245 }
246 
247 void
250  const StripGeomDetUnit& det,
251  std::vector<float>& localAmplitudes,
252  size_t& recordMinAffectedStrip,
253  size_t& recordMaxAffectedStrip,
254  const TrackerTopology *tTopo) const {
255 
256 
257  auto const & coupling = signalCoupling[typeOf(det,tTopo)];
258  const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
259  size_t Nstrips = topology.nstrips();
260 
261  if (!collection_points.empty()) count.dep(collection_points.size());
262 
263  for (auto signalpoint = collection_points.begin(); signalpoint != collection_points.end(); signalpoint++ ) {
264 
265  //In strip coordinates:
266  double chargePosition = topology.strip(signalpoint->position());
267  double chargeSpread = signalpoint->sigma() / topology.localPitch(signalpoint->position());
268 
269  size_t fromStrip = size_t(std::max( 0, int(std::floor( chargePosition - Nsigma*chargeSpread))));
270  size_t untilStrip = size_t(std::min( Nstrips, size_t(std::ceil( chargePosition + Nsigma*chargeSpread) )));
271 
272  count.str(std::max(0,int(untilStrip)-int(fromStrip)));
273  for (size_t strip = fromStrip; strip < untilStrip; strip++) {
274 
275  double chargeDepositedOnStrip = chargeDeposited( strip, Nstrips, signalpoint->amplitude() / geVperElectron, chargeSpread, chargePosition);
276 
277  size_t affectedFromStrip = size_t(std::max( 0, int(strip - coupling.size() + 1)));
278  size_t affectedUntilStrip = size_t(std::min( Nstrips, strip + coupling.size()) );
279  for (size_t affectedStrip = affectedFromStrip; affectedStrip < affectedUntilStrip; affectedStrip++) {
280  localAmplitudes.at( affectedStrip ) += chargeDepositedOnStrip * coupling.at(abs( affectedStrip - strip )) ;
281  }
282 
283  if( affectedFromStrip < recordMinAffectedStrip ) recordMinAffectedStrip = affectedFromStrip;
284  if( affectedUntilStrip > recordMaxAffectedStrip ) recordMaxAffectedStrip = affectedUntilStrip;
285  }
286  }
287 
288 }
289 
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
unsigned int tecRing(const DetId &id) const
ring id
#define abs(x)
Definition: mlp_lapack.h:159
std::vector< SignalPoint > collection_type
#define min(a, b)
Definition: mlp_lapack.h:161
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.
double charge(const std::vector< uint8_t > &Ampls)
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
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
virtual float strip(const LocalPoint &) const =0
float approx_erf(float x)
Definition: approx_erf.h:6
virtual float localPitch(const LocalPoint &) const =0
const T & max(const T &a, const T &b)
T sqrt(T t)
Definition: SSEVec.h:48
int j
Definition: DBlmapReader.cc:9
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:72
double f[11][100]
virtual StripGeomDetType & specificType() const
tuple conf
Definition: dbtoconf.py:185
int k[5][pyjets_maxn]
Definition: DetId.h:20
#define N
Definition: blowfish.cc:9
tuple cout
Definition: gather_cfg.py:121
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
#define constexpr
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