CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
SiTrivialInduceChargeOnStrips Class Reference

#include <SiTrivialInduceChargeOnStrips.h>

Inheritance diagram for SiTrivialInduceChargeOnStrips:
SiInduceChargeOnStrips

Public Member Functions

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 override
 
 SiTrivialInduceChargeOnStrips (const edm::ParameterSet &conf, double g)
 
 ~SiTrivialInduceChargeOnStrips () override
 
- Public Member Functions inherited from SiInduceChargeOnStrips
virtual ~SiInduceChargeOnStrips ()
 

Private Member Functions

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
 
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
 

Private Attributes

const float geVperElectron
 
const float Nsigma
 
const std::vector< std::vector< float > > signalCoupling
 

Detailed Description

Definition at line 9 of file SiTrivialInduceChargeOnStrips.h.

Constructor & Destructor Documentation

SiTrivialInduceChargeOnStrips::SiTrivialInduceChargeOnStrips ( const edm::ParameterSet conf,
double  g 
)

Definition at line 120 of file SiTrivialInduceChargeOnStrips.cc.

121  : signalCoupling(fillSignalCoupling(conf, Ntypes, type)), Nsigma(3.), geVperElectron(g) {}
type
Definition: HCALResponse.h:21
const std::vector< std::vector< float > > signalCoupling
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
SiTrivialInduceChargeOnStrips::~SiTrivialInduceChargeOnStrips ( )
inlineoverride

Definition at line 12 of file SiTrivialInduceChargeOnStrips.h.

References induce(), induceOriginal(), and induceVector().

12 {}

Member Function Documentation

void SiTrivialInduceChargeOnStrips::induce ( const SiChargeCollectionDrifter::collection_type collection_points,
const StripGeomDetUnit det,
std::vector< float > &  localAmplitudes,
size_t &  recordMinAffectedStrip,
size_t &  recordMaxAffectedStrip,
const TrackerTopology tTopo 
) const
overridevirtual

Implements SiInduceChargeOnStrips.

Definition at line 123 of file SiTrivialInduceChargeOnStrips.cc.

References induceVector().

Referenced by ~SiTrivialInduceChargeOnStrips().

128  {
129  induceVector(collection_points, det, localAmplitudes, recordMinAffectedStrip, recordMaxAffectedStrip, tTopo);
130 
131  /*
132  auto ominA=recordMinAffectedStrip, omaxA=recordMaxAffectedStrip;
133  std::vector<float> oampl(localAmplitudes);
134  induceOriginal(collection_points, det, oampl, ominA, omaxA, tTopo);
135 
136  // std::cout << "orig " << ominA << " " << omaxA << " ";
137  //for (auto a : oampl) std::cout << a << ",";
138  //std::cout << std::endl;
139 
140  auto minA=recordMinAffectedStrip, maxA=recordMaxAffectedStrip;
141  std::vector<float> ampl(localAmplitudes);
142  induceVector(collection_points, det, ampl, minA, maxA, tTopo);
143 
144  // std::cout << "vect " << minA << " " << maxA << " ";
145  //for (auto a : ampl) std::cout << a << ",";
146  //std::cout << std::endl;
147 
148  float diff=0;
149  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);}
150  if (diff> 1.e-4) {
151  std::cout << diff << std::endl;
152  std::cout << "orig " << ominA << " " << omaxA << " ";
153 // for (auto a : oampl) std::cout << a << ",";
154  std::cout << std::endl;
155  std::cout << "vect " << minA << " " << maxA << " ";
156 // for (auto a : ampl) std::cout << a << ",";
157  std::cout << std::endl;
158  }
159 
160  localAmplitudes.swap(ampl);
161  recordMinAffectedStrip=minA;
162  recordMaxAffectedStrip=maxA;
163  */
164 }
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
void SiTrivialInduceChargeOnStrips::induceOriginal ( const SiChargeCollectionDrifter::collection_type collection_points,
const StripGeomDetUnit det,
std::vector< float > &  localAmplitudes,
size_t &  recordMinAffectedStrip,
size_t &  recordMaxAffectedStrip,
const TrackerTopology tTopo 
) const
private

Definition at line 287 of file SiTrivialInduceChargeOnStrips.cc.

References funct::abs(), geVperElectron, StripTopology::localPitch(), SiStripPI::max, min(), Nsigma, StripTopology::nstrips(), StripGeomDetUnit::specificTopology(), StripTopology::strip(), and ecaldqm::topology().

Referenced by ~SiTrivialInduceChargeOnStrips().

292  {
293  auto const& coupling = signalCoupling[typeOf(det, tTopo)];
294  const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
295  size_t Nstrips = topology.nstrips();
296 
297  if (!collection_points.empty())
298  count.dep(collection_points.size());
299 
300  for (auto signalpoint = collection_points.begin(); signalpoint != collection_points.end(); signalpoint++) {
301  //In strip coordinates:
302  double chargePosition = topology.strip(signalpoint->position());
303  double chargeSpread = signalpoint->sigma() / topology.localPitch(signalpoint->position());
304 
305  size_t fromStrip = size_t(std::max(0, int(std::floor(chargePosition - Nsigma * chargeSpread))));
306  size_t untilStrip = size_t(std::min(Nstrips, size_t(std::ceil(chargePosition + Nsigma * chargeSpread))));
307 
308  count.str(std::max(0, int(untilStrip) - int(fromStrip)));
309  for (size_t strip = fromStrip; strip < untilStrip; strip++) {
310  double chargeDepositedOnStrip =
311  chargeDeposited(strip, Nstrips, signalpoint->amplitude() / geVperElectron, chargeSpread, chargePosition);
312 
313  size_t affectedFromStrip = size_t(std::max(0, int(strip - coupling.size() + 1)));
314  size_t affectedUntilStrip = size_t(std::min(Nstrips, strip + coupling.size()));
315  for (size_t affectedStrip = affectedFromStrip; affectedStrip < affectedUntilStrip; affectedStrip++) {
316  localAmplitudes.at(affectedStrip) +=
317  chargeDepositedOnStrip * coupling.at((size_t)std::abs((int)affectedStrip - (int)strip));
318  }
319 
320  if (affectedFromStrip < recordMinAffectedStrip)
321  recordMinAffectedStrip = affectedFromStrip;
322  if (affectedUntilStrip > recordMaxAffectedStrip)
323  recordMaxAffectedStrip = affectedUntilStrip;
324  }
325  }
326 }
CaloTopology const * topology(0)
virtual float strip(const LocalPoint &) const =0
const std::vector< std::vector< float > > signalCoupling
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
virtual int nstrips() const =0
virtual float localPitch(const LocalPoint &) const =0
void SiTrivialInduceChargeOnStrips::induceVector ( const SiChargeCollectionDrifter::collection_type collection_points,
const StripGeomDetUnit det,
std::vector< float > &  localAmplitudes,
size_t &  recordMinAffectedStrip,
size_t &  recordMaxAffectedStrip,
const TrackerTopology tTopo 
) const
private

include last strip

do crosstalk... (can be done better, most probably not worth)

Definition at line 166 of file SiTrivialInduceChargeOnStrips.cc.

References funct::abs(), CustomPhysics_cfi::amplitude, approx_erf(), ALCARECOTkAlJpsiMuMu_cff::charge, constexpr, delta, f, objects.autophobj::float, geVperElectron, mps_fire::i, gen::k, GetRecoTauVFromDQM_MC_cff::kk, StripTopology::localPitch(), SiStripPI::max, PFRecoTauDiscriminationByNProngs_cfi::MaxN, min(), N, Nsigma, StripTopology::nstrips(), position, SimDataFormats::CaloAnalysis::sc, StripGeomDetUnit::specificTopology(), mathSSE::sqrt(), StripTopology::strip(), and ecaldqm::topology().

Referenced by induce(), and ~SiTrivialInduceChargeOnStrips().

171  {
172  auto const& coupling = signalCoupling[typeOf(det, tTopo)];
173  const StripTopology& topology = dynamic_cast<const StripTopology&>(det.specificTopology());
174  const int Nstrips = topology.nstrips();
175 
176  if (Nstrips == 0)
177  return;
178 
179  const int NP = collection_points.size();
180  if (0 == NP)
181  return;
182 
183  constexpr int MaxN = 512;
184  // if NP too large split...
185 
186  for (int ip = 0; ip < NP; ip += MaxN) {
187  auto N = std::min(NP - ip, MaxN);
188 
189  count.dep(N);
190  float amplitude[N];
191  float chargePosition[N];
192  float chargeSpread[N];
193  int fromStrip[N];
194  int nStrip[N];
195 
196  // load not vectorize
197  //In strip coordinates:
198  for (int i = 0; i != N; ++i) {
199  auto j = ip + i;
200  if (0 == collection_points[j].amplitude())
201  count.zero();
202  chargePosition[i] = topology.strip(collection_points[j].position());
203  chargeSpread[i] = collection_points[j].sigma() / topology.localPitch(collection_points[j].position());
204  amplitude[i] = 0.5f * collection_points[j].amplitude() / geVperElectron;
205  }
206 
207  // this vectorize
208  for (int i = 0; i != N; ++i) {
209  fromStrip[i] = std::max(0, int(std::floor(chargePosition[i] - Nsigma * chargeSpread[i])));
210  nStrip[i] = std::min(Nstrips, int(std::ceil(chargePosition[i] + Nsigma * chargeSpread[i]))) - fromStrip[i];
211  }
212  int tot = 0;
213  for (int i = 0; i != N; ++i)
214  tot += nStrip[i];
215  tot += N; // add last strip
216  count.val(tot);
217  float value[tot];
218 
219  // assign relative position (lower bound of strip) in value;
220  int kk = 0;
221  for (int i = 0; i != N; ++i) {
222  auto delta = 1.f / (std::sqrt(2.f) * chargeSpread[i]);
223  auto pos = delta * (float(fromStrip[i]) - chargePosition[i]);
224 
225  // VI: before value[0] was not defined and value[tot] was filled
226  // to fix this the loop below was changed
227  for (int j = 0; j <= nStrip[i]; ++j) {
228  value[kk] = pos + float(j) * delta;
229  ++kk;
230  }
231  }
232  assert(kk == tot);
233 
234  // main loop fully vectorized
235  for (int k = 0; k != tot; ++k)
236  value[k] = approx_erf(value[k]);
237 
238  // saturate 0 & NStrips strip to 0 and 1???
239  kk = 0;
240  for (int i = 0; i != N; ++i) {
241  if (0 == fromStrip[i])
242  value[kk] = 0;
243  kk += nStrip[i];
244  if (Nstrips == fromStrip[i] + nStrip[i])
245  value[kk] = 1.f;
246  ++kk;
247  }
248  assert(kk == tot);
249 
250  // compute integral over strip (lower bound becomes the value)
251  for (int k = 0; k != tot - 1; ++k)
252  value[k] -= value[k + 1]; // this is negative!
253 
254  float charge[Nstrips];
255  for (int i = 0; i != Nstrips; ++i)
256  charge[i] = 0;
257  kk = 0;
258  for (int i = 0; i != N; ++i) {
259  for (int j = 0; j != nStrip[i]; ++j)
260  charge[fromStrip[i] + j] -= amplitude[i] * value[kk++];
261  ++kk; // skip last "strip"
262  }
263  assert(kk == tot);
264 
266  int minA = recordMinAffectedStrip, maxA = recordMaxAffectedStrip;
267  int sc = coupling.size();
268  for (int i = 0; i != Nstrips; ++i) {
269  int strip = i;
270  if (0 == charge[i])
271  continue;
272  auto affectedFromStrip = std::max(0, strip - sc + 1);
273  auto affectedUntilStrip = std::min(Nstrips, strip + sc);
274  for (auto affectedStrip = affectedFromStrip; affectedStrip < affectedUntilStrip; ++affectedStrip)
275  localAmplitudes[affectedStrip] += charge[i] * coupling[std::abs(affectedStrip - strip)];
276 
277  if (affectedFromStrip < minA)
278  minA = affectedFromStrip;
279  if (affectedUntilStrip > maxA)
280  maxA = affectedUntilStrip;
281  }
282  recordMinAffectedStrip = minA;
283  recordMaxAffectedStrip = maxA;
284  } // end loop ip
285 }
dbl * delta
Definition: mlp_gen.cc:36
CaloTopology const * topology(0)
virtual float strip(const LocalPoint &) const =0
const std::vector< std::vector< float > > signalCoupling
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
Definition: value.py:1
T min(T a, T b)
Definition: MathUtil.h:58
int k[5][pyjets_maxn]
#define N
Definition: blowfish.cc:9
virtual int nstrips() const =0
constexpr float approx_erf(float x)
Definition: approx_erf.h:6
virtual float localPitch(const LocalPoint &) const =0
static int position[264][3]
Definition: ReadPGInfo.cc:509
#define constexpr

Member Data Documentation

const float SiTrivialInduceChargeOnStrips::geVperElectron
private

Definition at line 38 of file SiTrivialInduceChargeOnStrips.h.

Referenced by induceOriginal(), and induceVector().

const float SiTrivialInduceChargeOnStrips::Nsigma
private

Definition at line 37 of file SiTrivialInduceChargeOnStrips.h.

Referenced by induceOriginal(), and induceVector().

const std::vector<std::vector<float> > SiTrivialInduceChargeOnStrips::signalCoupling
private

Definition at line 35 of file SiTrivialInduceChargeOnStrips.h.