CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/CalibCalorimetry/HcalStandardModules/src/HFLightCalRand.cc

Go to the documentation of this file.
00001 // Analysis of HF LED/Laser run: 
00002 // Case when signal has random time position in TS
00003 // SPE calibration for low light intensity or raw SPE calibration for high light intensity
00004 // and HF performance based on this analysis
00005 //
00006 // Igor Vodopiyanov. Oct-2007 .... update Sept-2008
00007 // Thanks G.Safronov, M.Mohammadi, F.Ratnikov
00008 //
00009 #include <memory>
00010 #include <string>
00011 #include <iostream>
00012 
00013 #include "TH1F.h"
00014 #include "TH2F.h"
00015 #include "TFile.h"
00016 #include "math.h"
00017 #include "TMath.h"
00018 #include "TF1.h"
00019 
00020 #include "CalibCalorimetry/HcalStandardModules/interface/HFLightCalRand.h"
00021 
00022 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00023 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00024 
00025 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00026 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00027 #include "CondFormats/HcalObjects/interface/HcalQIEShape.h"
00028 #include "CondFormats/HcalObjects/interface/HcalQIECoder.h"
00029 
00030 #include "FWCore/Framework/interface/Event.h"
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 #include "FWCore/Framework/interface/EventSetup.h"
00033 #include "FWCore/Framework/interface/ESHandle.h"
00034 
00035 using namespace std;
00036 Int_t NEvents, runNumb=0, EventN=0;
00037 
00038 namespace {
00039   //bool verbose = true;
00040   bool verbose = false;
00041 }
00042 
00043 HFLightCalRand::HFLightCalRand (const edm::ParameterSet& fConfiguration) :
00044   hfDigiCollectionTag_(fConfiguration.getParameter<edm::InputTag>("hfDigiCollectionTag")),
00045   hcalCalibDigiCollectionTag_(fConfiguration.getParameter<edm::InputTag>("hcalCalibDigiCollectionTag")) {
00046 
00047   //std::string histfile = fConfiguration.getUntrackedParameter<string>("rootFile");
00048   histfile = fConfiguration.getUntrackedParameter<string>("rootFile");
00049   textfile = fConfiguration.getUntrackedParameter<string>("textFile");
00050 }
00051 
00052 HFLightCalRand::~HFLightCalRand () {
00053   //delete mFile;
00054 }
00055 
00056 void HFLightCalRand::beginJob() {
00057 
00058   char htit[64];
00059   mFile = new TFile (histfile.c_str(),"RECREATE");
00060   if ((tFile = fopen(textfile.c_str(),"w"))==NULL)printf("\nNo textfile open\n\n");
00061   // General Histos
00062   htmax = new TH1F("htmax","Max TS",10,-0.5,9.5);
00063   htmean = new TH1F("htmean","Mean signal TS",100,0,10);
00064   hsignalmean = new TH1F("hsignalmean","Mean ADC 4maxTS",1201,-25,30000);
00065   hsignalrms = new TH1F("hsignalrms","RMS ADC 4maxTS",500,0,500);
00066   hpedmean = new TH1F("hpedmean","Mean ADC 4lowTS",200,-10,90);
00067   hpedrms = new TH1F("hpedrms","RMS ADC 4lowTS",200,0,100);
00068   hspes = new TH1F("hspes","SPE if measured",200,0,40);
00069   hnpevar = new TH1F("hnpevar","~N PE input",500,0,500);
00070   hsignalmapP = new TH2F("hsignalmapP","Mean(Response) - Mean(Pedestal) HFP;#eta;#phi",26,28.5,41.5,36,0,72);
00071   hsignalRMSmapP = new TH2F("hsignalRMSmapP","RMS Response HFP;#eta;#phi",26,28.5,41.5,36,0,72);
00072   hnpemapP = new TH2F("hnpemapP","~N PE input HFP;#eta;#phi",26,28.5,41.5,36,0,72);
00073   hnpemapP->SetOption("COLZ");hsignalmapP->SetOption("COLZ");hsignalRMSmapP->SetOption("COLZ");
00074   hsignalmapM = new TH2F("hsignalmapM","Mean(Response) - Mean(Pedestal) HFM;#eta;#phi",26,-41.5,-28.5,36,0,72);
00075   hsignalRMSmapM = new TH2F("hsignalRMSmapM","RMS Response HFM;#eta;#phi",26,-41.5,-28.5,36,0,72);
00076   hnpemapM = new TH2F("hnpemapM","~N PE input HFM;#eta;#phi",26,-41.5,-28.5,36,0,72);
00077   hnpemapM->SetOption("COLZ");hsignalmapM->SetOption("COLZ");hsignalRMSmapM->SetOption("COLZ");
00078   // Channel-by-channel histos
00079   for (int i=0;i<13;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
00080     if (i>10 && j%2==0) continue;
00081     sprintf(htit,"ts_+%d_%d_%d",i+29,j*2+1,k+1);
00082     hts[i][j][k] = new TH1F(htit,htit,10,-0.5,9.5); // TimeSlices (pulse shape)
00083     sprintf(htit,"tsmean_+%d_%d_%d",i+29,j*2+1,k+1);
00084     htsm[i][j][k] = new TH1F(htit,htit,100,0,10);   // Mean signal time estimated from TS 
00085     sprintf(htit,"sp_+%d_%d_%d",i+29,j*2+1,k+1);
00086     hsp[i][j][k] = new TH1F(htit,htit,1201,-25,30000); // Big-scale spectrum (linear ADC)
00087     sprintf(htit,"spe_+%d_%d_%d",i+29,j*2+1,k+1);
00088     hspe[i][j][k] = new TH1F(htit,htit,200,-9.5,190.5); // Small-scale spectrum (linear ADC)
00089     sprintf(htit,"ped_+%d_%d_%d",i+29,j*2+1,k+1);
00090     hped[i][j][k] = new TH1F(htit,htit,200,-9.5,190.5); // Pedestal spectrum
00091     sprintf(htit,"ts_-%d_%d_%d",i+29,j*2+1,k+1);
00092     hts[i+13][j][k] = new TH1F(htit,htit,10,-0.5,9.5);
00093     sprintf(htit,"tsmean_-%d_%d_%d",i+29,j*2+1,k+1);
00094     htsm[i+13][j][k] = new TH1F(htit,htit,100,0,10);  
00095     sprintf(htit,"sp_-%d_%d_%d",i+29,j*2+1,k+1);
00096     hsp[i+13][j][k] = new TH1F(htit,htit,1201,-25,30000);
00097     sprintf(htit,"spe_-%d_%d_%d",i+29,j*2+1,k+1);
00098     hspe[i+13][j][k] = new TH1F(htit,htit,200,-9.5,190.5); 
00099     sprintf(htit,"ped_-%d_%d_%d",i+29,j*2+1,k+1);
00100     hped[i+13][j][k] = new TH1F(htit,htit,200,-9.5,190.5); 
00101   } 
00102   // PIN-diodes histos
00103   for (int i=0;i<4;i++) for (int j=0;j<3;j++) {
00104     sprintf(htit,"ts_PIN%d_+Q%d",j+1,i+1);
00105     htspin[i][j] = new TH1F(htit,htit,10,-0.5,9.5);
00106     sprintf(htit,"sp_PIN%d_+Q%d",j+1,i+1);
00107     hsppin[i][j] = new TH1F(htit,htit,1601,-25,40000);
00108     sprintf(htit,"spe_PIN%d_+Q%d",j+1,i+1);
00109     hspepin[i][j] = new TH1F(htit,htit,200,-9.5,190.5);
00110     sprintf(htit,"ped_PIN%d_+Q%d",j+1,i+1);
00111     hpedpin[i][j] = new TH1F(htit,htit,200,-9.5,190.5);
00112     sprintf(htit,"tsmean_PIN%d_+Q%d",j+1,i+1);
00113     htsmpin[i][j] = new TH1F(htit,htit,100,0,10);  
00114     sprintf(htit,"ts_PIN%d_-Q%d",j+1,i+1);
00115     htspin[i+4][j] = new TH1F(htit,htit,10,-0.5,9.5);
00116     sprintf(htit,"sp_PIN%d_-Q%d",j+1,i+1);
00117     hsppin[i+4][j] = new TH1F(htit,htit,1601,-25,40000);
00118     sprintf(htit,"spe_PIN%d_-Q%d",j+1,i+1);
00119     hspepin[i+4][j] = new TH1F(htit,htit,200,-9.5,190.5);
00120     sprintf(htit,"ped_PIN%d_-Q%d",j+1,i+1);
00121     hpedpin[i+4][j] = new TH1F(htit,htit,200,-9.5,190.5);
00122     sprintf(htit,"tsmean_PIN%d_-Q%d",j+1,i+1);
00123     htsmpin[i+4][j] = new TH1F(htit,htit,100,0,10);  
00124   }
00125   std::cout<<std::endl<<"beginJob: histfile="<<histfile.c_str()<<"  textfile="<<textfile.c_str()<<std::endl;
00126   return;
00127 }
00128 
00129 void HistSpec(TH1F* hist, Double_t &mean, Double_t &rms, Double_t range=4) {
00130   Double_t xmin,xmax;
00131   mean=hist->GetMean();
00132   rms=hist->GetRMS();
00133   xmin=hist->GetXaxis()->GetXmin();
00134   xmax=hist->GetXaxis()->GetXmax();
00135   hist->SetAxisRange(mean-range*rms-100,mean+range*rms+100);
00136   mean=hist->GetMean();
00137   rms=hist->GetRMS();
00138   hist->SetAxisRange(mean-range*rms-100,mean+range*rms+100);
00139   mean=hist->GetMean();
00140   rms=hist->GetRMS();
00141   hist->SetAxisRange(xmin,xmax);
00142   return;
00143 }
00144 
00145 Double_t Fit3Peak(Double_t *x, Double_t *par) { 
00146 // Spectra fit function: Pedestal Gaussian + asymmetric 1PE + 2PE +3PE peaks
00147 
00148   Double_t sum,xx,A0,C0,r0,sigma0,mean1,sigma1,A1,C1,r1,mean2,sigma2,A2,C2,r2,mean3,sigma3,A3,C3,r3;
00149 
00150   const Double_t k0=2.0,k1=1.6, k2=2.0;
00151 
00152   xx=x[0];
00153   sigma0 = par[2];
00154   A0 = 2*NEvents/(2+2*par[0]+par[0]*par[0]+par[0]*par[0]*par[0]/3);
00155   r0 = ((xx-par[1])/sigma0);
00156   C0 = 1/(sigma0* TMath::Exp(-k0*k0/2)/k0 +
00157           sigma0*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k0/1.41421)));
00158   //sum = 1/(sqrt(2*3.14159)*par[2])*A0*TMath::Exp(-0.5*r0*r0);
00159   if(r0 < k0) sum = C0*A0*TMath::Exp(-0.5*r0*r0);
00160   else sum = C0*A0*TMath::Exp(0.5*k0*k0-k0*r0);
00161 
00162   mean1 = par[1]+par[3];
00163   //sigma1 = par[4];
00164   sigma1 = 1.547+0.125*par[3]+0.004042*par[3]*par[3];
00165   sigma1 = (sigma1+(9.1347e-3+3.845e-2*par[3])*par[4]*2.0)*par[2];
00166   A1 = A0*par[0];
00167   C1 = 1/(sigma1* TMath::Exp(-k1*k1/2)/k1 +
00168           sigma1*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k1/1.41421)));
00169   r1 = ((xx-mean1)/sigma1);
00170   if(r1 < k1) sum += C1*A1*TMath::Exp(-0.5*r1*r1);
00171   else sum += C1*A1*TMath::Exp(0.5*k1*k1-k1*r1);
00172 
00173   mean2 = 2*par[3]+par[1];
00174   sigma2 = sqrt(2*sigma1*sigma1 - par[2]*par[2]);
00175   //A2 = A0*par[5]*par[0]*par[0]/2;
00176   A2 = A0*par[0]*par[0]/2;
00177   C2 = 1/(sigma2* TMath::Exp(-k2*k2/2)/k2 +
00178           sigma2*sqrt(2*3.14159)*0.5*(1+TMath::Erf(k2/1.41421)));
00179   r2 = ((xx-mean2)/sigma2);
00180   if(r2 < k2) sum += C2*A2*TMath::Exp(-0.5*r2*r2);
00181   else sum += C2*A2*TMath::Exp(0.5*k2*k2-k2*r2);
00182 
00183   mean3 = 3*par[3]+par[1];
00184   sigma3 = sqrt(3*sigma1*sigma1 - 2*par[2]*par[2]);
00185   A3 = A0*par[0]*par[0]*par[0]/6;
00186   C3 = 1/(sigma3*sqrt(2*3.14159));
00187   r3 = ((xx-mean3)/sigma3);
00188   sum += C3*A3*TMath::Exp(-0.5*r3*r3);
00189 
00190   return sum;
00191 }
00192 
00193 void HFLightCalRand::endJob(void)
00194 {
00195   Double_t mean,rms,meanped,rmsped,npevar;
00196   Double_t par[5],parm[5],dspe=0,dnpe;
00197   Int_t tsmax,intspe;
00198   fprintf(tFile,"#RunN %d   Events processed %d",runNumb,EventN);
00199   std::cout<<"endJob: histos processing..."<<std::endl;
00200   std::cout<<"RunN= "<<runNumb<<"  Events processed= "<<EventN<<std::endl;
00201 
00202   for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
00203     if (i>10 && i<13 && j%2==0) continue;
00204     if (i>23 && j%2==0) continue;
00205     meanped=rmsped=mean=rms=0;
00206     if (hsp[i][j][k]->Integral()>0) {
00207       HistSpec(hped[i][j][k],meanped,rmsped);
00208       HistSpec(hsp[i][j][k],mean,rms);
00209       if (hspe[i][j][k]->Integral()>hsp[i][j][k]->Integral()*0.9 || mean<100) {
00210         HistSpec(hspe[i][j][k],mean,rms);
00211       }
00212       hsignalmean->Fill(mean);
00213       hsignalrms->Fill(rms);
00214       hpedmean->Fill(meanped);
00215       hpedrms->Fill(rmsped);
00216     }
00217   }
00218 
00219   meanped=hpedmean->GetMean();
00220   rmsped=hpedrms->GetMean();
00221   mean=hsignalmean->GetMean();
00222   rms=hsignalrms->GetMean();
00223   fprintf(tFile,"   MeanInput=<%.2f> [linADCcount]   RMS=<%.2f>\n",mean,rms);
00224   fprintf(tFile,"#eta/phi/depth  sum4maxTS     RMS      ~N_PE  sum4lowTS     RMS  maxTS  SPE +/- Err   Comment\n");
00225   TF1* fPed = new TF1("fPed","gaus",0,120);
00226   fPed->SetNpx(200);
00227   TF1 *fTot = new TF1("fTot",Fit3Peak ,0,200,5);
00228   fTot->SetNpx(800);
00229   for (int i=0;i<26;i++) for (int j=0;j<36;j++) for (int k=0;k<2;k++) {
00230     if (i>10 && i<13 && j%2==0) continue;
00231     if (i>23 && j%2==0) continue;
00232     HistSpec(hped[i][j][k],meanped,rmsped);
00233     HistSpec(hsp[i][j][k],mean,rms);
00234     par[3]=0;
00235     if (hspe[i][j][k]->Integral()>hsp[i][j][k]->Integral()*0.9 || mean<100) {
00236       HistSpec(hspe[i][j][k],mean,rms);
00237       if (hspe[i][j][k]->Integral(1,(int) (meanped+3*rmsped+12))/NEvents>0.1) {
00238         //if (hspe[i][j][k]->Integral()>100 && mean-meanped<100) {
00239         if (mean+rms*3-meanped-rmsped*3>1 && rmsped>0) { // SPE fit if low intensity>0
00240           par[1] = meanped;
00241           par[2] = rmsped;
00242           par[0] = hped[i][j][k]->GetMaximum();
00243           fPed->SetParameters(par);
00244           hped[i][j][k]->Fit(fPed,"BQ0");
00245           fPed->GetParameters(&par[0]);
00246           hped[i][j][k]->Fit(fPed,"B0Q","",par[1]-par[2]*3,par[1]+par[2]*3);
00247           fPed->GetParameters(par);
00248           hped[i][j][k]->Fit(fPed,"BLIQ","",par[1]-par[2]*3,par[1]+par[2]*3);
00249           fPed->GetParameters(&par[0]);
00250           fPed->GetParameters(&parm[0]);
00251           NEvents = (int) hspe[i][j][k]->Integral();
00252           par[0]=0.1;
00253           par[3]=10;
00254           par[4]=6;
00255           par[2]=par[2]/0.93;
00256           fTot->SetParameters(par);
00257           fTot->SetParLimits(0,0,2);
00258           //fTot->FixParameter(1,par[1]);
00259           //fTot->FixParameter(2,par[2]);
00260           fTot->SetParLimits(1,par[1]-4,par[1]+4);
00261           fTot->SetParLimits(2,par[2]*0.9,par[2]);
00262           fTot->SetParLimits(3,1.2,50);
00263           //fTot->SetParLimits(4,par[2]*1.1,100);
00264           par[4]=0;
00265           fTot->SetParLimits(4,-1.64,1.64);
00266           hspe[i][j][k]->Fit(fTot,"BLEQ","");
00267           fTot->GetParameters(par);
00268           dspe=fTot->GetParError(3);
00269           dnpe=fTot->GetParError(0);
00270           if (par[3]<1.21 || dnpe>par[0]) par[3]=-1;
00271           else if (par[0]>1.96 || par[3]>49) par[3]=0;
00272           else {
00273             hspes->Fill(par[3]);
00274           }
00275         } 
00276       }
00277     }
00278 
00279     // NPE
00280     npevar=0;
00281     if (par[3]>0) npevar=par[0];                          // NPE from SPE fit
00282     else {                                                // NPE from high intensity signal
00283       if (hspe[i][j][k]->Integral()>hsp[i][j][k]->Integral()*0.98) {
00284         HistSpec(hspe[i][j][k],mean,rms,3);
00285       }
00286       else {
00287         HistSpec(hsp[i][j][k],mean,rms,3);
00288       }
00289       if (rmsped>0) {
00290         if ((rms*rms-rmsped*rmsped/0.8836)>1 && mean>meanped+2) {
00291           npevar=(mean-meanped-2)*(mean-meanped-2)/(rms*rms-rmsped*rmsped/0.8836);
00292         }
00293         else if (mean<100) {
00294           intspe=int(hspe[i][j][k]->Integral());
00295           hspe[i][j][k]->SetAxisRange(meanped+2+rmsped*4,300);
00296           npevar=hspe[i][j][k]->Integral()/intspe;
00297           if (npevar>0.01) npevar=-1;
00298           else npevar=0;
00299           hspe[i][j][k]->SetAxisRange(-20,300);
00300         }
00301       }
00302     }
00303     if (npevar>5.0e-5) hnpevar->Fill(npevar);
00304 
00305     if (i<13) {
00306       hsignalmapP->Fill(i+28.6+k/2.0,j*2+1,mean-meanped-1.8); 
00307       hsignalRMSmapP->Fill(i+28.6+k/2.0,j*2+1,rms);
00308       if (npevar>0) hnpemapP->Fill(i+28.6+k/2.0,j*2+1,npevar);
00309       fprintf(tFile,"%3d%4d%5d  %11.2f%8.2f",i+29,j*2+1,k+1,mean,rms);
00310     }
00311     else {
00312       fprintf(tFile,"%3d%4d%5d  %11.2f%8.2f",13-i-29,j*2+1,k+1,mean,rms);
00313       hsignalmapM->Fill(13-i-28.6-k/2.0,j*2+1,mean-meanped-1.8);
00314       hsignalRMSmapM->Fill(13-i-28.6-k/2.0,j*2+1,rms);
00315       if (npevar>0) hnpemapM->Fill(13-i-28.6-k/2.0,j*2+1,npevar);
00316     }
00317     if (npevar>0) fprintf(tFile,"  %9.4f",npevar);
00318     else  fprintf(tFile,"      0    ");
00319     fprintf(tFile,"   %8.2f%8.2f",meanped,rmsped);
00320     tsmax=hts[i][j][k]->GetMaximumBin()-1;
00321     fprintf(tFile," %4d",tsmax);
00322     if (par[3]>0 && par[3]<99)  fprintf(tFile,"%8.2f%7.2f",par[3]+1,dspe);
00323     else if (npevar>0) fprintf(tFile,"%8.2f    0  ",(mean-meanped-1)/npevar);
00324     else fprintf(tFile,"     0      0  ");
00325 
00326     // Diagnostics 
00327     fprintf(tFile,"    ");
00328     if (hsp[i][j][k]->GetEntries()<=0) fprintf(tFile,"NoSignal\n");
00329     else if (hsp[i][j][k]->GetEntries()<=10) fprintf(tFile,"Nev<10\n");
00330     else {
00331       if (hsp[i][j][k]->Integral()<=10 || mean>12000)  fprintf(tFile,"SignalOffRange\n");
00332       else {
00333         if (hsp[i][j][k]->Integral()<100)  fprintf(tFile,"Nev<100/");
00334         if (npevar>0 && par[3]>0 && (npevar*NEvents<10 || npevar<0.001)) 
00335           fprintf(tFile,"LowSignal/");
00336         else if (npevar==0 && fabs(mean-meanped)<3) fprintf(tFile,"LowSignal/");
00337         if (par[3]<0)  fprintf(tFile,"BadFit/");
00338         else if (par[3]==0)  fprintf(tFile,"NoSPEFit/");
00339         else if (par[3]>0 && npevar>1)   fprintf(tFile,"NPE>1/");
00340         if (npevar<0)   fprintf(tFile,"Problem/");
00341         if (mean<2) fprintf(tFile,"LowMean/");
00342         if (rms<0.5) fprintf(tFile,"LowRMS/"); 
00343         if (meanped<-1) fprintf(tFile,"LowPed/");
00344         else if (meanped>25) fprintf(tFile,"HighPed/"); 
00345         if (rmsped<0.5 && rmsped>0) fprintf(tFile,"NarrowPed/"); 
00346         else if (rmsped>10) fprintf(tFile,"WidePed/");
00347         if (hped[i][j][k]->GetBinContent(201)>10) fprintf(tFile,"PedOffRange"); 
00348         fprintf(tFile,"-\n");
00349       }
00350     }
00351   }
00352 
00353   for (int i=0;i<8;i++) for (int j=0;j<3;j++) {
00354     HistSpec(hpedpin[i][j],meanped,rmsped);
00355     HistSpec(hsppin[i][j],mean,rms);
00356     if (hspepin[i][j]->Integral()>hsppin[i][j]->Integral()*0.9 || mean<100) {
00357       HistSpec(hspepin[i][j],mean,rms);
00358     }
00359     if (i<4) fprintf(tFile," PIN%d  +Q%d  %12.2f  %6.2f",j+1,i+1,mean,rms);
00360     else fprintf(tFile," PIN%d  -Q%d  %12.2f  %6.2f",j+1,i-3,mean,rms);
00361     fprintf(tFile,"  %15.2f  %6.2f",meanped,rmsped);
00362     tsmax=htspin[i][j]->GetMaximumBin()-1;
00363     fprintf(tFile,"  %4d\n",tsmax);
00364   } 
00365 
00366   mFile->Write();
00367   mFile->Close();
00368   fclose(tFile);
00369   std::cout<<std::endl<<" --endJob-- done"<<std::endl;
00370   return;
00371 }
00372 
00373 void HFLightCalRand::analyze(const edm::Event& fEvent, const edm::EventSetup& fSetup) {
00374 
00375   // event ID
00376   edm::EventID eventId = fEvent.id();
00377   int runNumber = eventId.run ();
00378   int eventNumber = eventId.event ();
00379   if (runNumb==0) runNumb=runNumber;
00380   EventN++;
00381   if (verbose) std::cout << "========================================="<<std::endl
00382                          << "run/event: "<<runNumber<<'/'<<eventNumber<<std::endl;
00383 
00384   Double_t buf[20];
00385   Double_t maxADC,signal,ped=0,meant;
00386   Int_t maxisample=0,i1=3,i2=6;
00387 
00388   // HF PIN-diodes
00389   edm::Handle<HcalCalibDigiCollection> calib;  
00390   fEvent.getByLabel(hcalCalibDigiCollectionTag_, calib);
00391   if (verbose) std::cout<<"Analysis-> total CAL digis= "<<calib->size()<<std::endl;
00392   /* COMMENTED OUT by J. Mans (7-28-2008) as major changes needed with new Calib DetId 
00393    re-commented out by R.Ofierzynski (11.Nov.2008) - to be able to provide a consistent code for CMSSW_3_0_0_pre3:
00394    major changes are needed for the new Calib DetId which does not have the old methods any more
00395 
00396   for (unsigned j = 0; j < calib->size (); ++j) {
00397     const HcalCalibDataFrame digi = (*calib)[j];
00398     HcalElectronicsId elecId = digi.elecId();
00399     HcalCalibDetId calibId = digi.id();
00400     if (verbose) std::cout<<calibId.sectorString().c_str()<<" "<<calibId.rbx()<<" "<<elecId.fiberChanId()<<std::endl;
00401     int isector = calibId.rbx()-1;
00402     int ipin = elecId.fiberChanId();
00403     int iside = -1;
00404     if (calibId.sectorString() == "HFP") iside = 0; 
00405     else if (calibId.sectorString() == "HFM") iside = 4;
00406 
00407     if (iside != -1) {
00408       maxADC=-99;
00409       for (int isample = 0; isample < digi.size(); ++isample) {
00410         int adc = digi[isample].adc();
00411         int capid = digi[isample].capid ();
00412         double linear_ADC = digi[isample].nominal_fC();
00413         if (verbose) std::cout<<"PIN linear_ADC = "<<linear_ADC<<std::endl;
00414         htspin[isector+iside][ipin]->Fill(isample,linear_ADC);
00415         buf[isample]=linear_ADC;
00416         if (maxADC<linear_ADC) {
00417           maxADC=linear_ADC;
00418           maxisample=isample;
00419         }
00420       }
00421       i1=maxisample-1;
00422       i2=maxisample+2;
00423       if (i1<0) {i1=0;i2=3;}
00424       else if (i2>9) {i1=6;i2=9;} 
00425       if      (i1==0) ped=buf[8]+buf[9]+buf[6]+buf[7];
00426       else if (i1==1) ped=buf[8]+buf[9]+buf[6]+buf[7];
00427       else if (i1==2) ped=buf[0]+buf[1]+buf[6]+buf[7];
00428       else if (i1==3) ped=buf[0]+buf[1]+buf[2]+buf[7];
00429       else if (i1>=4) ped=buf[0]+buf[1]+buf[2]+buf[3];
00430       signal=0;
00431       for (ii=0;ii<4;ii++) signal+=TMath::Max(buf[ii+i1],ped/4.0);
00432       hsppin[isector+iside][ipin]->Fill(signal);
00433       hspepin[isector+iside][ipin]->Fill(signal);
00434       hpedpin[isector+iside][ipin]->Fill(ped);
00435 
00436       // Mean PIN signal time estimation
00437       ped=ped/4;
00438       meant=0;
00439       for (ii=0;ii<4;ii++) meant+=(TMath::Max(buf[ii+i1],ped)-ped)*(ii+i1);
00440       if (signal-ped*4>0) meant/=(signal-ped*4); 
00441       else meant=i1+1;
00442       htsmpin[isector+iside][ipin]->Fill(meant);
00443     }
00444   }
00445   */  
00446   // HF
00447   edm::Handle<HFDigiCollection> hf_digi;
00448   fEvent.getByLabel(hfDigiCollectionTag_, hf_digi);
00449   if (verbose) std::cout<<"Analysis-> total HF digis= "<<hf_digi->size()<<std::endl;
00450 
00451   for (unsigned ihit = 0; ihit < hf_digi->size (); ++ihit) {
00452     const HFDataFrame& frame = (*hf_digi)[ihit];
00453     HcalDetId detId = frame.id();
00454     int ieta = detId.ieta();
00455     int iphi = detId.iphi();
00456     int depth = detId.depth();
00457     if (verbose) std::cout <<"HF digi # " <<ihit<<": eta/phi/depth: "
00458                            <<ieta<<'/'<<iphi<<'/'<< depth << std::endl;
00459 
00460     if (ieta>0) ieta = ieta-29;
00461     else ieta = 13-ieta-29;
00462 
00463     for (int ii=0; ii<10; ii++) buf[ii]=0;
00464     maxADC=-99;
00465     for (int isample = 0; isample < frame.size(); ++isample) {
00466       int adc = frame[isample].adc();
00467       int capid = frame[isample].capid ();
00468       double linear_ADC = frame[isample].nominal_fC();
00469       double nominal_fC = detId.subdet () == HcalForward ? 2.6 *  linear_ADC : linear_ADC;
00470 
00471       if (verbose) std::cout << "Analysis->     HF sample # " << isample 
00472                              << ", capid=" << capid 
00473                              << ": ADC=" << adc 
00474                              << ", linearized ADC=" << linear_ADC
00475                              << ", nominal fC=" << nominal_fC << std::endl;
00476 
00477       hts[ieta][(iphi-1)/2][depth-1]->Fill(isample,linear_ADC);
00478       buf[isample]=linear_ADC;
00479       if (maxADC<linear_ADC) {
00480         maxADC=linear_ADC;
00481         maxisample=isample;
00482       }
00483     }
00484      
00485     // Signal is four capIDs around maxTS, Pedestal is four capID off the signal
00486     maxADC=-99;    
00487     for (int ibf=0; ibf<frame.size()-1; ibf++) {
00488       Double_t sumbuf=0;
00489       for (int jbf=0; jbf<2; jbf++) {    
00490         sumbuf += buf[ibf+jbf];
00491         if (ibf+jbf<2) sumbuf -= (buf[ibf+jbf+4]+buf[ibf+jbf+8])/2.0;
00492         else if (ibf+jbf<4) sumbuf -= buf[ibf+jbf+4];
00493         else if (ibf+jbf<6) sumbuf -= (buf[ibf+jbf-4]+buf[ibf+jbf+4])/2.0;
00494         else if (ibf+jbf<8) sumbuf -= buf[ibf+jbf-4];
00495         else if (ibf+jbf<10) sumbuf -= (buf[ibf+jbf-4]+buf[ibf+jbf-8])/2.0;
00496       }
00497       if (sumbuf>maxADC) {
00498         maxADC=sumbuf;
00499         if (buf[ibf]>buf[ibf+1]) maxisample=ibf;
00500         else maxisample=ibf+1;
00501         //maxisample=ibf;
00502       }
00503     }
00504     htmax->Fill(maxisample);
00505     i1=maxisample-1;
00506     i2=maxisample+2;
00507     if (i1<0) {i1=0;i2=3;}
00508     else if (i2>9) {i1=6;i2=9;} 
00509     else if (i2<9) {
00510       if (buf[i1]<buf[i2+1]) {i1=i1+1;i2=i2+1;}
00511     }
00512     signal=buf[i1]+buf[i1+1]+buf[i1+2]+buf[i1+3];
00513     hsp[ieta][(iphi-1)/2][depth-1]->Fill(signal);
00514     hspe[ieta][(iphi-1)/2][depth-1]->Fill(signal);
00515     
00516     if      (i1==0) ped=(buf[4]+buf[8])/2.0+(buf[5]+buf[9])/2.0+buf[6]+buf[7];
00517     else if (i1==1) ped=(buf[0]+buf[8])/2.0+(buf[5]+buf[9])/2.0+buf[6]+buf[7];
00518     else if (i1==2) ped=(buf[0]+buf[8])/2.0+(buf[1]+buf[9])/2.0+buf[6]+buf[7];
00519     else if (i1==3) ped=(buf[0]+buf[8])/2.0+(buf[1]+buf[9])/2.0+buf[2]+buf[7];
00520     else if (i1==4) ped=(buf[0]+buf[8])/2.0+(buf[1]+buf[9])/2.0+buf[2]+buf[3];
00521     else if (i1==5) ped=(buf[0]+buf[4])/2.0+(buf[1]+buf[9])/2.0+buf[2]+buf[3];
00522     else if (i1==6) ped=(buf[0]+buf[4])/2.0+(buf[1]+buf[5])/2.0+buf[2]+buf[3];
00523     /*
00524     if      (i1==0) ped=buf[8]+buf[9]+buf[6]+buf[7];
00525     else if (i1==1) ped=buf[0]+buf[9]+buf[6]+buf[7];
00526     else if (i1==2) ped=buf[0]+buf[1]+buf[6]+buf[7];
00527     else if (i1==3) ped=buf[0]+buf[1]+buf[2]+buf[7];
00528     else if (i1>=4) ped=buf[0]+buf[1]+buf[2]+buf[3];
00529     */
00530     hped[ieta][(iphi-1)/2][depth-1]->Fill(ped);
00531 
00532     // Mean signal time estimation
00533     ped=ped/4;
00534     meant=(buf[i1]-ped)*i1+(buf[i1+1]-ped)*(i1+1)+(buf[i1+2]-ped)*(i1+2)+(buf[i1+3]-ped)*(i1+3);
00535     meant /= (buf[i1]-ped)+(buf[i1+1]-ped)+(buf[i1+2]-ped)+(buf[i1+3]-ped);
00536     htmean->Fill(meant);
00537     htsm[ieta][(iphi-1)/2][depth-1]->Fill(meant);
00538   }
00539 }
00540