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 // $Id: SiStripMeanCMExtractor.cc,v 1.2 2010/11/04 14:12:31 edwenger Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <sstream>
23 #include <memory>
24 #include <list>
25 #include <algorithm>
26 #include <cassert>
27 
28 
29 // user include files
37 
40 
42 
50 
51 
52 typedef std::map<uint32_t, std::vector<float> > CMMap;
53 
55  public:
56  explicit SiStripMeanCMExtractor( const edm::ParameterSet&);
58 
59  private:
60  virtual void beginJob() ;
61  virtual void produce(edm::Event&, const edm::EventSetup&);
62  virtual void endJob() ;
63 
64  void init(const edm::EventSetup&);
65 
68 
73  std::string _Algorithm;
74  uint16_t _nEventsToUse;
75  uint16_t _actualEvent;
76 
77  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.
78 };
79 
80 
82  _inputTag(conf.getParameter<edm::InputTag> ("CMCollection")),
83  _Algorithm(conf.getParameter<std::string>("Algorithm")),
84  _nEventsToUse(conf.getParameter<uint32_t>("NEvents"))
85 {
86 
87 
88  if(_nEventsToUse < 1) _nEventsToUse=1;
89  produces< edm::DetSetVector<SiStripProcessedRawDigi> > ("MEANAPVCM");
90 }
91 
92 
94 {
95 
96 }
97 
99 
100  uint32_t p_cache_id = es.get<SiStripPedestalsRcd>().cacheIdentifier();
101 
102  if(p_cache_id != pedestal_cache_id_) {
104  pedestal_cache_id_ = p_cache_id;
105  }
106 }
107 
108 // ------------ method called to produce the data ------------
109 void
111 {
112  using namespace edm;
113 
114 
115  //if(_actualEvent > _nEventsToUse) return;
116 
117  std::vector<edm::DetSet<SiStripProcessedRawDigi> > meancm;
118 
119  if(_Algorithm == "StoredCM"){
121  iEvent.getByLabel(_inputTag,inputCM);
122 
123  this->StoreMean(*inputCM);
124  this->ConvertMeanMapToDetSetVector(meancm);
125 
126  } else if (_Algorithm == "Pedestals"){
127  this->init(iSetup);
128 
130  iEvent.getByLabel(_inputTag,input);
131 
132  this->CMExtractorFromPedestals(*input,meancm);
133  }
134 
135  ++_actualEvent;
136 
137 
138 
139 
140  std::auto_ptr< edm::DetSetVector<SiStripProcessedRawDigi> > outputMeanCM(new edm::DetSetVector<SiStripProcessedRawDigi>(meancm) );
141  iEvent.put( outputMeanCM,"MEANAPVCM");
142 
143 }
144 
146  meancm.clear();
147  meancm.reserve(15000);
148 
150  rawDigis = input.begin(); rawDigis != input.end(); rawDigis++) {
151  SiStripPedestals::Range detPedestalRange = pedestalHandle_->getRange(rawDigis->id);
152  edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(rawDigis->id);
153 
154  for(uint16_t APV = 0; APV < rawDigis->size()/128; ++APV){
155  uint16_t MinPed =0;
156  for(uint16_t strip = APV*128; strip< (APV+1)*128; ++strip){
157  uint16_t ped = (uint16_t)pedestalHandle_->getPed(strip,detPedestalRange);
158  if(ped < MinPed) MinPed = ped;
159  }
160  if(MinPed>128) MinPed=128;
161  MeanCMDetSet.push_back(MinPed);
162  }
163 
164  meancm.push_back(MeanCMDetSet);
165  }
166 }
167 
169 
170  uint32_t detId;
171  CMMap::iterator itMap;
173 
174  for(itInput = Input.begin(); itInput != Input.end(); ++itInput){
175  detId = itInput->id;
176  itMap = _CMMap.find(detId);
178  std::vector<float> MeanCMNValue;
179  MeanCMNValue.clear();
180  if(itMap!=_CMMap.end()){ //the detId was already found
181  std::vector< float >& MapContent = itMap->second;
182  std::vector<float>::iterator itMapVector = MapContent.begin();
183  for(itCM = itInput->begin(); itCM != itInput->end(); ++itCM, ++itMapVector){
184  MeanCMNValue.push_back(itCM->adc() + *itMapVector);
185  }
186  _CMMap.erase(itMap);
187  _CMMap.insert(itMap, std::pair<uint32_t, std::vector<float> >(detId,MeanCMNValue));
188  } else { //no detId found
189  for(itCM = itInput->begin(); itCM != itInput->end(); ++itCM) MeanCMNValue.push_back(itCM->adc());
190  _CMMap.insert(std::pair<uint32_t, std::vector<float> >(detId,MeanCMNValue));
191  }
192  }
193 
194 }
195 
196 void
198  CMMap::iterator itMap;
199  std::vector<float>::const_iterator itMapVector;
200 
201  meancm.clear();
202  meancm.reserve(15000);
203 
204  for(itMap = _CMMap.begin(); itMap != _CMMap.end(); ++itMap){
205  edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(itMap->first);
206  for(itMapVector = (itMap->second).begin(); itMapVector != (itMap->second).end(); ++itMapVector) MeanCMDetSet.push_back(*itMapVector/(float)_actualEvent);
207  meancm.push_back(MeanCMDetSet);
208  }
209 
210 }
211 void
213 {
214  _actualEvent =1;
215 
216  _CMMap.clear();
217 
218 }
219 
220 void
222 
223 
224 }
225 
227 
#define Input(cl)
Definition: vmac.h:189
void push_back(const T &t)
Definition: DetSet.h:69
void strip(std::string &input, const std::string &blanks=" \n\t")
Definition: stringTools.cc:16
#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_
void init(const edm::EventSetup &)
void StoreMean(const edm::DetSetVector< SiStripProcessedRawDigi > &)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
void CMExtractorFromPedestals(const edm::DetSetVector< SiStripRawDigi > &, std::vector< edm::DetSet< SiStripProcessedRawDigi > > &)
#define end
Definition: vmac.h:38
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:356
tuple conf
Definition: dbtoconf.py:185
const T & get() const
Definition: EventSetup.h:55
std::map< uint32_t, std::vector< float > > CMMap
SiStripMeanCMExtractor(const edm::ParameterSet &)
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:341
collection_type::const_iterator const_iterator
Definition: DetSet.h:34
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
virtual void produce(edm::Event &, const edm::EventSetup &)