CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Attributes
edm::DataMixingSiPixelWorker Class Reference

#include <DataMixingSiPixelWorker.h>

Public Member Functions

void addSiPixelPileups (const int bcr, const edm::EventPrincipal *, unsigned int EventId, ModuleCallingContext const *)
 
void addSiPixelSignals (const edm::Event &e)
 
 DataMixingSiPixelWorker ()
 
 DataMixingSiPixelWorker (const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
 
void putSiPixel (edm::Event &e)
 
virtual ~DataMixingSiPixelWorker ()
 

Private Types

typedef std::multimap< int,
PixelDigi
OneDetectorMap
 
typedef std::map< uint32_t,
OneDetectorMap
SiGlobalIndex
 

Private Attributes

std::string label_
 
edm::InputTag pixeldigi_collectionPile_
 
edm::InputTag pixeldigi_collectionSig_
 
std::string PixelDigiCollectionDM_
 
edm::EDGetTokenT
< edm::DetSetVector< PixelDigi > > 
PixelDigiPToken_
 
edm::EDGetTokenT
< edm::DetSetVector< PixelDigi > > 
PixelDigiToken_
 
SiGlobalIndex SiHitStorage_
 

Detailed Description

Definition at line 38 of file DataMixingSiPixelWorker.h.

Member Typedef Documentation

typedef std::multimap<int, PixelDigi> edm::DataMixingSiPixelWorker::OneDetectorMap
private

Definition at line 67 of file DataMixingSiPixelWorker.h.

typedef std::map<uint32_t, OneDetectorMap> edm::DataMixingSiPixelWorker::SiGlobalIndex
private

Definition at line 68 of file DataMixingSiPixelWorker.h.

Constructor & Destructor Documentation

DataMixingSiPixelWorker::DataMixingSiPixelWorker ( )

Definition at line 28 of file DataMixingSiPixelWorker.cc.

28 { }
DataMixingSiPixelWorker::DataMixingSiPixelWorker ( const edm::ParameterSet ps,
edm::ConsumesCollector &&  iC 
)
explicit

standard constructor

Definition at line 31 of file DataMixingSiPixelWorker.cc.

References edm::ParameterSet::getParameter(), pixeldigi_collectionPile_, pixeldigi_collectionSig_, PixelDigiCollectionDM_, PixelDigiPToken_, PixelDigiToken_, SiHitStorage_, and AlCaHLTBitMon_QueryRunRegistry::string.

31  :
32  label_(ps.getParameter<std::string>("Label"))
33 
34  {
35 
36  // get the subdetector names
37  // this->getSubdetectorNames(); //something like this may be useful to check what we are supposed to do...
38 
39  // declare the products to produce
40 
41  pixeldigi_collectionSig_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig");
42  pixeldigi_collectionPile_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile");
43  PixelDigiCollectionDM_ = ps.getParameter<std::string>("PixelDigiCollectionDM");
44 
47 
48 
49  // clear local storage for this event
50  SiHitStorage_.clear();
51 
52  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiToken_
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiPToken_
DataMixingSiPixelWorker::~DataMixingSiPixelWorker ( )
virtual

Default destructor

Definition at line 56 of file DataMixingSiPixelWorker.cc.

56  {
57  }

Member Function Documentation

void DataMixingSiPixelWorker::addSiPixelPileups ( const int  bcr,
const edm::EventPrincipal ep,
unsigned int  EventId,
ModuleCallingContext const *  mcc 
)

Definition at line 97 of file DataMixingSiPixelWorker.cc.

References begin, edm::DetSetVector< T >::begin(), end, edm::DetSetVector< T >::end(), edm::EventPrincipal::id(), input, LogDebug, pixeldigi_collectionPile_, and SiHitStorage_.

Referenced by edm::DataMixingModule::pileWorker().

98  {
99 
100  LogDebug("DataMixingSiPixelWorker") <<"\n===============> adding pileups from event "<<ep->id()<<" for bunchcrossing "<<bcr;
101 
102  // fill in maps of hits; same code as addSignals, except now applied to the pileup events
103 
104  std::shared_ptr<Wrapper<edm::DetSetVector<PixelDigi> > const> inputPTR =
105  getProductByTag<edm::DetSetVector<PixelDigi> >(*ep, pixeldigi_collectionPile_, mcc);
106 
107  if(inputPTR ) {
108 
109  const edm::DetSetVector<PixelDigi> *input = const_cast< edm::DetSetVector<PixelDigi> * >(inputPTR->product());
110 
111 
112 
113  // Handle< edm::DetSetVector<PixelDigi> > input;
114 
115  // if( e->getByLabel(pixeldigi_collectionPile_,input) ) {
116 
117  //loop on all detsets (detectorIDs) inside the input collection
119  for (; DSViter!=input->end();DSViter++){
120 
121 #ifdef DEBUG
122  LogDebug("DataMixingSiPixelWorker") << "Pileups: Processing DetID " << DSViter->id;
123 #endif
124 
125  uint32_t detID = DSViter->id;
129 
130  // find correct local map (or new one) for this detector ID
131 
132  SiGlobalIndex::const_iterator itest;
133 
134  itest = SiHitStorage_.find(detID);
135 
136  if(itest!=SiHitStorage_.end()) { // this detID already has hits, add to existing map
137 
138  OneDetectorMap LocalMap = itest->second;
139 
140  // fill in local map with extra channels
141  for (icopy=begin; icopy!=end; icopy++) {
142  LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
143  }
144 
145  SiHitStorage_[detID]=LocalMap;
146 
147  }
148  else{ // fill local storage with this information, put in global collection
149 
150  OneDetectorMap LocalMap;
151 
152  for (icopy=begin; icopy!=end; icopy++) {
153  LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
154  }
155 
156  SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
157  }
158 
159  }
160  }
161  }
#define LogDebug(id)
EventID const & id() const
static std::string const input
Definition: EdmProvDump.cc:43
#define end
Definition: vmac.h:37
Container::value_type value_type
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:365
std::multimap< int, PixelDigi > OneDetectorMap
#define begin
Definition: vmac.h:30
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
void DataMixingSiPixelWorker::addSiPixelSignals ( const edm::Event e)

Definition at line 61 of file DataMixingSiPixelWorker.cc.

References begin, end, edm::Event::getByToken(), edm::EventBase::id(), input, LogDebug, PixelDigiToken_, and SiHitStorage_.

Referenced by edm::DataMixingModule::addSignals().

61  {
62  // fill in maps of hits
63 
64  LogDebug("DataMixingSiPixelWorker")<<"===============> adding MC signals for "<<e.id();
65 
66  Handle< edm::DetSetVector<PixelDigi> > input;
67 
68  if( e.getByToken(PixelDigiToken_,input) ) {
69 
70  //loop on all detsets (detectorIDs) inside the input collection
71  edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
72  for (; DSViter!=input->end();DSViter++){
73 
74 #ifdef DEBUG
75  LogDebug("DataMixingSiPixelWorker") << "Processing DetID " << DSViter->id;
76 #endif
77 
78  uint32_t detID = DSViter->id;
82 
83  OneDetectorMap LocalMap;
84 
85  for (icopy=begin; icopy!=end; icopy++) {
86  LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
87  }
88 
89  SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
90  }
91 
92  }
93  } // end of addSiPixelSignals
#define LogDebug(id)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiToken_
static std::string const input
Definition: EdmProvDump.cc:43
#define end
Definition: vmac.h:37
Container::value_type value_type
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:365
std::multimap< int, PixelDigi > OneDetectorMap
edm::EventID id() const
Definition: EventBase.h:60
#define begin
Definition: vmac.h:30
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
void DataMixingSiPixelWorker::putSiPixel ( edm::Event e)

Definition at line 165 of file DataMixingSiPixelWorker.cc.

References ecalMGPA::adc(), PixelDigiCollectionDM_, edm::DetSet< T >::push_back(), edm::Event::put(), and SiHitStorage_.

Referenced by edm::DataMixingModule::put().

165  {
166 
167  // collection of Digis to put in the event
168 
169  std::vector< edm::DetSet<PixelDigi> > vPixelDigi;
170 
171  // loop through our collection of detectors, merging hits and putting new ones in the output
172 
173  // big loop over Detector IDs:
174 
175  for(SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin();
176  IDet != SiHitStorage_.end(); IDet++) {
177 
178  edm::DetSet<PixelDigi> SPD(IDet->first); // Make empty collection with this detector ID
179 
180  OneDetectorMap LocalMap = IDet->second;
181 
182  //counter variables
183  int formerPixel = -1;
184  int currentPixel;
185  int ADCSum = 0;
186 
187 
188  OneDetectorMap::const_iterator iLocalchk;
189 
190  for(OneDetectorMap::const_iterator iLocal = LocalMap.begin();
191  iLocal != LocalMap.end(); ++iLocal) {
192 
193  currentPixel = iLocal->first;
194 
195  if (currentPixel == formerPixel) { // we have to add these digis together
196  ADCSum+=(iLocal->second).adc();
197  }
198  else{
199  if(formerPixel!=-1){ // ADC info stolen from SiStrips...
200  if (ADCSum > 511) ADCSum = 255;
201  else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
202  PixelDigi aHit(formerPixel, ADCSum);
203  SPD.push_back( aHit );
204  }
205  // save pointers for next iteration
206  formerPixel = currentPixel;
207  ADCSum = (iLocal->second).adc();
208  }
209 
210  iLocalchk = iLocal;
211  if((++iLocalchk) == LocalMap.end()) { //make sure not to lose the last one
212  if (ADCSum > 511) ADCSum = 255;
213  else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
214  SPD.push_back( PixelDigi(formerPixel, ADCSum) );
215  }
216 
217  }// end of loop over one detector
218 
219  // stick this into the global vector of detector info
220  vPixelDigi.push_back(SPD);
221 
222  } // end of big loop over all detector IDs
223 
224  // put the collection of digis in the event
225  LogInfo("DataMixingSiPixelWorker") << "total # Merged Pixels: " << vPixelDigi.size() ;
226 
227  // make new digi collection
228 
229  std::auto_ptr< edm::DetSetVector<PixelDigi> > MyPixelDigis(new edm::DetSetVector<PixelDigi>(vPixelDigi) );
230 
231  // put collection
232 
233  e.put( MyPixelDigis, PixelDigiCollectionDM_ );
234 
235  // clear local storage for this event
236  SiHitStorage_.clear();
237  }
int adc(sample_type sample)
get the ADC sample (12 bits)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
std::multimap< int, PixelDigi > OneDetectorMap

Member Data Documentation

std::string edm::DataMixingSiPixelWorker::label_
private
edm::InputTag edm::DataMixingSiPixelWorker::pixeldigi_collectionPile_
private

Definition at line 59 of file DataMixingSiPixelWorker.h.

Referenced by addSiPixelPileups(), and DataMixingSiPixelWorker().

edm::InputTag edm::DataMixingSiPixelWorker::pixeldigi_collectionSig_
private

Definition at line 58 of file DataMixingSiPixelWorker.h.

Referenced by DataMixingSiPixelWorker().

std::string edm::DataMixingSiPixelWorker::PixelDigiCollectionDM_
private

Definition at line 60 of file DataMixingSiPixelWorker.h.

Referenced by DataMixingSiPixelWorker(), and putSiPixel().

edm::EDGetTokenT<edm::DetSetVector<PixelDigi> > edm::DataMixingSiPixelWorker::PixelDigiPToken_
private

Definition at line 63 of file DataMixingSiPixelWorker.h.

Referenced by DataMixingSiPixelWorker().

edm::EDGetTokenT<edm::DetSetVector<PixelDigi> > edm::DataMixingSiPixelWorker::PixelDigiToken_
private

Definition at line 62 of file DataMixingSiPixelWorker.h.

Referenced by addSiPixelSignals(), and DataMixingSiPixelWorker().

SiGlobalIndex edm::DataMixingSiPixelWorker::SiHitStorage_
private