CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoLocalTracker/SiStripZeroSuppression/plugins/SiStripMeanCMExtractor.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripMeanCMExtractor
00004 // Class:      SiStripMeanCMExtractor
00005 // 
00013 //
00014 // Original Author:  Ivan Amos Cali,32 4-A08,+41227673039,
00015 //         Created:  Wed Oct 13 11:50:47 CEST 2010
00016 // $Id: SiStripMeanCMExtractor.cc,v 1.2 2010/11/04 14:12:31 edwenger Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <sstream>
00023 #include <memory>
00024 #include <list>
00025 #include <algorithm>
00026 #include <cassert>
00027 
00028 
00029 // user include files
00030 #include "FWCore/Framework/interface/Event.h"
00031 #include "FWCore/Framework/interface/EventSetup.h"
00032 #include "DataFormats/Common/interface/Handle.h"
00033 #include "FWCore/Framework/interface/EventSetup.h"
00034 #include "FWCore/Framework/interface/ESHandle.h"
00035 #include "FWCore/Framework/interface/Frameworkfwd.h"
00036 #include "FWCore/Framework/interface/EDProducer.h"
00037 
00038 #include "FWCore/Framework/interface/Event.h"
00039 #include "FWCore/Framework/interface/MakerMacros.h"
00040 
00041 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00042 
00043 #include "DataFormats/Common/interface/DetSet.h"
00044 #include "DataFormats/Common/interface/DetSetVector.h"
00045 #include "DataFormats/DetId/interface/DetId.h"
00046 #include "DataFormats/SiStripDigi/interface/SiStripProcessedRawDigi.h"
00047 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
00048 #include "CondFormats/SiStripObjects/interface/SiStripPedestals.h"
00049 #include "CondFormats/DataRecord/interface/SiStripPedestalsRcd.h"
00050 
00051 
00052 typedef std::map<uint32_t, std::vector<float> > CMMap;
00053 
00054 class SiStripMeanCMExtractor : public edm::EDProducer {
00055    public:
00056       explicit SiStripMeanCMExtractor( const edm::ParameterSet&);
00057       ~SiStripMeanCMExtractor();
00058 
00059    private:
00060       virtual void beginJob() ;
00061       virtual void produce(edm::Event&, const edm::EventSetup&);
00062       virtual void endJob() ;
00063       
00064           void init(const edm::EventSetup&);
00065           
00066           edm::ESHandle<SiStripPedestals> pedestalHandle_;
00067       uint32_t pedestal_cache_id_;
00068   
00069           void StoreMean(const edm::DetSetVector<SiStripProcessedRawDigi>& );
00070           void ConvertMeanMapToDetSetVector(std::vector<edm::DetSet<SiStripProcessedRawDigi> >&);
00071           void CMExtractorFromPedestals(const edm::DetSetVector<SiStripRawDigi>&,std::vector<edm::DetSet<SiStripProcessedRawDigi> >&);
00072           edm::InputTag _inputTag;
00073           std::string _Algorithm;
00074           uint16_t _nEventsToUse;
00075           uint16_t _actualEvent;
00076       
00077           CMMap _CMMap; //it contains the sum of the CM calculated before. The normalization for the number of events it is done at the end when it is written in the DetSetVector.
00078 };
00079 
00080 
00081 SiStripMeanCMExtractor::SiStripMeanCMExtractor(const edm::ParameterSet& conf) : 
00082      _inputTag(conf.getParameter<edm::InputTag> ("CMCollection")),
00083         _Algorithm(conf.getParameter<std::string>("Algorithm")),
00084         _nEventsToUse(conf.getParameter<uint32_t>("NEvents"))
00085 {
00086 
00087          
00088     if(_nEventsToUse < 1) _nEventsToUse=1;
00089         produces< edm::DetSetVector<SiStripProcessedRawDigi> > ("MEANAPVCM"); 
00090 }
00091 
00092 
00093 SiStripMeanCMExtractor::~SiStripMeanCMExtractor()
00094 {
00095  
00096 }
00097 
00098 void SiStripMeanCMExtractor::init(const edm::EventSetup& es){
00099      
00100         uint32_t p_cache_id = es.get<SiStripPedestalsRcd>().cacheIdentifier();
00101         
00102         if(p_cache_id != pedestal_cache_id_) {
00103                 es.get<SiStripPedestalsRcd>().get( pedestalHandle_ );
00104                 pedestal_cache_id_ = p_cache_id;
00105         }
00106 }
00107 
00108 // ------------ method called to produce the data  ------------
00109 void
00110 SiStripMeanCMExtractor::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00111 {
00112    using namespace edm;
00113    
00114      
00115    //if(_actualEvent > _nEventsToUse) return;
00116    
00117    std::vector<edm::DetSet<SiStripProcessedRawDigi> > meancm;
00118    
00119    if(_Algorithm == "StoredCM"){
00120         edm::Handle< edm::DetSetVector<SiStripProcessedRawDigi> > inputCM;
00121         iEvent.getByLabel(_inputTag,inputCM);
00122         
00123         this->StoreMean(*inputCM);
00124         this->ConvertMeanMapToDetSetVector(meancm);
00125         
00126    } else if (_Algorithm == "Pedestals"){
00127      this->init(iSetup);
00128          
00129      edm::Handle< edm::DetSetVector<SiStripRawDigi> > input;
00130      iEvent.getByLabel(_inputTag,input);
00131          
00132      this->CMExtractorFromPedestals(*input,meancm);
00133    }
00134    
00135    ++_actualEvent;
00136    
00137     
00138         
00139         
00140         std::auto_ptr< edm::DetSetVector<SiStripProcessedRawDigi> > outputMeanCM(new edm::DetSetVector<SiStripProcessedRawDigi>(meancm) );
00141     iEvent.put( outputMeanCM,"MEANAPVCM");
00142    
00143 }
00144 
00145 void SiStripMeanCMExtractor::CMExtractorFromPedestals(const edm::DetSetVector<SiStripRawDigi>& input, std::vector<edm::DetSet<SiStripProcessedRawDigi> >& meancm){
00146         meancm.clear();
00147         meancm.reserve(15000);    
00148         
00149          for ( edm::DetSetVector<SiStripRawDigi>::const_iterator 
00150           rawDigis = input.begin(); rawDigis != input.end(); rawDigis++) {
00151          SiStripPedestals::Range detPedestalRange = pedestalHandle_->getRange(rawDigis->id);
00152                  edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(rawDigis->id);
00153                 
00154                 for(uint16_t APV = 0; APV < rawDigis->size()/128; ++APV){
00155                         uint16_t MinPed =0;
00156                         for(uint16_t strip = APV*128; strip< (APV+1)*128; ++strip){
00157                           uint16_t ped =  (uint16_t)pedestalHandle_->getPed(strip,detPedestalRange);
00158                           if(ped < MinPed) MinPed = ped;
00159                         }
00160                         if(MinPed>128) MinPed=128;
00161                         MeanCMDetSet.push_back(MinPed);
00162                 }
00163                 
00164                 meancm.push_back(MeanCMDetSet); 
00165         }
00166 }
00167 
00168 void SiStripMeanCMExtractor::StoreMean(const edm::DetSetVector<SiStripProcessedRawDigi>& Input){
00169         
00170         uint32_t detId;
00171         CMMap::iterator itMap;
00172         edm::DetSetVector<SiStripProcessedRawDigi>::const_iterator itInput;
00173         
00174         for(itInput = Input.begin(); itInput != Input.end(); ++itInput){
00175                 detId = itInput->id;
00176                 itMap = _CMMap.find(detId);
00177                 edm::DetSet<SiStripProcessedRawDigi>::const_iterator itCM;
00178                 std::vector<float> MeanCMNValue;
00179                 MeanCMNValue.clear();
00180                 if(itMap!=_CMMap.end()){   //the detId was already found
00181                     std::vector< float >& MapContent = itMap->second;
00182                         std::vector<float>::iterator itMapVector = MapContent.begin();
00183                         for(itCM = itInput->begin(); itCM != itInput->end(); ++itCM, ++itMapVector){
00184                                 MeanCMNValue.push_back(itCM->adc() + *itMapVector); 
00185             }
00186                         _CMMap.erase(itMap);
00187             _CMMap.insert(itMap, std::pair<uint32_t, std::vector<float> >(detId,MeanCMNValue));                 
00188                 } else {                 //no detId found
00189                         for(itCM = itInput->begin(); itCM != itInput->end(); ++itCM) MeanCMNValue.push_back(itCM->adc());                       
00190                         _CMMap.insert(std::pair<uint32_t, std::vector<float> >(detId,MeanCMNValue));
00191                 }
00192         }
00193   
00194 }
00195 
00196 void 
00197 SiStripMeanCMExtractor::ConvertMeanMapToDetSetVector(std::vector<edm::DetSet<SiStripProcessedRawDigi> >& meancm){
00198         CMMap::iterator itMap;
00199         std::vector<float>::const_iterator itMapVector;
00200         
00201         meancm.clear();
00202         meancm.reserve(15000);    
00203         
00204     for(itMap = _CMMap.begin(); itMap != _CMMap.end(); ++itMap){
00205        edm::DetSet<SiStripProcessedRawDigi> MeanCMDetSet(itMap->first);
00206            for(itMapVector = (itMap->second).begin(); itMapVector != (itMap->second).end(); ++itMapVector) MeanCMDetSet.push_back(*itMapVector/(float)_actualEvent);
00207            meancm.push_back(MeanCMDetSet);      
00208         }       
00209     
00210 }
00211 void 
00212 SiStripMeanCMExtractor::beginJob()
00213 {
00214         _actualEvent =1;
00215                 
00216         _CMMap.clear();
00217         
00218 }
00219 
00220 void 
00221 SiStripMeanCMExtractor::endJob() {
00222     
00223            
00224 }
00225 
00226 DEFINE_FWK_MODULE(SiStripMeanCMExtractor);
00227