CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripMeanCMExtractor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripMeanCMExtractor
4 // Class: SiStripMeanCMExtractor
5 //
13 //
14 // Original Author: Ivan Amos Cali,32 4-A08,+41227673039,
15 // Created: Wed Oct 13 11:50:47 CEST 2010
16 //
17 //
18 
19 
20 // system include files
21 #include <sstream>
22 #include <memory>
23 #include <list>
24 #include <algorithm>
25 #include <cassert>
26 
27 
28 // user include files
36 
39 
41 
49 
50 
51 typedef std::map<uint32_t, std::vector<float> > CMMap;
52 
54  public:
55  explicit SiStripMeanCMExtractor( const edm::ParameterSet&);
57 
58  private:
59  virtual void beginJob() override ;
60  virtual void produce(edm::Event&, const edm::EventSetup&) override;
61  virtual void endJob() override ;
62 
63  void init(const edm::EventSetup&);
64 
67 
73  uint16_t _nEventsToUse;
74  uint16_t _actualEvent;
75 
76  CMMap _CMMap; //it contains the sum of the CM calculated before. The normalization for the number of events it is done at the end when it is written in the DetSetVector.
77 };
78 
79 
81  _inputTag(conf.getParameter<edm::InputTag> ("CMCollection")),
82  _Algorithm(conf.getParameter<std::string>("Algorithm")),
83  _nEventsToUse(conf.getParameter<uint32_t>("NEvents"))
84 {
85 
86 
87  if(_nEventsToUse < 1) _nEventsToUse=1;
88  produces< edm::DetSetVector<SiStripProcessedRawDigi> > ("MEANAPVCM");
89 }
90 
91 
93 {
94 
95 }
96 
98 
99  uint32_t p_cache_id = es.get<SiStripPedestalsRcd>().cacheIdentifier();
100 
101  if(p_cache_id != pedestal_cache_id_) {
103  pedestal_cache_id_ = p_cache_id;
104  }
105 }
106 
107 // ------------ method called to produce the data ------------
108 void
110 {
111  using namespace edm;
112 
113 
114  //if(_actualEvent > _nEventsToUse) return;
115 
116  std::vector<edm::DetSet<SiStripProcessedRawDigi> > meancm;
117 
118  if(_Algorithm == "StoredCM"){
120  iEvent.getByLabel(_inputTag,inputCM);
121 
122  this->StoreMean(*inputCM);
123  this->ConvertMeanMapToDetSetVector(meancm);
124 
125  } else if (_Algorithm == "Pedestals"){
126  this->init(iSetup);
127 
129  iEvent.getByLabel(_inputTag,input);
130 
131  this->CMExtractorFromPedestals(*input,meancm);
132  }
133 
134  ++_actualEvent;
135 
136 
137 
138 
139  std::auto_ptr< edm::DetSetVector<SiStripProcessedRawDigi> > outputMeanCM(new edm::DetSetVector<SiStripProcessedRawDigi>(meancm) );
140  iEvent.put( outputMeanCM,"MEANAPVCM");
141 
142 }
143 
145  meancm.clear();
146  meancm.reserve(15000);
147 
149  rawDigis = input.begin(); rawDigis != input.end(); rawDigis++) {
150  SiStripPedestals::Range detPedestalRange = pedestalHandle_->getRange(rawDigis->id);
151  edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(rawDigis->id);
152 
153  for(uint16_t APV = 0; APV < rawDigis->size()/128; ++APV){
154  uint16_t MinPed =0;
155  for(uint16_t strip = APV*128; strip< (APV+1)*128; ++strip){
156  uint16_t ped = (uint16_t)pedestalHandle_->getPed(strip,detPedestalRange);
157  if(ped < MinPed) MinPed = ped;
158  }
159  if(MinPed>128) MinPed=128;
160  MeanCMDetSet.push_back(MinPed);
161  }
162 
163  meancm.push_back(MeanCMDetSet);
164  }
165 }
166 
168 
169  uint32_t detId;
170  CMMap::iterator itMap;
172 
173  for(itInput = Input.begin(); itInput != Input.end(); ++itInput){
174  detId = itInput->id;
175  itMap = _CMMap.find(detId);
177  std::vector<float> MeanCMNValue;
178  MeanCMNValue.clear();
179  if(itMap!=_CMMap.end()){ //the detId was already found
180  std::vector< float >& MapContent = itMap->second;
181  std::vector<float>::iterator itMapVector = MapContent.begin();
182  for(itCM = itInput->begin(); itCM != itInput->end(); ++itCM, ++itMapVector){
183  MeanCMNValue.push_back(itCM->adc() + *itMapVector);
184  }
185  _CMMap.erase(itMap);
186  _CMMap.insert(itMap, std::pair<uint32_t, std::vector<float> >(detId,MeanCMNValue));
187  } else { //no detId found
188  for(itCM = itInput->begin(); itCM != itInput->end(); ++itCM) MeanCMNValue.push_back(itCM->adc());
189  _CMMap.insert(std::pair<uint32_t, std::vector<float> >(detId,MeanCMNValue));
190  }
191  }
192 
193 }
194 
195 void
197  CMMap::iterator itMap;
198  std::vector<float>::const_iterator itMapVector;
199 
200  meancm.clear();
201  meancm.reserve(15000);
202 
203  for(itMap = _CMMap.begin(); itMap != _CMMap.end(); ++itMap){
204  edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(itMap->first);
205  for(itMapVector = (itMap->second).begin(); itMapVector != (itMap->second).end(); ++itMapVector) MeanCMDetSet.push_back(*itMapVector/(float)_actualEvent);
206  meancm.push_back(MeanCMDetSet);
207  }
208 
209 }
210 void
212 {
213  _actualEvent =1;
214 
215  _CMMap.clear();
216 
217 }
218 
219 void
221 
222 
223 }
224 
226 
virtual void produce(edm::Event &, const edm::EventSetup &) override
#define Input(cl)
Definition: vmac.h:188
void push_back(const T &t)
Definition: DetSet.h:68
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void ConvertMeanMapToDetSetVector(std::vector< edm::DetSet< SiStripProcessedRawDigi > > &)
std::pair< ContainerIterator, ContainerIterator > Range
edm::ESHandle< SiStripPedestals > pedestalHandle_
static std::string const input
Definition: EdmProvDump.cc:43
void init(const edm::EventSetup &)
void StoreMean(const edm::DetSetVector< SiStripProcessedRawDigi > &)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
void CMExtractorFromPedestals(const edm::DetSetVector< SiStripRawDigi > &, std::vector< edm::DetSet< SiStripProcessedRawDigi > > &)
#define end
Definition: vmac.h:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:365
tuple conf
Definition: dbtoconf.py:185
const T & get() const
Definition: EventSetup.h:56
std::map< uint32_t, std::vector< float > > CMMap
SiStripMeanCMExtractor(const edm::ParameterSet &)
virtual void endJob() override
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:350
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
virtual void beginJob() override