CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
SiStripDigitizerAlgorithm Class Reference

#include <SiStripDigitizerAlgorithm.h>

Public Types

typedef float Amplitude
 
typedef
SiDigitalConverter::DigitalRawVecType 
DigitalRawVecType
 
typedef
SiDigitalConverter::DigitalVecType 
DigitalVecType
 
typedef std::map< int, float,
std::less< int > > 
hit_map_type
 
typedef
SiPileUpSignals::HitCounterToDigisMapType 
HitCounterToDigisMapType
 
typedef
SiPileUpSignals::HitToDigisMapType 
HitToDigisMapType
 

Public Member Functions

std::vector< StripDigiSimLinkmake_link ()
 
void run (edm::DetSet< SiStripDigi > &, edm::DetSet< SiStripRawDigi > &, const std::vector< std::pair< const PSimHit *, int > > &, StripGeomDetUnit *, GlobalVector, float, edm::ESHandle< SiStripGain > &, edm::ESHandle< SiStripThreshold > &, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripPedestals > &, edm::ESHandle< SiStripBadStrip > &)
 
void setParticleDataTable (const ParticleDataTable *pardt)
 
 SiStripDigitizerAlgorithm (const edm::ParameterSet &conf, CLHEP::HepRandomEngine &)
 
 ~SiStripDigitizerAlgorithm ()
 

Private Member Functions

void push_link (const DigitalVecType &, const HitToDigisMapType &, const HitCounterToDigisMapType &, const std::vector< double > &, unsigned int)
 
void push_link_raw (const DigitalRawVecType &, const HitToDigisMapType &, const HitCounterToDigisMapType &, const std::vector< double > &, unsigned int)
 

Private Attributes

bool APVSaturationFromHIP
 
double APVSaturationProb
 
bool BaselineShift
 
double cmnRMStec
 
double cmnRMStib
 
double cmnRMStid
 
double cmnRMStob
 
bool CommonModeNoise
 
edm::ParameterSet conf_
 
double cosmicShift
 
std::vector< double > detAmpl
 
DigitalVecType digis
 
size_t firstChannelWithSignal
 
double inefficiency
 
size_t lastChannelWithSignal
 
std::vector< StripDigiSimLinklink_coll
 
size_t localFirstChannel
 
size_t localLastChannel
 
std::vector< double > locAmpl
 
bool noise
 
int numStrips
 
const ParticleDataparticle
 
const ParticleDataTablepdt
 
bool peakMode
 
double pedOffset
 
DigitalRawVecType rawdigis
 
bool RealPedestals
 
CLHEP::HepRandomEngine & rndEngine
 
bool SingleStripNoise
 
int strip
 
double theElectronPerADC
 
int theFedAlgo
 
CLHEP::RandFlat * theFlatDistribution
 
SiTrivialDigitalConvertertheSiDigitalConverter
 
SiHitDigitizertheSiHitDigitizer
 
SiGaussianTailNoiseAddertheSiNoiseAdder
 
SiPileUpSignalstheSiPileUpSignals
 
SiStripFedZeroSuppressiontheSiZeroSuppress
 
double theThreshold
 
double theTOFCutForDeconvolution
 
double theTOFCutForPeak
 
double tofCut
 
bool zeroSuppression
 

Detailed Description

SiStripDigitizerAlgorithm converts hits to digis

Definition at line 39 of file SiStripDigitizerAlgorithm.h.

Member Typedef Documentation

Definition at line 46 of file SiStripDigitizerAlgorithm.h.

Definition at line 42 of file SiStripDigitizerAlgorithm.h.

Definition at line 41 of file SiStripDigitizerAlgorithm.h.

typedef std::map< int, float, std::less<int> > SiStripDigitizerAlgorithm::hit_map_type

Definition at line 45 of file SiStripDigitizerAlgorithm.h.

Definition at line 44 of file SiStripDigitizerAlgorithm.h.

Definition at line 43 of file SiStripDigitizerAlgorithm.h.

Constructor & Destructor Documentation

SiStripDigitizerAlgorithm::SiStripDigitizerAlgorithm ( const edm::ParameterSet conf,
CLHEP::HepRandomEngine &  eng 
)

Definition at line 18 of file SiStripDigitizerAlgorithm.cc.

References APVSaturationFromHIP, APVSaturationProb, BaselineShift, cmnRMStec, cmnRMStib, cmnRMStid, cmnRMStob, CommonModeNoise, conf_, cosmicShift, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), inefficiency, LogDebug, noise, peakMode, pedOffset, RealPedestals, rndEngine, SingleStripNoise, theElectronPerADC, theFedAlgo, theFlatDistribution, theSiDigitalConverter, theSiHitDigitizer, theSiNoiseAdder, theSiPileUpSignals, theSiZeroSuppress, theThreshold, theTOFCutForDeconvolution, theTOFCutForPeak, tofCut, and zeroSuppression.

18  :
19  conf_(conf),rndEngine(eng){
20  theThreshold = conf_.getParameter<double>("NoiseSigmaThreshold");
21  theFedAlgo = conf_.getParameter<int>("FedAlgorithm");
22  peakMode = conf_.getParameter<bool>("APVpeakmode");
23  theElectronPerADC = conf_.getParameter<double>( peakMode ? "electronPerAdcPeak" : "electronPerAdcDec" );
24  noise = conf_.getParameter<bool>("Noise");
25  zeroSuppression = conf_.getParameter<bool>("ZeroSuppression");
26  theTOFCutForPeak = conf_.getParameter<double>("TOFCutForPeak");
27  theTOFCutForDeconvolution = conf_.getParameter<double>("TOFCutForDeconvolution");
28  cosmicShift = conf_.getUntrackedParameter<double>("CosmicDelayShift");
29  inefficiency = conf_.getParameter<double>("Inefficiency");
30  RealPedestals = conf_.getParameter<bool>("RealPedestals");
31  SingleStripNoise = conf_.getParameter<bool>("SingleStripNoise");
32  CommonModeNoise = conf_.getParameter<bool>("CommonModeNoise");
33  BaselineShift = conf_.getParameter<bool>("BaselineShift");
34  APVSaturationFromHIP = conf_.getParameter<bool>("APVSaturationFromHIP");
35  APVSaturationProb = conf_.getParameter<double>("APVSaturationProb");
36  cmnRMStib = conf_.getParameter<double>("cmnRMStib");
37  cmnRMStob = conf_.getParameter<double>("cmnRMStob");
38  cmnRMStid = conf_.getParameter<double>("cmnRMStid");
39  cmnRMStec = conf_.getParameter<double>("cmnRMStec");
40  pedOffset = (unsigned int)conf_.getParameter<double>("PedestalsOffset");
41  if (peakMode) {
43  LogDebug("StripDigiInfo")<<"APVs running in peak mode (poor time resolution)";
44  } else {
46  LogDebug("StripDigiInfo")<<"APVs running in deconvolution mode (good time resolution)";
47  };
48 
52  theSiDigitalConverter = new SiTrivialDigitalConverter(theElectronPerADC);
54  theFlatDistribution = new CLHEP::RandFlat(rndEngine, 0., 1.);
55 
56 
57 }
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
CLHEP::HepRandomEngine & rndEngine
SiTrivialDigitalConverter * theSiDigitalConverter
SiStripFedZeroSuppression * theSiZeroSuppress
SiGaussianTailNoiseAdder * theSiNoiseAdder
SiStripDigitizerAlgorithm::~SiStripDigitizerAlgorithm ( )

Definition at line 59 of file SiStripDigitizerAlgorithm.cc.

References theFlatDistribution, theSiDigitalConverter, theSiHitDigitizer, theSiNoiseAdder, theSiPileUpSignals, and theSiZeroSuppress.

59  {
60  delete theSiHitDigitizer;
61  delete theSiPileUpSignals;
62  delete theSiNoiseAdder;
63  delete theSiDigitalConverter;
64  delete theSiZeroSuppress;
65  delete theFlatDistribution;
66  //delete particle;
67  //delete pdt;
68 }
SiTrivialDigitalConverter * theSiDigitalConverter
SiStripFedZeroSuppression * theSiZeroSuppress
SiGaussianTailNoiseAdder * theSiNoiseAdder

Member Function Documentation

std::vector<StripDigiSimLink> SiStripDigitizerAlgorithm::make_link ( )
inline

Definition at line 62 of file SiStripDigitizerAlgorithm.h.

References link_coll.

Referenced by SiStripDigitizer::produce().

62 { return link_coll; }
std::vector< StripDigiSimLink > link_coll
void SiStripDigitizerAlgorithm::push_link ( const DigitalVecType digis,
const HitToDigisMapType htd,
const HitCounterToDigisMapType hctd,
const std::vector< double > &  afterNoise,
unsigned int  detID 
)
private

Definition at line 297 of file SiStripDigitizerAlgorithm.cc.

References i, link_coll, and dtDQMClient_cfg::threshold.

Referenced by run().

301  {
302  link_coll.clear();
303  for ( DigitalVecType::const_iterator i=digis.begin(); i!=digis.end(); i++) {
304  // Instead of checking the validity of the links against the digis,
305  // let's loop over digis and push the corresponding link
306  HitToDigisMapType::const_iterator mi(htd.find(i->strip()));
307  if (mi == htd.end()) continue;
308  HitCounterToDigisMapType::const_iterator cmi(hctd.find(i->strip()));
309  std::map<const PSimHit *, Amplitude> totalAmplitudePerSimHit;
310  for (std::vector < std::pair < const PSimHit*, Amplitude > >::const_iterator simul =
311  (*mi).second.begin() ; simul != (*mi).second.end(); simul ++){
312  totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
313  }
314 
315  //--- include the noise as well
316  double totalAmplitude1 = afterNoise[(*mi).first];
317 
318  //--- digisimlink
319  int sim_counter=0;
320  for (std::map<const PSimHit *, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin();
321  iter != totalAmplitudePerSimHit.end(); iter++){
322  float threshold = 0.;
323  float fraction = (*iter).second/totalAmplitude1;
324  if ( fraction >= threshold) {
325  // Noise fluctuation could make fraction>1. Unphysical, set it by hand = 1.
326  if(fraction > 1.) fraction = 1.;
327  for (std::vector < std::pair < const PSimHit*, int > >::const_iterator
328  simcount = (*cmi).second.begin() ; simcount != (*cmi).second.end(); simcount ++){
329  if((*iter).first == (*simcount).first) sim_counter = (*simcount).second;
330  }
331  link_coll.push_back(StripDigiSimLink( (*mi).first, //channel
332  ((*iter).first)->trackId(), //simhit trackId
333  sim_counter, //simhit counter
334  ((*iter).first)->eventId(), //simhit eventId
335  fraction)); //fraction
336  }
337  }
338  }
339 }
int i
Definition: DBlmapReader.cc:9
std::vector< StripDigiSimLink > link_coll
void SiStripDigitizerAlgorithm::push_link_raw ( const DigitalRawVecType digis,
const HitToDigisMapType htd,
const HitCounterToDigisMapType hctd,
const std::vector< double > &  afterNoise,
unsigned int  detID 
)
private

Definition at line 341 of file SiStripDigitizerAlgorithm.cc.

References i, link_coll, and dtDQMClient_cfg::threshold.

Referenced by run().

345  {
346  link_coll.clear();
347  int nstrip = -1;
348  for ( DigitalRawVecType::const_iterator i=digis.begin(); i!=digis.end(); i++) {
349  nstrip++;
350  // Instead of checking the validity of the links against the digis,
351  // let's loop over digis and push the corresponding link
352  HitToDigisMapType::const_iterator mi(htd.find(nstrip));
353  HitCounterToDigisMapType::const_iterator cmi(hctd.find(nstrip));
354  if (mi == htd.end()) continue;
355  std::map<const PSimHit *, Amplitude> totalAmplitudePerSimHit;
356  for (std::vector < std::pair < const PSimHit*, Amplitude > >::const_iterator simul =
357  (*mi).second.begin() ; simul != (*mi).second.end(); simul ++){
358  totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
359  }
360 
361  //--- include the noise as well
362  double totalAmplitude1 = afterNoise[(*mi).first];
363 
364  //--- digisimlink
365  int sim_counter_raw=0;
366  for (std::map<const PSimHit *, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin();
367  iter != totalAmplitudePerSimHit.end(); iter++){
368  float threshold = 0.;
369  float fraction = (*iter).second/totalAmplitude1;
370  if (fraction >= threshold) {
371  //Noise fluctuation could make fraction>1. Unphysical, set it by hand.
372  if(fraction >1.) fraction = 1.;
373  //add counter information
374  for (std::vector < std::pair < const PSimHit*, int > >::const_iterator
375  simcount = (*cmi).second.begin() ; simcount != (*cmi).second.end(); simcount ++){
376  if((*iter).first == (*simcount).first) sim_counter_raw = (*simcount).second;
377  }
378  link_coll.push_back(StripDigiSimLink( (*mi).first, //channel
379  ((*iter).first)->trackId(), //simhit trackId
380  sim_counter_raw, //simhit counter
381  ((*iter).first)->eventId(), //simhit eventId
382  fraction)); //fraction
383  }
384  }
385  }
386 }
int i
Definition: DBlmapReader.cc:9
std::vector< StripDigiSimLink > link_coll
void SiStripDigitizerAlgorithm::run ( edm::DetSet< SiStripDigi > &  outdigi,
edm::DetSet< SiStripRawDigi > &  outrawdigi,
const std::vector< std::pair< const PSimHit *, int > > &  input,
StripGeomDetUnit det,
GlobalVector  bfield,
float  langle,
edm::ESHandle< SiStripGain > &  gainHandle,
edm::ESHandle< SiStripThreshold > &  thresholdHandle,
edm::ESHandle< SiStripNoises > &  noiseHandle,
edm::ESHandle< SiStripPedestals > &  pedestalHandle,
edm::ESHandle< SiStripBadStrip > &  deadChannelHandle 
)

Definition at line 73 of file SiStripDigitizerAlgorithm.cc.

References SiPileUpSignals::add(), SiGaussianTailNoiseAdder::addBaselineShift(), SiGaussianTailNoiseAdder::addCMNoise(), SiGaussianTailNoiseAdder::addNoise(), SiGaussianTailNoiseAdder::addNoiseVR(), SiGaussianTailNoiseAdder::addPedestals(), APVSaturationFromHIP, APVSaturationProb, BaselineShift, DeDxDiscriminatorTools::charge(), cmnRMStec, cmnRMStib, cmnRMStid, cmnRMStob, CommonModeNoise, SiTrivialDigitalConverter::convert(), SiTrivialDigitalConverter::convertRaw(), cosmicShift, gather_cfg::cout, edm::DetSet< T >::data, detAmpl, digis, SiPileUpSignals::dumpCounterLink(), SiPileUpSignals::dumpLink(), firstChannelWithSignal, SiStripBadStrip::data::firstStrip, GeomDet::geographicalId(), inefficiency, LaserDQM_cfg::input, lastChannelWithSignal, localFirstChannel, localLastChannel, locAmpl, mag(), noise, NULL, numStrips, particle, pdt, pedOffset, SiHitDigitizer::processHit(), push_link(), push_link_raw(), SiStripBadStrip::data::range, rawdigis, DetId::rawId(), RealPedestals, SiPileUpSignals::reset(), SingleStripNoise, StripGeomDetUnit::specificTopology(), strip, DetId::subdetId(), SiStripFedZeroSuppression::suppress(), GeomDet::surface(), theElectronPerADC, theFlatDistribution, theSiDigitalConverter, theSiHitDigitizer, theSiNoiseAdder, theSiPileUpSignals, theSiZeroSuppress, tofCut, Surface::toGlobal(), and zeroSuppression.

Referenced by SiStripDigitizer::produce().

83  {
85  unsigned int detID = det->geographicalId().rawId();
86  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
87  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
88  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
89  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detID);
90  numStrips = (det->specificTopology()).nstrips();
91 
92  //stroing the bad stip of the the module. the module is not removed but just signal put to 0
93  std::vector<bool> badChannels;
94  badChannels.clear();
95  badChannels.insert(badChannels.begin(),numStrips,false);
97  for(SiStripBadStrip::ContainerIterator it=detBadStripRange.first;it!=detBadStripRange.second;++it){
98  fs=deadChannelHandle->decode(*it);
99  for(int strip = fs.firstStrip; strip <fs.firstStrip+fs.range; ++strip )badChannels[strip] = true;
100  }
101 
102 
103  // local amplitude of detector channels (from processed PSimHit)
104 // locAmpl.clear();
105  detAmpl.clear();
106 // locAmpl.insert(locAmpl.begin(),numStrips,0.);
107  // total amplitude of detector channels
108  detAmpl.insert(detAmpl.begin(),numStrips,0.);
109 
112 
113  // First: loop on the SimHits
114  std::vector<std::pair<const PSimHit*, int > >::const_iterator simHitIter = input.begin();
115  std::vector<std::pair<const PSimHit*, int > >::const_iterator simHitIterEnd = input.end();
116  if(theFlatDistribution->fire()>inefficiency) {
117  for (;simHitIter != simHitIterEnd; ++simHitIter) {
118  locAmpl.clear();
119  locAmpl.insert(locAmpl.begin(),numStrips,0.);
120  // check TOF
121  if ( std::fabs( ((*simHitIter).first)->tof() - cosmicShift - det->surface().toGlobal(((*simHitIter).first)->localPosition()).mag()/30.) < tofCut && ((*simHitIter).first)->energyLoss()>0) {
123  localLastChannel = 0;
124  // process the hit
125  theSiHitDigitizer->processHit(((*simHitIter).first),*det,bfield,langle, locAmpl, localFirstChannel, localLastChannel);
126 
127  //APV Killer to simulate HIP effect
128  //------------------------------------------------------
129 
131  int pdg_id = ((*simHitIter).first)->particleType();
132  particle = pdt->particle(pdg_id);
133  if(particle != NULL){
134  float charge = particle->charge();
135  bool isHadron = particle->isHadron();
136  if(charge!=0 && isHadron){
138  int FirstAPV = localFirstChannel/128;
139  int LastAPV = localLastChannel/128;
140  std::cout << "-------------------HIP--------------" << std::endl;
141  std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
142  for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) badChannels[strip] = true; //doing like that I remove the signal information only after the
143  //stip that got the HIP but it remains the signal of the previous
144  //one. I'll make a further loop to remove all signal
145 
146  }
147  }
148  }
149  }
150 
151 
152  theSiPileUpSignals->add(locAmpl, localFirstChannel, localLastChannel, ((*simHitIter).first), (*simHitIter).second);
153 
154  // sum signal on strips
155  for (size_t iChannel=localFirstChannel; iChannel<localLastChannel; iChannel++) {
156  if(locAmpl[iChannel]>0.) {
157  //if(!badChannels[iChannel]) detAmpl[iChannel]+=locAmpl[iChannel];
158  //locAmpl[iChannel]=0;
159  detAmpl[iChannel]+=locAmpl[iChannel];
160  }
161  }
164  }
165  }
166  }
167 
168  //removing signal from the dead (and HIP effected) strips
169  for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] =0;
170 
173 
174  if(zeroSuppression){
175  if(noise){
176  int RefStrip = int(numStrips/2.);
177  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
178  RefStrip++;
179  }
180  if(RefStrip<numStrips){
181  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
182  float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
184  }
185  }
186  digis.clear();
187  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
188  push_link(digis, theLink, theCounterLink, detAmpl,detID);
189  outdigi.data = digis;
190  }
191 
192 
193  if(!zeroSuppression){
194  //if(noise){
195  // the constant pedestal offset is needed because
196  // negative adc counts are not allowed in case
197  // Pedestal and CMN subtraction is performed.
198  // The pedestal value read from the conditions
199  // is pedValue and after the pedestal subtraction
200  // the baseline is zero. The Common Mode Noise
201  // is not subtracted from the negative adc counts
202  // channels. Adding pedOffset the baseline is set
203  // to pedOffset after pedestal subtraction and CMN
204  // is subtracted to all the channels since none of
205  // them has negative adc value. The pedOffset is
206  // treated as a constant component in the CMN
207  // estimation and subtracted as CMN.
208 
209 
210  //calculating the charge deposited on each APV and subtracting the shift
211  //------------------------------------------------------
212  if(BaselineShift){
214  }
215 
216  //Adding the strip noise
217  //------------------------------------------------------
218  if(noise){
219  std::vector<float> noiseRMSv;
220  noiseRMSv.clear();
221  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
222 
223  if(SingleStripNoise){
224  for(int strip=0; strip< numStrips; ++strip){
225  if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
226  }
227 
228  } else {
229  int RefStrip = 0; //int(numStrips/2.);
230  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
231  RefStrip++;
232  }
233  if(RefStrip<numStrips){
234  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
235  for(int strip=0; strip< numStrips; ++strip){
236  if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
237  }
238  }
239  }
240 
241  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv);
242  }
243 
244  //adding the CMN
245  //------------------------------------------------------
246  if(CommonModeNoise){
247  float cmnRMS = 0.;
248  DetId detId(detID);
249  uint32_t SubDet = detId.subdetId();
250  if(SubDet==3){
251  cmnRMS = cmnRMStib;
252  }else if(SubDet==4){
253  cmnRMS = cmnRMStid;
254  }else if(SubDet==5){
255  cmnRMS = cmnRMStob;
256  }else if(SubDet==6){
257  cmnRMS = cmnRMStec;
258  }
259  cmnRMS *= theElectronPerADC;
260  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels);
261  }
262 
263 
264  //Adding the pedestals
265  //------------------------------------------------------
266 
267  std::vector<float> vPeds;
268  vPeds.clear();
269  vPeds.insert(vPeds.begin(),numStrips,0.);
270 
271  if(RealPedestals){
272  for(int strip=0; strip< numStrips; ++strip){
273  if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
274  }
275  } else {
276  for(int strip=0; strip< numStrips; ++strip){
277  if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
278  }
279  }
280 
282 
283 
284  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
285  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
286  //}else{
287 
288  rawdigis.clear();
289  rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
290  push_link_raw(rawdigis, theLink, theCounterLink, detAmpl,detID);
291  outrawdigi.data = rawdigis;
292 
293  //}
294  }
295 }
unsigned short range
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:78
DigitalVecType convert(const std::vector< double > &, edm::ESHandle< SiStripGain > &, unsigned int detid)
void addCMNoise(std::vector< double > &, float, std::vector< bool > &)
const HitCounterToDigisMapType & dumpCounterLink() const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< unsigned int >::const_iterator ContainerIterator
#define NULL
Definition: scimark2.h:8
std::pair< ContainerIterator, ContainerIterator > Range
void addBaselineShift(std::vector< double > &, std::vector< bool > &)
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
double charge(const std::vector< uint8_t > &Ampls)
void addNoiseVR(std::vector< double > &, std::vector< float > &)
void push_link_raw(const DigitalRawVecType &, const HitToDigisMapType &, const HitCounterToDigisMapType &, const std::vector< double > &, unsigned int)
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
void push_link(const DigitalVecType &, const HitToDigisMapType &, const HitCounterToDigisMapType &, const std::vector< double > &, unsigned int)
void suppress(const std::vector< SiStripDigi > &, std::vector< SiStripDigi > &, const uint32_t &, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripThreshold > &)
std::pair< ContainerIterator, ContainerIterator > Range
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:72
std::map< int, std::vector< std::pair< const PSimHit *, int > >, std::less< int > > HitCounterToDigisMapType
std::map< int, std::vector< std::pair< const PSimHit *, Amplitude > >, std::less< int > > HitToDigisMapType
void processHit(const PSimHit *, const StripGeomDetUnit &, GlobalVector, float, std::vector< double > &, size_t &, size_t &)
Definition: DetId.h:20
unsigned short firstStrip
SiTrivialDigitalConverter * theSiDigitalConverter
void addNoise(std::vector< double > &, size_t &, size_t &, int, float)
collection_type data
Definition: DetSet.h:79
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
std::pair< ContainerIterator, ContainerIterator > Range
SiStripFedZeroSuppression * theSiZeroSuppress
SiGaussianTailNoiseAdder * theSiNoiseAdder
tuple cout
Definition: gather_cfg.py:121
virtual void add(const std::vector< double > &locAmpl, const size_t &firstChannelWithSignal, const size_t &lastChannelWithSignal, const PSimHit *hit, const int &counter)
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:41
DigitalRawVecType convertRaw(const std::vector< double > &, edm::ESHandle< SiStripGain > &, unsigned int detid)
void addPedestals(std::vector< double > &, std::vector< float > &)
const ParticleDataTable * pdt
const HitToDigisMapType & dumpLink() const
void SiStripDigitizerAlgorithm::setParticleDataTable ( const ParticleDataTable pardt)
inline

Definition at line 65 of file SiStripDigitizerAlgorithm.h.

References pdt, SiHitDigitizer::setParticleDataTable(), and theSiHitDigitizer.

Referenced by SiStripDigitizer::produce().

65  {
67  pdt= pardt;
68  }
void setParticleDataTable(const ParticleDataTable *pdt)
const ParticleDataTable * pdt

Member Data Documentation

bool SiStripDigitizerAlgorithm::APVSaturationFromHIP
private

Definition at line 85 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::APVSaturationProb
private

Definition at line 78 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

bool SiStripDigitizerAlgorithm::BaselineShift
private

Definition at line 84 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::cmnRMStec
private

Definition at line 77 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::cmnRMStib
private

Definition at line 74 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::cmnRMStid
private

Definition at line 76 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::cmnRMStob
private

Definition at line 75 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

bool SiStripDigitizerAlgorithm::CommonModeNoise
private

Definition at line 83 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

edm::ParameterSet SiStripDigitizerAlgorithm::conf_
private

Definition at line 71 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::cosmicShift
private

Definition at line 96 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

std::vector<double> SiStripDigitizerAlgorithm::detAmpl
private

Definition at line 108 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

DigitalVecType SiStripDigitizerAlgorithm::digis
private

Definition at line 120 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

size_t SiStripDigitizerAlgorithm::firstChannelWithSignal
private

Definition at line 100 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

double SiStripDigitizerAlgorithm::inefficiency
private

Definition at line 97 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

size_t SiStripDigitizerAlgorithm::lastChannelWithSignal
private

Definition at line 101 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

std::vector<StripDigiSimLink> SiStripDigitizerAlgorithm::link_coll
private

Definition at line 122 of file SiStripDigitizerAlgorithm.h.

Referenced by make_link(), push_link(), and push_link_raw().

size_t SiStripDigitizerAlgorithm::localFirstChannel
private

Definition at line 102 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

size_t SiStripDigitizerAlgorithm::localLastChannel
private

Definition at line 103 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

std::vector<double> SiStripDigitizerAlgorithm::locAmpl
private

Definition at line 106 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

bool SiStripDigitizerAlgorithm::noise
private

Definition at line 80 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

int SiStripDigitizerAlgorithm::numStrips
private

Definition at line 92 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

const ParticleData* SiStripDigitizerAlgorithm::particle
private

Definition at line 111 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

const ParticleDataTable* SiStripDigitizerAlgorithm::pdt
private

Definition at line 110 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and setParticleDataTable().

bool SiStripDigitizerAlgorithm::peakMode
private

Definition at line 79 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::pedOffset
private

Definition at line 98 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

DigitalRawVecType SiStripDigitizerAlgorithm::rawdigis
private

Definition at line 121 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

bool SiStripDigitizerAlgorithm::RealPedestals
private

Definition at line 81 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

CLHEP::HepRandomEngine& SiStripDigitizerAlgorithm::rndEngine
private

Definition at line 118 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

bool SiStripDigitizerAlgorithm::SingleStripNoise
private

Definition at line 82 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

int SiStripDigitizerAlgorithm::strip
private

Definition at line 93 of file SiStripDigitizerAlgorithm.h.

Referenced by run().

double SiStripDigitizerAlgorithm::theElectronPerADC
private

Definition at line 72 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

int SiStripDigitizerAlgorithm::theFedAlgo
private

Definition at line 87 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

CLHEP::RandFlat* SiStripDigitizerAlgorithm::theFlatDistribution
private
SiTrivialDigitalConverter* SiStripDigitizerAlgorithm::theSiDigitalConverter
private
SiHitDigitizer* SiStripDigitizerAlgorithm::theSiHitDigitizer
private
SiGaussianTailNoiseAdder* SiStripDigitizerAlgorithm::theSiNoiseAdder
private
SiPileUpSignals* SiStripDigitizerAlgorithm::theSiPileUpSignals
private
SiStripFedZeroSuppression* SiStripDigitizerAlgorithm::theSiZeroSuppress
private
double SiStripDigitizerAlgorithm::theThreshold
private

Definition at line 73 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::theTOFCutForDeconvolution
private

Definition at line 90 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::theTOFCutForPeak
private

Definition at line 89 of file SiStripDigitizerAlgorithm.h.

Referenced by SiStripDigitizerAlgorithm().

double SiStripDigitizerAlgorithm::tofCut
private

Definition at line 91 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().

bool SiStripDigitizerAlgorithm::zeroSuppression
private

Definition at line 88 of file SiStripDigitizerAlgorithm.h.

Referenced by run(), and SiStripDigitizerAlgorithm().