CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/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() { } 
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     // declare the products to produce
00039 
00040     pixeldigi_collectionSig_   = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig");
00041     pixeldigi_collectionPile_   = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile");
00042     PixelDigiCollectionDM_  = ps.getParameter<std::string>("PixelDigiCollectionDM");
00043 
00044     // clear local storage for this event                                                                     
00045     SiHitStorage_.clear();
00046 
00047   }
00048                
00049 
00050   // Virtual destructor needed.
00051   DataMixingSiPixelWorker::~DataMixingSiPixelWorker() { 
00052   }  
00053 
00054 
00055 
00056   void DataMixingSiPixelWorker::addSiPixelSignals(const edm::Event &e) { 
00057     // fill in maps of hits
00058 
00059     LogDebug("DataMixingSiPixelWorker")<<"===============> adding MC signals for "<<e.id();
00060 
00061     Handle< edm::DetSetVector<PixelDigi> >  input;
00062 
00063     if( e.getByLabel(pixeldigi_collectionSig_,input) ) {
00064 
00065       //loop on all detsets (detectorIDs) inside the input collection
00066       edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
00067       for (; DSViter!=input->end();DSViter++){
00068 
00069 #ifdef DEBUG
00070         LogDebug("DataMixingSiPixelWorker")  << "Processing DetID " << DSViter->id;
00071 #endif
00072 
00073         uint32_t detID = DSViter->id;
00074         edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin();
00075         edm::DetSet<PixelDigi>::const_iterator end   =(DSViter->data).end();
00076         edm::DetSet<PixelDigi>::const_iterator icopy;
00077   
00078         OneDetectorMap LocalMap;
00079 
00080         for (icopy=begin; icopy!=end; icopy++) {
00081           LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00082         }
00083 
00084         SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
00085       }
00086  
00087     }    
00088   } // end of addSiPixelSignals
00089 
00090 
00091 
00092   void DataMixingSiPixelWorker::addSiPixelPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr) {
00093   
00094     LogDebug("DataMixingSiPixelWorker") <<"\n===============> adding pileups from event  "<<ep->id()<<" for bunchcrossing "<<bcr;
00095 
00096     // fill in maps of hits; same code as addSignals, except now applied to the pileup events
00097 
00098     boost::shared_ptr<Wrapper<edm::DetSetVector<PixelDigi> >  const> inputPTR =
00099       getProductByTag<edm::DetSetVector<PixelDigi> >(*ep, pixeldigi_collectionPile_ );
00100 
00101     if(inputPTR ) {
00102 
00103       const edm::DetSetVector<PixelDigi>  *input = const_cast< edm::DetSetVector<PixelDigi> * >(inputPTR->product());
00104 
00105 
00106 
00107       //   Handle< edm::DetSetVector<PixelDigi> >  input;
00108 
00109       //   if( e->getByLabel(pixeldigi_collectionPile_,input) ) {
00110 
00111       //loop on all detsets (detectorIDs) inside the input collection
00112       edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
00113       for (; DSViter!=input->end();DSViter++){
00114 
00115 #ifdef DEBUG
00116         LogDebug("DataMixingSiPixelWorker")  << "Pileups: Processing DetID " << DSViter->id;
00117 #endif
00118 
00119         uint32_t detID = DSViter->id;
00120         edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin();
00121         edm::DetSet<PixelDigi>::const_iterator end   =(DSViter->data).end();
00122         edm::DetSet<PixelDigi>::const_iterator icopy;
00123 
00124         // find correct local map (or new one) for this detector ID
00125 
00126         SiGlobalIndex::const_iterator itest;
00127 
00128         itest = SiHitStorage_.find(detID);
00129 
00130         if(itest!=SiHitStorage_.end()) {  // this detID already has hits, add to existing map
00131 
00132           OneDetectorMap LocalMap = itest->second;
00133 
00134           // fill in local map with extra channels
00135           for (icopy=begin; icopy!=end; icopy++) {
00136             LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00137           }
00138 
00139           SiHitStorage_[detID]=LocalMap;
00140           
00141         }
00142         else{ // fill local storage with this information, put in global collection
00143 
00144           OneDetectorMap LocalMap;
00145 
00146           for (icopy=begin; icopy!=end; icopy++) {
00147             LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00148           }
00149 
00150           SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
00151         }
00152 
00153       }
00154     }
00155   }
00156 
00157 
00158  
00159   void DataMixingSiPixelWorker::putSiPixel(edm::Event &e) {
00160 
00161     // collection of Digis to put in the event
00162 
00163     std::vector< edm::DetSet<PixelDigi> > vPixelDigi;
00164 
00165     // loop through our collection of detectors, merging hits and putting new ones in the output
00166 
00167     // big loop over Detector IDs:
00168 
00169     for(SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin();
00170         IDet != SiHitStorage_.end(); IDet++) {
00171 
00172       edm::DetSet<PixelDigi> SPD(IDet->first); // Make empty collection with this detector ID
00173         
00174       OneDetectorMap LocalMap = IDet->second;
00175 
00176       //counter variables
00177       int formerPixel = -1;
00178       int currentPixel;
00179       int ADCSum = 0;
00180 
00181 
00182       OneDetectorMap::const_iterator iLocalchk;
00183 
00184       for(OneDetectorMap::const_iterator iLocal  = LocalMap.begin();
00185           iLocal != LocalMap.end(); ++iLocal) {
00186 
00187         currentPixel = iLocal->first; 
00188 
00189         if (currentPixel == formerPixel) { // we have to add these digis together
00190           ADCSum+=(iLocal->second).adc();
00191         }
00192         else{
00193           if(formerPixel!=-1){             // ADC info stolen from SiStrips...
00194             if (ADCSum > 511) ADCSum = 255;
00195             else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
00196             PixelDigi aHit(formerPixel, ADCSum);
00197             SPD.push_back( aHit );        
00198           }
00199           // save pointers for next iteration
00200           formerPixel = currentPixel;
00201           ADCSum = (iLocal->second).adc();
00202         }
00203 
00204         iLocalchk = iLocal;
00205         if((++iLocalchk) == LocalMap.end()) {  //make sure not to lose the last one
00206           if (ADCSum > 511) ADCSum = 255;
00207           else if (ADCSum > 253 && ADCSum < 512) ADCSum = 254;
00208           SPD.push_back( PixelDigi(formerPixel, ADCSum) );        
00209         } 
00210 
00211       }// end of loop over one detector
00212 
00213       // stick this into the global vector of detector info
00214       vPixelDigi.push_back(SPD);
00215 
00216     } // end of big loop over all detector IDs
00217 
00218     // put the collection of digis in the event   
00219     LogInfo("DataMixingSiPixelWorker") << "total # Merged Pixels: " << vPixelDigi.size() ;
00220 
00221     // make new digi collection
00222     
00223     std::auto_ptr< edm::DetSetVector<PixelDigi> > MyPixelDigis(new edm::DetSetVector<PixelDigi>(vPixelDigi) );
00224 
00225     // put collection
00226 
00227     e.put( MyPixelDigis, PixelDigiCollectionDM_ );
00228 
00229     // clear local storage for this event
00230     SiHitStorage_.clear();
00231   }
00232 
00233 } //edm