CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

SiStripBaselineAnalyzer Class Reference

#include <Validation/SiStripAnalyzer/src/SiStripBaselineAnalyzer.cc>

Inheritance diagram for SiStripBaselineAnalyzer:
edm::EDAnalyzer

List of all members.

Public Member Functions

 SiStripBaselineAnalyzer (const edm::ParameterSet &)
 ~SiStripBaselineAnalyzer ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void endJob ()

Private Attributes

uint16_t actualModule_
TCanvas * Canvas_
TH1F * h1BadAPVperEvent_
TH1F * h1Baseline_
TH1F * h1Clusters_
TH1F * h1ProcessedRawDigis_
uint16_t nModuletoDisplay_
TFile * oFile_
std::string outputFile_
edm::InputTag srcBaseline_
edm::InputTag srcProcessedRawDigi_
std::auto_ptr
< SiStripPedestalsSubtractor
subtractorPed_
std::vector< TH1F > vBaselineHisto_
std::vector< TH1F > vClusterHisto_
std::vector< TH1F > vProcessedRawDigiHisto_

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 68 of file SiStripBaselineAnalyzer.cc.


Constructor & Destructor Documentation

SiStripBaselineAnalyzer::SiStripBaselineAnalyzer ( const edm::ParameterSet conf) [explicit]

Definition at line 103 of file SiStripBaselineAnalyzer.cc.

References SiStripRawProcessingFactory::create_SubtractorPed(), edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), h1BadAPVperEvent_, nModuletoDisplay_, outputFile_, srcBaseline_, srcProcessedRawDigi_, subtractorPed_, vBaselineHisto_, vClusterHisto_, and vProcessedRawDigiHisto_.

                                                                           {
   
  outputFile_ = conf.getUntrackedParameter<std::string>("outputFile", "StripHistos.root");
  srcBaseline_ =  conf.getParameter<edm::InputTag>( "srcBaseline" );
  srcProcessedRawDigi_ =  conf.getParameter<edm::InputTag>( "srcProcessedRawDigi" );
  subtractorPed_ = SiStripRawProcessingFactory::create_SubtractorPed(conf.getParameter<edm::ParameterSet>("Algorithms"));
  nModuletoDisplay_ = conf.getParameter<uint32_t>( "nModuletoDisplay_" );
  
  vProcessedRawDigiHisto_.clear();
  vProcessedRawDigiHisto_.reserve(10000);
  
  vBaselineHisto_.clear();
  vBaselineHisto_.reserve(10000);
  
  vClusterHisto_.clear();
  vClusterHisto_.reserve(10000);
 
  
  h1BadAPVperEvent_ = new TH1F("BadAPV/Event","BadAPV/Event", 2001, -0.5, 2000.5);
  h1BadAPVperEvent_->SetXTitle("# Modules with Bad APVs");
  h1BadAPVperEvent_->SetYTitle("Entries");
  h1BadAPVperEvent_->SetLineWidth(2);
  h1BadAPVperEvent_->SetLineStyle(2);
 
}
SiStripBaselineAnalyzer::~SiStripBaselineAnalyzer ( )

Definition at line 130 of file SiStripBaselineAnalyzer.cc.

{
 
   

}

Member Function Documentation

void SiStripBaselineAnalyzer::analyze ( const edm::Event e,
const edm::EventSetup es 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 138 of file SiStripBaselineAnalyzer.cc.

References actualModule_, edmNew::DetSetVector< T >::begin(), edm::DetSet< T >::begin(), gather_cfg::cout, edmNew::DetSetVector< T >::end(), edm::EventID::event(), event(), edm::Event::getByLabel(), h1BadAPVperEvent_, h1Baseline_, h1Clusters_, h1ProcessedRawDigis_, edm::EventBase::id(), nModuletoDisplay_, edm::EventID::run(), CrabTask::run, getRunRegistry::runs, srcBaseline_, srcProcessedRawDigi_, strip(), subtractorPed_, vBaselineHisto_, vClusterHisto_, and vProcessedRawDigiHisto_.

{
   using namespace edm;
   
   
   subtractorPed_->init(es);
   
   edm::Handle<edm::DetSetVector<SiStripProcessedRawDigi> > moduleBaseline;
   e.getByLabel(srcBaseline_, moduleBaseline);
   
   edm::Handle< edm::DetSetVector<SiStripRawDigi> > moduleRawDigi;
   e.getByLabel(srcProcessedRawDigi_,moduleRawDigi);

   edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters;
   edm::InputTag clusLabel("siStripClusters");
   e.getByLabel(clusLabel, clusters);

   
   char detIds[20];
   char evs[20];
   char runs[20];    
    
   edm::DetSetVector<SiStripProcessedRawDigi>::const_iterator itDSBaseline = moduleBaseline->begin();
   edm::DetSetVector<SiStripRawDigi>::const_iterator itRawDigis = moduleRawDigi->begin();
   
   uint32_t NBabAPVs = moduleRawDigi->size();     
   std::cout<< "Number of module with HIP in this event: " << NBabAPVs << std::endl;
   h1BadAPVperEvent_->Fill(NBabAPVs);
   
   for (; itRawDigis != moduleRawDigi->end(); ++itRawDigis, ++itDSBaseline) {
      if(actualModule_ > nModuletoDisplay_) return;
      uint32_t detId = itRawDigis->id;
          
          
      if(itDSBaseline->id != detId){
                std::cout << "Collections out of Synch. Something of fishy is going on ;-)" << std::endl;
                return;
      }   
      
          actualModule_++;
          uint32_t event = e.id().event();
          uint32_t run = e.id().run();
          //std::cout << "processing module N: " << actualModule_<< " detId: " << detId << " event: "<< event << std::endl; 
          
          sprintf(detIds,"%ul", detId);
          sprintf(evs,"%ul", event);
          sprintf(runs,"%ul", run);
          char* dHistoName = Form("Id:%s_run:%s_ev:%s",detIds, runs, evs);
          h1ProcessedRawDigis_ = new TH1F(dHistoName,dHistoName, 768, -0.5, 767.5); 
          h1Baseline_ = new TH1F(dHistoName,dHistoName, 768, -0.5, 767.5); 
      h1Clusters_ = new TH1F(dHistoName,dHistoName, 768, -0.5, 767.5);
          
          
          h1ProcessedRawDigis_->SetXTitle("strip#");  
          h1ProcessedRawDigis_->SetYTitle("ADC");
          h1ProcessedRawDigis_->SetMaximum(1024.);
      h1ProcessedRawDigis_->SetMinimum(-300.);
          h1ProcessedRawDigis_->SetLineWidth(2);

   
      h1Baseline_->SetXTitle("strip#");
      h1Baseline_->SetYTitle("ADC");
      h1Baseline_->SetMaximum(1024.);
      h1Baseline_->SetMinimum(-300.);
      h1Baseline_->SetLineWidth(2);
          h1Baseline_->SetLineStyle(2);
          h1Baseline_->SetLineColor(2);
         
          h1Clusters_->SetXTitle("strip#");
      h1Clusters_->SetYTitle("ADC");
      h1Clusters_->SetMaximum(1024.);
      h1Clusters_->SetMinimum(-300.);
      h1Clusters_->SetLineWidth(2);
          h1Clusters_->SetLineStyle(2);
          h1Clusters_->SetLineColor(3);
          
          edm::DetSet<SiStripProcessedRawDigi>::const_iterator  itBaseline; 
          std::vector<int16_t>::const_iterator itProcessedRawDigis;
          
          
          std::vector<int16_t> ProcessedRawDigis(itRawDigis->size());
          subtractorPed_->subtract( *itRawDigis, ProcessedRawDigis);
          
          
          int strip =0;
      for(itProcessedRawDigis = ProcessedRawDigis.begin(), itBaseline = itDSBaseline->begin();itProcessedRawDigis != ProcessedRawDigis.end(); ++itProcessedRawDigis, ++itBaseline){
                h1ProcessedRawDigis_->Fill(strip, *itProcessedRawDigis);
                h1Baseline_->Fill(strip, itBaseline->adc()); 
                ++strip;
      }   
          
          edmNew::DetSetVector<SiStripCluster>::const_iterator itClusters = clusters->begin();
          for ( ; itClusters != clusters->end(); ++itClusters ){
                for ( edmNew::DetSet<SiStripCluster>::const_iterator clus =     itClusters->begin(); clus != itClusters->end(); ++clus){
            if(clus->geographicalId() == detId){
                                int firststrip = clus->firstStrip();
                    //std::cout << "Found cluster in detId " << detId << " " << firststrip << " " << clus->amplitudes().size() << " -----------------------------------------------" << std::endl;              
                        strip=0;
                                for( std::vector<uint8_t>::const_iterator itAmpl = clus->amplitudes().begin(); itAmpl != clus->amplitudes().end(); ++itAmpl){
                   h1Clusters_->Fill(firststrip+strip, *itAmpl);
                                   ++strip;
                                }
                        }              
        }
          }
          
          
         
         vProcessedRawDigiHisto_.push_back(*h1ProcessedRawDigis_);
         vBaselineHisto_.push_back(*h1Baseline_);
         vClusterHisto_.push_back(*h1Clusters_);
         
        }
                
    
}
void SiStripBaselineAnalyzer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 257 of file SiStripBaselineAnalyzer.cc.

References actualModule_.

{
  
  
  actualModule_ =0;
   
 
}
void SiStripBaselineAnalyzer::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 268 of file SiStripBaselineAnalyzer.cc.

References h1BadAPVperEvent_, oFile_, outputFile_, vBaselineHisto_, vClusterHisto_, and vProcessedRawDigiHisto_.

                                {
    oFile_ = new TFile((const char*)outputFile_.c_str(), "RECREATE");
    oFile_->mkdir("ProcessedRawDigis");
    oFile_->mkdir("Baseline");
        oFile_->mkdir("Clusters");
        
        
    std::vector<TH1F>::iterator itvProcessedRawDigis, itvBaseline, itvClusters; 
    itvProcessedRawDigis = vProcessedRawDigiHisto_.begin();
    itvBaseline = vBaselineHisto_.begin();
    itvClusters = vClusterHisto_.begin();
        
        oFile_->cd();
        h1BadAPVperEvent_->Write();
        for(; itvProcessedRawDigis != vProcessedRawDigiHisto_.end(); ++itvProcessedRawDigis, ++itvBaseline, ++itvClusters){
            oFile_->cd("ProcessedRawDigis");
                itvProcessedRawDigis->Write();
                oFile_->cd("Baseline");
            itvBaseline->Write();
                oFile_->cd("Clusters");
                itvClusters->Write();
        }
        oFile_->Write();
    oFile_->Close();
  
}

Member Data Documentation

Definition at line 99 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and beginJob().

Definition at line 93 of file SiStripBaselineAnalyzer.cc.

Definition at line 85 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().

Definition at line 88 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

Definition at line 89 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

Definition at line 87 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze().

Definition at line 98 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

Definition at line 91 of file SiStripBaselineAnalyzer.cc.

Referenced by endJob().

std::string SiStripBaselineAnalyzer::outputFile_ [private]

Definition at line 81 of file SiStripBaselineAnalyzer.cc.

Referenced by endJob(), and SiStripBaselineAnalyzer().

Definition at line 82 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

Definition at line 83 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

Definition at line 79 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), and SiStripBaselineAnalyzer().

std::vector<TH1F> SiStripBaselineAnalyzer::vBaselineHisto_ [private]

Definition at line 95 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().

std::vector<TH1F> SiStripBaselineAnalyzer::vClusterHisto_ [private]

Definition at line 96 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().

Definition at line 94 of file SiStripBaselineAnalyzer.cc.

Referenced by analyze(), endJob(), and SiStripBaselineAnalyzer().