CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

EcalPulseShapeGrapher Class Reference

#include <EcalPulseShapeGrapher.h>

Inheritance diagram for EcalPulseShapeGrapher:
edm::EDAnalyzer

List of all members.

Public Member Functions

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

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void endJob ()
std::string intToString (int)

Private Attributes

int abscissa [10]
int ampCut_
std::map< int, TH1F * > ampHistMap_
std::map< int, TH1F * > cutAmpHistMap_
edm::InputTag EBDigis_
edm::InputTag EBUncalibratedRecHitCollection_
edm::InputTag EEDigis_
edm::InputTag EEUncalibratedRecHitCollection_
EcalFedMapfedMap_
TFile * file_
std::map< int, TH1F * > firstSampleHistMap_
std::vector< int > listChannels_
int ordinate [10]
std::map< int, TH2F * > pulseShapeHistMap_
std::map< int, TH2F * > rawPulseShapeHistMap_
std::string rootFilename_

Detailed Description

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

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

Definition at line 47 of file EcalPulseShapeGrapher.h.


Constructor & Destructor Documentation

EcalPulseShapeGrapher::EcalPulseShapeGrapher ( const edm::ParameterSet iConfig) [explicit]

Definition at line 33 of file EcalPulseShapeGrapher.cc.

References abscissa, ampHistMap_, cutAmpHistMap_, fedMap_, firstSampleHistMap_, edm::ParameterSet::getUntrackedParameter(), i, intToString(), listChannels_, mergeVDriftHistosByStation::name, pulseShapeHistMap_, rawPulseShapeHistMap_, and indexGen::title.

                                                                           :
EBUncalibratedRecHitCollection_ (iConfig.getParameter<edm::InputTag>("EBUncalibratedRecHitCollection")),
EBDigis_ (iConfig.getParameter<edm::InputTag>("EBDigiCollection")),
EEUncalibratedRecHitCollection_ (iConfig.getParameter<edm::InputTag>("EEUncalibratedRecHitCollection")),
EEDigis_ (iConfig.getParameter<edm::InputTag>("EEDigiCollection")),
ampCut_ (iConfig.getUntrackedParameter<int>("AmplitudeCutADC", 13)),
rootFilename_ (iConfig.getUntrackedParameter<std::string>("rootFilename","pulseShapeGrapher"))
{
   //now do what ever initialization is needed

   std::vector<int> listDefaults;
   listDefaults.push_back(-1);
   listChannels_ = iConfig.getUntrackedParameter<std::vector<int> >("listChannels", listDefaults);

   for(std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr)
   {
     std::string title = "Amplitude of cry "+intToString(*itr);
     std::string name = "ampOfCry"+intToString(*itr);
     ampHistMap_[*itr] = new TH1F(name.c_str(),title.c_str(),100,0,100);
     ampHistMap_[*itr]->GetXaxis()->SetTitle("ADC");

     title = "Amplitude (over 13 ADC) of cry "+intToString(*itr);
     name = "cutAmpOfCry"+intToString(*itr);
     cutAmpHistMap_[*itr] = new TH1F(name.c_str(),title.c_str(),100,0,100);
     cutAmpHistMap_[*itr]->GetXaxis()->SetTitle("ADC");

     title = "Pulse shape of cry "+intToString(*itr);
     name = "PulseShapeCry"+intToString(*itr);
     pulseShapeHistMap_[*itr] = new TH2F(name.c_str(),title.c_str(),10,0,10,220,-20,2);
     pulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");

     title = "Raw Pulse shape of cry "+intToString(*itr);
     name = "RawPulseShapeCry"+intToString(*itr);
     rawPulseShapeHistMap_[*itr] = new TH2F(name.c_str(),title.c_str(),10,0,10,500,0,500);
     rawPulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
     
     title = "Amplitude of first sample, cry "+intToString(*itr);
     name = "AmpOfFirstSampleCry"+intToString(*itr);
     firstSampleHistMap_[*itr] = new TH1F(name.c_str(),title.c_str(),300,100,400);
     firstSampleHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
   }
   
   fedMap_ = new EcalFedMap();
   
   for (int i=0; i<10; i++)
     abscissa[i] = i;
}
EcalPulseShapeGrapher::~EcalPulseShapeGrapher ( )

Definition at line 82 of file EcalPulseShapeGrapher.cc.

{
 
   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

void EcalPulseShapeGrapher::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 97 of file EcalPulseShapeGrapher.cc.

References EcalMGPASample::adc(), ampCut_, ampHistMap_, EcalUncalibratedRecHit::amplitude(), cutAmpHistMap_, EcalElectronicsId::dccId(), EBDigis_, EBUncalibratedRecHitCollection_, EEDigis_, EEUncalibratedRecHitCollection_, L1Comparator_cfi::FEDid, fedMap_, firstSampleHistMap_, EcalMGPASample::gainId(), edm::Event::getByLabel(), EcalFedMap::getSliceFromFed(), EEDetId::hashedIndex(), EBDetId::hashedIndex(), i, EcalUncalibratedRecHit::id(), listChannels_, pulseShapeHistMap_, rawPulseShapeHistMap_, EcalDataFrame::sample(), and findQualityFiles::size.

{
   using namespace edm;
   using namespace std;

   int numHitsWithActivity = 0;
   
   //int eventNum = iEvent.id().event();
   //vector<EcalUncalibratedRecHit> sampleHitsPastCut;
   
   Handle<EcalUncalibratedRecHitCollection> EBHits;
   iEvent.getByLabel(EBUncalibratedRecHitCollection_, EBHits);
   Handle<EcalUncalibratedRecHitCollection> EEHits;
   iEvent.getByLabel(EEUncalibratedRecHitCollection_, EEHits);
   //cout << "event: " << eventNum << " sample hits collection size: " << sampleHits->size() << endl;
   //Handle<EcalUncalibratedRecHitCollection> fittedHits;
   //iEvent.getByLabel(FittedUncalibratedRecHitCollection_, fittedHits);
   Handle<EBDigiCollection>  EBdigis;
   iEvent.getByLabel(EBDigis_, EBdigis);
   Handle<EEDigiCollection>  EEdigis;
   iEvent.getByLabel(EEDigis_, EEdigis);
   //cout << "event: " << eventNum << " digi collection size: " << digis->size() << endl;

   auto_ptr<EcalElectronicsMapping> ecalElectronicsMap(new EcalElectronicsMapping);

   // Loop over the hits
   for(EcalUncalibratedRecHitCollection::const_iterator hitItr = EBHits->begin(); hitItr != EBHits->end(); ++hitItr)
   {
     EcalUncalibratedRecHit hit = (*hitItr);
     float amplitude = hit.amplitude();
     EBDetId hitDetId = hit.id();

     // Get the Fedid
     EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
     int FEDid = 600+elecId.dccId();
     string SMname = fedMap_->getSliceFromFed(FEDid);
     
     vector<int>::const_iterator itr = listChannels_.begin();
     while(itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex())
     {
       itr++;
     }
     if(itr==listChannels_.end())
       continue;

     ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
     //cout << "Cry hash:" << hitDetId.hashedIndex() << " amplitude: " << amplitude << endl;
     if(amplitude < ampCut_)
       continue;

     cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
     numHitsWithActivity++;
     EBDigiCollection::const_iterator digiItr= EBdigis->begin();
     while(digiItr != EBdigis->end() && digiItr->id() != hitItr->id())
     {
       digiItr++;
     }
     if(digiItr == EBdigis->end())
       continue;

     double sampleADC[10];
     EBDataFrame df(*digiItr);
     double pedestal = 200;

     if(df.sample(0).gainId()!=1 || df.sample(1).gainId()!=1) continue; //goes to the next digi
     else {
       sampleADC[0] = df.sample(0).adc();
       sampleADC[1] = df.sample(1).adc();
       pedestal = (double)(sampleADC[0]+sampleADC[1])/(double)2;
     } 

     for (int i=0; (unsigned int)i< digiItr->size(); ++i ) {
       EBDataFrame df(*digiItr);
       double gain = 12.;
       if(df.sample(i).gainId()==1)
         gain = 1.;
       else if(df.sample(i).gainId()==2)
         gain = 2.;
       sampleADC[i] = pedestal+(df.sample(i).adc()-pedestal)*gain;
     }

     //cout << "1) maxsample amp:" << maxSampleAmp << " maxSampleIndex:" << maxSampleIndex << endl;
     for(int i=0; i<10;++i)
     {
       //cout << "Filling hist for:" << hitDetId.hashedIndex() << " with sample:" << i 
         //<< " amp:" <<(float)(df.sample(i).adc()-baseline)/(maxSampleAmp-baseline) << endl;
       //cout << "ADC of sample:" << df.sample(i).adc() << " baseline:" << baseline << " maxSampleAmp:" 
         //<< maxSampleAmp << endl << endl;
       pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]-pedestal)/amplitude);
       rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
     }
     firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
   }

   // Now do the same for the EE hits
   for(EcalUncalibratedRecHitCollection::const_iterator hitItr = EEHits->begin(); hitItr != EEHits->end(); ++hitItr)
   {
     EcalUncalibratedRecHit hit = (*hitItr);
     float amplitude = hit.amplitude();
     EEDetId hitDetId = hit.id();

     // Get the Fedid
     EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
     int FEDid = 600+elecId.dccId();
     string SMname = fedMap_->getSliceFromFed(FEDid);
     
     vector<int>::const_iterator itr = listChannels_.begin();
     while(itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex())
     {
       itr++;
     }
     if(itr==listChannels_.end())
       continue;

     ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
     //cout << "Cry hash:" << hitDetId.hashedIndex() << " amplitude: " << amplitude << endl;
     if(amplitude < ampCut_)
       continue;

     cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
     numHitsWithActivity++;
     EEDigiCollection::const_iterator digiItr= EEdigis->begin();
     while(digiItr != EEdigis->end() && digiItr->id() != hitItr->id())
     {
       digiItr++;
     }
     if(digiItr == EEdigis->end())
       continue;

     double sampleADC[10];
     EEDataFrame df(*digiItr);
     double pedestal = 200;

     if(df.sample(0).gainId()!=1 || df.sample(1).gainId()!=1) continue; //goes to the next digi
     else {
       sampleADC[0] = df.sample(0).adc();
       sampleADC[1] = df.sample(1).adc();
       pedestal = (double)(sampleADC[0]+sampleADC[1])/(double)2;
     } 

     for (int i=0; (unsigned int)i< digiItr->size(); ++i ) {
       EEDataFrame df(*digiItr);
       double gain = 12.;
       if(df.sample(i).gainId()==1)
         gain = 1.;
       else if(df.sample(i).gainId()==2)
         gain = 2.;
       sampleADC[i] = pedestal+(df.sample(i).adc()-pedestal)*gain;
     }

     //cout << "1) maxsample amp:" << maxSampleAmp << " maxSampleIndex:" << maxSampleIndex << endl;
     for(int i=0; i<10;++i)
     {
       //cout << "Filling hist for:" << hitDetId.hashedIndex() << " with sample:" << i 
         //<< " amp:" <<(float)(df.sample(i).adc()-baseline)/(maxSampleAmp-baseline) << endl;
       //cout << "ADC of sample:" << df.sample(i).adc() << " baseline:" << baseline << " maxSampleAmp:" 
         //<< maxSampleAmp << endl << endl;
       pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]-pedestal)/amplitude);
       rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
     }
     firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
   }
}
void EcalPulseShapeGrapher::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 264 of file EcalPulseShapeGrapher.cc.

References ampHistMap_, cutAmpHistMap_, file_, firstSampleHistMap_, listChannels_, rawPulseShapeHistMap_, and rootFilename_.

                              {
  
   rootFilename_+=".root";
   file_ = new TFile(rootFilename_.c_str(),"RECREATE");
   TH1::AddDirectory(false);
  
   for(std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr)
   {
     ampHistMap_[*itr]->Write();
     cutAmpHistMap_[*itr]->Write();
     firstSampleHistMap_[*itr]->Write();

     rawPulseShapeHistMap_[*itr]->Write();
     TProfile* t2 = (TProfile*) (rawPulseShapeHistMap_[*itr]->ProfileX());
     t2->Write();
     //TODO: fix the normalization so these are correct
     //pulseShapeHistMap_[*itr]->Write();
     //TProfile* t1 = (TProfile*) (pulseShapeHistMap_[*itr]->ProfileX());
     //t1->Write();
   }

   
  file_->Write();
  file_->Close();
}
std::string EcalPulseShapeGrapher::intToString ( int  num) [private]

Definition at line 290 of file EcalPulseShapeGrapher.cc.

Referenced by EcalPulseShapeGrapher().

{
  using namespace std;
  ostringstream myStream;
  myStream << num << flush;
  return(myStream.str()); //returns the string form of the stringstream object
}

Member Data Documentation

Definition at line 64 of file EcalPulseShapeGrapher.h.

Referenced by EcalPulseShapeGrapher().

Definition at line 73 of file EcalPulseShapeGrapher.h.

Referenced by analyze().

std::map<int,TH1F*> EcalPulseShapeGrapher::ampHistMap_ [private]

Definition at line 67 of file EcalPulseShapeGrapher.h.

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

std::map<int,TH1F*> EcalPulseShapeGrapher::cutAmpHistMap_ [private]

Definition at line 71 of file EcalPulseShapeGrapher.h.

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

Definition at line 60 of file EcalPulseShapeGrapher.h.

Referenced by analyze().

Definition at line 59 of file EcalPulseShapeGrapher.h.

Referenced by analyze().

Definition at line 62 of file EcalPulseShapeGrapher.h.

Referenced by analyze().

Definition at line 61 of file EcalPulseShapeGrapher.h.

Referenced by analyze().

Definition at line 78 of file EcalPulseShapeGrapher.h.

Referenced by analyze(), and EcalPulseShapeGrapher().

TFile* EcalPulseShapeGrapher::file_ [private]

Definition at line 76 of file EcalPulseShapeGrapher.h.

Referenced by endJob().

std::map<int,TH1F*> EcalPulseShapeGrapher::firstSampleHistMap_ [private]

Definition at line 69 of file EcalPulseShapeGrapher.h.

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

std::vector<int> EcalPulseShapeGrapher::listChannels_ [private]

Definition at line 66 of file EcalPulseShapeGrapher.h.

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

Definition at line 65 of file EcalPulseShapeGrapher.h.

std::map<int,TH2F*> EcalPulseShapeGrapher::pulseShapeHistMap_ [private]

Definition at line 68 of file EcalPulseShapeGrapher.h.

Referenced by analyze(), and EcalPulseShapeGrapher().

std::map<int,TH2F*> EcalPulseShapeGrapher::rawPulseShapeHistMap_ [private]

Definition at line 70 of file EcalPulseShapeGrapher.h.

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

std::string EcalPulseShapeGrapher::rootFilename_ [private]

Definition at line 74 of file EcalPulseShapeGrapher.h.

Referenced by endJob().