CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTBCalo/EcalSimpleTBAnalysis/src/EcalSimpleTBAnalyzer.cc

Go to the documentation of this file.
00001 
00008 //
00009 // $Id: EcalSimpleTBAnalyzer.cc,v 1.10 2012/02/01 19:41:58 vskarupe Exp $
00010 //
00011 //
00012 
00013 #include "RecoTBCalo/EcalSimpleTBAnalysis/interface/EcalSimpleTBAnalyzer.h"
00014 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00015 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00016 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00017 #include "TBDataFormats/EcalTBObjects/interface/EcalTBHodoscopeRecInfo.h"
00018 #include "TBDataFormats/EcalTBObjects/interface/EcalTBTDCRecInfo.h"
00019 #include "TBDataFormats/EcalTBObjects/interface/EcalTBEventHeader.h"
00020 
00021 
00022 //#include<fstream>
00023 
00024 #include "TFile.h"
00025 #include "TH1.h"
00026 #include "TH2.h"
00027 #include "TF1.h"
00028 
00029 #include <iostream>
00030 #include <string>
00031 #include <stdexcept>
00032 //
00033 // constants, enums and typedefs
00034 //
00035 
00036 //
00037 // static data member definitions
00038 //
00039 
00040 //
00041 // constructors and destructor
00042 //
00043 
00044 
00045 //========================================================================
00046 EcalSimpleTBAnalyzer::EcalSimpleTBAnalyzer( const edm::ParameterSet& iConfig ) : xtalInBeam_(0)
00047 //========================================================================
00048 {
00049    //now do what ever initialization is needed
00050    rootfile_          = iConfig.getUntrackedParameter<std::string>("rootfile","ecalSimpleTBanalysis.root");
00051    digiCollection_     = iConfig.getParameter<std::string>("digiCollection");
00052    digiProducer_       = iConfig.getParameter<std::string>("digiProducer");
00053    hitCollection_     = iConfig.getParameter<std::string>("hitCollection");
00054    hitProducer_       = iConfig.getParameter<std::string>("hitProducer");
00055    hodoRecInfoCollection_     = iConfig.getParameter<std::string>("hodoRecInfoCollection");
00056    hodoRecInfoProducer_       = iConfig.getParameter<std::string>("hodoRecInfoProducer");
00057    tdcRecInfoCollection_     = iConfig.getParameter<std::string>("tdcRecInfoCollection");
00058    tdcRecInfoProducer_       = iConfig.getParameter<std::string>("tdcRecInfoProducer");
00059    eventHeaderCollection_     = iConfig.getParameter<std::string>("eventHeaderCollection");
00060    eventHeaderProducer_       = iConfig.getParameter<std::string>("eventHeaderProducer");
00061 
00062 
00063    std::cout << "EcalSimpleTBAnalyzer: fetching hitCollection: " << hitCollection_.c_str()
00064         << " produced by " << hitProducer_.c_str() << std::endl;
00065 
00066 }
00067 
00068 
00069 //========================================================================
00070 EcalSimpleTBAnalyzer::~EcalSimpleTBAnalyzer()
00071 //========================================================================
00072 {
00073   // do anything here that needs to be done at desctruction time
00074   // (e.g. close files, deallocate resources etc.)
00075   // Amplitude vs TDC offset
00076 //   if (h_ampltdc)
00077 //   delete h_ampltdc;
00078   
00079 //   // Reconstructed energies
00080 //   delete h_e1x1;
00081 //   delete h_e3x3; 
00082 //   delete h_e5x5; 
00083   
00084 //   delete h_bprofx; 
00085 //   delete h_bprofy; 
00086   
00087 //   delete h_qualx; 
00088 //   delete h_qualy; 
00089   
00090 //   delete h_slopex; 
00091 //   delete h_slopey; 
00092   
00093 //   delete h_mapx; 
00094 //   delete h_mapy; 
00095 
00096 }
00097 
00098 //========================================================================
00099 void
00100 EcalSimpleTBAnalyzer::beginJob() {
00101 //========================================================================
00102 
00103   // Amplitude vs TDC offset
00104   h_ampltdc = new TH2F("h_ampltdc","Max Amplitude vs TDC offset", 100,0.,1.,1000, 0., 4000.);
00105 
00106   // Reconstructed energies
00107   h_tableIsMoving = new TH1F("h_tableIsMoving","TableIsMoving", 100000, 0., 100000.);
00108 
00109   h_e1x1 = new TH1F("h_e1x1","E1x1 energy", 1000, 0., 4000.);
00110   h_e3x3 = new TH1F("h_e3x3","E3x3 energy", 1000, 0., 4000.);
00111   h_e5x5 = new TH1F("h_e5x5","E5x5 energy", 1000, 0., 4000.);
00112 
00113   h_e1x1_center = new TH1F("h_e1x1_center","E1x1 energy", 1000, 0., 4000.);
00114   h_e3x3_center = new TH1F("h_e3x3_center","E3x3 energy", 1000, 0., 4000.);
00115   h_e5x5_center = new TH1F("h_e5x5_center","E5x5 energy", 1000, 0., 4000.);
00116 
00117   h_e1e9 = new TH1F("h_e1e9","E1/E9 ratio", 600, 0., 1.2);
00118   h_e1e25 = new TH1F("h_e1e25","E1/E25 ratio", 600, 0., 1.2);
00119   h_e9e25 = new TH1F("h_e9e25","E9/E25 ratio", 600, 0., 1.2);
00120 
00121   h_bprofx = new TH1F("h_bprofx","Beam Profile X",100,-20.,20.);
00122   h_bprofy = new TH1F("h_bprofy","Beam Profile Y",100,-20.,20.);
00123 
00124   h_qualx = new TH1F("h_qualx","Beam Quality X",5000,0.,5.);
00125   h_qualy = new TH1F("h_qualy","Beam Quality X",5000,0.,5.);
00126 
00127   h_slopex = new TH1F("h_slopex","Beam Slope X",500, -5e-4 , 5e-4 );
00128   h_slopey = new TH1F("h_slopey","Beam Slope Y",500, -5e-4 , 5e-4 );
00129 
00130   char hname[50];
00131   char htitle[50];
00132   for (unsigned int icry=0;icry<25;icry++)
00133     {       
00134       sprintf(hname,"h_mapx_%d",icry);
00135       sprintf(htitle,"Max Amplitude vs X %d",icry);
00136       h_mapx[icry] = new TH2F(hname,htitle,80,-20,20,1000,0.,4000.);
00137       sprintf(hname,"h_mapy_%d",icry);
00138       sprintf(htitle,"Max Amplitude vs Y %d",icry);
00139       h_mapy[icry] = new TH2F(hname,htitle,80,-20,20,1000,0.,4000.);
00140     }
00141   
00142   h_e1e9_mapx = new TH2F("h_e1e9_mapx","E1/E9 vs X",80,-20,20,600,0.,1.2);
00143   h_e1e9_mapy = new TH2F("h_e1e9_mapy","E1/E9 vs Y",80,-20,20,600,0.,1.2);
00144 
00145   h_e1e25_mapx = new TH2F("h_e1e25_mapx","E1/E25 vs X",80,-20,20,600,0.,1.2);
00146   h_e1e25_mapy = new TH2F("h_e1e25_mapy","E1/E25 vs Y",80,-20,20,600,0.,1.2);
00147 
00148   h_e9e25_mapx = new TH2F("h_e9e25_mapx","E9/E25 vs X",80,-20,20,600,0.,1.2);
00149   h_e9e25_mapy = new TH2F("h_e9e25_mapy","E9/E25 vs Y",80,-20,20,600,0.,1.2);
00150 
00151   h_Shape_ = new TH2F("h_Shape_","Xtal in Beam Shape",250,0,10,350,0,3500);
00152 
00153 }
00154 
00155 //========================================================================
00156 void
00157 EcalSimpleTBAnalyzer::endJob() {
00158 //========================================================================
00159 
00160   TFile f(rootfile_.c_str(),"RECREATE");
00161 
00162   // Amplitude vs TDC offset
00163   h_ampltdc->Write(); 
00164 
00165   // Reconstructed energies
00166   h_e1x1->Write(); 
00167   h_e3x3->Write(); 
00168   h_e5x5->Write(); 
00169 
00170   h_e1x1_center->Write(); 
00171   h_e3x3_center->Write(); 
00172   h_e5x5_center->Write(); 
00173 
00174   h_e1e9->Write(); 
00175   h_e1e25->Write(); 
00176   h_e9e25->Write(); 
00177 
00178   h_bprofx->Write(); 
00179   h_bprofy->Write(); 
00180 
00181   h_qualx->Write(); 
00182   h_qualy->Write(); 
00183 
00184   h_slopex->Write(); 
00185   h_slopey->Write(); 
00186   
00187   h_Shape_->Write();
00188 
00189   for (unsigned int icry=0;icry<25;icry++)
00190     {       
00191       h_mapx[icry]->Write(); 
00192       h_mapy[icry]->Write(); 
00193     }
00194 
00195   h_e1e9_mapx->Write(); 
00196   h_e1e9_mapy->Write(); 
00197 
00198   h_e1e25_mapx->Write(); 
00199   h_e1e25_mapy->Write(); 
00200 
00201   h_e9e25_mapx->Write(); 
00202   h_e9e25_mapy->Write(); 
00203 
00204   h_tableIsMoving->Write();
00205 
00206   f.Close();
00207 }
00208 
00209 //
00210 // member functions
00211 //
00212 
00213 //========================================================================
00214 void
00215 EcalSimpleTBAnalyzer::analyze( const edm::Event& iEvent, const edm::EventSetup& iSetup ) {
00216 //========================================================================
00217 
00218    using namespace edm;
00219    using namespace cms;
00220 
00221 
00222 
00223    Handle<EBDigiCollection> pdigis;
00224    const EBDigiCollection* digis=0;
00225    //std::cout << "EcalSimpleTBAnalyzer::analyze getting product with label: " << digiProducer_.c_str()<< " prodname: " << digiCollection_.c_str() << endl;
00226    iEvent.getByLabel( digiProducer_, digiCollection_,pdigis);
00227    if ( pdigis.isValid() ) {
00228      digis = pdigis.product(); // get a ptr to the product
00229      //iEvent.getByLabel( hitProducer_, phits);
00230    } else {
00231            edm::LogError("EcalSimpleTBAnalyzerError") << "Error! can't get the product " << digiCollection_;
00232    }
00233 
00234    // fetch the digis and compute signal amplitude
00235    Handle<EBUncalibratedRecHitCollection> phits;
00236    const EBUncalibratedRecHitCollection* hits=0;
00237    //std::cout << "EcalSimpleTBAnalyzer::analyze getting product with label: " << digiProducer_.c_str()<< " prodname: " << digiCollection_.c_str() << endl;
00238    iEvent.getByLabel( hitProducer_, hitCollection_,phits);
00239    if (phits.isValid()) {
00240      hits = phits.product(); // get a ptr to the product
00241      //iEvent.getByLabel( hitProducer_, phits);
00242    } else {
00243            edm::LogError("EcalSimpleTBAnalyzerError") << "Error! can't get the product " << hitCollection_;
00244    }
00245 
00246    Handle<EcalTBHodoscopeRecInfo> pHodo;
00247    const EcalTBHodoscopeRecInfo* recHodo=0;
00248    //std::cout << "EcalSimpleTBAnalyzer::analyze getting product with label: " << digiProducer_.c_str()<< " prodname: " << digiCollection_.c_str() << endl;
00249    iEvent.getByLabel( hodoRecInfoProducer_, hodoRecInfoCollection_, pHodo);
00250    if ( pHodo.isValid() ) {
00251      recHodo = pHodo.product(); // get a ptr to the product
00252    } else {
00253            edm::LogError("EcalSimpleTBAnalyzerError") << "Error! can't get the product " << hodoRecInfoCollection_;
00254    }
00255 
00256    Handle<EcalTBTDCRecInfo> pTDC;
00257    const EcalTBTDCRecInfo* recTDC=0;
00258    //std::cout << "EcalSimpleTBAnalyzer::analyze getting product with label: " << digiProducer_.c_str()<< " prodname: " << digiCollection_.c_str() << endl;
00259    iEvent.getByLabel( tdcRecInfoProducer_, tdcRecInfoCollection_, pTDC);
00260    if ( pTDC.isValid() ) {
00261      recTDC = pTDC.product(); // get a ptr to the product
00262    } else {
00263            edm::LogError("EcalSimpleTBAnalyzerError") << "Error! can't get the product " << tdcRecInfoCollection_;
00264    }
00265 
00266    Handle<EcalTBEventHeader> pEventHeader;
00267    const EcalTBEventHeader* evtHeader=0;
00268    //std::cout << "EcalSimpleTBAnalyzer::analyze getting product with label: " << digiProducer_.c_str()<< " prodname: " << digiCollection_.c_str() << endl;
00269    iEvent.getByLabel( eventHeaderProducer_ , pEventHeader );
00270    if ( pEventHeader.isValid() ) {
00271      evtHeader = pEventHeader.product(); // get a ptr to the product
00272    } else {
00273            edm::LogError("EcalSimpleTBAnalyzerError") << "Error! can't get the product " << eventHeaderProducer_;
00274    }
00275    
00276    if (!hits)
00277      return;
00278 
00279    if (!recTDC)
00280      return;
00281 
00282    if (!recHodo)
00283      return;
00284 
00285    if (!evtHeader)
00286      return;
00287 
00288    if (hits->size() == 0)
00289      return;
00290 
00291    if (evtHeader->tableIsMoving())
00292      h_tableIsMoving->Fill(evtHeader->eventNumber());
00293 
00294    // Crystal hit by beam
00295    if (xtalInBeam_.null())
00296      {
00297        xtalInBeam_ = EBDetId(1,evtHeader->crystalInBeam(),EBDetId::SMCRYSTALMODE);
00298        std::cout<< "Xtal In Beam is " << xtalInBeam_.ic() << std::endl;
00299      }
00300    else if (xtalInBeam_ != EBDetId(1,evtHeader->crystalInBeam(),EBDetId::SMCRYSTALMODE))
00301      return;
00302 
00303    if (evtHeader->tableIsMoving())
00304      return;
00305    
00306    
00307 //    EBDetId maxHitId(0); 
00308 //    float maxHit= -999999.;
00309    
00310 //    for(EBUncalibratedRecHitCollection::const_iterator ithit = hits->begin(); ithit != hits->end(); ++ithit) 
00311 //      {
00312 //        if (ithit->amplitude()>=maxHit)
00313 //       {
00314 //         maxHit=ithit->amplitude();
00315 //         maxHitId=ithit->id();
00316 //       }
00317        
00318 //      }   
00319 
00320 //    if (maxHitId==EBDetId(0))
00321 //      return;
00322 
00323 //   EBDetId maxHitId(1,704,EBDetId::SMCRYSTALMODE); 
00324 
00325    //Find EBDetId in a 5x5 Matrix (to be substituted by the Selector code)
00326    // Something like 
00327    // EBFixedWindowSelector<EcalUncalibratedRecHit> Simple5x5Matrix(hits,maxHitId,5,5);
00328    // std::vector<EcalUncalibratedRecHit> Energies5x5 = Simple5x5Matrix.getHits();
00329 
00330 
00331    EBDetId Xtals5x5[25];
00332    for (unsigned int icry=0;icry<25;icry++)
00333      {
00334        unsigned int row = icry / 5;
00335        unsigned int column= icry %5;
00336        int ieta=xtalInBeam_.ieta()+column-2;
00337        int iphi=xtalInBeam_.iphi()+row-2;
00338        EBDetId tempId(ieta, iphi,EBDetId::ETAPHIMODE);
00339        if (tempId.ism()==1) 
00340                Xtals5x5[icry]=tempId;
00341        else
00342                Xtals5x5[icry]=EBDetId(0);
00344      } 
00345 
00346 
00347    
00348    double samples_save[10]; for(int i=0; i < 10; ++i) samples_save[i]=0.0;
00349    
00350    // find the rechit corresponding digi and the max sample
00351    EBDigiCollection::const_iterator myDg = digis->find(xtalInBeam_);
00352    double eMax = 0.;
00353    if (myDg != digis->end())
00354      {
00355        EBDataFrame myDigi = (*myDg);
00356        for (int sample = 0; sample < myDigi.size(); ++sample)
00357          {
00358            double analogSample = myDigi.sample(sample).adc();
00359            samples_save[sample] = analogSample;
00360            //  std::cout << analogSample << " ";
00361            if ( eMax < analogSample )
00362              {
00363                eMax = analogSample;
00364              }
00365          }
00366        // std::cout << std::endl;
00367      }
00368 
00369    for(int i =0; i < 10; ++i) {
00370      h_Shape_->Fill(double(i)+recTDC->offset(),samples_save[i]);
00371    }
00372 
00373    double amplitude[25];
00374    
00375    double amplitude3x3=0;  
00376    double amplitude5x5=0;  
00377    
00378    for (unsigned int icry=0;icry<25;icry++)
00379      {
00380        if (!Xtals5x5[icry].null())
00381          {
00382            amplitude[icry]=(hits->find(Xtals5x5[icry]))->amplitude();
00383            amplitude5x5 += amplitude[icry];
00384            // Is in 3x3?
00385            if ( icry == 6  || icry == 7  || icry == 8 ||
00386                 icry == 11 || icry == 12 || icry ==13 ||
00387                 icry == 16 || icry == 17 || icry ==18   )
00388              {
00389                amplitude3x3+=amplitude[icry];
00390              }
00391          }
00392      }
00393 
00394 
00395    h_e1x1->Fill(amplitude[12]);
00396    h_e3x3->Fill(amplitude3x3);
00397    h_e5x5->Fill(amplitude5x5);
00398 
00399    h_e1e9->Fill(amplitude[12]/amplitude3x3);
00400    h_e1e25->Fill(amplitude[12]/amplitude5x5);
00401    h_e9e25->Fill(amplitude3x3/amplitude5x5);
00402 
00403    if (recTDC)
00404      h_ampltdc->Fill(recTDC->offset(),amplitude[12]);
00405 
00406    if (recHodo)
00407      {
00408        float x=recHodo->posX();
00409        float y=recHodo->posY();
00410        float xslope=recHodo->slopeX();
00411        float yslope=recHodo->slopeY();
00412        float xqual=recHodo->qualX();
00413        float yqual=recHodo->qualY();
00414        
00415        //Filling beam profiles
00416        h_bprofx->Fill(x);
00417        h_bprofy->Fill(y);
00418        h_qualx->Fill(xqual);
00419        h_qualy->Fill(yqual);
00420        h_slopex->Fill(xslope);
00421        h_slopey->Fill(yslope);
00422 
00423        //Fill central events
00424 
00425        
00426        if ( fabs(x + 2.5) < 2.5 && fabs(y + 0.5) < 2.5)
00427          {
00428            h_e1x1_center->Fill(amplitude[12]);
00429            h_e3x3_center->Fill(amplitude3x3);
00430            h_e5x5_center->Fill(amplitude5x5);
00431            
00432            h_e1e9->Fill(amplitude[12]/amplitude3x3);
00433            h_e1e25->Fill(amplitude[12]/amplitude5x5);
00434            h_e9e25->Fill(amplitude3x3/amplitude5x5);
00435          }
00436 
00437        for (unsigned int icry=0;icry<25;icry++)
00438          {       
00439            h_mapx[icry]->Fill(x,amplitude[icry]);
00440            h_mapy[icry]->Fill(y,amplitude[icry]);
00441          }
00442 
00443        h_e1e9_mapx->Fill(x,amplitude[12]/amplitude3x3);
00444        h_e1e9_mapy->Fill(y,amplitude[12]/amplitude3x3);
00445 
00446        h_e1e25_mapx->Fill(x,amplitude[12]/amplitude5x5);
00447        h_e1e25_mapy->Fill(y,amplitude[12]/amplitude5x5);
00448 
00449        h_e9e25_mapx->Fill(x,amplitude3x3/amplitude5x5);
00450        h_e9e25_mapy->Fill(y,amplitude3x3/amplitude5x5);
00451      }
00452 
00453 }
00454 
00455