CMS 3D CMS Logo

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_collection_   = ps.getParameter<edm::InputTag>("pixeldigiCollection");
00049     PixelDigiCollectionDM_  = ps.getParameter<std::string>("PixelDigiCollectionDM");
00050 
00051     // clear local storage for this event                                                                     
00052     SiHitStorage_.clear();
00053 
00054   }
00055                
00056 
00057   // Virtual destructor needed.
00058   DataMixingSiPixelWorker::~DataMixingSiPixelWorker() { 
00059     delete sel_;
00060     sel_=0;
00061   }  
00062 
00063 
00064 
00065   void DataMixingSiPixelWorker::addSiPixelSignals(const edm::Event &e) { 
00066     // fill in maps of hits
00067 
00068     LogDebug("DataMixingSiPixelWorker")<<"===============> adding MC signals for "<<e.id();
00069 
00070     Handle< edm::DetSetVector<PixelDigi> >  input;
00071 
00072     if( e.getByLabel(pixeldigi_collection_,input) ) {
00073 
00074       //loop on all detsets (detectorIDs) inside the input collection
00075       edm::DetSetVector<PixelDigi>::const_iterator DSViter=input->begin();
00076       for (; DSViter!=input->end();DSViter++){
00077 
00078 #ifdef DEBUG
00079         LogDebug("DataMixingSiPixelWorker")  << "Processing DetID " << DSViter->id;
00080 #endif
00081 
00082         uint32_t detID = DSViter->id;
00083         edm::DetSet<PixelDigi>::const_iterator begin =(DSViter->data).begin();
00084         edm::DetSet<PixelDigi>::const_iterator end   =(DSViter->data).end();
00085         edm::DetSet<PixelDigi>::const_iterator icopy;
00086   
00087         OneDetectorMap LocalMap;
00088 
00089         for (icopy=begin; icopy!=end; icopy++) {
00090           LocalMap.insert(OneDetectorMap::value_type( (icopy->channel()), *icopy ));
00091         }
00092 
00093         SiHitStorage_.insert( SiGlobalIndex::value_type( detID, LocalMap ) );
00094       }
00095  
00096     }    
00097   } // end of addSiPixelSignals
00098 
00099 
00100 
00101   void DataMixingSiPixelWorker::addSiPixelPileups(const int bcr, Event *e, unsigned int eventNr) {
00102   
00103     LogDebug("DataMixingSiPixelWorker") <<"\n===============> adding pileups from event  "<<e->id()<<" for bunchcrossing "<<bcr;
00104 
00105     // fill in maps of hits; same code as addSignals, except now applied to the pileup events
00106 
00107     Handle< edm::DetSetVector<PixelDigi> >  input;
00108 
00109     if( e->getByLabel(pixeldigi_collection_,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 = 0;
00178       int currentPixel;
00179       int ADCSum = 0;
00180       PixelDigi OldHit;
00181       int nmatch=0;
00182 
00183       OneDetectorMap::const_iterator iLocalchk;
00184 
00185       for(OneDetectorMap::const_iterator iLocal  = LocalMap.begin();
00186           iLocal != LocalMap.end(); ++iLocal) {
00187 
00188         currentPixel = iLocal->first; 
00189 
00190         if (currentPixel == formerPixel) { // we have to add these digis together
00191           nmatch++;                  // use this to avoid using the "count" function
00192           ADCSum+=(iLocal->second).adc();          // on every element...
00193 
00194           iLocalchk = iLocal;
00195           if((iLocalchk++) == LocalMap.end()) {  //make sure not to lose the last one
00196             PixelDigi aHit(formerPixel, ADCSum);
00197             SPD.push_back( aHit );        
00198             // reset adc sum, nmatch
00199             ADCSum = 0 ;
00200             nmatch=0;     
00201           }
00202         }
00203         else {
00204           if(nmatch>0) {
00205             PixelDigi aHit(formerPixel, ADCSum);
00206             SPD.push_back( aHit );        
00207             // reset adc sum, nmatch
00208             ADCSum = 0 ;
00209             nmatch=0;     
00210           }
00211           else {
00212             SPD.push_back( OldHit );
00213           }
00214         
00215           iLocalchk = iLocal;
00216           if((iLocalchk++) == LocalMap.end()) {  //make sure not to lose the last one
00217             SPD.push_back( iLocal->second );
00218           }
00219 
00220           // save pointers for next iteration
00221           OldHit = iLocal->second;
00222           formerPixel = currentPixel;
00223           ADCSum = (iLocal->second).adc();
00224         }
00225       }  // end of loop over one detector
00226 
00227       // stick this into the global vector of detector info
00228       vPixelDigi.push_back(SPD);
00229 
00230     } // end of big loop over all detector IDs
00231 
00232     // put the collection of digis in the event   
00233     LogInfo("DataMixingSiPixelWorker") << "total # Merged Pixels: " << vPixelDigi.size() ;
00234 
00235     // make new digi collection
00236     
00237     std::auto_ptr< edm::DetSetVector<PixelDigi> > MyPixelDigis(new edm::DetSetVector<PixelDigi>(vPixelDigi) );
00238 
00239     // put collection
00240 
00241     e.put( MyPixelDigis, PixelDigiCollectionDM_ );
00242 
00243     // clear local storage for this event
00244     SiHitStorage_.clear();
00245   }
00246 
00247 } //edm

Generated on Tue Jun 9 17:47:22 2009 for CMSSW by  doxygen 1.5.4