CMS 3D CMS Logo

SiStripZeroSuppression.cc
Go to the documentation of this file.
2 
12 #include <memory>
13 
15  : algorithms(SiStripRawProcessingFactory::create(conf.getParameter<edm::ParameterSet>("Algorithms"))),
16  produceRawDigis(conf.getParameter<bool>("produceRawDigis")),
17  storeCM(conf.getParameter<bool>("storeCM")),
18  fixCM(conf.getParameter<bool>("fixCM")),
19  produceCalculatedBaseline(conf.getParameter<bool>("produceCalculatedBaseline")),
20  produceBaselinePoints(conf.getParameter<bool>("produceBaselinePoints")),
21  storeInZScollBadAPV(conf.getParameter<bool>("storeInZScollBadAPV")),
22  produceHybridFormat(conf.getParameter<bool>("produceHybridFormat"))
23 {
24  for ( const auto& inputTag : conf.getParameter<std::vector<edm::InputTag> >("RawDigiProducersList") )
25  {
26  const auto& tagName = inputTag.instance();
27  produces<edm::DetSetVector<SiStripDigi>>(tagName);
28  if (produceRawDigis)
29  produces<edm::DetSetVector<SiStripRawDigi>>(tagName);
30  if (storeCM)
31  produces<edm::DetSetVector<SiStripProcessedRawDigi>>("APVCM"+tagName);
33  produces<edm::DetSetVector<SiStripProcessedRawDigi>>("BADAPVBASELINE"+tagName);
35  produces<edm::DetSetVector<SiStripDigi>>("BADAPVBASELINEPOINTS"+tagName);
36 
37  RawType inputType = RawType::Unknown;
38  if ( tagName == "ProcessedRaw" ) {
39  inputType = RawType::ProcessedRaw;
40  if (produceHybridFormat) throw cms::Exception("Processed Raw Cannot be converted in hybrid Format");
41  } else if ( tagName == "VirginRaw" ) {
42  inputType = RawType::VirginRaw;
43  } else if ( tagName == "ScopeMode" ) {
44  inputType = RawType::ScopeMode;
45  if (produceHybridFormat) throw cms::Exception("Scope Mode cannot be converted in hybrid Format");
46  }
47  if ( RawType::Unknown != inputType ) {
48  rawInputs.emplace_back(tagName, inputType, consumes<edm::DetSetVector<SiStripRawDigi>>(inputTag));
49  } else if ( tagName == "ZeroSuppressed" ) {
51  } else {
52  throw cms::Exception("Unknown input type") << tagName << " unknown. "
53  << "SiStripZeroZuppression can only process types \"VirginRaw\", \"ProcessedRaw\" and \"ZeroSuppressed\"";
54  }
55  }
56 
57  if ( produceHybridFormat && ( "HybridEmulation" != conf.getParameter<edm::ParameterSet>("Algorithms").getParameter<std::string>("APVInspectMode") ) )
58  throw cms::Exception("Invalid option") << "When producing data in the hybrid format, the APV restorer must be configured with APVInspectMode='HybridEmulation'";
59 
60  if ( ! ( rawInputs.empty() && hybridInputs.empty() ) ) {
61  output_base.reserve(16000);
62  if (produceRawDigis) output_base_raw.reserve(16000);
63  if (storeCM) output_apvcm.reserve(16000);
64  if (produceCalculatedBaseline) output_baseline.reserve(16000);
66  }
67 }
68 
70 {
71  algorithms->initialize(es, e);
72 
73  for ( const auto& input : rawInputs ) {
74  clearOutputs();
76  e.getByToken(std::get<rawtoken_t>(input), inDigis);
77  if ( ! inDigis->empty() )
78  processRaw(*inDigis, std::get<RawType>(input));
79  putOutputs(e, std::get<std::string>(input));
80  }
81  for ( const auto& input : hybridInputs ) {
82  clearOutputs();
84  e.getByToken(std::get<zstoken_t>(input), inDigis);
85  if ( ! inDigis->empty() ) {
86  processHybrid(*inDigis);
87  }
88  putOutputs(e, std::get<std::string>(input));
89  }
90 }
91 
93 {
94  output_base.clear();
95  output_base_raw.clear();
96  output_baseline.clear();
97  output_baseline_points.clear();
98  output_apvcm.clear();
99 }
101 {
102  evt.put(std::make_unique<edm::DetSetVector<SiStripDigi>>(output_base), tagName);
103  if (produceRawDigis)
104  evt.put(std::make_unique<edm::DetSetVector<SiStripRawDigi>>(output_base_raw), tagName);
106  evt.put(std::make_unique<edm::DetSetVector<SiStripProcessedRawDigi>>(output_baseline), "BADAPVBASELINE"+tagName);
108  evt.put(std::make_unique<edm::DetSetVector<SiStripDigi>>(output_baseline_points), "BADAPVBASELINEPOINTS"+tagName);
109  if (storeCM)
110  evt.put(std::make_unique<edm::DetSetVector<SiStripProcessedRawDigi>>(output_apvcm), "APVCM"+tagName);
111 }
112 
114 {
115  for ( const auto& rawDigis : input ) {
116  edm::DetSet<SiStripDigi> suppressedDigis(rawDigis.id);
117 
118  int16_t nAPVflagged = 0;
119  if ( RawType::ProcessedRaw == inType ) {
120  nAPVflagged = algorithms->suppressProcessedRawData(rawDigis, suppressedDigis);
121  } else if ( RawType::ScopeMode == inType) {
122  nAPVflagged = algorithms->suppressVirginRawData(rawDigis, suppressedDigis);
123  } else if ( RawType::VirginRaw == inType ) {
124  if ( produceHybridFormat ) {
125  nAPVflagged = algorithms->convertVirginRawToHybrid(rawDigis, suppressedDigis);
126  } else {
127  nAPVflagged = algorithms->suppressVirginRawData(rawDigis, suppressedDigis);
128  }
129  }
130 
131  storeExtraOutput(rawDigis.id, nAPVflagged);
132  if (!suppressedDigis.empty() && (storeInZScollBadAPV || nAPVflagged ==0))
133  output_base.push_back(std::move(suppressedDigis));
134 
135  if (produceRawDigis && nAPVflagged > 0) {
136  output_base_raw.push_back(formatRawDigis(rawDigis));
137  }
138  }
139 }
140 
142 {
143  for ( const auto& inDigis : input ) {
144  edm::DetSet<SiStripDigi> suppressedDigis(inDigis.id);
145 
146  std::vector<int16_t> rawDigis;
147  const auto nAPVflagged = algorithms->suppressHybridData(inDigis, suppressedDigis, rawDigis);
148 
149  storeExtraOutput(inDigis.id, nAPVflagged);
150  if (!suppressedDigis.empty() && (storeInZScollBadAPV || nAPVflagged ==0))
151  output_base.push_back(std::move(suppressedDigis));
152 
153  if (produceRawDigis && nAPVflagged > 0) {
154  output_base_raw.push_back(formatRawDigis(inDigis.id, rawDigis));
155  }
156  }
157 }
158 
160 {
161  edm::DetSet<SiStripRawDigi> outRawDigis(rawDigis.id);
162  outRawDigis.reserve(rawDigis.size());
163  const std::vector<bool>& apvf = algorithms->getAPVFlags();
164  uint32_t strip=0;
165  for ( const auto rawDigi : rawDigis ) {
166  int16_t apvN = strip/128;
167  if (apvf[apvN]) outRawDigis.push_back(rawDigi);
168  else outRawDigis.push_back(SiStripRawDigi(0));
169  ++strip;
170  }
171  return outRawDigis;
172 }
173 
174 inline edm::DetSet<SiStripRawDigi> SiStripZeroSuppression::formatRawDigis(uint32_t detId, const std::vector<int16_t>& rawDigis)
175 {
176  edm::DetSet<SiStripRawDigi> outRawDigis(detId);
177  outRawDigis.reserve(rawDigis.size());
178  const std::vector<bool>& apvf = algorithms->getAPVFlags();
179  uint32_t strip=0;
180  for ( const auto rawDigi : rawDigis ) {
181  int16_t apvN = strip/128;
182  if (apvf[apvN]) outRawDigis.push_back(SiStripRawDigi(rawDigi));
183  else outRawDigis.push_back(SiStripRawDigi(0));
184  ++strip;
185  }
186  return outRawDigis;
187 }
188 
189 inline void SiStripZeroSuppression::storeExtraOutput(uint32_t id, int16_t nAPVflagged)
190 {
191  const auto& vmedians = algorithms->getAPVsCM();
192  if (storeCM) storeCMN(id, vmedians);
193  if (nAPVflagged > 0){
194  if (produceCalculatedBaseline) storeBaseline(id, vmedians);
196  }
197 }
198 
199 inline void SiStripZeroSuppression::storeBaseline(uint32_t id, const medians_t& vmedians)
200 {
201  const auto& baselinemap = algorithms->getBaselineMap();
202 
203  edm::DetSet<SiStripProcessedRawDigi> baselineDetSet(id);
204  baselineDetSet.reserve(vmedians.size()*128);
205  for ( const auto& vmed : vmedians ) {
206  const uint16_t apvN = vmed.first;
207  const float median = vmed.second;
208  auto itBaselineMap = baselinemap.find(apvN);
209  if ( baselinemap.end() == itBaselineMap ) {
210  for (size_t strip=0; strip < 128; ++strip)
211  baselineDetSet.push_back(SiStripProcessedRawDigi(median));
212  } else {
213  for (size_t strip=0; strip < 128; ++strip)
214  baselineDetSet.push_back(SiStripProcessedRawDigi((itBaselineMap->second)[strip]));
215  }
216  }
217 
218  if ( ! baselineDetSet.empty() )
219  output_baseline.push_back(baselineDetSet);
220 }
221 
223 {
224  edm::DetSet<SiStripDigi> baspointDetSet(id);
225  for ( const auto& itBaselinePointVect : algorithms->getSmoothedPoints() ) {
226  const uint16_t apvN = itBaselinePointVect.first;
227  for ( const auto& itBaselinePointMap : itBaselinePointVect.second ) {
228  const uint16_t bpstrip = itBaselinePointMap.first + apvN*128;
229  const int16_t bp = itBaselinePointMap.second;
230  baspointDetSet.push_back(SiStripDigi(bpstrip, bp+128));
231  }
232  }
233 
234  if ( ! baspointDetSet.empty() )
235  output_baseline_points.push_back(std::move(baspointDetSet));
236 }
237 
238 inline void SiStripZeroSuppression::storeCMN(uint32_t id, const medians_t& vmedians)
239 {
240  std::vector<bool> apvf(6, false);
241  if (fixCM) {
242  const auto& apvFlagged = algorithms->getAPVFlags();
243  std::copy(std::begin(apvFlagged), std::end(apvFlagged), std::begin(apvf));
244  }
245 
247  short apvNb=0;
248  for ( const auto& vmed : vmedians ) {
249  if ( vmed.first > apvNb ) {
250  for ( int i{0}; i < vmed.first-apvNb; ++i ) {
251  apvDetSet.push_back(SiStripProcessedRawDigi(-999.));
252  apvNb++;
253  }
254  }
255 
256  if ( fixCM && apvf[vmed.first] ) {
257  apvDetSet.push_back(SiStripProcessedRawDigi(-999.));
258  } else {
259  apvDetSet.push_back(SiStripProcessedRawDigi(vmed.second));
260  }
261  apvNb++;
262  }
263 
264  if(!apvDetSet.empty())
265  output_apvcm.push_back(std::move(apvDetSet));
266 }
T getParameter(std::string const &) const
std::vector< std::pair< short, float >> medians_t
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
void push_back(const T &t)
Definition: DetSet.h:68
std::vector< edm::DetSet< SiStripDigi > > output_base
def copy(args, dbName)
def create(alignables, pedeDump, additionalData, outputFile, config)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< std::tuple< std::string, zstoken_t > > hybridInputs
A signed Digi for the silicon strip detector, containing only adc information, and suitable for stori...
std::vector< edm::DetSet< SiStripDigi > > output_baseline_points
void processRaw(const edm::DetSetVector< SiStripRawDigi > &input, RawType inType)
static std::string const input
Definition: EdmProvDump.cc:45
size_type size() const
Definition: DetSet.h:63
void reserve(size_t s)
Definition: DetSet.h:67
void storeExtraOutput(uint32_t, int16_t)
std::vector< edm::DetSet< SiStripProcessedRawDigi > > output_apvcm
std::unique_ptr< SiStripRawProcessingAlgorithms > algorithms
std::vector< edm::DetSet< SiStripProcessedRawDigi > > output_baseline
#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
std::vector< std::tuple< std::string, RawType, rawtoken_t > > rawInputs
void storeBaseline(uint32_t, const medians_t &)
std::vector< edm::DetSet< SiStripRawDigi > > output_base_raw
void storeCMN(uint32_t, const medians_t &)
void processHybrid(const edm::DetSetVector< SiStripDigi > &input)
SiStripZeroSuppression(const edm::ParameterSet &)
bool empty() const
Definition: DetSet.h:64
edm::DetSet< SiStripRawDigi > formatRawDigis(const edm::DetSet< SiStripRawDigi > &rawDigis)
#define begin
Definition: vmac.h:32
HLT enums.
det_id_type id
Definition: DetSet.h:77
void putOutputs(edm::Event &evt, const std::string &tagName)
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
def move(src, dest)
Definition: eostools.py:511