Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:00:04

0001 // -*- C++ -*-
0002 //
0003 // Package:    EcalPulseShapeGrapher
0004 // Class:      EcalPulseShapeGrapher
0005 //
0006 /**\class EcalPulseShapeGrapher EcalPulseShapeGrapher.cc Analyzers/EcalPulseShapeGrapher/src/EcalPulseShapeGrapher.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Seth Cooper
0015 //         Created:  Tue Feb  5 11:35:45 CST 2008
0016 //
0017 //
0018 
0019 #include "CaloOnlineTools/EcalTools/plugins/EcalPulseShapeGrapher.h"
0020 
0021 //
0022 // constants, enums and typedefs
0023 //
0024 
0025 //
0026 // static data member definitions
0027 //
0028 
0029 //
0030 // constructors and destructor
0031 //
0032 EcalPulseShapeGrapher::EcalPulseShapeGrapher(const edm::ParameterSet& iConfig)
0033     : EBUncalibratedRecHitCollection_(consumes<EcalUncalibratedRecHitCollection>(
0034           iConfig.getParameter<edm::InputTag>("EBUncalibratedRecHitCollection"))),
0035       EBDigis_(consumes<EBDigiCollection>(iConfig.getParameter<edm::InputTag>("EBDigiCollection"))),
0036       EEUncalibratedRecHitCollection_(consumes<EcalUncalibratedRecHitCollection>(
0037           iConfig.getParameter<edm::InputTag>("EEUncalibratedRecHitCollection"))),
0038       EEDigis_(consumes<EEDigiCollection>(iConfig.getParameter<edm::InputTag>("EEDigiCollection"))),
0039       ampCut_(iConfig.getUntrackedParameter<int>("AmplitudeCutADC", 13)),
0040       rootFilename_(iConfig.getUntrackedParameter<std::string>("rootFilename", "pulseShapeGrapher")) {
0041   //now do what ever initialization is needed
0042 
0043   std::vector<int> listDefaults;
0044   listDefaults.push_back(-1);
0045   listChannels_ = iConfig.getUntrackedParameter<std::vector<int> >("listChannels", listDefaults);
0046 
0047   for (std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr) {
0048     std::string title = "Amplitude of cry " + intToString(*itr);
0049     std::string name = "ampOfCry" + intToString(*itr);
0050     ampHistMap_[*itr] = new TH1F(name.c_str(), title.c_str(), 100, 0, 100);
0051     ampHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
0052 
0053     title = "Amplitude (over 13 ADC) of cry " + intToString(*itr);
0054     name = "cutAmpOfCry" + intToString(*itr);
0055     cutAmpHistMap_[*itr] = new TH1F(name.c_str(), title.c_str(), 100, 0, 100);
0056     cutAmpHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
0057 
0058     title = "Pulse shape of cry " + intToString(*itr);
0059     name = "PulseShapeCry" + intToString(*itr);
0060     pulseShapeHistMap_[*itr] = new TH2F(name.c_str(), title.c_str(), 10, 0, 10, 220, -20, 2);
0061     pulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
0062 
0063     title = "Raw Pulse shape of cry " + intToString(*itr);
0064     name = "RawPulseShapeCry" + intToString(*itr);
0065     rawPulseShapeHistMap_[*itr] = new TH2F(name.c_str(), title.c_str(), 10, 0, 10, 500, 0, 500);
0066     rawPulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
0067 
0068     title = "Amplitude of first sample, cry " + intToString(*itr);
0069     name = "AmpOfFirstSampleCry" + intToString(*itr);
0070     firstSampleHistMap_[*itr] = new TH1F(name.c_str(), title.c_str(), 300, 100, 400);
0071     firstSampleHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
0072   }
0073 
0074   fedMap_ = new EcalFedMap();
0075 
0076   for (int i = 0; i < 10; i++)
0077     abscissa[i] = i;
0078 }
0079 
0080 //
0081 // member functions
0082 //
0083 
0084 // ------------ method called to for each event  ------------
0085 void EcalPulseShapeGrapher::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0086   using namespace edm;
0087   using namespace std;
0088 
0089   const Handle<EcalUncalibratedRecHitCollection>& EBHits = iEvent.getHandle(EBUncalibratedRecHitCollection_);
0090   const Handle<EcalUncalibratedRecHitCollection>& EEHits = iEvent.getHandle(EEUncalibratedRecHitCollection_);
0091   const Handle<EBDigiCollection>& EBdigis = iEvent.getHandle(EBDigis_);
0092   const Handle<EEDigiCollection>& EEdigis = iEvent.getHandle(EEDigis_);
0093 
0094   unique_ptr<EcalElectronicsMapping> ecalElectronicsMap(new EcalElectronicsMapping);
0095 
0096   // Loop over the hits
0097   for (EcalUncalibratedRecHitCollection::const_iterator hitItr = EBHits->begin(); hitItr != EBHits->end(); ++hitItr) {
0098     EcalUncalibratedRecHit hit = (*hitItr);
0099     float amplitude = hit.amplitude();
0100     EBDetId hitDetId = hit.id();
0101 
0102     // Get the Fedid
0103     EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
0104     int FEDid = 600 + elecId.dccId();
0105     string SMname = fedMap_->getSliceFromFed(FEDid);
0106 
0107     vector<int>::const_iterator itr = listChannels_.begin();
0108     while (itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex()) {
0109       itr++;
0110     }
0111     if (itr == listChannels_.end())
0112       continue;
0113 
0114     ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
0115     if (amplitude < ampCut_)
0116       continue;
0117 
0118     cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
0119     EBDigiCollection::const_iterator digiItr = EBdigis->begin();
0120     while (digiItr != EBdigis->end() && digiItr->id() != hitItr->id()) {
0121       digiItr++;
0122     }
0123     if (digiItr == EBdigis->end())
0124       continue;
0125 
0126     double sampleADC[10];
0127     EBDataFrame df(*digiItr);
0128     double pedestal = 200;
0129 
0130     if (df.sample(0).gainId() != 1 || df.sample(1).gainId() != 1)
0131       continue;  //goes to the next digi
0132     else {
0133       sampleADC[0] = df.sample(0).adc();
0134       sampleADC[1] = df.sample(1).adc();
0135       pedestal = (double)(sampleADC[0] + sampleADC[1]) / (double)2;
0136     }
0137 
0138     for (int i = 0; (unsigned int)i < digiItr->size(); ++i) {
0139       EBDataFrame df(*digiItr);
0140       double gain = 12.;
0141       if (df.sample(i).gainId() == 1)
0142         gain = 1.;
0143       else if (df.sample(i).gainId() == 2)
0144         gain = 2.;
0145       sampleADC[i] = pedestal + (df.sample(i).adc() - pedestal) * gain;
0146     }
0147 
0148     for (int i = 0; i < 10; ++i) {
0149       pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i] - pedestal) / amplitude);
0150       rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
0151     }
0152     firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
0153   }
0154 
0155   // Now do the same for the EE hits
0156   for (EcalUncalibratedRecHitCollection::const_iterator hitItr = EEHits->begin(); hitItr != EEHits->end(); ++hitItr) {
0157     EcalUncalibratedRecHit hit = (*hitItr);
0158     float amplitude = hit.amplitude();
0159     EEDetId hitDetId = hit.id();
0160 
0161     // Get the Fedid
0162     EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
0163     int FEDid = 600 + elecId.dccId();
0164     string SMname = fedMap_->getSliceFromFed(FEDid);
0165 
0166     vector<int>::const_iterator itr = listChannels_.begin();
0167     while (itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex()) {
0168       itr++;
0169     }
0170     if (itr == listChannels_.end())
0171       continue;
0172 
0173     ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
0174     if (amplitude < ampCut_)
0175       continue;
0176 
0177     cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
0178     EEDigiCollection::const_iterator digiItr = EEdigis->begin();
0179     while (digiItr != EEdigis->end() && digiItr->id() != hitItr->id()) {
0180       digiItr++;
0181     }
0182     if (digiItr == EEdigis->end())
0183       continue;
0184 
0185     double sampleADC[10];
0186     EEDataFrame df(*digiItr);
0187     double pedestal = 200;
0188 
0189     if (df.sample(0).gainId() != 1 || df.sample(1).gainId() != 1)
0190       continue;  //goes to the next digi
0191     else {
0192       sampleADC[0] = df.sample(0).adc();
0193       sampleADC[1] = df.sample(1).adc();
0194       pedestal = (double)(sampleADC[0] + sampleADC[1]) / (double)2;
0195     }
0196 
0197     for (int i = 0; (unsigned int)i < digiItr->size(); ++i) {
0198       EEDataFrame df(*digiItr);
0199       double gain = 12.;
0200       if (df.sample(i).gainId() == 1)
0201         gain = 1.;
0202       else if (df.sample(i).gainId() == 2)
0203         gain = 2.;
0204       sampleADC[i] = pedestal + (df.sample(i).adc() - pedestal) * gain;
0205     }
0206 
0207     for (int i = 0; i < 10; ++i) {
0208       pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i] - pedestal) / amplitude);
0209       rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
0210     }
0211     firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
0212   }
0213 }
0214 
0215 // ------------ method called once each job just after ending the event loop  ------------
0216 void EcalPulseShapeGrapher::endJob() {
0217   rootFilename_ += ".root";
0218   file_ = new TFile(rootFilename_.c_str(), "RECREATE");
0219   TH1::AddDirectory(false);
0220 
0221   for (std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr) {
0222     ampHistMap_[*itr]->Write();
0223     cutAmpHistMap_[*itr]->Write();
0224     firstSampleHistMap_[*itr]->Write();
0225 
0226     rawPulseShapeHistMap_[*itr]->Write();
0227     TProfile* t2 = (TProfile*)(rawPulseShapeHistMap_[*itr]->ProfileX());
0228     t2->Write();
0229     //TODO: fix the normalization so these are correct
0230   }
0231 
0232   file_->Write();
0233   file_->Close();
0234 }
0235 
0236 std::string EcalPulseShapeGrapher::intToString(int num) {
0237   using namespace std;
0238   ostringstream myStream;
0239   myStream << num << flush;
0240   return (myStream.str());  //returns the string form of the stringstream object
0241 }