CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/CaloOnlineTools/EcalTools/plugins/EcalPulseShapeGrapher.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    EcalPulseShapeGrapher
00004 // Class:      EcalPulseShapeGrapher
00005 // 
00013 //
00014 // Original Author:  Seth Cooper
00015 //         Created:  Tue Feb  5 11:35:45 CST 2008
00016 // $Id: EcalPulseShapeGrapher.cc,v 1.2 2010/01/04 15:07:40 ferriff Exp $
00017 //
00018 //
00019 
00020 #include "CaloOnlineTools/EcalTools/plugins/EcalPulseShapeGrapher.h"
00021 
00022 //
00023 // constants, enums and typedefs
00024 //
00025 
00026 //
00027 // static data member definitions
00028 //
00029 
00030 //
00031 // constructors and destructor
00032 //
00033 EcalPulseShapeGrapher::EcalPulseShapeGrapher(const edm::ParameterSet& iConfig) :
00034 EBUncalibratedRecHitCollection_ (iConfig.getParameter<edm::InputTag>("EBUncalibratedRecHitCollection")),
00035 EBDigis_ (iConfig.getParameter<edm::InputTag>("EBDigiCollection")),
00036 EEUncalibratedRecHitCollection_ (iConfig.getParameter<edm::InputTag>("EEUncalibratedRecHitCollection")),
00037 EEDigis_ (iConfig.getParameter<edm::InputTag>("EEDigiCollection")),
00038 ampCut_ (iConfig.getUntrackedParameter<int>("AmplitudeCutADC", 13)),
00039 rootFilename_ (iConfig.getUntrackedParameter<std::string>("rootFilename","pulseShapeGrapher"))
00040 {
00041    //now do what ever initialization is needed
00042 
00043    std::vector<int> listDefaults;
00044    listDefaults.push_back(-1);
00045    listChannels_ = iConfig.getUntrackedParameter<std::vector<int> >("listChannels", listDefaults);
00046 
00047    for(std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr)
00048    {
00049      std::string title = "Amplitude of cry "+intToString(*itr);
00050      std::string name = "ampOfCry"+intToString(*itr);
00051      ampHistMap_[*itr] = new TH1F(name.c_str(),title.c_str(),100,0,100);
00052      ampHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
00053 
00054      title = "Amplitude (over 13 ADC) of cry "+intToString(*itr);
00055      name = "cutAmpOfCry"+intToString(*itr);
00056      cutAmpHistMap_[*itr] = new TH1F(name.c_str(),title.c_str(),100,0,100);
00057      cutAmpHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
00058 
00059      title = "Pulse shape of cry "+intToString(*itr);
00060      name = "PulseShapeCry"+intToString(*itr);
00061      pulseShapeHistMap_[*itr] = new TH2F(name.c_str(),title.c_str(),10,0,10,220,-20,2);
00062      pulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
00063 
00064      title = "Raw Pulse shape of cry "+intToString(*itr);
00065      name = "RawPulseShapeCry"+intToString(*itr);
00066      rawPulseShapeHistMap_[*itr] = new TH2F(name.c_str(),title.c_str(),10,0,10,500,0,500);
00067      rawPulseShapeHistMap_[*itr]->GetXaxis()->SetTitle("sample");
00068      
00069      title = "Amplitude of first sample, cry "+intToString(*itr);
00070      name = "AmpOfFirstSampleCry"+intToString(*itr);
00071      firstSampleHistMap_[*itr] = new TH1F(name.c_str(),title.c_str(),300,100,400);
00072      firstSampleHistMap_[*itr]->GetXaxis()->SetTitle("ADC");
00073    }
00074    
00075    fedMap_ = new EcalFedMap();
00076    
00077    for (int i=0; i<10; i++)
00078      abscissa[i] = i;
00079 }
00080 
00081 
00082 EcalPulseShapeGrapher::~EcalPulseShapeGrapher()
00083 {
00084  
00085    // do anything here that needs to be done at desctruction time
00086    // (e.g. close files, deallocate resources etc.)
00087 
00088 }
00089 
00090 
00091 //
00092 // member functions
00093 //
00094 
00095 // ------------ method called to for each event  ------------
00096 void
00097 EcalPulseShapeGrapher::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00098 {
00099    using namespace edm;
00100    using namespace std;
00101 
00102    int numHitsWithActivity = 0;
00103    
00104    //int eventNum = iEvent.id().event();
00105    //vector<EcalUncalibratedRecHit> sampleHitsPastCut;
00106    
00107    Handle<EcalUncalibratedRecHitCollection> EBHits;
00108    iEvent.getByLabel(EBUncalibratedRecHitCollection_, EBHits);
00109    Handle<EcalUncalibratedRecHitCollection> EEHits;
00110    iEvent.getByLabel(EEUncalibratedRecHitCollection_, EEHits);
00111    //cout << "event: " << eventNum << " sample hits collection size: " << sampleHits->size() << endl;
00112    //Handle<EcalUncalibratedRecHitCollection> fittedHits;
00113    //iEvent.getByLabel(FittedUncalibratedRecHitCollection_, fittedHits);
00114    Handle<EBDigiCollection>  EBdigis;
00115    iEvent.getByLabel(EBDigis_, EBdigis);
00116    Handle<EEDigiCollection>  EEdigis;
00117    iEvent.getByLabel(EEDigis_, EEdigis);
00118    //cout << "event: " << eventNum << " digi collection size: " << digis->size() << endl;
00119 
00120    auto_ptr<EcalElectronicsMapping> ecalElectronicsMap(new EcalElectronicsMapping);
00121 
00122    // Loop over the hits
00123    for(EcalUncalibratedRecHitCollection::const_iterator hitItr = EBHits->begin(); hitItr != EBHits->end(); ++hitItr)
00124    {
00125      EcalUncalibratedRecHit hit = (*hitItr);
00126      float amplitude = hit.amplitude();
00127      EBDetId hitDetId = hit.id();
00128 
00129      // Get the Fedid
00130      EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
00131      int FEDid = 600+elecId.dccId();
00132      string SMname = fedMap_->getSliceFromFed(FEDid);
00133      
00134      vector<int>::const_iterator itr = listChannels_.begin();
00135      while(itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex())
00136      {
00137        itr++;
00138      }
00139      if(itr==listChannels_.end())
00140        continue;
00141 
00142      ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
00143      //cout << "Cry hash:" << hitDetId.hashedIndex() << " amplitude: " << amplitude << endl;
00144      if(amplitude < ampCut_)
00145        continue;
00146 
00147      cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
00148      numHitsWithActivity++;
00149      EBDigiCollection::const_iterator digiItr= EBdigis->begin();
00150      while(digiItr != EBdigis->end() && digiItr->id() != hitItr->id())
00151      {
00152        digiItr++;
00153      }
00154      if(digiItr == EBdigis->end())
00155        continue;
00156 
00157      double sampleADC[10];
00158      EBDataFrame df(*digiItr);
00159      double pedestal = 200;
00160 
00161      if(df.sample(0).gainId()!=1 || df.sample(1).gainId()!=1) continue; //goes to the next digi
00162      else {
00163        sampleADC[0] = df.sample(0).adc();
00164        sampleADC[1] = df.sample(1).adc();
00165        pedestal = (double)(sampleADC[0]+sampleADC[1])/(double)2;
00166      } 
00167 
00168      for (int i=0; (unsigned int)i< digiItr->size(); ++i ) {
00169        EBDataFrame df(*digiItr);
00170        double gain = 12.;
00171        if(df.sample(i).gainId()==1)
00172          gain = 1.;
00173        else if(df.sample(i).gainId()==2)
00174          gain = 2.;
00175        sampleADC[i] = pedestal+(df.sample(i).adc()-pedestal)*gain;
00176      }
00177 
00178      //cout << "1) maxsample amp:" << maxSampleAmp << " maxSampleIndex:" << maxSampleIndex << endl;
00179      for(int i=0; i<10;++i)
00180      {
00181        //cout << "Filling hist for:" << hitDetId.hashedIndex() << " with sample:" << i 
00182          //<< " amp:" <<(float)(df.sample(i).adc()-baseline)/(maxSampleAmp-baseline) << endl;
00183        //cout << "ADC of sample:" << df.sample(i).adc() << " baseline:" << baseline << " maxSampleAmp:" 
00184          //<< maxSampleAmp << endl << endl;
00185        pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]-pedestal)/amplitude);
00186        rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
00187      }
00188      firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
00189    }
00190 
00191    // Now do the same for the EE hits
00192    for(EcalUncalibratedRecHitCollection::const_iterator hitItr = EEHits->begin(); hitItr != EEHits->end(); ++hitItr)
00193    {
00194      EcalUncalibratedRecHit hit = (*hitItr);
00195      float amplitude = hit.amplitude();
00196      EEDetId hitDetId = hit.id();
00197 
00198      // Get the Fedid
00199      EcalElectronicsId elecId = ecalElectronicsMap->getElectronicsId(hitDetId);
00200      int FEDid = 600+elecId.dccId();
00201      string SMname = fedMap_->getSliceFromFed(FEDid);
00202      
00203      vector<int>::const_iterator itr = listChannels_.begin();
00204      while(itr != listChannels_.end() && (*itr) != hitDetId.hashedIndex())
00205      {
00206        itr++;
00207      }
00208      if(itr==listChannels_.end())
00209        continue;
00210 
00211      ampHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
00212      //cout << "Cry hash:" << hitDetId.hashedIndex() << " amplitude: " << amplitude << endl;
00213      if(amplitude < ampCut_)
00214        continue;
00215 
00216      cutAmpHistMap_[hitDetId.hashedIndex()]->Fill(amplitude);
00217      numHitsWithActivity++;
00218      EEDigiCollection::const_iterator digiItr= EEdigis->begin();
00219      while(digiItr != EEdigis->end() && digiItr->id() != hitItr->id())
00220      {
00221        digiItr++;
00222      }
00223      if(digiItr == EEdigis->end())
00224        continue;
00225 
00226      double sampleADC[10];
00227      EEDataFrame df(*digiItr);
00228      double pedestal = 200;
00229 
00230      if(df.sample(0).gainId()!=1 || df.sample(1).gainId()!=1) continue; //goes to the next digi
00231      else {
00232        sampleADC[0] = df.sample(0).adc();
00233        sampleADC[1] = df.sample(1).adc();
00234        pedestal = (double)(sampleADC[0]+sampleADC[1])/(double)2;
00235      } 
00236 
00237      for (int i=0; (unsigned int)i< digiItr->size(); ++i ) {
00238        EEDataFrame df(*digiItr);
00239        double gain = 12.;
00240        if(df.sample(i).gainId()==1)
00241          gain = 1.;
00242        else if(df.sample(i).gainId()==2)
00243          gain = 2.;
00244        sampleADC[i] = pedestal+(df.sample(i).adc()-pedestal)*gain;
00245      }
00246 
00247      //cout << "1) maxsample amp:" << maxSampleAmp << " maxSampleIndex:" << maxSampleIndex << endl;
00248      for(int i=0; i<10;++i)
00249      {
00250        //cout << "Filling hist for:" << hitDetId.hashedIndex() << " with sample:" << i 
00251          //<< " amp:" <<(float)(df.sample(i).adc()-baseline)/(maxSampleAmp-baseline) << endl;
00252        //cout << "ADC of sample:" << df.sample(i).adc() << " baseline:" << baseline << " maxSampleAmp:" 
00253          //<< maxSampleAmp << endl << endl;
00254        pulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]-pedestal)/amplitude);
00255        rawPulseShapeHistMap_[hitDetId.hashedIndex()]->Fill(i, (float)(sampleADC[i]));
00256      }
00257      firstSampleHistMap_[hitDetId.hashedIndex()]->Fill(sampleADC[0]);
00258    }
00259 }
00260 
00261 
00262 // ------------ method called once each job just after ending the event loop  ------------
00263 void 
00264 EcalPulseShapeGrapher::endJob() {
00265   
00266    rootFilename_+=".root";
00267    file_ = new TFile(rootFilename_.c_str(),"RECREATE");
00268    TH1::AddDirectory(false);
00269   
00270    for(std::vector<int>::const_iterator itr = listChannels_.begin(); itr != listChannels_.end(); ++itr)
00271    {
00272      ampHistMap_[*itr]->Write();
00273      cutAmpHistMap_[*itr]->Write();
00274      firstSampleHistMap_[*itr]->Write();
00275 
00276      rawPulseShapeHistMap_[*itr]->Write();
00277      TProfile* t2 = (TProfile*) (rawPulseShapeHistMap_[*itr]->ProfileX());
00278      t2->Write();
00279      //TODO: fix the normalization so these are correct
00280      //pulseShapeHistMap_[*itr]->Write();
00281      //TProfile* t1 = (TProfile*) (pulseShapeHistMap_[*itr]->ProfileX());
00282      //t1->Write();
00283    }
00284 
00285    
00286   file_->Write();
00287   file_->Close();
00288 }
00289 
00290 std::string EcalPulseShapeGrapher::intToString(int num)
00291 {
00292   using namespace std;
00293   ostringstream myStream;
00294   myStream << num << flush;
00295   return(myStream.str()); //returns the string form of the stringstream object
00296 }