CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/CalibCalorimetry/EcalLaserAnalyzer/plugins/EcalLaserAnalyzerYousi.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    EcalLaserAnalyzerYousi
00004 // Class:      EcalLaserAnalyzerYousi
00005 // 
00013 //
00014 // Original Author:  Yousi Ma
00015 //         Created:  Tue Jun 19 23:06:36 CEST 2007
00016 // $Id: EcalLaserAnalyzerYousi.cc,v 1.4 2012/02/09 10:07:36 eulisse Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 #include <iostream>
00024 #include <fstream>
00025 
00026 
00027 // user include files
00028 #include "FWCore/Framework/interface/Frameworkfwd.h"
00029 #include "FWCore/Framework/interface/EDAnalyzer.h"
00030 
00031 #include "FWCore/Framework/interface/Event.h"
00032 #include "FWCore/Framework/interface/MakerMacros.h"
00033 
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00036 
00037 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
00038 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00039 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00040 #include "DataFormats/EcalDigi/interface/EBDataFrame.h"
00041 #include "DataFormats/EcalRawData/interface/EcalRawDataCollections.h"
00042 
00043 #include "TFile.h"
00044 #include "TH1F.h"
00045 #include "TH2F.h"
00046 #include "TROOT.h"
00047 #include "TF1.h"
00048 #include "TNtuple.h"
00049 #include "TDirectory.h"
00050 
00051 
00052 
00053 //
00054 // class decleration
00055 //
00056 
00057 class EcalLaserAnalyzerYousi : public edm::EDAnalyzer {
00058    public:
00059       explicit EcalLaserAnalyzerYousi(const edm::ParameterSet&);
00060       ~EcalLaserAnalyzerYousi();
00061 
00062 
00063    private:
00064       virtual void beginJob();
00065       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00066       virtual void endJob();
00067 
00068       // ----------member data ---------------------------
00069   // Declare histograms and ROOT trees, etc.
00070   TH1F *C_APD[1700];
00071   TH1F *C_APDPN[1700];
00072   TH1F *C_PN[1700];
00073   TH1F *C_J[1700];
00074   TH2F *C_APDPN_J[1700];
00075   
00076   TH1F *peakAPD[2];
00077   TH1F *peakAPDPN[2];
00078   TH1F *APD_LM[9];
00079   TH1F *APDPN_LM[9];
00080   TH2F *APDPN_J_LM[9];
00081   TH2F *APDPN_J_H[2];
00082   
00083   
00084   //fixme make declare and init separate
00085   TH2F *APD; 
00086   TH2F *APD_RMS; 
00087   TH2F *APDPN; 
00088   TH2F *APDPN_RMS; 
00089   TH2F *PN; 
00090   TH2F *APDPN_J; 
00091   TH2F *APDPN_C; 
00092   
00093   TH1F *FitHist;
00094   TH2F *Count; 
00095   
00096   
00097   TFile *fPN;
00098   TFile *fAPD;
00099   TFile *fROOT;
00100 
00101   TNtuple *C_Tree[1700];
00102 
00103   //parameters
00104 
00105   std::string hitCollection_ ;
00106   std::string hitProducer_ ;
00107   //  std::string PNFileName_ ;
00108   //  std::string ABFileName_ ;
00109   std::string outFileName_ ;
00110   std::string SM_ ;
00111   std::string Run_ ;
00112   std::string digiProducer_ ;
00113   std::string PNdigiCollection_ ;
00114 
00115 
00116 };
00117 
00118 //
00119 // constants, enums and typedefs
00120 //
00121 
00122 //
00123 // static data member definitions
00124 //
00125 
00126 //
00127 // constructors and destructor
00128 //
00129 EcalLaserAnalyzerYousi::EcalLaserAnalyzerYousi(const edm::ParameterSet& iConfig)
00130 {
00131    //now do what ever initialization is needed
00132   //get the PN and AB file names
00133   //get the output file names, digi producers, etc
00134 
00135   hitCollection_ = iConfig.getUntrackedParameter<std::string>("hitCollection");
00136   hitProducer_ = iConfig.getUntrackedParameter<std::string>("hitProducer");
00137   //  PNFileName_ = iConfig.getUntrackedParameter<std::string>("PNFileName");
00138   //  ABFileName_ = iConfig.getUntrackedParameter<std::string>("ABFileName");
00139   outFileName_ = iConfig.getUntrackedParameter<std::string>("outFileName");
00140   SM_ = iConfig.getUntrackedParameter<std::string>("SM");
00141   Run_ = iConfig.getUntrackedParameter<std::string>("Run");
00142   digiProducer_ = iConfig.getUntrackedParameter<std::string>("digiProducer");
00143   PNdigiCollection_ = iConfig.getUntrackedParameter<std::string>("PNdigiCollection");
00144 
00145 
00146 
00147 
00148 }
00149 
00150 
00151 EcalLaserAnalyzerYousi::~EcalLaserAnalyzerYousi()
00152 {
00153  
00154    // do anything here that needs to be done at desctruction time
00155    // (e.g. close files, deallocate resources etc.)
00156 
00157 
00158 
00159 }
00160 
00161 
00162 //
00163 // member functions
00164 //
00165 
00166 // ------------ method called to for each event  ------------
00167 void
00168 EcalLaserAnalyzerYousi::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00169 {
00170    using namespace edm;
00171   
00172    //  if ( fPN->IsOpen() ) { edm::LogInfo("EcalLaserAnalyzerYousi") <<"fPN is open in analyze OKAAAAAAAAYYY \n\n"; }
00173 
00174 
00175 
00176    edm::Handle<EcalRawDataCollection> DCCHeaders;
00177    iEvent.getByLabel(digiProducer_, DCCHeaders);
00178  
00179    EcalDCCHeaderBlock::EcalDCCEventSettings settings = DCCHeaders->begin()->getEventSettings();
00180 
00181    int wavelength = settings.wavelength;
00182 
00183    //   std::cout<<"wavelength: "<<wavelength<<"\n\n";
00184 
00185    if ( wavelength !=0 ) { return; }//only process blue laser
00186 
00187    edm::Handle<EBUncalibratedRecHitCollection> hits;
00188 
00189    try{ 
00190      iEvent.getByLabel(hitProducer_ , hitCollection_ , hits);
00191      //     iEvent.getByType(hits);
00192    } catch ( std::exception& ex ) {
00193     LogError("EcalLaserAnalyzerYousi") << "Cannot get product:  EBRecHitCollection from: " 
00194                                     << hitCollection_ << " - returning.\n\n";
00195     //    return;
00196    }
00197    
00198    edm::Handle<EcalPnDiodeDigiCollection> pndigis;
00199 
00200    try{ 
00201      //     iEvent.getByLabel(hitProducer_, hits);
00202           iEvent.getByLabel(digiProducer_, PNdigiCollection_, pndigis);
00203      //iEvent.getByType( pndigis );
00204    } catch ( std::exception& ex ) {
00205     LogError("EcalLaserAnalyzerYousi") << "Cannot get product:  EBdigiCollection from: " 
00206                                     << "getbytype" << " - returning.\n\n";
00207     //    return;
00208    }
00209 
00210 
00211    Float_t PN_amp[5];
00212 
00213    //do some averaging over each pair of PNs
00214    for (int j = 0; j < 5; ++j) {
00215      PN_amp[j] = 0;
00216      for (int z = 0; z<2; ++z) {
00217        FitHist->Reset();
00218        TF1 peakFit("peakFit", "[0] +[1]*x +[2]*x^2", 30, 50);
00219        TF1 pedFit("pedFit","[0]",0,5);
00220        
00221        for(int k = 0; k < 50; k++){
00222          FitHist->SetBinContent(k, (*pndigis)[j+z*5].sample(k).adc() );
00223        }
00224        pedFit.SetParameter(0,750);
00225        FitHist->Fit(&pedFit,"RQI");
00226        Float_t  ped = pedFit.GetParameter(0);
00227        
00228        Int_t maxbin = FitHist->GetMaximumBin();
00229        peakFit.SetRange(FitHist->GetBinCenter(maxbin)-4*FitHist->GetBinWidth(maxbin),
00230                         FitHist->GetBinCenter(maxbin)+4*FitHist->GetBinWidth(maxbin));
00231        peakFit.SetParameters(750,4,-.05);
00232        FitHist->Fit(&peakFit,"RQI");
00233        Float_t max = peakFit.Eval(-peakFit.GetParameter(1)/(2*peakFit.GetParameter(2)));
00234        if(ped != max){
00235          PN_amp[j] = PN_amp[j] + max - ped;
00236        } else {
00237          PN_amp[j] = PN_amp[j] + max;
00238        }
00239        
00240      }//end of z loop
00241      PN_amp[j] = PN_amp[j]/2.0;
00242    }//end of j loop
00243 
00244 
00245 
00246 
00247    //do some real PN, APD calculations
00248 
00249 
00250 
00251    //FIXME. previously used .info files to get time, what to do now?
00252    
00253    //   TNtuple *Time = new TNtuple("Time", "Time", "Time");
00254    //   Int_t iTime = Get_Time(Input_File);
00255    //   Time->Fill(iTime);
00256    
00257 
00258    Float_t fTree[7];
00259    
00260 
00261    //   b->GetEntry(EVT);
00262    EBDetId ID;
00263    Float_t theAPD;
00264    Float_t thePN;
00265    Float_t Jitter;
00266    Float_t Chi2;
00267    Int_t CN = hits->size();
00268    //                   cout<<"num CN: "<<CN<<endl;
00269    for(int j = 0; j < CN; j++){
00270      ID = (*hits)[j].id();
00271      theAPD = (*hits)[j].amplitude();
00272      Jitter = (*hits)[j].jitter();
00273      Chi2 = (*hits)[j].chi2();
00274      thePN = PN_amp[ (ID.ic() + 299)/400 ];
00275      
00276      //                                 cout<<"THE APD: "<<theAPD<<endl;
00277      //                                 cout<<"THE PN: "<<thePN<<endl;
00278      
00279      C_APD[ID.ic()-1]->Fill(theAPD);
00280      C_APDPN[ID.ic()-1]->Fill(theAPD/thePN);
00281      C_PN[ID.ic()-1]->Fill(thePN);
00282      C_J[ID.ic()-1]->Fill(Jitter);
00283      C_APDPN_J[ID.ic()-1]->Fill(Jitter, theAPD/thePN);
00284      APDPN_J->Fill(Jitter, theAPD/thePN);
00285      APDPN_C->Fill(Chi2, theAPD/thePN);
00286      fTree[0] = theAPD;
00287      fTree[1] = thePN;
00288      fTree[2] = theAPD/thePN;
00289      fTree[3] = Jitter;
00290      fTree[4] = Chi2;
00291      fTree[5] = (*hits)[j].pedestal();
00292      fTree[6] = iEvent.id().event();
00293      C_Tree[ID.ic()-1]->Fill(fTree);
00294      if(((ID.ic()-1)%20 > 9)||((ID.ic()-1)<100)){
00295        peakAPD[0]->Fill(theAPD);
00296        peakAPDPN[0]->Fill(theAPD/thePN);
00297        APDPN_J_H[0]->Fill(Jitter, theAPD/thePN);
00298      } else {
00299        peakAPD[1]->Fill(theAPD);
00300        peakAPDPN[1]->Fill(theAPD/thePN);
00301        APDPN_J_H[1]->Fill(Jitter, theAPD/thePN);
00302      }
00303      if((ID.ic()-1) < 100){
00304        APD_LM[0]->Fill(theAPD);
00305        APDPN_LM[0]->Fill(theAPD/thePN);
00306        APDPN_J_LM[0]->Fill(Jitter, theAPD/thePN);
00307      } 
00308      else {
00309        Int_t index;
00310        if(((ID.ic()-1)%20) < 10){ 
00311          index = ((ID.ic()-101)/400)*2 + 1;
00312          APD_LM[index]->Fill(theAPD);
00313          APDPN_LM[index]->Fill(theAPD/thePN);
00314          APDPN_J_LM[index]->Fill(Jitter, theAPD/thePN);
00315        }
00316        else {
00317          index = ((ID.ic()-101)/400)*2 + 2;
00318          APD_LM[index]->Fill(theAPD);
00319          APDPN_LM[index]->Fill(theAPD/thePN);
00320          APDPN_J_LM[index]->Fill(Jitter, theAPD/thePN);
00321        }
00322      }
00323    }//end of CN loop
00324    
00325 
00326 
00327 
00328    //now that you got the PN and APD's, make the ntuples. done
00329 
00330    //vec from ROOT version should correspond to hits_itr or something similar. done
00331 
00332    //check WL from PNdiodedigi, should be ==0, o.w (blue data only). don't process. done
00333 
00334    //get PN pulse, and do fitting of pulse. i.e. fill hist with PN.apd() or equivalent. done
00335 
00336    //fit to first 5 for PED, and 30-50 bins for pulse (poly2 for the moment). done
00337 
00338 
00339 
00340 
00341 }
00342 
00343 
00344 // ------------ method called once each job just before starting event loop  ------------
00345 void 
00346 EcalLaserAnalyzerYousi::beginJob()
00347 {
00348 
00349   edm::LogInfo("EcalLaserAnalyzerYousi") << "running laser analyzer \n\n";
00350 
00351   fROOT = new TFile(outFileName_.c_str(), "RECREATE");
00352   fROOT->cd();
00353 
00354   //init all the histos and files?
00355   APD = new TH2F("APD", "APD", 85, 0., 85., 20, 0., 20.);
00356   APD_RMS = new TH2F("APD_RMS", "APD_RMS", 85, 0., 85., 20, 0., 20.);
00357   APDPN = new TH2F("APDPN", "APDPN", 85, 0., 85., 20, 0., 20.);
00358   APDPN_RMS = new TH2F("APDPN_RMS", "APDPN_RMS", 85, 0., 85., 20, 0., 20.);
00359   PN = new TH2F("PN", "PN", 85, 0., 85., 20, 0., 20.);
00360   APDPN_J = new TH2F("JittervAPDPN", "JittervAPDPN",250, 3., 7., 250, 1., 2.);
00361   APDPN_C = new TH2F("Chi2vAPDPN", "Chi2vAPDPN", 250, 0., 50., 250, 0., 5.0); 
00362   FitHist = new TH1F("FitHist", "FitHist", 50, 0, 50);
00363   Count = new TH2F("Count", "Count", 85, 0., 1., 20, 0., 1.);
00364 
00365 
00366   for(int i = 0; i < 1700; i++){
00367     std::ostringstream name_1;
00368     std::ostringstream name_2;
00369     std::ostringstream name_3;
00370     std::ostringstream name_4;
00371     std::ostringstream name_5;
00372     name_1 << "C_APD_" << i+1;
00373     name_2 << "C_APDPN_" << i+1;
00374     name_3 << "C_PN_" << i+1;
00375     name_4 << "C_J_" << i+1;
00376     name_5 << "C_APDPN_J_" << i+1;
00377     C_APD[i] = new TH1F(name_1.str().c_str(), name_1.str().c_str(), 2500, 0., 5000.);
00378     C_APDPN[i] = new TH1F(name_2.str().c_str(), name_2.str().c_str(), 20000, 0., 25.);
00379     C_PN[i] = new TH1F(name_3.str().c_str(), name_3.str().c_str(), 1000, 0., 4000.);
00380     C_J[i] = new TH1F(name_4.str().c_str(), name_4.str().c_str(), 250, 3.0, 7.);
00381     C_APDPN_J[i] = new TH2F(name_5.str().c_str(), name_5.str().c_str(), 250, 3.0, 6., 250, 1., 2.2);
00382   }
00383 
00384   for(int i = 0; i < 2; i++){
00385     std::ostringstream aname_1;
00386     std::ostringstream aname_2;
00387     std::ostringstream aname_3;
00388     aname_1 << "peakAPD_" << i;
00389     aname_2 << "peakAPDPN_" << i;
00390     aname_3 << "JittervAPDPN_Half_" << i;
00391     peakAPD[i] = new TH1F(aname_1.str().c_str(), aname_1.str().c_str(), 1000, 0., 5000.);
00392     peakAPDPN[i] = new TH1F(aname_2.str().c_str(), aname_2.str().c_str(), 1000, 0., 8.);
00393     APDPN_J_H[i] = new TH2F(aname_3.str().c_str(), aname_3.str().c_str(), 250, 3., 7., 250, 1., 2.2);
00394   }
00395   
00396   for(int i = 0; i < 9; i++){
00397     std::ostringstream bname_1;
00398     std::ostringstream bname_2;
00399     std::ostringstream bname_3;
00400     bname_1 << "APD_LM_" << i;
00401     bname_2 << "APDPN_LM_" << i;
00402     bname_3 << "APDPN_J_LM_" << i;
00403     APD_LM[i] = new TH1F(bname_1.str().c_str(), bname_1.str().c_str(), 500, 0., 5000.);
00404     APDPN_LM[i] = new TH1F(bname_2.str().c_str(), bname_2.str().c_str(), 500, 0., 8.);
00405     APDPN_J_LM[i] = new TH2F(bname_3.str().c_str(), bname_3.str().c_str(), 250, 3., 7., 250, 1., 2.2);
00406   }
00407 
00408   //get the PN file. or don't get, and read from event.
00409 
00410 
00411   //don't need to get AB, it will be read in via framework poolsource = ???
00412 
00413   //configure the final NTuple
00414    std::ostringstream varlist;
00415    varlist << "APD:PN:APDPN:Jitter:Chi2:ped:EVT";
00416    for(int i = 0; i < 1700; i++){
00417      std::ostringstream name;
00418      name << "C_Tree_" << i+1;
00419      C_Tree[i] = (TNtuple*) new TNtuple(name.str().c_str(), name.str().c_str(),
00420                                         varlist.str().c_str());
00421    }
00422 
00423 }
00424 
00425 // ------------ method called once each job just after ending the event loop  ------------
00426 void 
00427 EcalLaserAnalyzerYousi::endJob() {
00428 
00429   //write the file (get ouput file name first).
00430   TFile *fROOT = (TFile*) new TFile(outFileName_.c_str(), "RECREATE");
00431   
00432   //  TDirectory *DIR = fROOT->Get(Run_.c_str());
00433   TDirectory *DIR;
00434   //  if(DIR == NULL){
00435   DIR = fROOT->mkdir(Run_.c_str());
00436     //  }
00437   DIR->cd();
00438   for(int j = 0; j < 1700; j++){
00439     Float_t min_r, max_r;
00440     Float_t RMS, Sigma, K;
00441     Int_t iCount;
00442     TF1 *gs1;
00443     TF1 *gs2;
00444     TF1 *gs3;
00445     
00446     RMS = C_APD[j]->GetRMS();
00447     APD_RMS->SetBinContent(85-(j/20), 20-(j%20), RMS);
00448     Sigma = 999999;
00449     K = 2.5;
00450     iCount = 0;
00451     while(Sigma > RMS){
00452       min_r = C_APD[j]->GetBinCenter(C_APD[j]->GetMaximumBin()) - K*RMS;
00453       max_r = C_APD[j]->GetBinCenter(C_APD[j]->GetMaximumBin()) + K*RMS;
00454       gs1 = new TF1("gs1", "gaus", min_r, max_r);
00455       C_APD[j]->Fit("gs1", "RQI");
00456       Sigma = gs1->GetParameter(2);
00457       K = K*1.5;
00458       iCount++;
00459       if(iCount > 2){
00460         C_APD[j]->Fit("gaus", "QI");
00461         gs1 = C_APD[j]->GetFunction("gaus");
00462         break;
00463       }
00464     }
00465     
00466     RMS = C_APDPN[j]->GetRMS();
00467     APDPN_RMS->SetBinContent(85-(j/20), 20-(j%20), RMS);
00468     Sigma = 999999;
00469     K = 2.5;
00470     iCount = 0;
00471     while(Sigma > RMS){
00472       min_r = C_APDPN[j]->GetBinCenter(C_APDPN[j]->GetMaximumBin()) - K*RMS;
00473       max_r = C_APDPN[j]->GetBinCenter(C_APDPN[j]->GetMaximumBin()) + K*RMS;
00474       gs2 = new TF1("gs2", "gaus", min_r, max_r);
00475       C_APDPN[j]->Fit("gs2", "RQI");
00476       Sigma = gs2->GetParameter(2);
00477       K = K*1.5;
00478       iCount++;
00479       if(iCount > 2){
00480         C_APDPN[j]->Fit("gaus", "QI");
00481         gs2 = C_APDPN[j]->GetFunction("gaus");
00482         break;
00483       }
00484     }
00485 
00486     TF1 *newgs1;
00487     TF1 *newgs2;
00488     
00489     C_PN[j]->Fit("gaus", "Q");
00490     C_APD[j]->Fit("gaus","QI");
00491     C_APDPN[j]->Fit("gaus","QI");
00492     C_APD[j]->Write("", TObject::kOverwrite);
00493     C_APDPN[j]->Write("", TObject::kOverwrite);
00494     C_PN[j]->Write("", TObject::kOverwrite);
00495     C_J[j]->Write("", TObject::kOverwrite);
00496     C_APDPN_J[j]->Write("", TObject::kOverwrite);
00497     newgs1 = C_APD[j]->GetFunction("gaus");
00498     newgs2 = C_APDPN[j]->GetFunction("gaus");
00499     gs3 = C_PN[j]->GetFunction("gaus");
00500     Float_t theAPD = newgs1->GetParameter(1);
00501     APD->SetBinContent(85-(j/20), 20-(j%20), theAPD);
00502     Float_t theAPDPN = newgs2->GetParameter(1);
00503     APDPN->SetBinContent(85-(j/20), 20-(j%20), theAPDPN);
00504     Float_t thePN = gs3->GetParameter(1);
00505     //          cout<<"LOOK HERE thePN = "<< thePN<<endl;
00506     PN->SetBinContent(85-(j/20), 20-(j%20), thePN);
00507     C_Tree[j]->Write("", TObject::kOverwrite);
00508   }
00509   
00510   for(int i = 0; i < 9; i++){
00511     APD_LM[i]->Write("", TObject::kOverwrite);
00512     APDPN_LM[i]->Write("", TObject::kOverwrite);
00513     APDPN_J_LM[i]->Write("", TObject::kOverwrite);
00514   }
00515   
00516   
00517   //  Time->Write("", TObject::kOverwrite);
00518   APD->Write("", TObject::kOverwrite);
00519   APD_RMS->Write("", TObject::kOverwrite);
00520   APDPN_RMS->Write("", TObject::kOverwrite);
00521   APDPN->Write("", TObject::kOverwrite);
00522   APDPN_J->Write("", TObject::kOverwrite);
00523   APDPN_C->Write("", TObject::kOverwrite);
00524   PN->Write("", TObject::kOverwrite);
00525   peakAPD[0]->Write("", TObject::kOverwrite);
00526   peakAPD[1]->Write("", TObject::kOverwrite);
00527   peakAPDPN[0]->Write("", TObject::kOverwrite);
00528   peakAPDPN[1]->Write("", TObject::kOverwrite);
00529   APDPN_J_H[0]->Write("", TObject::kOverwrite);
00530   APDPN_J_H[1]->Write("", TObject::kOverwrite);
00531   
00532   
00533   
00534   // don't Make plots
00535   //  fROOT->Close();
00536   
00537 //   fPN->Close();
00538 //   fAPD->Close();
00539   
00540   
00541     fROOT->Write();
00542 //   fROOT->Close();
00543 
00544 }
00545 
00546 //define this as a plug-in
00547 DEFINE_FWK_MODULE(EcalLaserAnalyzerYousi);