CMS 3D CMS Logo

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