#include <DataMixingSiPixelWorker.h>
Public Member Functions | |
void | addSiPixelPileups (const int bcr, edm::EventPrincipal *, unsigned int EventId) |
void | addSiPixelSignals (const edm::Event &e) |
DataMixingSiPixelWorker (const edm::ParameterSet &ps) | |
DataMixingSiPixelWorker () | |
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_ |
Selector * | sel_ |
SiGlobalIndex | SiHitStorage_ |
Definition at line 36 of file DataMixingSiPixelWorker.h.
typedef std::multimap<int, PixelDigi> edm::DataMixingSiPixelWorker::OneDetectorMap [private] |
Definition at line 62 of file DataMixingSiPixelWorker.h.
typedef std::map<uint32_t, OneDetectorMap> edm::DataMixingSiPixelWorker::SiGlobalIndex [private] |
Definition at line 63 of file DataMixingSiPixelWorker.h.
DataMixingSiPixelWorker::DataMixingSiPixelWorker | ( | ) |
Definition at line 27 of file DataMixingSiPixelWorker.cc.
{ sel_=0;}
DataMixingSiPixelWorker::DataMixingSiPixelWorker | ( | const edm::ParameterSet & | ps | ) | [explicit] |
standard constructor
Definition at line 30 of file DataMixingSiPixelWorker.cc.
References edm::ParameterSet::getParameter(), label_, pixeldigi_collectionPile_, pixeldigi_collectionSig_, PixelDigiCollectionDM_, sel_, and SiHitStorage_.
: label_(ps.getParameter<std::string>("Label")) { // get the subdetector names // this->getSubdetectorNames(); //something like this may be useful to check what we are supposed to do... // create input selector if (label_.size()>0){ sel_=new Selector( ModuleLabelSelector(label_)); } else { sel_=new Selector( MatchAllSelector()); } // declare the products to produce pixeldigi_collectionSig_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig"); pixeldigi_collectionPile_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile"); PixelDigiCollectionDM_ = ps.getParameter<std::string>("PixelDigiCollectionDM"); // clear local storage for this event SiHitStorage_.clear(); }
DataMixingSiPixelWorker::~DataMixingSiPixelWorker | ( | ) | [virtual] |
void DataMixingSiPixelWorker::addSiPixelPileups | ( | const int | bcr, |
edm::EventPrincipal * | ep, | ||
unsigned int | EventId | ||
) |
Definition at line 102 of file DataMixingSiPixelWorker.cc.
References edm::DetSetVector< T >::begin(), begin, end, edm::DetSetVector< T >::end(), edm::EventPrincipal::id(), collect_tpl::input, LogDebug, pixeldigi_collectionPile_, and SiHitStorage_.
Referenced by edm::DataMixingModule::addPileups().
{ LogDebug("DataMixingSiPixelWorker") <<"\n===============> adding pileups from event "<<ep->id()<<" for bunchcrossing "<<bcr; // fill in maps of hits; same code as addSignals, except now applied to the pileup events boost::shared_ptr<Wrapper<edm::DetSetVector<PixelDigi> > const> inputPTR = getProductByTag<edm::DetSetVector<PixelDigi> >(*ep, pixeldigi_collectionPile_ ); if(inputPTR ) { const edm::DetSetVector<PixelDigi> *input = const_cast< edm::DetSetVector<PixelDigi> * >(inputPTR->product()); // Handle< edm::DetSetVector<PixelDigi> > input; // if( e->getByLabel(pixeldigi_collectionPile_,input) ) { //loop on all detsets (detectorIDs) inside the input collection edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin(); for (; DSViter!=input->end();DSViter++){ #ifdef DEBUG LogDebug("DataMixingSiPixelWorker") << "Pileups: Processing DetID " << DSViter->id; #endif uint32_t detID = DSViter->id; edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin(); edm::DetSet<PixelDigi>::const_iterator end =(DSViter->data).end(); edm::DetSet<PixelDigi>::const_iterator icopy; // find correct local map (or new one) for this detector ID SiGlobalIndex::const_iterator itest; itest = SiHitStorage_.find(detID); if(itest!=SiHitStorage_.end()) { // this detID already has hits, add to existing map OneDetectorMap LocalMap = itest->second; // fill in local map with extra channels for (icopy=begin; icopy!=end; icopy++) { LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy )); } SiHitStorage_[detID]=LocalMap; } else{ // fill local storage with this information, put in global collection OneDetectorMap LocalMap; for (icopy=begin; icopy!=end; icopy++) { LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy )); } SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) ); } } } }
void DataMixingSiPixelWorker::addSiPixelSignals | ( | const edm::Event & | e | ) |
Definition at line 66 of file DataMixingSiPixelWorker.cc.
References begin, end, edm::Event::getByLabel(), edm::EventBase::id(), collect_tpl::input, LogDebug, pixeldigi_collectionSig_, and SiHitStorage_.
Referenced by edm::DataMixingModule::addSignals().
{ // fill in maps of hits LogDebug("DataMixingSiPixelWorker")<<"===============> adding MC signals for "<<e.id(); Handle< edm::DetSetVector<PixelDigi> > input; if( e.getByLabel(pixeldigi_collectionSig_,input) ) { //loop on all detsets (detectorIDs) inside the input collection edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin(); for (; DSViter!=input->end();DSViter++){ #ifdef DEBUG LogDebug("DataMixingSiPixelWorker") << "Processing DetID " << DSViter->id; #endif uint32_t detID = DSViter->id; edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin(); edm::DetSet<PixelDigi>::const_iterator end =(DSViter->data).end(); edm::DetSet<PixelDigi>::const_iterator icopy; OneDetectorMap LocalMap; for (icopy=begin; icopy!=end; icopy++) { LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy )); } SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) ); } } } // end of addSiPixelSignals
void DataMixingSiPixelWorker::putSiPixel | ( | edm::Event & | e | ) |
Definition at line 169 of file DataMixingSiPixelWorker.cc.
References ecalMGPA::adc(), PixelDigiCollectionDM_, edm::DetSet< T >::push_back(), edm::Event::put(), and SiHitStorage_.
Referenced by edm::DataMixingModule::put().
{ // collection of Digis to put in the event std::vector< edm::DetSet<PixelDigi> > vPixelDigi; // loop through our collection of detectors, merging hits and putting new ones in the output // big loop over Detector IDs: for(SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin(); IDet != SiHitStorage_.end(); IDet++) { edm::DetSet<PixelDigi> SPD(IDet->first); // Make empty collection with this detector ID OneDetectorMap LocalMap = IDet->second; //counter variables int formerPixel = -1; int currentPixel; int ADCSum = 0; OneDetectorMap::const_iterator iLocalchk; for(OneDetectorMap::const_iterator iLocal = LocalMap.begin(); iLocal != LocalMap.end(); ++iLocal) { currentPixel = iLocal->first; if (currentPixel == formerPixel) { // we have to add these digis together ADCSum+=(iLocal->second).adc(); } else{ if(formerPixel!=-1){ // ADC info stolen from SiStrips... if (ADCSum > 511) ADCSum = 255; else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254; PixelDigi aHit(formerPixel, ADCSum); SPD.push_back( aHit ); } // save pointers for next iteration formerPixel = currentPixel; ADCSum = (iLocal->second).adc(); } iLocalchk = iLocal; if((++iLocalchk) == LocalMap.end()) { //make sure not to lose the last one if (ADCSum > 511) ADCSum = 255; else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254; SPD.push_back( PixelDigi(formerPixel, ADCSum) ); } }// end of loop over one detector // stick this into the global vector of detector info vPixelDigi.push_back(SPD); } // end of big loop over all detector IDs // put the collection of digis in the event LogInfo("DataMixingSiPixelWorker") << "total # Merged Pixels: " << vPixelDigi.size() ; // make new digi collection std::auto_ptr< edm::DetSetVector<PixelDigi> > MyPixelDigis(new edm::DetSetVector<PixelDigi>(vPixelDigi) ); // put collection e.put( MyPixelDigis, PixelDigiCollectionDM_ ); // clear local storage for this event SiHitStorage_.clear(); }
std::string edm::DataMixingSiPixelWorker::label_ [private] |
Definition at line 71 of file DataMixingSiPixelWorker.h.
Referenced by DataMixingSiPixelWorker().
Definition at line 57 of file DataMixingSiPixelWorker.h.
Referenced by addSiPixelPileups(), and DataMixingSiPixelWorker().
Definition at line 56 of file DataMixingSiPixelWorker.h.
Referenced by addSiPixelSignals(), and DataMixingSiPixelWorker().
std::string edm::DataMixingSiPixelWorker::PixelDigiCollectionDM_ [private] |
Definition at line 58 of file DataMixingSiPixelWorker.h.
Referenced by DataMixingSiPixelWorker(), and putSiPixel().
Selector* edm::DataMixingSiPixelWorker::sel_ [private] |
Definition at line 70 of file DataMixingSiPixelWorker.h.
Referenced by DataMixingSiPixelWorker(), and ~DataMixingSiPixelWorker().
Definition at line 65 of file DataMixingSiPixelWorker.h.
Referenced by addSiPixelPileups(), addSiPixelSignals(), DataMixingSiPixelWorker(), and putSiPixel().