CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Attributes
PreMixingSiPixelWorker Class Reference
Inheritance diagram for PreMixingSiPixelWorker:
PreMixingWorker

Public Member Functions

void addPileups (PileUpEventPrincipal const &, edm::EventSetup const &es) override
 
void addSignals (edm::Event const &e, edm::EventSetup const &es) override
 
void initializeEvent (edm::Event const &e, edm::EventSetup const &c) override
 
 PreMixingSiPixelWorker (const edm::ParameterSet &ps, edm::ProducesCollector, edm::ConsumesCollector &&iC)
 
void put (edm::Event &e, edm::EventSetup const &iSetup, std::vector< PileupSummaryInfo > const &ps, int bs) override
 
 ~PreMixingSiPixelWorker () override=default
 
- Public Member Functions inherited from PreMixingWorker
virtual void beginLuminosityBlock (edm::LuminosityBlock const &iLumi, edm::EventSetup const &iSetup)
 
virtual void beginRun (edm::Run const &iRun, edm::EventSetup const &iSetup)
 
virtual void endRun ()
 
virtual void finalizeBunchCrossing (edm::Event &iEvent, edm::EventSetup const &iSetup, int bunchCrossing)
 
virtual void initializeBunchCrossing (edm::Event const &iEvent, edm::EventSetup const &iSetup, int bunchCrossing)
 
 PreMixingWorker ()=default
 
virtual ~PreMixingWorker ()=default
 

Private Types

typedef int Amplitude
 
typedef std::multimap< int, PixelDigiOneDetectorMap
 
typedef std::map< uint32_t, OneDetectorMapSiGlobalIndex
 
typedef std::map< int, Amplitude, std::less< int > > signal_map_type
 
typedef std::map< uint32_t, signal_map_typesignalMaps
 

Private Attributes

SiPixelDigitizerAlgorithm digitizer_
 
bool firstFinalizeEvent_ = true
 
bool firstInitializeEvent_ = true
 
const std::string geometryType_
 
edm::ESHandle< TrackerGeometrypDD
 
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 41 of file PreMixingSiPixelWorker.cc.

Member Typedef Documentation

typedef int PreMixingSiPixelWorker::Amplitude
private

Definition at line 66 of file PreMixingSiPixelWorker.cc.

typedef std::multimap<int, PixelDigi> PreMixingSiPixelWorker::OneDetectorMap
private

Definition at line 71 of file PreMixingSiPixelWorker.cc.

typedef std::map<uint32_t, OneDetectorMap> PreMixingSiPixelWorker::SiGlobalIndex
private

Definition at line 72 of file PreMixingSiPixelWorker.cc.

typedef std::map<int, Amplitude, std::less<int> > PreMixingSiPixelWorker::signal_map_type
private

Definition at line 67 of file PreMixingSiPixelWorker.cc.

typedef std::map<uint32_t, signal_map_type> PreMixingSiPixelWorker::signalMaps
private

Definition at line 68 of file PreMixingSiPixelWorker.cc.

Constructor & Destructor Documentation

PreMixingSiPixelWorker::PreMixingSiPixelWorker ( const edm::ParameterSet ps,
edm::ProducesCollector  producesCollector,
edm::ConsumesCollector &&  iC 
)

Definition at line 83 of file PreMixingSiPixelWorker.cc.

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

86  : digitizer_(ps), geometryType_(ps.getParameter<std::string>("PixGeometryType")) {
87  // declare the products to produce
88 
89  pixeldigi_collectionSig_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig");
90  pixeldigi_collectionPile_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile");
91  PixelDigiCollectionDM_ = ps.getParameter<std::string>("PixelDigiCollectionDM");
92 
95 
98 
99  // clear local storage for this event
100  SiHitStorage_.clear();
101 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiPToken_
ProductRegistryHelper::BranchAliasSetterT< ProductType > produces()
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiToken_
SiPixelDigitizerAlgorithm digitizer_
PreMixingSiPixelWorker::~PreMixingSiPixelWorker ( )
overridedefault

Member Function Documentation

void PreMixingSiPixelWorker::addPileups ( PileUpEventPrincipal const &  pep,
edm::EventSetup const &  es 
)
overridevirtual

Implements PreMixingWorker.

Definition at line 145 of file PreMixingSiPixelWorker.cc.

References begin, PileUpEventPrincipal::bunchCrossing(), end, PileUpEventPrincipal::getByLabel(), edm::EventPrincipal::id(), input, edm::HandleBase::isValid(), LogDebug, pixeldigi_collectionPile_, PileUpEventPrincipal::principal(), and SiHitStorage_.

145  {
146  LogDebug("PreMixingSiPixelWorker") << "\n===============> adding pileups from event " << pep.principal().id()
147  << " for bunchcrossing " << pep.bunchCrossing();
148 
149  // fill in maps of hits; same code as addSignals, except now applied to the pileup events
150 
152  pep.getByLabel(pixeldigi_collectionPile_, inputHandle);
153 
154  if (inputHandle.isValid()) {
155  const auto& input = *inputHandle;
156 
157  //loop on all detsets (detectorIDs) inside the input collection
159  for (; DSViter != input.end(); DSViter++) {
160 #ifdef DEBUG
161  LogDebug("PreMixingSiPixelWorker") << "Pileups: Processing DetID " << DSViter->id;
162 #endif
163 
164  uint32_t detID = DSViter->id;
166  edm::DetSet<PixelDigi>::const_iterator end = (DSViter->data).end();
168 
169  // find correct local map (or new one) for this detector ID
170 
171  SiGlobalIndex::const_iterator itest;
172 
173  itest = SiHitStorage_.find(detID);
174 
175  if (itest != SiHitStorage_.end()) { // this detID already has hits, add to existing map
176 
177  OneDetectorMap LocalMap = itest->second;
178 
179  // fill in local map with extra channels
180  for (icopy = begin; icopy != end; icopy++) {
181  LocalMap.insert(OneDetectorMap::value_type((icopy->channel()), *icopy));
182  }
183 
184  SiHitStorage_[detID] = LocalMap;
185 
186  } else { // fill local storage with this information, put in global collection
187 
188  OneDetectorMap LocalMap;
189 
190  for (icopy = begin; icopy != end; icopy++) {
191  LocalMap.insert(OneDetectorMap::value_type((icopy->channel()), *icopy));
192  }
193 
194  SiHitStorage_.insert(SiGlobalIndex::value_type(detID, LocalMap));
195  }
196  }
197  }
198 }
#define LogDebug(id)
std::multimap< int, PixelDigi > OneDetectorMap
static std::string const input
Definition: EdmProvDump.cc:48
#define end
Definition: vmac.h:39
bool isValid() const
Definition: HandleBase.h:70
#define begin
Definition: vmac.h:32
collection_type::const_iterator const_iterator
Definition: DetSet.h:32
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
void PreMixingSiPixelWorker::addSignals ( edm::Event const &  e,
edm::EventSetup const &  es 
)
overridevirtual

Implements PreMixingWorker.

Definition at line 114 of file PreMixingSiPixelWorker.cc.

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

114  {
115  // fill in maps of hits
116 
117  LogDebug("PreMixingSiPixelWorker") << "===============> adding MC signals for " << e.id();
118 
120 
121  if (e.getByToken(PixelDigiToken_, input)) {
122  //loop on all detsets (detectorIDs) inside the input collection
123  edm::DetSetVector<PixelDigi>::const_iterator DSViter = input->begin();
124  for (; DSViter != input->end(); DSViter++) {
125 #ifdef DEBUG
126  LogDebug("PreMixingSiPixelWorker") << "Processing DetID " << DSViter->id;
127 #endif
128 
129  uint32_t detID = DSViter->id;
131  edm::DetSet<PixelDigi>::const_iterator end = (DSViter->data).end();
133 
134  OneDetectorMap LocalMap;
135 
136  for (icopy = begin; icopy != end; icopy++) {
137  LocalMap.insert(OneDetectorMap::value_type((icopy->channel()), *icopy));
138  }
139 
140  SiHitStorage_.insert(SiGlobalIndex::value_type(detID, LocalMap));
141  }
142  }
143 } // end of addSiPixelSignals
#define LogDebug(id)
std::multimap< int, PixelDigi > OneDetectorMap
static std::string const input
Definition: EdmProvDump.cc:48
edm::EDGetTokenT< edm::DetSetVector< PixelDigi > > PixelDigiToken_
#define end
Definition: vmac.h:39
#define begin
Definition: vmac.h:32
collection_type::const_iterator const_iterator
Definition: DetSet.h:32
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
void PreMixingSiPixelWorker::initializeEvent ( edm::Event const &  e,
edm::EventSetup const &  c 
)
overridevirtual
void PreMixingSiPixelWorker::put ( edm::Event e,
edm::EventSetup const &  iSetup,
std::vector< PileupSummaryInfo > const &  ps,
int  bs 
)
overridevirtual

Implements PreMixingWorker.

Definition at line 200 of file PreMixingSiPixelWorker.cc.

References ecalMGPA::adc(), SiPixelDigitizerAlgorithm::calculateInstlumiFactor(), SiPixelDigitizerAlgorithm::chooseScenario(), DEFINE_PREMIXING_WORKER, TrackerGeometry::detUnits(), SiPixelDigitizerAlgorithm::digitize(), digitizer_, Exception, firstFinalizeEvent_, edm::EventSetup::get(), edm::RandomNumberGenerator::getEngine(), SiPixelDigitizerAlgorithm::init_DynIneffDB(), SiPixelDigitizerAlgorithm::killBadFEDChannels(), eostools::move(), pDD, PixelDigiCollectionDM_, edm::ESHandle< T >::product(), edm::Event::put(), SiPixelDigitizerAlgorithm::setSimAccumulator(), SiHitStorage_, and edm::Event::streamID().

203  {
204  // collection of Digis to put in the event
205 
206  std::vector<edm::DetSet<PixelDigi>> vPixelDigi;
207 
208  // loop through our collection of detectors, merging hits and putting new ones in the output
209  signalMaps signal;
210 
211  // big loop over Detector IDs:
212 
213  for (SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin(); IDet != SiHitStorage_.end(); IDet++) {
214  uint32_t detID = IDet->first;
215 
216  OneDetectorMap LocalMap = IDet->second;
217 
218  signal_map_type Signals;
219  Signals.clear();
220 
221  //counter variables
222  int formerPixel = -1;
223  int currentPixel;
224  int ADCSum = 0;
225 
226  OneDetectorMap::const_iterator iLocalchk;
227 
228  for (OneDetectorMap::const_iterator iLocal = LocalMap.begin(); iLocal != LocalMap.end(); ++iLocal) {
229  currentPixel = iLocal->first;
230 
231  if (currentPixel == formerPixel) { // we have to add these digis together
232  ADCSum += (iLocal->second).adc();
233  } else {
234  if (formerPixel != -1) { // ADC info stolen from SiStrips...
235  if (ADCSum > 511)
236  ADCSum = 255;
237  else if (ADCSum > 253 && ADCSum < 512)
238  ADCSum = 254;
239 
240  Signals.insert(std::make_pair(formerPixel, ADCSum));
241  }
242  // save pointers for next iteration
243  formerPixel = currentPixel;
244  ADCSum = (iLocal->second).adc();
245  }
246 
247  iLocalchk = iLocal;
248  if ((++iLocalchk) == LocalMap.end()) { //make sure not to lose the last one
249  if (ADCSum > 511)
250  ADCSum = 255;
251  else if (ADCSum > 253 && ADCSum < 512)
252  ADCSum = 254;
253  Signals.insert(std::make_pair(formerPixel, ADCSum));
254  }
255 
256  } // end of loop over one detector
257 
258  // stick this into the global vector of detector info
259  signal.insert(std::make_pair(detID, Signals));
260 
261  } // end of big loop over all detector IDs
262 
263  // put the collection of digis in the event
264  edm::LogInfo("PreMixingSiPixelWorker") << "total # Merged Pixels: " << signal.size();
265 
266  std::vector<edm::DetSet<PixelDigi>> theDigiVector;
267 
268  // Load inefficiency constants (1st pass), set pileup information.
269  if (firstFinalizeEvent_) {
270  digitizer_.init_DynIneffDB(iSetup, bs);
271  firstFinalizeEvent_ = false;
272  }
273 
276 
278  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
279 
281  iSetup.get<TrackerTopologyRcd>().get(tTopoHand);
282  const TrackerTopology* tTopo = tTopoHand.product();
283 
285  std::unique_ptr<PixelFEDChannelCollection> PixelFEDChannelCollection_ = digitizer_.chooseScenario(ps, engine);
286  if (PixelFEDChannelCollection_ == nullptr) {
287  throw cms::Exception("NullPointerError") << "PixelFEDChannelCollection not set in chooseScenario function.\n";
288  }
289  e.put(std::move(PixelFEDChannelCollection_), PixelDigiCollectionDM_);
290  }
291 
292  for (const auto& iu : pDD->detUnits()) {
293  if (iu->type().isTrackerPixel()) {
294  edm::DetSet<PixelDigi> collector(iu->geographicalId().rawId());
295  edm::DetSet<PixelDigiSimLink> linkcollector(
296  iu->geographicalId().rawId()); // ignored as DigiSimLinks are combined separately
297 
298  digitizer_.digitize(dynamic_cast<const PixelGeomDetUnit*>(iu), collector.data, linkcollector.data, tTopo, engine);
299  if (!collector.data.empty()) {
300  theDigiVector.push_back(std::move(collector));
301  }
302  }
303  }
304 
305  e.put(std::make_unique<edm::DetSetVector<PixelDigi>>(theDigiVector), PixelDigiCollectionDM_);
306 
307  // clear local storage for this event
308  SiHitStorage_.clear();
309 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
void init_DynIneffDB(const edm::EventSetup &, const unsigned int &)
std::multimap< int, PixelDigi > OneDetectorMap
const DetContainer & detUnits() const override
Returm a vector of all GeomDet.
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
SiPixelDigitizerAlgorithm digitizer_
void digitize(const PixelGeomDetUnit *pixdet, std::vector< PixelDigi > &digis, std::vector< PixelDigiSimLink > &simlinks, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
std::unique_ptr< PixelFEDChannelCollection > chooseScenario(PileupMixingContent *puInfo, CLHEP::HepRandomEngine *)
void setSimAccumulator(const std::map< uint32_t, std::map< int, int > > &signalMap)
void calculateInstlumiFactor(PileupMixingContent *puInfo)
std::map< uint32_t, signal_map_type > signalMaps
edm::ESHandle< TrackerGeometry > pDD
StreamID streamID() const
Definition: Event.h:96
std::map< int, Amplitude, std::less< int > > signal_map_type
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

SiPixelDigitizerAlgorithm PreMixingSiPixelWorker::digitizer_
private

Definition at line 61 of file PreMixingSiPixelWorker.cc.

Referenced by initializeEvent(), and put().

bool PreMixingSiPixelWorker::firstFinalizeEvent_ = true
private

Definition at line 79 of file PreMixingSiPixelWorker.cc.

Referenced by put().

bool PreMixingSiPixelWorker::firstInitializeEvent_ = true
private

Definition at line 78 of file PreMixingSiPixelWorker.cc.

Referenced by initializeEvent().

const std::string PreMixingSiPixelWorker::geometryType_
private

Definition at line 76 of file PreMixingSiPixelWorker.cc.

Referenced by initializeEvent().

edm::ESHandle<TrackerGeometry> PreMixingSiPixelWorker::pDD
private

Definition at line 59 of file PreMixingSiPixelWorker.cc.

Referenced by initializeEvent(), and put().

edm::InputTag PreMixingSiPixelWorker::pixeldigi_collectionPile_
private

Definition at line 53 of file PreMixingSiPixelWorker.cc.

Referenced by addPileups(), and PreMixingSiPixelWorker().

edm::InputTag PreMixingSiPixelWorker::pixeldigi_collectionSig_
private

Definition at line 52 of file PreMixingSiPixelWorker.cc.

Referenced by PreMixingSiPixelWorker().

std::string PreMixingSiPixelWorker::PixelDigiCollectionDM_
private

Definition at line 54 of file PreMixingSiPixelWorker.cc.

Referenced by PreMixingSiPixelWorker(), and put().

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

Definition at line 57 of file PreMixingSiPixelWorker.cc.

Referenced by PreMixingSiPixelWorker().

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

Definition at line 56 of file PreMixingSiPixelWorker.cc.

Referenced by addSignals(), and PreMixingSiPixelWorker().

SiGlobalIndex PreMixingSiPixelWorker::SiHitStorage_
private

Definition at line 74 of file PreMixingSiPixelWorker.cc.

Referenced by addPileups(), addSignals(), PreMixingSiPixelWorker(), and put().