CMS 3D CMS Logo

SiStripRawProcessingAlgorithms.cc
Go to the documentation of this file.
2 
12 #include <memory>
13 
14 SiStripRawProcessingAlgorithms::SiStripRawProcessingAlgorithms(std::unique_ptr<SiStripPedestalsSubtractor> ped,
15  std::unique_ptr<SiStripCommonModeNoiseSubtractor> cmn,
16  std::unique_ptr<SiStripFedZeroSuppression> zs,
17  std::unique_ptr<SiStripAPVRestorer> res,
18  bool doAPVRest,
19  bool useCMMap)
20  : subtractorPed(std::move(ped)),
21  subtractorCMN(std::move(cmn)),
22  suppressor(std::move(zs)),
23  restorer(std::move(res)),
24  doAPVRestore(doAPVRest),
25  useCMMeanMap(useCMMap) {}
26 
28  subtractorPed->init(es);
29  subtractorCMN->init(es);
30  suppressor->init(es);
31  if (restorer.get())
32  restorer->init(es);
33 
36  trGeo = tracker.product();
37 }
38 
40  initialize(es);
41  if (restorer.get() && doAPVRestore && useCMMeanMap)
42  restorer->loadMeanCMMap(e);
43 }
44 
59  edm::DetSet<SiStripDigi>& suppressedDigis,
60  uint16_t firstAPV) {
61  uint16_t nAPVFlagged = 0;
62  auto beginAPV = hybridDigis.begin();
63  const auto indigis_end = hybridDigis.end();
64  auto iAPV = firstAPV;
65  while (beginAPV != indigis_end) {
66  const auto endAPV = std::lower_bound(beginAPV, indigis_end, SiStripDigi((iAPV + 1) * 128, 0));
67  const auto nDigisInAPV = std::distance(beginAPV, endAPV);
68  if (nDigisInAPV > 64) {
69  digivector_t workDigis(128, -1024);
70  for (auto it = beginAPV; it != endAPV; ++it) {
71  workDigis[it->strip() - 128 * iAPV] = it->adc() * 2 - 1024;
72  }
73  digivector_t workDigisPedSubtracted(workDigis);
74  subtractorCMN->subtract(hybridDigis.id, iAPV, workDigis);
75  const auto apvFlagged = restorer->inspectAndRestore(
76  hybridDigis.id, iAPV, workDigisPedSubtracted, workDigis, subtractorCMN->getAPVsCM());
77  nAPVFlagged += apvFlagged;
78  if (getAPVFlags()[iAPV]) {
79  suppressor->suppress(workDigis, iAPV, suppressedDigis);
80  } else { // bad APV: more than 64 but not flagged
81  for (uint16_t i = 0; i != 128; ++i) {
82  const auto digi = workDigisPedSubtracted[i];
83  if (digi > 0) {
84  suppressedDigis.push_back(SiStripDigi(iAPV * 128 + i, suppressor->truncate(digi)));
85  }
86  }
87  }
88  } else { // already zero-suppressed, copy and truncate
89  std::transform(beginAPV, endAPV, std::back_inserter(suppressedDigis), [this](const SiStripDigi inDigi) {
90  return SiStripDigi(inDigi.strip(), suppressor->truncate(inDigi.adc()));
91  });
92  }
93  beginAPV = endAPV;
94  ++iAPV;
95  }
96  return nAPVFlagged;
97 }
98 
112  uint16_t firstAPV,
113  digivector_t& procRawDigis,
115  subtractorPed->subtract(id, firstAPV * 128, procRawDigis);
116  return suppressProcessedRawData(id, firstAPV, procRawDigis, output);
117 }
118 
131  digivector_t rawdigis;
132  rawdigis.reserve(rawDigis.size());
133  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
134  return digi.adc();
135  });
136  return suppressVirginRawData(rawDigis.id, 0, rawdigis, output);
137 }
138 
152  uint16_t firstAPV,
153  digivector_t& procRawDigis,
155  digivector_t procRawDigisPedSubtracted;
156 
157  int16_t nAPVFlagged = 0;
158  if (doAPVRestore)
159  procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
160  subtractorCMN->subtract(id, firstAPV, procRawDigis);
161  if (doAPVRestore)
162  nAPVFlagged =
163  restorer->inspectAndRestore(id, firstAPV, procRawDigisPedSubtracted, procRawDigis, subtractorCMN->getAPVsCM());
164  suppressor->suppress(procRawDigis, firstAPV, output);
165  return nAPVFlagged;
166 }
167 
180  digivector_t rawdigis;
181  rawdigis.reserve(rawDigis.size());
182  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
183  return digi.adc();
184  });
185  return suppressProcessedRawData(rawDigis.id, 0, rawdigis, output);
186 }
187 
204  uint16_t firstAPV,
205  digivector_t& procRawDigis,
207  digivector_t procRawDigisPedSubtracted;
208 
209  for (auto& digi : procRawDigis) {
210  digi += 1024;
211  } // adding one MSB
212 
213  subtractorPed->subtract(id, firstAPV * 128, procRawDigis); // all strips are pedestals subtracted
214 
215  for (auto& digi : procRawDigis) {
216  digi /= 2;
217  }
218 
219  procRawDigisPedSubtracted.assign(procRawDigis.begin(), procRawDigis.end());
220 
221  subtractorCMN->subtract(id, firstAPV, procRawDigis);
222 
223  const auto nAPVFlagged = restorer->inspect(id, firstAPV, procRawDigis, subtractorCMN->getAPVsCM());
224 
225  for (auto& digi : procRawDigis) {
226  digi *= 2;
227  }
228 
229  const std::vector<bool>& apvf = getAPVFlags();
230  const std::size_t nAPVs = procRawDigis.size() / 128;
231  for (uint16_t iAPV = firstAPV; iAPV < nAPVs + firstAPV; ++iAPV) {
232  if (apvf[iAPV]) {
233  //GB 23/6/08: truncation should be done at the very beginning
234  for (uint16_t i = 0; i < 128; ++i) {
235  const int16_t digi = procRawDigisPedSubtracted[128 * (iAPV - firstAPV) + i];
236  output.push_back(SiStripDigi(128 * iAPV + i, (digi < 0 ? 0 : suppressor->truncate(digi))));
237  }
238  } else {
239  const auto firstDigiIt = std::begin(procRawDigis) + 128 * (iAPV - firstAPV);
240  std::vector<int16_t> singleAPVdigi(firstDigiIt, firstDigiIt + 128);
241  suppressor->suppress(singleAPVdigi, iAPV, output);
242  }
243  }
244 
245  return nAPVFlagged;
246 }
247 
261  edm::DetSet<SiStripDigi>& suppressedDigis) {
262  digivector_t rawdigis;
263  rawdigis.reserve(rawDigis.size());
264  std::transform(std::begin(rawDigis), std::end(rawDigis), std::back_inserter(rawdigis), [](SiStripRawDigi digi) {
265  return digi.adc();
266  });
267  return convertVirginRawToHybrid(rawDigis.id, 0, rawdigis, suppressedDigis);
268 }
edm::DetSet::push_back
void push_back(const T &t)
Definition: DetSet.h:66
SiStripRawProcessingAlgorithms::subtractorPed
const std::unique_ptr< SiStripPedestalsSubtractor > subtractorPed
Definition: SiStripRawProcessingAlgorithms.h:54
SiStripRawProcessingFactory.h
Handle.h
mps_fire.i
i
Definition: mps_fire.py:355
SiStripRawProcessingAlgorithms::initialize
void initialize(const edm::EventSetup &)
Definition: SiStripRawProcessingAlgorithms.cc:27
edm::DetSet::size
size_type size() const
Definition: DetSet.h:61
edm::DetSet< SiStripDigi >
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:32
edm::DetSet::id
det_id_type id
Definition: DetSet.h:79
SiStripRawDigi.h
edm::DetSet::begin
iterator begin()
Definition: DetSet.h:57
SiStripRawProcessingAlgorithms::convertVirginRawToHybrid
uint16_t convertVirginRawToHybrid(uint32_t detId, uint16_t firstAPV, digivector_t &inDigis, edm::DetSet< SiStripDigi > &rawDigis)
Definition: SiStripRawProcessingAlgorithms.cc:203
HLT_2018_cff.distance
distance
Definition: HLT_2018_cff.py:6417
SiStripDetId.h
SiStripRawProcessingAlgorithms::trGeo
const TrackerGeometry * trGeo
Definition: SiStripRawProcessingAlgorithms.h:63
SiStripRawProcessingAlgorithms::suppressProcessedRawData
uint16_t suppressProcessedRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
Definition: SiStripRawProcessingAlgorithms.cc:151
HLT_2018_cff.doAPVRestore
doAPVRestore
Definition: HLT_2018_cff.py:8174
end
#define end
Definition: vmac.h:39
SiStripRawProcessingAlgorithms::suppressHybridData
uint16_t suppressHybridData(const edm::DetSet< SiStripDigi > &inDigis, edm::DetSet< SiStripDigi > &suppressedDigis, uint16_t firstAPV=0)
Definition: SiStripRawProcessingAlgorithms.cc:58
HLT_2018_cff.useCMMeanMap
useCMMeanMap
Definition: HLT_2018_cff.py:8172
SiStripRawDigi
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
Definition: SiStripRawDigi.h:15
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
SiStripRawDigi::adc
uint16_t adc() const
Definition: SiStripRawDigi.h:22
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
SiStripDigi.h
edm::ESHandle< TrackerGeometry >
SiStripDigi::adc
const uint16_t & adc() const
Definition: SiStripDigi.h:34
HcalDetIdTransform::transform
unsigned transform(const HcalDetId &id, unsigned transformCode)
Definition: HcalDetIdTransform.cc:7
SiStripRawProcessingAlgorithms::suppressor
const std::unique_ptr< SiStripFedZeroSuppression > suppressor
Definition: SiStripRawProcessingAlgorithms.h:56
SiStripRawProcessingAlgorithms::SiStripRawProcessingAlgorithms
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)
Definition: SiStripRawProcessingAlgorithms.cc:14
cuda_std::lower_bound
__host__ constexpr __device__ RandomIt lower_bound(RandomIt first, RandomIt last, const T &value, Compare comp={})
Definition: cudastdAlgorithm.h:27
PbPb_ZMuSkimMuonDPG_cff.tracker
tracker
Definition: PbPb_ZMuSkimMuonDPG_cff.py:60
SiStripRawProcessingAlgorithms::suppressVirginRawData
uint16_t suppressVirginRawData(uint32_t detId, uint16_t firstAPV, digivector_t &procRawDigis, edm::DetSet< SiStripDigi > &output)
Definition: SiStripRawProcessingAlgorithms.cc:111
SiStripDigi::strip
const uint16_t & strip() const
Definition: SiStripDigi.h:33
Event.h
SiStripRawProcessingAlgorithms::getAPVFlags
const std::vector< bool > & getAPVFlags() const
Definition: SiStripRawProcessingAlgorithms.h:47
SiStripRawProcessingAlgorithms::digivector_t
SiStripAPVRestorer::digivector_t digivector_t
Definition: SiStripRawProcessingAlgorithms.h:19
SiStripRawProcessingAlgorithms.h
SiStripRawProcessingAlgorithms::useCMMeanMap
const bool useCMMeanMap
Definition: SiStripRawProcessingAlgorithms.h:61
edm::EventSetup
Definition: EventSetup.h:57
DetSetVector.h
get
#define get
res
Definition: Electron.h:6
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
SiStripRawProcessingAlgorithms::restorer
const std::unique_ptr< SiStripAPVRestorer > restorer
Definition: SiStripRawProcessingAlgorithms.h:57
EventSetup.h
SiStripDigi
A Digi for the silicon strip detector, containing both strip and adc information, and suitable for st...
Definition: SiStripDigi.h:12
edm::DetSet::end
iterator end()
Definition: DetSet.h:58
ParameterSet.h
SiStripRawProcessingAlgorithms::subtractorCMN
const std::unique_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
Definition: SiStripRawProcessingAlgorithms.h:55
edm::Event
Definition: Event.h:73
SiStripRawProcessingAlgorithms::doAPVRestore
const bool doAPVRestore
Definition: SiStripRawProcessingAlgorithms.h:60
begin
#define begin
Definition: vmac.h:32
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37