CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | Friends
SiStripRawProcessingAlgorithms Class Reference

#include <SiStripRawProcessingAlgorithms.h>

Public Types

using digivector_t = SiStripAPVRestorer::digivector_t
 

Public Member Functions

uint16_t convertVirginRawToHybrid (uint32_t detId, uint16_t firstAPV, digivector_t &inDigis, edm::DetSet< SiStripDigi > &rawDigis)
 
uint16_t convertVirginRawToHybrid (const edm::DetSet< SiStripRawDigi > &rawDigis, edm::DetSet< SiStripDigi > &suppressedDigis)
 
const std::vector< bool > & getAPVFlags () const
 
const SiStripAPVRestorer::medians_tgetAPVsCM () const
 
const SiStripAPVRestorer::baselinemap_tgetBaselineMap () const
 
const std::map< uint16_t, SiStripAPVRestorer::digimap_t > & getSmoothedPoints () const
 
void initialize (const edm::EventSetup &)
 
void initialize (const edm::EventSetup &, const edm::Event &)
 
uint16_t suppressHybridData (const uint16_t maxStrip, const edm::DetSet< SiStripDigi > &inDigis, edm::DetSet< SiStripDigi > &suppressedDigis, uint16_t firstAPV=0)
 
uint16_t suppressProcessedRawData (uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
 
uint16_t suppressProcessedRawData (const edm::DetSet< SiStripRawDigi > &rawDigis, edm::DetSet< SiStripDigi > &output)
 
uint16_t suppressVirginRawData (uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
 
uint16_t suppressVirginRawData (const edm::DetSet< SiStripRawDigi > &rawDigis, edm::DetSet< SiStripDigi > &output)
 

Public Attributes

const std::unique_ptr< SiStripAPVRestorerrestorer
 
const std::unique_ptr< SiStripCommonModeNoiseSubtractorsubtractorCMN
 
const std::unique_ptr< SiStripPedestalsSubtractorsubtractorPed
 
const std::unique_ptr< SiStripFedZeroSuppressionsuppressor
 

Private Member Functions

 SiStripRawProcessingAlgorithms (edm::ConsumesCollector iC, std::unique_ptr< SiStripPedestalsSubtractor > ped, std::unique_ptr< SiStripCommonModeNoiseSubtractor > cmn, std::unique_ptr< SiStripFedZeroSuppression > zs, std::unique_ptr< SiStripAPVRestorer > res, bool doAPVRest, bool useCMMap)
 

Private Attributes

const bool doAPVRestore
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtkGeomToken_
 
const TrackerGeometrytrGeo
 
const bool useCMMeanMap
 

Friends

class SiStripRawProcessingFactory
 

Detailed Description

Definition at line 15 of file SiStripRawProcessingAlgorithms.h.

Member Typedef Documentation

◆ digivector_t

Definition at line 19 of file SiStripRawProcessingAlgorithms.h.

Constructor & Destructor Documentation

◆ SiStripRawProcessingAlgorithms()

SiStripRawProcessingAlgorithms::SiStripRawProcessingAlgorithms ( edm::ConsumesCollector  iC,
std::unique_ptr< SiStripPedestalsSubtractor ped,
std::unique_ptr< SiStripCommonModeNoiseSubtractor cmn,
std::unique_ptr< SiStripFedZeroSuppression zs,
std::unique_ptr< SiStripAPVRestorer res,
bool  doAPVRest,
bool  useCMMap 
)
private

Definition at line 14 of file SiStripRawProcessingAlgorithms.cc.

23  suppressor(std::move(zs)),
25  doAPVRestore(doAPVRest),
26  useCMMeanMap(useCMMap),
const std::unique_ptr< SiStripPedestalsSubtractor > subtractorPed
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const std::unique_ptr< SiStripAPVRestorer > restorer
Definition: Electron.h:6
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
const std::unique_ptr< SiStripFedZeroSuppression > suppressor
def move(src, dest)
Definition: eostools.py:511

Member Function Documentation

◆ convertVirginRawToHybrid() [1/2]

uint16_t SiStripRawProcessingAlgorithms::convertVirginRawToHybrid ( uint32_t  id,
uint16_t  firstAPV,
digivector_t procRawDigis,
edm::DetSet< SiStripDigi > &  output 
)

Zero-suppress virgin raw data in "hybrid" mode

Subtracts pedestals (in 11bit mode, x->(x+1024-ped)/2) and common-mode noise, and inspects the digis then. If not flagged by the hybrid APV inspector, the zero-suppression is performed as usual (evaluation and subtraction of the baseline, truncation). Otherwise, the pedestal-subtracted digis (as above) are saved in one 128-strip cluster. Note: the APV restorer is used, it must be configured with APVInspectMode='HybridEmulation' if this method is called.

Parameters
idmodule DetId
firstAPVindex of the first APV considered
procRawDigisinput (virgin raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
outputzero-suppressed digis (or pedestal-subtracted digis, see above)
Returns
number of restored APVs

Definition at line 235 of file SiStripRawProcessingAlgorithms.cc.

References getAPVFlags(), mps_fire::i, restorer, subtractorCMN, subtractorPed, suppressor, and HLT_2024v14_cff::truncate.

Referenced by convertVirginRawToHybrid().

238  {
239  digivector_t procRawDigisPedSubtracted;
240 
241  for (auto& digi : procRawDigis) {
242  digi += 1024;
243  } // adding one MSB
244 
245  subtractorPed->subtract(id, firstAPV * 128, procRawDigis); // all strips are pedestals subtracted
246 
247  for (auto& digi : procRawDigis) {
248  digi /= 2;
249  }
250 
251  procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
252 
253  subtractorCMN->subtract(id, firstAPV, procRawDigis);
254 
255  const auto nAPVFlagged = restorer->inspect(id, firstAPV, procRawDigis, subtractorCMN->getAPVsCM());
256 
257  for (auto& digi : procRawDigis) {
258  digi *= 2;
259  }
260 
261  const std::vector<bool>& apvf = getAPVFlags();
262  const std::size_t nAPVs = procRawDigis.size() / 128;
263  for (uint16_t iAPV = firstAPV; iAPV < nAPVs + firstAPV; ++iAPV) {
264  if (apvf[iAPV]) {
265  //GB 23/6/08: truncation should be done at the very beginning
266  for (uint16_t i = 0; i < 128; ++i) {
267  const int16_t digi = procRawDigisPedSubtracted[128 * (iAPV - firstAPV) + i];
268  output.push_back(SiStripDigi(128 * iAPV + i, (digi < 0 ? 0 : suppressor->truncate(digi))));
269  }
270  } else {
271  const auto firstDigiIt = std::begin(procRawDigis) + 128 * (iAPV - firstAPV);
272  std::vector<int16_t> singleAPVdigi(firstDigiIt, firstDigiIt + 128);
273  suppressor->suppress(singleAPVdigi, iAPV, output);
274  }
275  }
276 
277  return nAPVFlagged;
278 }
const std::unique_ptr< SiStripPedestalsSubtractor > subtractorPed
SiStripAPVRestorer::digivector_t digivector_t
const std::unique_ptr< SiStripAPVRestorer > restorer
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
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
const std::vector< bool > & getAPVFlags() const
Definition: output.py:1

◆ convertVirginRawToHybrid() [2/2]

uint16_t SiStripRawProcessingAlgorithms::convertVirginRawToHybrid ( const edm::DetSet< SiStripRawDigi > &  rawDigis,
edm::DetSet< SiStripDigi > &  suppressedDigis 
)

Zero-suppress virgin raw data in "hybrid" mode

Subtracts pedestals (in 11bit mode, x->(x+1024-ped)/2) and common-mode noise, and inspects the digis then. If flagged by the hybrid APV inspector, the zero-suppression is performed as usual (evaluation and subtraction of the baseline, truncation). Otherwise, the pedestal-subtracted digis are saved in one 128-strip cluster.

Parameters
rawDigisinput (virgin) raw digis
outputzero-suppressed digis (or pedestal-subtracted digis, see above)
Returns
number of restored APVs

Definition at line 292 of file SiStripRawProcessingAlgorithms.cc.

References SiStripRawDigi::adc(), convertVirginRawToHybrid(), edm::DetSet< T >::id, edm::DetSet< T >::size(), and HcalDetIdTransform::transform().

293  {
294  digivector_t rawdigis;
295  rawdigis.reserve(rawDigis.size());
296  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
297  return digi.adc();
298  });
299  return convertVirginRawToHybrid(rawDigis.id, 0, rawdigis, suppressedDigis);
300 }
SiStripAPVRestorer::digivector_t digivector_t
uint16_t convertVirginRawToHybrid(uint32_t detId, uint16_t firstAPV, digivector_t &inDigis, edm::DetSet< SiStripDigi > &rawDigis)
uint16_t adc() const
det_id_type id
Definition: DetSet.h:79
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
size_type size() const
Definition: DetSet.h:61
unsigned transform(const HcalDetId &id, unsigned transformCode)

◆ getAPVFlags()

const std::vector<bool>& SiStripRawProcessingAlgorithms::getAPVFlags ( ) const
inline

Definition at line 48 of file SiStripRawProcessingAlgorithms.h.

References restorer.

Referenced by convertVirginRawToHybrid(), and suppressHybridData().

48 { return restorer->getAPVFlags(); }
const std::unique_ptr< SiStripAPVRestorer > restorer

◆ getAPVsCM()

const SiStripAPVRestorer::medians_t& SiStripRawProcessingAlgorithms::getAPVsCM ( ) const
inline

Definition at line 53 of file SiStripRawProcessingAlgorithms.h.

References subtractorCMN.

53 { return subtractorCMN->getAPVsCM(); }
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN

◆ getBaselineMap()

const SiStripAPVRestorer::baselinemap_t& SiStripRawProcessingAlgorithms::getBaselineMap ( ) const
inline

Definition at line 49 of file SiStripRawProcessingAlgorithms.h.

References restorer.

49 { return restorer->getBaselineMap(); }
const std::unique_ptr< SiStripAPVRestorer > restorer

◆ getSmoothedPoints()

const std::map<uint16_t, SiStripAPVRestorer::digimap_t>& SiStripRawProcessingAlgorithms::getSmoothedPoints ( ) const
inline

Definition at line 50 of file SiStripRawProcessingAlgorithms.h.

References restorer.

50  {
51  return restorer->getSmoothedPoints();
52  }
const std::unique_ptr< SiStripAPVRestorer > restorer

◆ initialize() [1/2]

void SiStripRawProcessingAlgorithms::initialize ( const edm::EventSetup es)

Definition at line 29 of file SiStripRawProcessingAlgorithms.cc.

References edm::EventSetup::getData(), restorer, subtractorCMN, subtractorPed, suppressor, tkGeomToken_, and trGeo.

Referenced by initialize().

29  {
30  subtractorPed->init(es);
31  subtractorCMN->init(es);
32  suppressor->init(es);
33  if (restorer.get())
34  restorer->init(es);
35 
36  trGeo = &es.getData(tkGeomToken_);
37 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const std::unique_ptr< SiStripPedestalsSubtractor > subtractorPed
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const std::unique_ptr< SiStripAPVRestorer > restorer
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
const std::unique_ptr< SiStripFedZeroSuppression > suppressor

◆ initialize() [2/2]

void SiStripRawProcessingAlgorithms::initialize ( const edm::EventSetup es,
const edm::Event e 
)

◆ suppressHybridData()

uint16_t SiStripRawProcessingAlgorithms::suppressHybridData ( const uint16_t  maxNStrips,
const edm::DetSet< SiStripDigi > &  hybridDigis,
edm::DetSet< SiStripDigi > &  suppressedDigis,
uint16_t  firstAPV = 0 
)

Zero-suppress "hybrid" raw data

Subtracts common-mode noise, and inspects the digis then. If flagged by the APV inspector, the zero-suppression is performed as usual. Otherwise, the positive inputs are copied.

Parameters
hybridDigisinput ADCs in ZS format (regular ZS or "hybrid", i.e. processed as x->(x+1024-ped)/2)
suppressedDigiszero-suppressed digis
firstAPV(optional) number of the first APV for which digis are should be handled (otherwise all present)
Returns
number of restored APVs

Definition at line 58 of file SiStripRawProcessingAlgorithms.cc.

References SiStripDigi::adc(), edm::DetSet< T >::begin(), HLT_2024v14_cff::distance, edm::DetSet< T >::end(), getAPVFlags(), mps_fire::i, edm::DetSet< T >::id, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, pfDeepBoostedJetPreprocessParams_cfi::lower_bound, edm::DetSet< T >::push_back(), restorer, SiStripDigi::strip(), subtractorCMN, suppressor, and HcalDetIdTransform::transform().

61  {
62  uint16_t nAPVFlagged = 0; // Count of flagged APVs
63  auto currentDigi = hybridDigis.begin();
64  const auto endDigi = hybridDigis.end();
65  auto currentAPV = firstAPV;
66 
67  // Loop through the APVs in the range
68  while (currentDigi != endDigi) {
69  // Determine the range of digis belonging to the current APV
70  const auto nextAPVBoundary = SiStripDigi((currentAPV + 1) * 128, 0);
71 
72  // Reject any APV larger than the max possible
73  if (nextAPVBoundary.strip() > maxNStrips) {
74  edm::LogError("SiStripRawProcessingAlgorithms")
75  << "In DetId " << suppressedDigis.id << " encountered APV boundary with strip number "
76  << nextAPVBoundary.strip() << ", which exceeds the maximum allowed value for this module (" << maxNStrips
77  << "). Exiting loop.";
78  break;
79  }
80 
81  const auto nextAPVDigi = std::lower_bound(currentDigi, endDigi, nextAPVBoundary);
82  const auto nDigisInAPV = std::distance(currentDigi, nextAPVDigi);
83 
84  // Handle based on the number of digis in the current APV
85  if (nDigisInAPV > 64) {
86  // Process hybrid data for noisy APV
87  digivector_t workDigis(128, -1024);
88 
89  // Populate `workDigis` with values from `currentDigi`
90  for (auto it = currentDigi; it != nextAPVDigi; ++it) {
91  workDigis[it->strip() - 128 * currentAPV] = it->adc() * 2 - 1024;
92  }
93 
94  // Perform pedestal subtraction
95  digivector_t workDigisPedSubtracted(workDigis);
96  subtractorCMN->subtract(hybridDigis.id, currentAPV, workDigis);
97 
98  // Inspect and restore digis
99  const auto apvFlagged = restorer->inspectAndRestore(
100  hybridDigis.id, currentAPV, workDigisPedSubtracted, workDigis, subtractorCMN->getAPVsCM());
101  nAPVFlagged += apvFlagged;
102 
103  // Process based on the APV flag
104  if (getAPVFlags()[currentAPV]) {
105  // Suppress flagged APVs
106  suppressor->suppress(workDigis, currentAPV, suppressedDigis);
107  } else {
108  // Handle bad APV: not flagged but exceeds threshold
109  for (uint16_t i = 0; i < 128; ++i) {
110  const auto digi = workDigisPedSubtracted[i];
111  if (digi > 0) {
112  suppressedDigis.push_back(SiStripDigi(currentAPV * 128 + i, suppressor->truncate(digi)));
113  }
114  }
115  }
116  } else {
117  // Already zero-suppressed: copy and truncate
118  std::transform(currentDigi, nextAPVDigi, std::back_inserter(suppressedDigis), [this](const SiStripDigi& inDigi) {
119  return SiStripDigi(inDigi.strip(), suppressor->truncate(inDigi.adc()));
120  });
121  }
122 
123  // Move to the next APV
124  currentDigi = nextAPVDigi;
125  ++currentAPV;
126  }
127 
128  return nAPVFlagged;
129 }
iterator end()
Definition: DetSet.h:58
void push_back(const T &t)
Definition: DetSet.h:66
SiStripAPVRestorer::digivector_t digivector_t
Log< level::Error, false > LogError
const std::unique_ptr< SiStripAPVRestorer > restorer
const uint16_t & strip() const
Definition: SiStripDigi.h:33
const uint16_t & adc() const
Definition: SiStripDigi.h:34
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
A Digi for the silicon strip detector, containing both strip and adc information, and suitable for st...
Definition: SiStripDigi.h:12
iterator begin()
Definition: DetSet.h:57
const std::unique_ptr< SiStripFedZeroSuppression > suppressor
const std::vector< bool > & getAPVFlags() const
det_id_type id
Definition: DetSet.h:79
unsigned transform(const HcalDetId &id, unsigned transformCode)

◆ suppressProcessedRawData() [1/2]

uint16_t SiStripRawProcessingAlgorithms::suppressProcessedRawData ( uint32_t  id,
uint16_t  firstAPV,
digivector_t procRawDigis,
edm::DetSet< SiStripDigi > &  output 
)

Zero-suppress processed (pedestals-subtracted) raw data.

Subtracts common-mode noise and (optionally, if doAPVRestore) re-evaluates and subtracts the baseline.

Parameters
idmodule DetId
firstAPVindex of the first APV to consider
procRawDigisinput (processed raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
outputzero-suppressed digis
Returns
number of restored APVs

Definition at line 183 of file SiStripRawProcessingAlgorithms.cc.

References doAPVRestore, restorer, subtractorCMN, and suppressor.

Referenced by suppressProcessedRawData(), and suppressVirginRawData().

186  {
187  digivector_t procRawDigisPedSubtracted;
188 
189  int16_t nAPVFlagged = 0;
190  if (doAPVRestore)
191  procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
192  subtractorCMN->subtract(id, firstAPV, procRawDigis);
193  if (doAPVRestore)
194  nAPVFlagged =
195  restorer->inspectAndRestore(id, firstAPV, procRawDigisPedSubtracted, procRawDigis, subtractorCMN->getAPVsCM());
196  suppressor->suppress(procRawDigis, firstAPV, output);
197  return nAPVFlagged;
198 }
SiStripAPVRestorer::digivector_t digivector_t
const std::unique_ptr< SiStripAPVRestorer > restorer
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
const std::unique_ptr< SiStripFedZeroSuppression > suppressor
Definition: output.py:1

◆ suppressProcessedRawData() [2/2]

uint16_t SiStripRawProcessingAlgorithms::suppressProcessedRawData ( const edm::DetSet< SiStripRawDigi > &  rawDigis,
edm::DetSet< SiStripDigi > &  output 
)

Zero-suppress processed (pedestals-subtracted) raw data.

Subtracts common-mode noise and (optionally, if doAPVRestore) re-evaluates and subtracts the baseline.

Parameters
rawDigisinput (processed) raw digis
outputzero-suppressed digis
Returns
number of restored APVs

Definition at line 210 of file SiStripRawProcessingAlgorithms.cc.

References SiStripRawDigi::adc(), edm::DetSet< T >::id, edm::DetSet< T >::size(), suppressProcessedRawData(), and HcalDetIdTransform::transform().

211  {
212  digivector_t rawdigis;
213  rawdigis.reserve(rawDigis.size());
214  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
215  return digi.adc();
216  });
217  return suppressProcessedRawData(rawDigis.id, 0, rawdigis, output);
218 }
uint16_t suppressProcessedRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
SiStripAPVRestorer::digivector_t digivector_t
uint16_t adc() const
det_id_type id
Definition: DetSet.h:79
Definition: output.py:1
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
size_type size() const
Definition: DetSet.h:61
unsigned transform(const HcalDetId &id, unsigned transformCode)

◆ suppressVirginRawData() [1/2]

uint16_t SiStripRawProcessingAlgorithms::suppressVirginRawData ( uint32_t  id,
uint16_t  firstAPV,
digivector_t procRawDigis,
edm::DetSet< SiStripDigi > &  output 
)

Zero-suppress virgin raw data.

Subtracts pedestals and common-mode noise, and (optionally, if doAPVRestore) re-evaluates and subtracts the baseline.

Parameters
idmodule DetId
firstAPVindex of the first APV to consider
procRawDigisinput (virgin raw) ADCs. Output: the ADCs after all subtractions, but before zero-suppression
outputzero-suppressed digis
Returns
number of restored APVs

Definition at line 143 of file SiStripRawProcessingAlgorithms.cc.

References subtractorPed, and suppressProcessedRawData().

Referenced by suppressVirginRawData().

146  {
147  subtractorPed->subtract(id, firstAPV * 128, procRawDigis);
148  return suppressProcessedRawData(id, firstAPV, procRawDigis, output);
149 }
const std::unique_ptr< SiStripPedestalsSubtractor > subtractorPed
uint16_t suppressProcessedRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
Definition: output.py:1

◆ suppressVirginRawData() [2/2]

uint16_t SiStripRawProcessingAlgorithms::suppressVirginRawData ( const edm::DetSet< SiStripRawDigi > &  rawDigis,
edm::DetSet< SiStripDigi > &  output 
)

Zero-suppress virgin raw data.

Subtracts pedestals and common-mode noise, and (optionally, if doAPVRestore) re-evaluates and subtracts the baseline.

Parameters
rawDigisinput (virgin) raw digis
outputzero-suppressed digis
Returns
number of restored APVs

Definition at line 161 of file SiStripRawProcessingAlgorithms.cc.

References SiStripRawDigi::adc(), edm::DetSet< T >::id, edm::DetSet< T >::size(), suppressVirginRawData(), and HcalDetIdTransform::transform().

162  {
163  digivector_t rawdigis;
164  rawdigis.reserve(rawDigis.size());
165  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
166  return digi.adc();
167  });
168  return suppressVirginRawData(rawDigis.id, 0, rawdigis, output);
169 }
uint16_t suppressVirginRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
SiStripAPVRestorer::digivector_t digivector_t
uint16_t adc() const
det_id_type id
Definition: DetSet.h:79
Definition: output.py:1
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
size_type size() const
Definition: DetSet.h:61
unsigned transform(const HcalDetId &id, unsigned transformCode)

Friends And Related Function Documentation

◆ SiStripRawProcessingFactory

friend class SiStripRawProcessingFactory
friend

Definition at line 16 of file SiStripRawProcessingAlgorithms.h.

Member Data Documentation

◆ doAPVRestore

const bool SiStripRawProcessingAlgorithms::doAPVRestore
private

Definition at line 61 of file SiStripRawProcessingAlgorithms.h.

Referenced by initialize(), and suppressProcessedRawData().

◆ restorer

const std::unique_ptr<SiStripAPVRestorer> SiStripRawProcessingAlgorithms::restorer

◆ subtractorCMN

const std::unique_ptr<SiStripCommonModeNoiseSubtractor> SiStripRawProcessingAlgorithms::subtractorCMN

◆ subtractorPed

const std::unique_ptr<SiStripPedestalsSubtractor> SiStripRawProcessingAlgorithms::subtractorPed

◆ suppressor

const std::unique_ptr<SiStripFedZeroSuppression> SiStripRawProcessingAlgorithms::suppressor

◆ tkGeomToken_

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiStripRawProcessingAlgorithms::tkGeomToken_
private

Definition at line 64 of file SiStripRawProcessingAlgorithms.h.

Referenced by initialize().

◆ trGeo

const TrackerGeometry* SiStripRawProcessingAlgorithms::trGeo
private

Definition at line 65 of file SiStripRawProcessingAlgorithms.h.

Referenced by initialize().

◆ useCMMeanMap

const bool SiStripRawProcessingAlgorithms::useCMMeanMap
private

Definition at line 62 of file SiStripRawProcessingAlgorithms.h.

Referenced by initialize().