CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DataMixingSiPixelWorker.cc
Go to the documentation of this file.
1 // File: DataMixingSiPixelWorker.cc
2 // Description: see DataMixingSiPixelWorker.h
3 // Author: Mike Hildreth, University of Notre Dame
4 //
5 //--------------------------------------------
6 
7 #include <map>
15 //
16 //
18 
19 
20 using namespace std;
21 
22 namespace edm
23 {
24 
25  // Virtual constructor
26 
27  DataMixingSiPixelWorker::DataMixingSiPixelWorker() { }
28 
29  // Constructor
30  DataMixingSiPixelWorker::DataMixingSiPixelWorker(const edm::ParameterSet& ps, edm::ConsumesCollector && iC) :
31  label_(ps.getParameter<std::string>("Label"))
32 
33  {
34 
35  // get the subdetector names
36  // this->getSubdetectorNames(); //something like this may be useful to check what we are supposed to do...
37 
38  // declare the products to produce
39 
40  pixeldigi_collectionSig_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig");
41  pixeldigi_collectionPile_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile");
42  PixelDigiCollectionDM_ = ps.getParameter<std::string>("PixelDigiCollectionDM");
43 
46 
47 
48  // clear local storage for this event
49  SiHitStorage_.clear();
50 
51  }
52 
53 
54  // Virtual destructor needed.
56  }
57 
58 
59 
61  // fill in maps of hits
62 
63  LogDebug("DataMixingSiPixelWorker")<<"===============> adding MC signals for "<<e.id();
64 
66 
67  if( e.getByToken(PixelDigiToken_,input) ) {
68 
69  //loop on all detsets (detectorIDs) inside the input collection
70  edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
71  for (; DSViter!=input->end();DSViter++){
72 
73 #ifdef DEBUG
74  LogDebug("DataMixingSiPixelWorker") << "Processing DetID " << DSViter->id;
75 #endif
76 
77  uint32_t detID = DSViter->id;
81 
82  OneDetectorMap LocalMap;
83 
84  for (icopy=begin; icopy!=end; icopy++) {
85  LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
86  }
87 
88  SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
89  }
90 
91  }
92  } // end of addSiPixelSignals
93 
94 
95 
96  void DataMixingSiPixelWorker::addSiPixelPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr,
97  ModuleCallingContext const* mcc) {
98 
99  LogDebug("DataMixingSiPixelWorker") <<"\n===============> adding pileups from event "<<ep->id()<<" for bunchcrossing "<<bcr;
100 
101  // fill in maps of hits; same code as addSignals, except now applied to the pileup events
102 
103  boost::shared_ptr<Wrapper<edm::DetSetVector<PixelDigi> > const> inputPTR =
104  getProductByTag<edm::DetSetVector<PixelDigi> >(*ep, pixeldigi_collectionPile_, mcc);
105 
106  if(inputPTR ) {
107 
108  const edm::DetSetVector<PixelDigi> *input = const_cast< edm::DetSetVector<PixelDigi> * >(inputPTR->product());
109 
110 
111 
112  // Handle< edm::DetSetVector<PixelDigi> > input;
113 
114  // if( e->getByLabel(pixeldigi_collectionPile_,input) ) {
115 
116  //loop on all detsets (detectorIDs) inside the input collection
118  for (; DSViter!=input->end();DSViter++){
119 
120 #ifdef DEBUG
121  LogDebug("DataMixingSiPixelWorker") << "Pileups: Processing DetID " << DSViter->id;
122 #endif
123 
124  uint32_t detID = DSViter->id;
128 
129  // find correct local map (or new one) for this detector ID
130 
131  SiGlobalIndex::const_iterator itest;
132 
133  itest = SiHitStorage_.find(detID);
134 
135  if(itest!=SiHitStorage_.end()) { // this detID already has hits, add to existing map
136 
137  OneDetectorMap LocalMap = itest->second;
138 
139  // fill in local map with extra channels
140  for (icopy=begin; icopy!=end; icopy++) {
141  LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
142  }
143 
144  SiHitStorage_[detID]=LocalMap;
145 
146  }
147  else{ // fill local storage with this information, put in global collection
148 
149  OneDetectorMap LocalMap;
150 
151  for (icopy=begin; icopy!=end; icopy++) {
152  LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
153  }
154 
155  SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
156  }
157 
158  }
159  }
160  }
161 
162 
163 
165 
166  // collection of Digis to put in the event
167 
168  std::vector< edm::DetSet<PixelDigi> > vPixelDigi;
169 
170  // loop through our collection of detectors, merging hits and putting new ones in the output
171 
172  // big loop over Detector IDs:
173 
174  for(SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin();
175  IDet != SiHitStorage_.end(); IDet++) {
176 
177  edm::DetSet<PixelDigi> SPD(IDet->first); // Make empty collection with this detector ID
178 
179  OneDetectorMap LocalMap = IDet->second;
180 
181  //counter variables
182  int formerPixel = -1;
183  int currentPixel;
184  int ADCSum = 0;
185 
186 
187  OneDetectorMap::const_iterator iLocalchk;
188 
189  for(OneDetectorMap::const_iterator iLocal = LocalMap.begin();
190  iLocal != LocalMap.end(); ++iLocal) {
191 
192  currentPixel = iLocal->first;
193 
194  if (currentPixel == formerPixel) { // we have to add these digis together
195  ADCSum+=(iLocal->second).adc();
196  }
197  else{
198  if(formerPixel!=-1){ // ADC info stolen from SiStrips...
199  if (ADCSum > 511) ADCSum = 255;
200  else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
201  PixelDigi aHit(formerPixel, ADCSum);
202  SPD.push_back( aHit );
203  }
204  // save pointers for next iteration
205  formerPixel = currentPixel;
206  ADCSum = (iLocal->second).adc();
207  }
208 
209  iLocalchk = iLocal;
210  if((++iLocalchk) == LocalMap.end()) { //make sure not to lose the last one
211  if (ADCSum > 511) ADCSum = 255;
212  else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
213  SPD.push_back( PixelDigi(formerPixel, ADCSum) );
214  }
215 
216  }// end of loop over one detector
217 
218  // stick this into the global vector of detector info
219  vPixelDigi.push_back(SPD);
220 
221  } // end of big loop over all detector IDs
222 
223  // put the collection of digis in the event
224  LogInfo("DataMixingSiPixelWorker") << "total # Merged Pixels: " << vPixelDigi.size() ;
225 
226  // make new digi collection
227 
228  std::auto_ptr< edm::DetSetVector<PixelDigi> > MyPixelDigis(new edm::DetSetVector<PixelDigi>(vPixelDigi) );
229 
230  // put collection
231 
232  e.put( MyPixelDigis, PixelDigiCollectionDM_ );
233 
234  // clear local storage for this event
235  SiHitStorage_.clear();
236  }
237 
238 } //edm
int adc(sample_type sample)
get the ADC sample (12 bits)
#define LogDebug(id)
T getParameter(std::string const &) const
void push_back(const T &t)
Definition: DetSet.h:68
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
EventID const & id() const
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiToken_
void addSiPixelPileups(const int bcr, const edm::EventPrincipal *, unsigned int EventId, ModuleCallingContext const *)
static std::string const input
Definition: EdmProvDump.cc:44
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiPToken_
#define end
Definition: vmac.h:37
Container::value_type value_type
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:363
std::multimap< int, PixelDigi > OneDetectorMap
edm::EventID id() const
Definition: EventBase.h:56
#define begin
Definition: vmac.h:30
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
void addSiPixelSignals(const edm::Event &e)