CMS 3D CMS Logo

SiStripRawProcessingAlgorithms.cc
Go to the documentation of this file.
2 
12 #include <memory>
13 
14 
16  std::unique_ptr<SiStripPedestalsSubtractor> ped,
17  std::unique_ptr<SiStripCommonModeNoiseSubtractor> cmn,
18  std::unique_ptr<SiStripFedZeroSuppression> zs,
19  std::unique_ptr<SiStripAPVRestorer> res,
20  bool doAPVRest, bool useCMMap)
21 : subtractorPed(std::move(ped)),
22  subtractorCMN(std::move(cmn)),
23  suppressor(std::move(zs)),
24  restorer(std::move(res)),
25  doAPVRestore(doAPVRest),
26  useCMMeanMap(useCMMap)
27 {}
28 
29 
31 {
32  subtractorPed->init(es);
33  subtractorCMN->init(es);
34  suppressor->init(es);
35  if ( restorer.get() ) restorer->init(es);
36 
38  es.get<TrackerDigiGeometryRecord>().get(tracker);
39  trGeo = tracker.product();
40 }
41 
43 {
44  initialize(es);
45  if ( restorer.get() && doAPVRestore&&useCMMeanMap ) restorer->loadMeanCMMap(e);
46 }
47 
63 //IMPORTANT: don't forget the conversion from hybrids on the bad APVs (*2 -1024)
64 uint16_t SiStripRawProcessingAlgorithms::suppressHybridData(uint32_t id, uint16_t firstAPV, digivector_t& procRawDigis, edm::DetSet<SiStripDigi>& suppressedDigis)
65 {
66  digivector_t procRawDigisPedSubtracted(procRawDigis);
67 
68  subtractorCMN->subtract(id, firstAPV, procRawDigis);
69 
70  const auto nAPVFlagged = restorer->inspectAndRestore(id, firstAPV, procRawDigisPedSubtracted, procRawDigis, subtractorCMN->getAPVsCM());
71 
72  const std::vector<bool>& apvf = getAPVFlags();
73  const std::size_t nAPVs = procRawDigis.size()/128;
74  for ( uint16_t iAPV = firstAPV; iAPV < firstAPV+nAPVs; ++iAPV ) {
75  if ( apvf[iAPV] ) {
76  const auto firstDigiIt = std::begin(procRawDigis)+128*(iAPV-firstAPV);
77  std::vector<int16_t> singleAPVdigi(firstDigiIt, firstDigiIt+128);
78  suppressor->suppress(singleAPVdigi, iAPV, suppressedDigis);
79  } else {
80  for ( uint16_t i = 0; i < 128; ++i ) {
81  const int16_t digi = procRawDigisPedSubtracted[128*(iAPV-firstAPV)+i];
82  if ( digi > 0 ) { suppressedDigis.push_back(SiStripDigi(iAPV*128+i, suppressor->truncate(digi))); }
83  }
84  }
85  }
86 
87  return nAPVFlagged;
88 }
89 
103 {
104  convertHybridDigiToRawDigiVector(hybridDigis, rawDigis);
105  return suppressHybridData(hybridDigis.id, 0, rawDigis, suppressedDigis);
106 }
107 
119 {
120  const auto stripModuleGeom = dynamic_cast<const StripGeomDetUnit*>(trGeo->idToDetUnit(inDigis.id));
121  const std::size_t nStrips = stripModuleGeom->specificTopology().nstrips();
122  const std::size_t nAPVs = nStrips/128;
123 
124  rawDigis.assign(nStrips, 0);
125  std::vector<uint16_t> stripsPerAPV(nAPVs, 0);
126 
127  for ( SiStripDigi digi : inDigis ) {
128  rawDigis[digi.strip()] = digi.adc();
129  ++stripsPerAPV[digi.strip()/128];
130  }
131 
132  for ( uint16_t iAPV = 0; iAPV < nAPVs; ++iAPV ) {
133  if ( stripsPerAPV[iAPV] > 64 ) {
134  for ( uint16_t strip = iAPV*128; strip < (iAPV+1)*128; ++strip )
135  rawDigis[strip] = rawDigis[strip] * 2 - 1024;
136  }
137  }
138 }
139 
153 {
154  subtractorPed->subtract(id, firstAPV*128, procRawDigis);
155  return suppressProcessedRawData(id, firstAPV, procRawDigis, output);
156 }
157 
169 {
170  digivector_t rawdigis; rawdigis.reserve(rawDigis.size());
171  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis),
172  [] ( SiStripRawDigi digi ) { return digi.adc(); } );
173  return suppressVirginRawData(rawDigis.id, 0, rawdigis, output);
174 }
175 
189 {
190  digivector_t procRawDigisPedSubtracted ;
191 
192  int16_t nAPVFlagged =0;
193  if ( doAPVRestore )
194  procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
195  subtractorCMN->subtract(id, firstAPV, procRawDigis);
196  if ( doAPVRestore )
197  nAPVFlagged = restorer->inspectAndRestore(id, firstAPV, procRawDigisPedSubtracted, procRawDigis, subtractorCMN->getAPVsCM());
198  suppressor->suppress(procRawDigis, firstAPV, output);
199  return nAPVFlagged;
200 }
201 
213 {
214  digivector_t rawdigis; rawdigis.reserve(rawDigis.size());
215  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis),
216  [] ( SiStripRawDigi digi ) { return digi.adc(); } );
217  return suppressProcessedRawData(rawDigis.id, 0, rawdigis, output);
218 }
219 
236 {
237  digivector_t procRawDigisPedSubtracted;
238 
239  for ( auto& digi : procRawDigis ) { digi += 1024; } // adding one MSB
240 
241  subtractorPed->subtract(id, firstAPV*128, procRawDigis); // all strips are pedestals subtracted
242 
243  for ( auto& digi : procRawDigis ) { digi /= 2; }
244 
245  procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
246 
247  subtractorCMN->subtract(id, firstAPV, procRawDigis);
248 
249  const auto nAPVFlagged = restorer->inspect(id, firstAPV, procRawDigis, subtractorCMN->getAPVsCM());
250 
251  for ( auto& digi : procRawDigis ) { digi *= 2; }
252 
253  const std::vector<bool>& apvf = getAPVFlags();
254  const std::size_t nAPVs = procRawDigis.size()/128;
255  for ( uint16_t iAPV = firstAPV; iAPV < nAPVs+firstAPV; ++iAPV ) {
256  if ( apvf[iAPV] ) {
257  //GB 23/6/08: truncation should be done at the very beginning
258  for ( uint16_t i = 0; i < 128; ++i ) {
259  const int16_t digi = procRawDigisPedSubtracted[128*(iAPV-firstAPV)+i];
260  output.push_back(SiStripDigi(128*iAPV+i, ( digi < 0 ? 0 : suppressor->truncate(digi) )));
261  }
262  } else {
263  const auto firstDigiIt = std::begin(procRawDigis)+128*(iAPV-firstAPV);
264  std::vector<int16_t> singleAPVdigi(firstDigiIt, firstDigiIt+128);
265  suppressor->suppress(singleAPVdigi, iAPV, output);
266  }
267  }
268 
269  return nAPVFlagged;
270 }
271 
285 {
286  digivector_t rawdigis; rawdigis.reserve(rawDigis.size());
287  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis),
288  [] ( SiStripRawDigi digi ) { return digi.adc(); } );
289  return convertVirginRawToHybrid(rawDigis.id, 0, rawdigis, suppressedDigis);
290 }
const uint16_t & adc() const
uint16_t suppressVirginRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
void push_back(const T &t)
Definition: DetSet.h:68
const std::unique_ptr< SiStripPedestalsSubtractor > subtractorPed
uint16_t suppressProcessedRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
SiStripAPVRestorer::digivector_t digivector_t
uint16_t convertVirginRawToHybrid(uint32_t detId, uint16_t firstAPV, digivector_t &inDigis, edm::DetSet< SiStripDigi > &rawDigis)
SiStripRawProcessingAlgorithms(std::unique_ptr< SiStripPedestalsSubtractor > ped, std::unique_ptr< SiStripCommonModeNoiseSubtractor > cmn, std::unique_ptr< SiStripFedZeroSuppression > zs, std::unique_ptr< SiStripAPVRestorer > res, bool doAPVRest, bool useCMMap)
void convertHybridDigiToRawDigiVector(const edm::DetSet< SiStripDigi > &inDigis, digivector_t &rawDigis)
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
const std::unique_ptr< SiStripAPVRestorer > restorer
Definition: Electron.h:6
size_type size() const
Definition: DetSet.h:63
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
#define end
Definition: vmac.h:39
A Digi for the silicon strip detector, containing both strip and adc information, and suitable for st...
Definition: SiStripDigi.h:12
const std::unique_ptr< SiStripFedZeroSuppression > suppressor
virtual int nstrips() const =0
const std::vector< bool > & getAPVFlags() const
#define begin
Definition: vmac.h:32
det_id_type id
Definition: DetSet.h:77
T get() const
Definition: EventSetup.h:62
uint16_t suppressHybridData(const edm::DetSet< SiStripDigi > &inDigis, edm::DetSet< SiStripDigi > &suppressedDigis, digivector_t &rawDigis)
void initialize(const edm::EventSetup &)
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511