CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/SimGeneral/DataMixingModule/plugins/DataMixingSiPixelWorker.cc

Go to the documentation of this file.
00001 // File: DataMixingSiPixelWorker.cc
00002 // Description:  see DataMixingSiPixelWorker.h
00003 // Author:  Mike Hildreth, University of Notre Dame
00004 //
00005 //--------------------------------------------
00006 
00007 #include <map>
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "FWCore/Utilities/interface/EDMException.h"
00010 #include "FWCore/Framework/interface/ConstProductRegistry.h"
00011 #include "FWCore/ServiceRegistry/interface/Service.h"
00012 #include "DataFormats/Common/interface/Handle.h"
00013 #include "DataFormats/Provenance/interface/Provenance.h"
00014 #include "DataFormats/Provenance/interface/BranchDescription.h"
00015 //
00016 //
00017 #include "DataMixingSiPixelWorker.h"
00018 
00019 
00020 using namespace std;
00021 
00022 namespace edm
00023 {
00024 
00025   // Virtual constructor
00026 
00027   DataMixingSiPixelWorker::DataMixingSiPixelWorker() { sel_=0;} 
00028 
00029   // Constructor 
00030   DataMixingSiPixelWorker::DataMixingSiPixelWorker(const edm::ParameterSet& ps) : 
00031                                                             label_(ps.getParameter<std::string>("Label"))
00032 
00033   {                                                         
00034 
00035     // get the subdetector names
00036     //    this->getSubdetectorNames();  //something like this may be useful to check what we are supposed to do...
00037 
00038     // create input selector
00039     if (label_.size()>0){
00040       sel_=new Selector( ModuleLabelSelector(label_));
00041     }
00042     else {
00043       sel_=new Selector( MatchAllSelector());
00044     }
00045 
00046     // declare the products to produce
00047 
00048     pixeldigi_collectionSig_   = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig");
00049     pixeldigi_collectionPile_   = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile");
00050     PixelDigiCollectionDM_  = ps.getParameter<std::string>("PixelDigiCollectionDM");
00051 
00052     // clear local storage for this event                                                                     
00053     SiHitStorage_.clear();
00054 
00055   }
00056                
00057 
00058   // Virtual destructor needed.
00059   DataMixingSiPixelWorker::~DataMixingSiPixelWorker() { 
00060     delete sel_;
00061     sel_=0;
00062   }  
00063 
00064 
00065 
00066   void DataMixingSiPixelWorker::addSiPixelSignals(const edm::Event &e) { 
00067     // fill in maps of hits
00068 
00069     LogDebug("DataMixingSiPixelWorker")<<"===============> adding MC signals for "<<e.id();
00070 
00071     Handle< edm::DetSetVector<PixelDigi> >  input;
00072 
00073     if( e.getByLabel(pixeldigi_collectionSig_,input) ) {
00074 
00075       //loop on all detsets (detectorIDs) inside the input collection
00076       edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
00077       for (; DSViter!=input->end();DSViter++){
00078 
00079 #ifdef DEBUG
00080         LogDebug("DataMixingSiPixelWorker")  << "Processing DetID " << DSViter->id;
00081 #endif
00082 
00083         uint32_t detID = DSViter->id;
00084         edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin();
00085         edm::DetSet<PixelDigi>::const_iterator end   =(DSViter->data).end();
00086         edm::DetSet<PixelDigi>::const_iterator icopy;
00087   
00088         OneDetectorMap LocalMap;
00089 
00090         for (icopy=begin; icopy!=end; icopy++) {
00091           LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00092         }
00093 
00094         SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
00095       }
00096  
00097     }    
00098   } // end of addSiPixelSignals
00099 
00100 
00101 
00102   void DataMixingSiPixelWorker::addSiPixelPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr) {
00103   
00104     LogDebug("DataMixingSiPixelWorker") <<"\n===============> adding pileups from event  "<<ep->id()<<" for bunchcrossing "<<bcr;
00105 
00106     // fill in maps of hits; same code as addSignals, except now applied to the pileup events
00107 
00108     boost::shared_ptr<Wrapper<edm::DetSetVector<PixelDigi> >  const> inputPTR =
00109       getProductByTag<edm::DetSetVector<PixelDigi> >(*ep, pixeldigi_collectionPile_ );
00110 
00111     if(inputPTR ) {
00112 
00113       const edm::DetSetVector<PixelDigi>  *input = const_cast< edm::DetSetVector<PixelDigi> * >(inputPTR->product());
00114 
00115 
00116 
00117       //   Handle< edm::DetSetVector<PixelDigi> >  input;
00118 
00119       //   if( e->getByLabel(pixeldigi_collectionPile_,input) ) {
00120 
00121       //loop on all detsets (detectorIDs) inside the input collection
00122       edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
00123       for (; DSViter!=input->end();DSViter++){
00124 
00125 #ifdef DEBUG
00126         LogDebug("DataMixingSiPixelWorker")  << "Pileups: Processing DetID " << DSViter->id;
00127 #endif
00128 
00129         uint32_t detID = DSViter->id;
00130         edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin();
00131         edm::DetSet<PixelDigi>::const_iterator end   =(DSViter->data).end();
00132         edm::DetSet<PixelDigi>::const_iterator icopy;
00133 
00134         // find correct local map (or new one) for this detector ID
00135 
00136         SiGlobalIndex::const_iterator itest;
00137 
00138         itest = SiHitStorage_.find(detID);
00139 
00140         if(itest!=SiHitStorage_.end()) {  // this detID already has hits, add to existing map
00141 
00142           OneDetectorMap LocalMap = itest->second;
00143 
00144           // fill in local map with extra channels
00145           for (icopy=begin; icopy!=end; icopy++) {
00146             LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00147           }
00148 
00149           SiHitStorage_[detID]=LocalMap;
00150           
00151         }
00152         else{ // fill local storage with this information, put in global collection
00153 
00154           OneDetectorMap LocalMap;
00155 
00156           for (icopy=begin; icopy!=end; icopy++) {
00157             LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00158           }
00159 
00160           SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
00161         }
00162 
00163       }
00164     }
00165   }
00166 
00167 
00168  
00169   void DataMixingSiPixelWorker::putSiPixel(edm::Event &e) {
00170 
00171     // collection of Digis to put in the event
00172 
00173     std::vector< edm::DetSet<PixelDigi> > vPixelDigi;
00174 
00175     // loop through our collection of detectors, merging hits and putting new ones in the output
00176 
00177     // big loop over Detector IDs:
00178 
00179     for(SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin();
00180         IDet != SiHitStorage_.end(); IDet++) {
00181 
00182       edm::DetSet<PixelDigi> SPD(IDet->first); // Make empty collection with this detector ID
00183         
00184       OneDetectorMap LocalMap = IDet->second;
00185 
00186       //counter variables
00187       int formerPixel = -1;
00188       int currentPixel;
00189       int ADCSum = 0;
00190 
00191 
00192       OneDetectorMap::const_iterator iLocalchk;
00193 
00194       for(OneDetectorMap::const_iterator iLocal  = LocalMap.begin();
00195           iLocal != LocalMap.end(); ++iLocal) {
00196 
00197         currentPixel = iLocal->first; 
00198 
00199         if (currentPixel == formerPixel) { // we have to add these digis together
00200           ADCSum+=(iLocal->second).adc();
00201         }
00202         else{
00203           if(formerPixel!=-1){             // ADC info stolen from SiStrips...
00204             if (ADCSum > 511) ADCSum = 255;
00205             else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
00206             PixelDigi aHit(formerPixel, ADCSum);
00207             SPD.push_back( aHit );        
00208           }
00209           // save pointers for next iteration
00210           formerPixel = currentPixel;
00211           ADCSum = (iLocal->second).adc();
00212         }
00213 
00214         iLocalchk = iLocal;
00215         if((++iLocalchk) == LocalMap.end()) {  //make sure not to lose the last one
00216           if (ADCSum > 511) ADCSum = 255;
00217           else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
00218           SPD.push_back( PixelDigi(formerPixel, ADCSum) );        
00219         } 
00220 
00221       }// end of loop over one detector
00222 
00223       // stick this into the global vector of detector info
00224       vPixelDigi.push_back(SPD);
00225 
00226     } // end of big loop over all detector IDs
00227 
00228     // put the collection of digis in the event   
00229     LogInfo("DataMixingSiPixelWorker") << "total # Merged Pixels: " << vPixelDigi.size() ;
00230 
00231     // make new digi collection
00232     
00233     std::auto_ptr< edm::DetSetVector<PixelDigi> > MyPixelDigis(new edm::DetSetVector<PixelDigi>(vPixelDigi) );
00234 
00235     // put collection
00236 
00237     e.put( MyPixelDigis, PixelDigiCollectionDM_ );
00238 
00239     // clear local storage for this event
00240     SiHitStorage_.clear();
00241   }
00242 
00243 } //edm