CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripProcessedRawDigiProducer.cc
Go to the documentation of this file.
2 
8 
11 
15 
16 #include <functional>
17 
19  : inputTags(conf.getParameter<std::vector<edm::InputTag> >("DigiProducersList")),
20  inputTokensDigi(edm::vector_transform(inputTags, [this](edm::InputTag const & tag){return consumes<edm::DetSetVector<SiStripDigi> >(tag);})),
21  inputTokensRawDigi(edm::vector_transform(inputTags, [this](edm::InputTag const & tag){return consumes<edm::DetSetVector<SiStripRawDigi> >(tag);})),
24 
25  produces< edm::DetSetVector<SiStripProcessedRawDigi> >("");
26 }
27 
30 
31  std::auto_ptr< edm::DetSetVector<SiStripProcessedRawDigi> > output(new edm::DetSetVector<SiStripProcessedRawDigi>());
34 
35  es.get<SiStripGainRcd>().get(gainHandle);
36  subtractorPed->init(es);
37  subtractorCMN->init(es);
38 
39  std::string label = findInput(inputRawdigis, inputTokensRawDigi, e);
40  if( "VirginRaw" == label ) vr_process(*inputRawdigis, *output);
41  else if( "ProcessedRaw" == label ) pr_process(*inputRawdigis, *output);
42  else if( "ZeroSuppressed" == findInput(inputDigis, inputTokensDigi, e) ) zs_process(*inputDigis, *output);
43  else
44  edm::LogError("Input Not Found");
45 
46  e.put(output);
47 }
48 
49 template<class T>
50 inline
52 findInput(edm::Handle<T>& handle, const std::vector<edm::EDGetTokenT<T> >& tokens, const edm::Event& e ) {
53 
54  for( typename std::vector<edm::EDGetTokenT<T> >::const_iterator
55  token = tokens.begin(); token != tokens.end(); ++token ) {
56  unsigned index(token - tokens.begin());
57  e.getByToken(*token, handle);
58  if( handle.isValid() && !handle->empty() ) {
59  edm::LogInfo("Input") << inputTags.at(index);
60  return inputTags.at(index).instance();
61  }
62  }
63  return "Input Not Found";
64 }
65 
68  std::vector<float> digis;
69  for(edm::DetSetVector<SiStripDigi>::const_iterator detset = input.begin(); detset != input.end(); detset++ ) {
70  digis.clear();
71  for(edm::DetSet<SiStripDigi>::const_iterator digi = detset->begin(); digi != detset->end(); digi++) {
72  digis.resize( digi->strip(), 0);
73  digis.push_back( digi->adc() );
74  }
75  common_process( detset->id, digis, output);
76  }
77 }
78 
81  for(edm::DetSetVector<SiStripRawDigi>::const_iterator detset=input.begin(); detset!=input.end(); detset++) {
82  std::vector<float> digis;
83  transform(detset->begin(), detset->end(), back_inserter(digis), boost::bind(&SiStripRawDigi::adc , _1));
84  subtractorCMN->subtract(detset->id, 0, digis);
85  common_process( detset->id, digis, output);
86  }
87 }
88 
91  for(edm::DetSetVector<SiStripRawDigi>::const_iterator detset=input.begin(); detset!=input.end(); detset++) {
92  std::vector<int16_t> int_digis(detset->size());
93  subtractorPed->subtract(*detset,int_digis);
94  std::vector<float> digis(int_digis.begin(), int_digis.end());
95  subtractorCMN->subtract(detset->id, 0, digis);
96  common_process( detset->id, digis, output);
97  }
98 }
99 
101 common_process(const uint32_t detId, std::vector<float> & digis, edm::DetSetVector<SiStripProcessedRawDigi>& output) {
102 
103  //Apply Gains
104  SiStripApvGain::Range detGainRange = gainHandle->getRange(detId);
105  for(std::vector<float>::iterator it=digis.begin(); it<digis.end(); it++)
106  (*it)/= (gainHandle->getStripGain(it-digis.begin(), detGainRange));
107 
108  //Insert as DetSet
110  copy(digis.begin(), digis.end(), back_inserter(ds.data) );
111  output.insert(ds);
112 }
const uint16_t & adc() const
std::vector< edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > > inputTokensRawDigi
std::auto_ptr< SiStripPedestalsSubtractor > subtractorPed
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
std::string findInput(edm::Handle< T > &handle, const std::vector< edm::EDGetTokenT< T > > &tokens, const edm::Event &e)
subtractorPed(SiStripRawProcessingFactory::create_SubtractorPed(conf))
inputTokensRawDigi(edm::vector_transform(inputTags, [this](edm::InputTag const &tag){return consumes< edm::DetSetVector< SiStripRawDigi > >(tag);}))
static std::auto_ptr< SiStripCommonModeNoiseSubtractor > create_SubtractorCMN(const edm::ParameterSet &)
static std::string const input
Definition: EdmProvDump.cc:44
void pr_process(const edm::DetSetVector< SiStripRawDigi > &, edm::DetSetVector< SiStripProcessedRawDigi > &)
SiStripProcessedRawDigiProducer(edm::ParameterSet const &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
subtractorCMN(SiStripRawProcessingFactory::create_SubtractorCMN(conf))
tuple handle
Definition: patZpeak.py:22
std::pair< ContainerIterator, ContainerIterator > Range
void common_process(const uint32_t, std::vector< float > &, edm::DetSetVector< SiStripProcessedRawDigi > &)
bool isValid() const
Definition: HandleBase.h:76
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:363
std::vector< edm::EDGetTokenT< edm::DetSetVector< SiStripDigi > > > inputTokensDigi
void produce(edm::Event &e, const edm::EventSetup &es)
tuple conf
Definition: dbtoconf.py:185
std::auto_ptr< SiStripCommonModeNoiseSubtractor > subtractorCMN
const T & get() const
Definition: EventSetup.h:55
string const
Definition: compareJSON.py:14
void zs_process(const edm::DetSetVector< SiStripDigi > &, edm::DetSetVector< SiStripProcessedRawDigi > &)
void insert(detset const &s)
Insert the given DetSet.
Definition: DetSetVector.h:237
collection_type data
Definition: DetSet.h:78
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:348
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
static std::auto_ptr< SiStripPedestalsSubtractor > create_SubtractorPed(const edm::ParameterSet &)
void vr_process(const edm::DetSetVector< SiStripRawDigi > &, edm::DetSetVector< SiStripProcessedRawDigi > &)