CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/CalibCalorimetry/CastorCalib/src/CastorPedestalAnalysis.cc

Go to the documentation of this file.
00001 #include "CalibFormats/CastorObjects/interface/CastorCoderDb.h"
00002 #include "CalibFormats/CastorObjects/interface/CastorCalibrations.h"
00003 #include "CalibFormats/CastorObjects/interface/CastorDbService.h"
00004 #include "CalibFormats/CastorObjects/interface/CastorDbRecord.h"
00005 #include "CondFormats/CastorObjects/interface/CastorQIECoder.h"
00006 #include "CondFormats/CastorObjects/interface/CastorPedestals.h"
00007 #include "CondFormats/CastorObjects/interface/CastorPedestalWidths.h"
00008 
00009 #include "CalibCalorimetry/CastorCalib/interface/CastorPedestalAnalysis.h"
00010 #include <TFile.h>
00011 #include <math.h>
00012 
00013 using namespace std;
00014 
00015 CastorPedestalAnalysis::CastorPedestalAnalysis(const edm::ParameterSet& ps)
00016   : fRefPedestals (0),
00017     fRefPedestalWidths (0),
00018     fRawPedestals (0),
00019     fRawPedestalWidths (0),
00020     fValPedestals (0),
00021     fValPedestalWidths (0)
00022 {
00023   evt=0;
00024   sample=0;
00025   m_file=0;
00026   m_AllPedsOK=0;
00027   for(int i=0; i<4; i++) m_stat[i]=0;
00028   for(int k=0;k<4;k++) state.push_back(true);
00029 
00030 // user cfg parameters
00031   m_outputFileMean = ps.getUntrackedParameter<string>("outputFileMeans", "");
00032   if ( m_outputFileMean.size() != 0 ) {
00033     cout << "Castor pedestal means will be saved to " << m_outputFileMean.c_str() << endl;
00034   } 
00035   m_outputFileWidth = ps.getUntrackedParameter<string>("outputFileWidths", "");
00036   if ( m_outputFileWidth.size() != 0 ) {
00037     cout << "Castor pedestal widths will be saved to " << m_outputFileWidth.c_str() << endl;
00038   } 
00039   m_outputFileROOT = ps.getUntrackedParameter<string>("outputFileHist", "");
00040   if ( m_outputFileROOT.size() != 0 ) {
00041     cout << "Castor pedestal histograms will be saved to " << m_outputFileROOT.c_str() << endl;
00042   } 
00043   m_nevtsample = ps.getUntrackedParameter<int>("nevtsample",0);
00044 // for compatibility with previous versions
00045   if(m_nevtsample==9999999) m_nevtsample=0;
00046   m_pedsinADC = ps.getUntrackedParameter<int>("pedsinADC",0);
00047   m_hiSaveflag = ps.getUntrackedParameter<int>("hiSaveflag",0);
00048   m_pedValflag = ps.getUntrackedParameter<int>("pedValflag",0);
00049   if(m_pedValflag<0) m_pedValflag=0;
00050   if (m_nevtsample>0 && m_pedValflag>0) {
00051     cout<<"WARNING - incompatible cfg options: nevtsample = "<<m_nevtsample<<", pedValflag = "<<m_pedValflag<<endl;
00052     cout<<"Setting pedValflag = 0"<<endl;
00053     m_pedValflag=0;
00054   }
00055   if(m_pedValflag>1) m_pedValflag=1;
00056   m_startTS = ps.getUntrackedParameter<int>("firstTS", 0);
00057   if(m_startTS<0) m_startTS=0;
00058   m_endTS = ps.getUntrackedParameter<int>("lastTS", 9);
00059 
00060 //  m_logFile.open("CastorPedestalAnalysis.log");
00061 
00062   castorHists.ALLPEDS = new TH1F("Castor All Pedestals","HF All Peds",10,0,9);
00063   castorHists.PEDRMS= new TH1F("Castor All Pedestal Widths","HF All Pedestal RMS",100,0,3);
00064   castorHists.PEDMEAN= new TH1F("Castor All Pedestal Means","HF All Pedestal Means",100,0,9);
00065   castorHists.CHI2= new TH1F("Castor Chi2/ndf for whole range Gauss fit","HF Chi2/ndf Gauss",200,0.,50.);
00066 }
00067 
00068 //-----------------------------------------------------------------------------
00069 CastorPedestalAnalysis::~CastorPedestalAnalysis(){
00070 
00071   for(_meot=castorHists.PEDTRENDS.begin(); _meot!=castorHists.PEDTRENDS.end(); _meot++){
00072     for(int i=0; i<16; i++) _meot->second[i].first->Delete();
00073   }
00074 
00075   castorHists.ALLPEDS->Delete();
00076   castorHists.PEDRMS->Delete();
00077   castorHists.PEDMEAN->Delete();
00078   castorHists.CHI2->Delete();
00079 }
00080 
00081 //-----------------------------------------------------------------------------
00082 void CastorPedestalAnalysis::setup(const std::string& m_outputFileROOT) {
00083   // open the histogram file, create directories within
00084   m_file=new TFile(m_outputFileROOT.c_str(),"RECREATE");
00085   m_file->mkdir("Castor");
00086   m_file->cd();
00087 }
00088 
00089 //-----------------------------------------------------------------------------
00090 void CastorPedestalAnalysis::processEvent(const CastorDigiCollection& castor,
00091                                         const CastorDbService& cond)
00092 {
00093   evt++;
00094   sample=1;
00095   evt_curr=evt;
00096   if(m_nevtsample>0) {
00097     sample = (evt-1)/m_nevtsample +1;
00098     evt_curr = evt%m_nevtsample;
00099     if(evt_curr==0)evt_curr=m_nevtsample;
00100   }
00101 
00102   m_shape = cond.getCastorShape();
00103   // HF
00104   try{
00105     if(!castor.size()) throw (int)castor.size();
00106     for (CastorDigiCollection::const_iterator j=castor.begin(); j!=castor.end(); j++){
00107       const CastorDataFrame digi = (const CastorDataFrame)(*j);
00108       m_coder = cond.getCastorCoder(digi.id());
00109       for (int i=m_startTS; i<digi.size() && i<=m_endTS; i++) {
00110         for(int flag=0; flag<4; flag++){
00111           if(i+flag<digi.size() && i+flag<=m_endTS){
00112             per2CapsHists(flag,2,digi.id(),digi.sample(i),digi.sample(i+flag),castorHists.PEDTRENDS,cond);
00113           }
00114         }
00115       }
00116       if(m_startTS==0 && m_endTS>4){
00117         AllChanHists(digi.id(),digi.sample(0),digi.sample(1),digi.sample(2),digi.sample(3),digi.sample(4),digi.sample(5),castorHists.PEDTRENDS);
00118       }
00119     }
00120   } 
00121   catch (int i ) {
00122 //    m_logFile << "Event with " << i<<" Castor Digis passed." << std::endl;
00123   } 
00124   // Call the function every m_nevtsample events
00125   if(m_nevtsample>0) {
00126     if(evt%m_nevtsample==0) SampleAnalysis();
00127   }
00128 }
00129 
00130 //-----------------------------------------------------------------------------
00131 void CastorPedestalAnalysis::per2CapsHists(int flag, int id, const HcalDetId detid, const HcalQIESample& qie1, const HcalQIESample& qie2, map<HcalDetId, map<int,PEDBUNCH> > &toolT, const CastorDbService& cond) {
00132 
00133 // this function is due to be called for every time slice, it fills either a charge
00134 // histo for a single capID (flag=0) or a product histo for two capIDs (flag>0)
00135 
00136   static const int bins=10;
00137   static const int bins2=100;
00138   float lo=-0.5; float hi=9.5;
00139   map<int,PEDBUNCH> _mei;
00140   static map<HcalDetId, map<int,float> > QieCalibMap;
00141   string type = "Castor";
00142 
00143   /*
00144   if(id==0){
00145     if(detid.ieta()<16) type = "HB";
00146     if(detid.ieta()>16) type = "HE";
00147     if(detid.ieta()==16){
00148       if(detid.depth()<3) type = "HB";
00149       if(detid.depth()==3) type = "HE";
00150     }
00151   } 
00152   else if(id==1) type = "HO";
00153   else if(id==2) type = "HF"; 
00154   */
00155 
00156   _meot = toolT.find(detid);
00157 
00158 // if histos for the current channel do not exist, first create them,
00159   if (_meot==toolT.end()){
00160     map<int,PEDBUNCH> insert;
00161     map<int,float> qiecalib;
00162     char name[1024];
00163     for(int i=0; i<4; i++){
00164       lo=-0.5;
00165       // fix from Andy: if you convert to fC and then bin in units of 1, you may 'skip' a bin while
00166       // filling, since the ADCs are quantized
00167       if (m_pedsinADC) hi=9.5;
00168       else hi = 11.5;
00169       sprintf(name,"%s Pedestal, eta=%d phi=%d d=%d cap=%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth(),i);  
00170       insert[i].first =  new TH1F(name,name,bins,lo,hi);
00171       sprintf(name,"%s Product, eta=%d phi=%d d=%d caps=%d*%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth(),i,(i+1)%4);  
00172       insert[4+i].first = new TH1F(name,name,bins2,0.,100.);
00173       sprintf(name,"%s Product, eta=%d phi=%d d=%d caps=%d*%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth(),i,(i+2)%4);  
00174       insert[8+i].first = new TH1F(name,name,bins2,0.,100.);
00175       sprintf(name,"%s Product, eta=%d phi=%d d=%d caps=%d*%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth(),i,(i+3)%4);  
00176       insert[12+i].first = new TH1F(name,name,bins2,0.,100.);
00177     }
00178     sprintf(name,"%s Signal in TS 4+5, eta=%d phi=%d d=%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth());  
00179     insert[16].first = new TH1F(name,name,21,-0.5,20.5);
00180     sprintf(name,"%s Signal in TS 4+5-2-3, eta=%d phi=%d d=%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth());  
00181     insert[17].first = new TH1F(name,name,21,-10.5,10.5);
00182     sprintf(name,"%s Signal in TS 4+5-(0+1+2+3)/2., eta=%d phi=%d d=%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth());  
00183     insert[18].first = new TH1F(name,name,21,-10.5,10.5);
00184     toolT[detid] = insert;
00185     _meot = toolT.find(detid);
00186 // store QIE calibrations in a map for later reuse
00187     QieCalibMap[detid]=qiecalib;
00188   }
00189 
00190   _mei = _meot->second;
00191 
00192   const CastorQIECoder* coder = cond.getCastorCoder(detid);
00193   const CastorQIEShape* shape = cond.getCastorShape();
00194   float charge1 = coder->charge(*shape,qie1.adc(),qie1.capid());
00195   float charge2 = coder->charge(*shape,qie2.adc(),qie2.capid());
00196 
00197 // fill single capID histo
00198   if(flag==0){
00199     if(m_nevtsample>0) {
00200       if((evt-1)%m_nevtsample==0 && state[qie1.capid()]){
00201         state[qie1.capid()]=false; 
00202         _mei[qie1.capid()].first->Reset();
00203         _mei[qie1.capid()+4].first->Reset();
00204         _mei[qie1.capid()+8].first->Reset();
00205         _mei[qie1.capid()+12].first->Reset();
00206       }
00207     }
00208     if (qie1.adc()<bins){
00209       if (m_pedsinADC) _mei[qie1.capid()].first->Fill(qie1.adc());
00210       else _mei[qie1.capid()].first->Fill(charge1); 
00211     }
00212     else if(qie1.adc()>=bins){
00213       _mei[qie1.capid()].first->AddBinContent(bins+1,1);
00214     }
00215   }
00216 
00217 // fill 2 capID histo
00218   if(flag>0){
00219     map<int,float> qiecalib = QieCalibMap[detid];
00220     //float charge1=(qie1.adc()-qiecalib[qie1.capid()+4])/qiecalib[qie1.capid()];
00221     //float charge2=(qie2.adc()-qiecalib[qie2.capid()+4])/qiecalib[qie2.capid()];
00222     if (charge1*charge2<bins2){
00223       _mei[qie1.capid()+4*flag].first->Fill(charge1*charge2);
00224     }
00225     else{
00226       _mei[qie1.capid()+4*flag].first->Fill(bins2);
00227     }
00228   }
00229 
00230   if(flag==0){
00231     //    if(id==0) hbHists.ALLPEDS->Fill(qie1.adc());
00232     //   else if(id==1) hoHists.ALLPEDS->Fill(qie1.adc());
00233     //   else if(id==2) castorHists.ALLPEDS->Fill(qie1.adc());
00234     castorHists.ALLPEDS->Fill(qie1.adc());
00235   }
00236 }
00237 
00238 //-----------------------------------------------------------------------------
00239 void CastorPedestalAnalysis::AllChanHists(const HcalDetId detid, const HcalQIESample& qie0, const HcalQIESample& qie1, const HcalQIESample& qie2, const HcalQIESample& qie3, const HcalQIESample& qie4, const HcalQIESample& qie5, map<HcalDetId, map<int,PEDBUNCH> > &toolT) { 
00240 
00241 // this function is due to be called for every channel
00242 
00243   _meot = toolT.find(detid);
00244   map<int,PEDBUNCH> _mei = _meot->second;
00245   _mei[16].first->Fill(qie4.adc()+qie5.adc()-1.);
00246   _mei[17].first->Fill(qie4.adc()+qie5.adc()-qie2.adc()-qie3.adc());
00247   _mei[18].first->Fill(qie4.adc()+qie5.adc()-(qie0.adc()+qie1.adc()+qie2.adc()+qie3.adc())/2.);
00248 }
00249 
00250 //-----------------------------------------------------------------------------
00251 void CastorPedestalAnalysis::SampleAnalysis(){
00252   // it is called every m_nevtsample events (a sample) and the end of run
00253   char PedSampleNum[20];
00254 
00255 // Compute pedestal constants for each HBHE, HO, HF
00256   sprintf(PedSampleNum,"Castor_Sample%d",sample);
00257   m_file->cd();
00258   m_file->mkdir(PedSampleNum);
00259   m_file->cd(PedSampleNum);
00260   GetPedConst(castorHists.PEDTRENDS,castorHists.PEDMEAN,castorHists.PEDRMS);
00261 }
00262 
00263 //-----------------------------------------------------------------------------
00264 void CastorPedestalAnalysis::GetPedConst(map<HcalDetId, map<int,PEDBUNCH> > &toolT, TH1F* PedMeans, TH1F* PedWidths)
00265 {
00266 // Completely rewritten version oct 2006
00267 // Compute pedestal constants and fill into CastorPedestals and CastorPedestalWidths objects
00268   float cap[4]; float sig[4][4]; float dcap[4]; float dsig[4][4]; float chi2[4];
00269 
00270   for(_meot=toolT.begin(); _meot!=toolT.end(); _meot++){
00271     HcalDetId detid = _meot->first;
00272 
00273 // take mean and width from a Gaussian fit or directly from the histo
00274     if(fitflag>0){
00275       for (int i=0; i<4; i++) {
00276         TF1 *fit = _meot->second[i].first->GetFunction("gaus");
00277         chi2[i]=0;
00278         if(fit->GetNDF()!=0) chi2[i]=fit->GetChisquare()/fit->GetNDF();
00279         cap[i]=fit->GetParameter(1);
00280         sig[i][i]=fit->GetParameter(2);
00281         dcap[i]=fit->GetParError(1);
00282         dsig[i][i]=fit->GetParError(2);
00283       }
00284     }
00285     else{
00286       for (int i=0; i<4; i++) {
00287         cap[i]=_meot->second[i].first->GetMean();
00288         sig[i][i]=_meot->second[i].first->GetRMS();
00289         m_stat[i]=0;
00290 
00291         for(int j=m_startTS; j<m_endTS+1; j++){
00292           m_stat[i]+=_meot->second[i].first->GetBinContent(j+1);
00293         }
00294         dcap[i] = sig[i][i]/sqrt(m_stat[i]);
00295 //        dsig[i][i] = dcap[i]*sig[i][i]/cap[i];
00296         dsig[i][i] = sig[i][i]/sqrt(2.*m_stat[i]);
00297         chi2[i]=0.;
00298       }
00299     }
00300 
00301     for (int i=0; i<4; i++) {
00302       if(m_hiSaveflag>0) {
00303         if (m_pedsinADC)
00304         _meot->second[i].first->GetXaxis()->SetTitle("ADC");
00305         else _meot->second[i].first->GetXaxis()->SetTitle("Charge, fC");
00306         _meot->second[i].first->GetYaxis()->SetTitle("CapID samplings");
00307         _meot->second[i].first->Write();
00308       }
00309       if(m_nevtsample>0) {
00310         _meot->second[i].second.first[0].push_back(cap[i]);
00311         _meot->second[i].second.first[1].push_back(dcap[i]);
00312         _meot->second[i].second.first[2].push_back(sig[i][i]);
00313         _meot->second[i].second.first[3].push_back(dsig[i][i]);
00314         _meot->second[i].second.first[4].push_back(chi2[i]);
00315       }
00316       PedMeans->Fill(cap[i]);
00317       PedWidths->Fill(sig[i][i]);
00318     }
00319 
00320 // special histos for Shuichi
00321     if(m_hiSaveflag==-100){
00322       for(int i=16; i<19; i++){
00323         if (m_pedsinADC)
00324         _meot->second[i].first->GetXaxis()->SetTitle("ADC");
00325         else _meot->second[i].first->GetXaxis()->SetTitle("Charge, fC");
00326         _meot->second[i].first->GetYaxis()->SetTitle("Events");
00327         _meot->second[i].first->Write();
00328       }
00329     }
00330 
00331 // diagonal sigma is width squared
00332     sig[0][0]=sig[0][0]*sig[0][0];
00333     sig[1][1]=sig[1][1]*sig[1][1];
00334     sig[2][2]=sig[2][2]*sig[2][2];
00335     sig[3][3]=sig[3][3]*sig[3][3];
00336 
00337 // off diagonal sigmas (correlations) are computed from 3 histograms
00338 // here we still have all 4*3=12 combinations
00339     sig[0][1]= _meot->second[4].first->GetMean()-cap[0]*cap[1];
00340     sig[0][2]= _meot->second[8].first->GetMean()-cap[0]*cap[2];
00341     sig[1][2]= _meot->second[5].first->GetMean()-cap[1]*cap[2];
00342     sig[1][3]= _meot->second[9].first->GetMean()-cap[1]*cap[3];
00343     sig[2][3]= _meot->second[6].first->GetMean()-cap[2]*cap[3];
00344     sig[0][3]= _meot->second[12].first->GetMean()-cap[0]*cap[3];
00345     sig[1][0]= _meot->second[13].first->GetMean()-cap[1]*cap[0];
00346     sig[2][0]= _meot->second[10].first->GetMean()-cap[2]*cap[0];
00347     sig[2][1]= _meot->second[14].first->GetMean()-cap[2]*cap[1];
00348     sig[3][1]= _meot->second[11].first->GetMean()-cap[3]*cap[1];
00349     sig[3][2]= _meot->second[15].first->GetMean()-cap[3]*cap[2];
00350     sig[3][0]= _meot->second[7].first->GetMean()-cap[3]*cap[0];
00351 
00352 // there is no proper error calculation for the correlation coefficients
00353     for(int i=0; i<4; i++){
00354       if(m_nevtsample>0) {
00355         _meot->second[i].second.first[5].push_back(sig[i][(i+1)%4]);
00356         _meot->second[i].second.first[6].push_back(2*sig[i][i]*dsig[i][i]);
00357         _meot->second[i].second.first[7].push_back(sig[i][(i+2)%4]);
00358         _meot->second[i].second.first[8].push_back(2*sig[i][i]*dsig[i][i]);
00359         _meot->second[i].second.first[9].push_back(sig[i][(i+3)%4]);
00360         _meot->second[i].second.first[10].push_back(2*sig[i][i]*dsig[i][i]);
00361       }
00362 // save product histos if desired
00363       if(m_hiSaveflag>10) {
00364         if (m_pedsinADC)
00365         _meot->second[i+4].first->GetXaxis()->SetTitle("ADC^2");
00366         else _meot->second[i+4].first->GetXaxis()->SetTitle("Charge^2, fC^2");
00367         _meot->second[i+4].first->GetYaxis()->SetTitle("2-CapID samplings");
00368         _meot->second[i+4].first->Write();
00369         if (m_pedsinADC)
00370         _meot->second[i+8].first->GetXaxis()->SetTitle("ADC^2");
00371         else _meot->second[i+8].first->GetXaxis()->SetTitle("Charge^2, fC^2");
00372         _meot->second[i+8].first->GetYaxis()->SetTitle("2-CapID samplings");
00373         _meot->second[i+8].first->Write();
00374         if (m_pedsinADC)
00375         _meot->second[i+12].first->GetXaxis()->SetTitle("ADC^2");
00376         else _meot->second[i+12].first->GetXaxis()->SetTitle("Charge^2, fC^2");
00377         _meot->second[i+12].first->GetYaxis()->SetTitle("2-CapID samplings");
00378         _meot->second[i+12].first->Write();
00379       }
00380     }
00381 
00382 // fill the objects - at this point only close and medium correlations are stored
00383 // and the matrix is assumed symmetric
00384     if (m_nevtsample<1) {
00385       sig[1][0]=sig[0][1];
00386       sig[2][0]=sig[0][2];
00387       sig[2][1]=sig[1][2];
00388       sig[3][1]=sig[1][3];
00389       sig[3][2]=sig[2][3];
00390       sig[0][3]=sig[3][0];
00391       if (fRawPedestals) {
00392           CastorPedestal item(detid,cap[0],cap[1],cap[2],cap[3]);
00393           fRawPedestals->addValues(item);
00394       }
00395       if (fRawPedestalWidths) {
00396         CastorPedestalWidth widthsp(detid);
00397         widthsp.setSigma(0,0,sig[0][0]);
00398         widthsp.setSigma(0,1,sig[0][1]);
00399         widthsp.setSigma(0,2,sig[0][2]);
00400         widthsp.setSigma(1,1,sig[1][1]);
00401         widthsp.setSigma(1,2,sig[1][2]);
00402         widthsp.setSigma(1,3,sig[1][3]);
00403         widthsp.setSigma(2,2,sig[2][2]);
00404         widthsp.setSigma(2,3,sig[2][3]);
00405         widthsp.setSigma(3,3,sig[3][3]);
00406         widthsp.setSigma(3,0,sig[0][3]);
00407         fRawPedestalWidths->addValues(widthsp);
00408       }
00409     }
00410   }
00411 }
00412 
00413 //-----------------------------------------------------------------------------
00414 int CastorPedestalAnalysis::done(const CastorPedestals* fInputPedestals, 
00415                                 const CastorPedestalWidths* fInputPedestalWidths,
00416                                 CastorPedestals* fOutputPedestals, 
00417                                 CastorPedestalWidths* fOutputPedestalWidths)
00418 {
00419    int nstat[4];
00420 
00421 // Pedestal objects
00422   // inputs...
00423   fRefPedestals = fInputPedestals;
00424   fRefPedestalWidths = fInputPedestalWidths;
00425   
00426   // outputs...
00427   if(m_pedValflag>0) {
00428     fValPedestals = fOutputPedestals;
00429     fValPedestalWidths = fOutputPedestalWidths;
00430     fRawPedestals = new CastorPedestals();
00431     fRawPedestalWidths = new CastorPedestalWidths();
00432   }
00433   else {
00434     fRawPedestals = fOutputPedestals;
00435     fRawPedestalWidths = fOutputPedestalWidths;
00436     fValPedestals = new CastorPedestals();
00437     fValPedestalWidths = new CastorPedestalWidths();
00438   }
00439 
00440 // compute pedestal constants
00441   if(m_nevtsample<1) SampleAnalysis();
00442   if(m_nevtsample>0) {
00443     if(evt%m_nevtsample!=0) SampleAnalysis();
00444   }
00445 
00446 // trending histos
00447   if(m_nevtsample>0){
00448     m_file->cd();
00449     m_file->cd("Castor");
00450     Trendings(castorHists.PEDTRENDS,castorHists.CHI2,castorHists.CAPID_AVERAGE,castorHists.CAPID_CHI2);
00451   }
00452 
00453   if (m_nevtsample<1) {
00454 
00455 // pedestal validation: m_AllPedsOK=-1 means not validated,
00456 //                                   0 everything OK,
00457 //                                   N>0 : mod(N,100000) drifts + width changes
00458 //                                         int(N/100000) missing channels
00459     m_AllPedsOK=-1;
00460     if(m_pedValflag>0) {
00461       for (int i=0; i<4; i++) nstat[i]=(int)m_stat[i];
00462       int NPedErrors=CastorPedVal(nstat,fRefPedestals,fRefPedestalWidths,
00463                             fRawPedestals,fRawPedestalWidths,
00464                             fValPedestals,fValPedestalWidths);
00465       m_AllPedsOK=NPedErrors;
00466     }
00467 // setting m_AllPedsOK=-2 will inhibit writing pedestals out
00468 //    if(m_pedValflag==1){
00469 //      if(evt<100)m_AllPedsOK=-2;
00470 //    }
00471 
00472   }
00473 
00474   // Write other histograms.
00475 
00476   // Castor
00477   m_file->cd();
00478   m_file->cd("Castor");
00479   castorHists.ALLPEDS->Write();
00480   castorHists.PEDRMS->Write();
00481   castorHists.PEDMEAN->Write();
00482 
00483   m_file->Close();
00484   cout << "Hcal/Castor histograms written to " << m_outputFileROOT.c_str() << endl;
00485   return (int)m_AllPedsOK;
00486 }
00487 
00488 //-----------------------------------------------------------------------------
00489 void CastorPedestalAnalysis::Trendings(map<HcalDetId, map<int,PEDBUNCH> > &toolT, TH1F* Chi2, TH1F* CapidAverage, TH1F* CapidChi2){
00490 
00491 // check stability of pedestal constants in a single long run
00492 
00493   map<int, std::vector<double> > AverageValues;
00494 
00495   for(_meot=toolT.begin(); _meot!=toolT.end(); _meot++){
00496     for(int i=0; i<4; i++){
00497       char name[1024];
00498       HcalDetId detid = _meot->first;
00499       sprintf(name,"Pedestal trend, eta=%d phi=%d d=%d cap=%d",detid.ieta(),detid.iphi(),detid.depth(),i);
00500       int bins = _meot->second[i].second.first[0].size();
00501       float lo =0.5;
00502       float hi = (float)bins+0.5;
00503       _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00504       sprintf(name,"Width trend, eta=%d phi=%d d=%d cap=%d",detid.ieta(),detid.iphi(),detid.depth(),i);
00505       bins = _meot->second[i].second.first[2].size();
00506       hi = (float)bins+0.5;
00507       _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00508       sprintf(name,"Correlation trend, eta=%d phi=%d d=%d caps=%d*%d",detid.ieta(),detid.iphi(),detid.depth(),i,(i+1)%4);
00509       bins = _meot->second[i].second.first[5].size();
00510       hi = (float)bins+0.5;
00511       _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00512 /*      sprintf(name,"Correlation trend, eta=%d phi=%d d=%d caps=%d*%d",detid.ieta(),detid.iphi(),detid.depth(),i,(i+2)%4);
00513       bins = _meot->second[i].second.first[7].size();
00514       hi = (float)bins+0.5;
00515       _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00516       sprintf(name,"Correlation trend, eta=%d phi=%d d=%d caps=%d*%d",detid.ieta(),detid.iphi(),detid.depth(),i,(i+3)%4);
00517       bins = _meot->second[i].second.first[9].size();
00518       hi = (float)bins+0.5;
00519       _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi)); */
00520 
00521       std::vector<double>::iterator sample_it;
00522       // Pedestal mean - put content and errors
00523       int j=0;
00524       for(sample_it=_meot->second[i].second.first[0].begin();
00525           sample_it!=_meot->second[i].second.first[0].end();sample_it++){
00526         _meot->second[i].second.second[0]->SetBinContent(++j,*sample_it);
00527       }
00528       j=0;
00529       for(sample_it=_meot->second[i].second.first[1].begin();
00530           sample_it!=_meot->second[i].second.first[1].end();sample_it++){
00531         _meot->second[i].second.second[0]->SetBinError(++j,*sample_it);
00532       }
00533       // fit with a constant - extract parameters
00534       _meot->second[i].second.second[0]->Fit("pol0","Q");
00535       TF1 *fit = _meot->second[i].second.second[0]->GetFunction("pol0");
00536       AverageValues[0].push_back(fit->GetParameter(0));
00537       AverageValues[1].push_back(fit->GetParError(0));
00538       if(sample>1)
00539       AverageValues[2].push_back(fit->GetChisquare()/fit->GetNDF());
00540       else
00541       AverageValues[2].push_back(fit->GetChisquare());
00542       sprintf(name,"Sample (%d events)",m_nevtsample);
00543       _meot->second[i].second.second[0]->GetXaxis()->SetTitle(name);
00544       _meot->second[i].second.second[0]->GetYaxis()->SetTitle("Pedestal value");
00545       _meot->second[i].second.second[0]->Write();
00546       // Pedestal width - put content and errors
00547       j=0;
00548       for(sample_it=_meot->second[i].second.first[2].begin();
00549           sample_it!=_meot->second[i].second.first[2].end();sample_it++){
00550         _meot->second[i].second.second[1]->SetBinContent(++j,*sample_it);
00551       }
00552       j=0;
00553       for(sample_it=_meot->second[i].second.first[3].begin();
00554           sample_it!=_meot->second[i].second.first[3].end();sample_it++){
00555         _meot->second[i].second.second[1]->SetBinError(++j,*sample_it);
00556       }
00557       _meot->second[i].second.second[1]->GetXaxis()->SetTitle(name);
00558       _meot->second[i].second.second[1]->GetYaxis()->SetTitle("Pedestal width");
00559       _meot->second[i].second.second[1]->Write();
00560       // Correlation coeffs - put contents and errors
00561       j=0;
00562       for(sample_it=_meot->second[i].second.first[5].begin();
00563           sample_it!=_meot->second[i].second.first[5].end();sample_it++){
00564         _meot->second[i].second.second[2]->SetBinContent(++j,*sample_it);
00565       }
00566       j=0;
00567       for(sample_it=_meot->second[i].second.first[6].begin();
00568           sample_it!=_meot->second[i].second.first[6].end();sample_it++){
00569         _meot->second[i].second.second[2]->SetBinError(++j,*sample_it);
00570       }
00571       _meot->second[i].second.second[2]->GetXaxis()->SetTitle(name);
00572       _meot->second[i].second.second[2]->GetYaxis()->SetTitle("Close correlation");
00573       _meot->second[i].second.second[2]->Write();
00574  /*     j=0;
00575       for(sample_it=_meot->second[i].second.first[7].begin();
00576           sample_it!=_meot->second[i].second.first[7].end();sample_it++){
00577         _meot->second[i].second.second[3]->SetBinContent(++j,*sample_it);
00578       }
00579       j=0;
00580       for(sample_it=_meot->second[i].second.first[8].begin();
00581           sample_it!=_meot->second[i].second.first[8].end();sample_it++){
00582         _meot->second[i].second.second[3]->SetBinError(++j,*sample_it);
00583       }
00584       _meot->second[i].second.second[3]->GetXaxis()->SetTitle(name);
00585       _meot->second[i].second.second[3]->GetYaxis()->SetTitle("Intermediate correlation");
00586       _meot->second[i].second.second[3]->Write();
00587       j=0;
00588       for(sample_it=_meot->second[i].second.first[9].begin();
00589           sample_it!=_meot->second[i].second.first[9].end();sample_it++){
00590         _meot->second[i].second.second[4]->SetBinContent(++j,*sample_it);
00591       }
00592       j=0;
00593       for(sample_it=_meot->second[i].second.first[10].begin();
00594           sample_it!=_meot->second[i].second.first[10].end();sample_it++){
00595         _meot->second[i].second.second[4]->SetBinError(++j,*sample_it);
00596       }
00597       _meot->second[i].second.second[4]->GetXaxis()->SetTitle(name);
00598       _meot->second[i].second.second[4]->GetYaxis()->SetTitle("Distant correlation");
00599       _meot->second[i].second.second[4]->Write(); */
00600       // chi2
00601       j=0;
00602       for(sample_it=_meot->second[i].second.first[4].begin();
00603           sample_it!=_meot->second[i].second.first[4].end();sample_it++){
00604         Chi2->Fill(*sample_it);
00605       }
00606     }
00607   }
00608   CapidAverage= new TH1F("Constant fit: Pedestal Values",
00609                          "Constant fit: Pedestal Values",
00610                          AverageValues[0].size(),0.,AverageValues[0].size());
00611   std::vector<double>::iterator sample_it;
00612   int j=0;
00613   for(sample_it=AverageValues[0].begin();
00614       sample_it!=AverageValues[0].end();sample_it++){
00615     CapidAverage->SetBinContent(++j,*sample_it);
00616   }
00617   j=0;
00618   for(sample_it=AverageValues[1].begin();
00619       sample_it!=AverageValues[1].end();sample_it++){
00620     CapidAverage->SetBinError(++j,*sample_it);
00621   }
00622   CapidChi2= new TH1F("Constant fit: Chi2/ndf",
00623                       "Constant fit: Chi2/ndf",
00624                       AverageValues[2].size(),0.,AverageValues[2].size());
00625   j=0;
00626   for(sample_it=AverageValues[2].begin();
00627       sample_it!=AverageValues[2].end();sample_it++){
00628     CapidChi2->SetBinContent(++j,*sample_it);
00629     //CapidChi2->SetBinError(++j,0);
00630   }
00631   Chi2->GetXaxis()->SetTitle("Chi2/ndf");
00632   Chi2->GetYaxis()->SetTitle("50 x [(16+2) x 4 x 4] `events`");
00633   Chi2->Write();
00634   CapidAverage->GetYaxis()->SetTitle("Pedestal value");
00635   CapidAverage->GetXaxis()->SetTitle("(16+2) x 4 x 4 `events`");
00636   CapidAverage->Write();
00637   CapidChi2->GetYaxis()->SetTitle("Chi2/ndf");
00638   CapidChi2->GetXaxis()->SetTitle("(16+2) x 4 x 4 `events`");
00639   CapidChi2->Write();
00640 
00641 }
00642 
00643 //-----------------------------------------------------------------------------
00644 int CastorPedestalAnalysis::CastorPedVal(int nstat[4], const CastorPedestals* fRefPedestals, 
00645                                 const CastorPedestalWidths* fRefPedestalWidths,
00646                                 CastorPedestals* fRawPedestals,
00647                                 CastorPedestalWidths* fRawPedestalWidths,
00648                                 CastorPedestals* fValPedestals, 
00649                                 CastorPedestalWidths* fValPedestalWidths)
00650 {
00651 // new version of pedestal validation - it is designed to be as independent of
00652 // all the rest as possible - you only need to provide valid pedestal objects
00653 // and a vector of statistics per capID to use this as standalone code
00654   HcalDetId detid;
00655   float RefPedVals[4]; float RefPedSigs[4][4];
00656   float RawPedVals[4]; float RawPedSigs[4][4];
00657   map<HcalDetId,bool> isinRaw;
00658   map<HcalDetId,bool> isinRef;
00659   std::vector<DetId> RefChanns=fRefPedestals->getAllChannels();
00660   std::vector<DetId> RawChanns=fRawPedestals->getAllChannels();
00661   std::ofstream PedValLog;
00662   PedValLog.open("CastorPedVal.log");
00663 
00664   if(nstat[0]+nstat[1]+nstat[2]+nstat[3]<2500) PedValLog<<"CastorPedVal: warning - low statistics"<<std::endl;
00665 // find complete list of channels in current data and reference
00666   for (int i=0; i<(int)RawChanns.size(); i++){
00667     isinRef[HcalDetId(RawChanns[i])]=false;
00668   }
00669   for (int i=0; i<(int)RefChanns.size(); i++){
00670     detid=HcalDetId(RefChanns[i]);
00671     isinRaw[detid]=false;
00672     isinRef[detid]=true;
00673   }
00674   for (int i=0; i<(int)RawChanns.size(); i++){
00675     detid=HcalDetId(RawChanns[i]);
00676     isinRaw[detid]=true;
00677     if (isinRef[detid]==false) {
00678       PedValLog<<"CastorPedVal: channel "<<detid<<" not found in reference set"<<std::endl;
00679       std::cerr<<"CastorPedVal: channel "<<detid<<" not found in reference set"<<std::endl;
00680     }
00681   }
00682 
00683 // main loop over channels
00684   int erflag=0;
00685   for (int i=0; i<(int)RefChanns.size(); i++){
00686     detid=HcalDetId(RefChanns[i]);
00687     for (int icap=0; icap<4; icap++) {
00688       RefPedVals[icap]=fRefPedestals->getValues(detid)->getValue(icap);
00689       for (int icap2=icap; icap2<4; icap2++) {
00690         RefPedSigs[icap][icap2]=fRefPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00691         if(icap2!=icap)RefPedSigs[icap2][icap]=RefPedSigs[icap][icap2];
00692       }
00693     }
00694 
00695 // read new raw values
00696     if(isinRaw[detid]) {
00697       for (int icap=0; icap<4; icap++) {
00698         RawPedVals[icap]=fRawPedestals->getValues(detid)->getValue(icap);
00699         for (int icap2=icap; icap2<4; icap2++) {
00700           RawPedSigs[icap][icap2]=fRawPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00701           if(icap2!=icap)RawPedSigs[icap2][icap]=RawPedSigs[icap][icap2];
00702         }
00703       }
00704 
00705 // first quick check if raw values make sense: if not, the channel is treated like absent
00706       for (int icap=0; icap<4; icap++) {
00707         if(RawPedVals[icap]<1. || RawPedSigs[icap][icap]<0.01) isinRaw[detid]=false;
00708         for (int icap2=icap; icap2<4; icap2++){
00709           if(fabs(RawPedSigs[icap][icap2]/sqrt(RawPedSigs[icap][icap]*RawPedSigs[icap2][icap2]))>1.) isinRaw[detid]=false;
00710         }
00711       }
00712     }
00713 
00714 // check raw values against reference
00715     if(isinRaw[detid]) {
00716       for (int icap=0; icap<4; icap++) {
00717         int icap2=(icap+1)%4;
00718         float width=sqrt(RawPedSigs[icap][icap]);
00719         float erof1=width/sqrt((float)nstat[icap]);
00720         float erof2=sqrt(erof1*erof1+RawPedSigs[icap][icap]/(float)nstat[icap]);
00721         float erofwidth=width/sqrt(2.*nstat[icap]);
00722         float diffof1=RawPedVals[icap]-RefPedVals[icap];
00723         float diffof2=RawPedVals[icap]+RawPedVals[icap2]-RefPedVals[icap]-RefPedVals[icap2];
00724         float diffofw=width-sqrt(RefPedSigs[icap][icap]);
00725 
00726 // validation in 2 TS for HB, HE, HO, in 1 TS for HF
00727         int nTS=2;
00728         if(detid.subdet()==HcalForward) nTS=1;
00729         if(nTS==1 && fabs(diffof1)>0.5+erof1) { 
00730           erflag+=1;
00731           PedValLog<<"HcalPedVal: drift in channel "<<detid<<" cap "<<icap<<": "<<RawPedVals[icap]<<" - "<<RefPedVals[icap]<<" = "<<diffof1<<std::endl;
00732         }
00733         if(nTS==2 && fabs(diffof2)>0.5+erof2) { 
00734           erflag+=1;
00735           PedValLog<<"HcalPedVal: drift in channel "<<detid<<" caps "<<icap<<"+"<<icap2<<": "<<RawPedVals[icap]<<"+"<<RawPedVals[icap2]<<" - "<<RefPedVals[icap]<<"+"<<RefPedVals[icap2]<<" = "<<diffof2<<std::endl;
00736         }
00737         if(fabs(diffofw)>0.15*width+erofwidth) {
00738           erflag+=1;
00739           PedValLog<<"HcalPedVal: width changed in channel "<<detid<<" cap "<<icap<<": "<<width<<" - "<<sqrt(RefPedSigs[icap][icap])<<" = "<<diffofw<<std::endl;
00740         }
00741       }
00742     }
00743 
00744 // for disconnected/bad channels restore reference values
00745     else {
00746       PedValLog<<"HcalPedVal: no valid data from channel "<<detid<<std::endl;
00747       erflag+=100000;
00748       CastorPedestal item(detid,RefPedVals[0],RefPedVals[1],RefPedVals[2],RefPedVals[3]);
00749       fValPedestals->addValues(item);
00750       CastorPedestalWidth widthsp(detid);
00751       for (int icap=0; icap<4; icap++) {
00752         for (int icap2=icap; icap2<4; icap2++) widthsp.setSigma(icap2,icap,RefPedSigs[icap2][icap]);
00753       }
00754       fValPedestalWidths->addValues(widthsp);
00755     }
00756 
00757 // end of channel loop
00758   }
00759 
00760   if(erflag==0) PedValLog<<"HcalPedVal: all pedestals checked OK"<<std::endl;
00761 
00762 // now construct the remaining part of the validated objects
00763 // if nothing changed outside tolerance, validated set = reference set
00764   if(erflag%100000 == 0) {
00765     for (int i=0; i<(int)RefChanns.size(); i++){
00766       detid=HcalDetId(RefChanns[i]);
00767       if (isinRaw[detid]) {
00768         CastorPedestalWidth widthsp(detid);
00769         for (int icap=0; icap<4; icap++) {
00770           RefPedVals[icap]=fRefPedestals->getValues(detid)->getValue(icap);
00771           for (int icap2=icap; icap2<4; icap2++) {
00772             RefPedSigs[icap][icap2]=fRefPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00773             if(icap2!=icap)RefPedSigs[icap2][icap]=RefPedSigs[icap][icap2];
00774             widthsp.setSigma(icap2,icap,RefPedSigs[icap2][icap]);
00775           }
00776         }
00777         fValPedestalWidths->addValues(widthsp);
00778         CastorPedestal item(detid,RefPedVals[0],RefPedVals[1],RefPedVals[2],RefPedVals[3]);
00779         fValPedestals->addValues(item);
00780       }
00781     }
00782   }
00783 
00784 // if anything changed, validated set = raw set + reference for missing/bad channels
00785   else {
00786     for (int i=0; i<(int)RawChanns.size(); i++){
00787       detid=HcalDetId(RawChanns[i]);
00788       if (isinRaw[detid]) {
00789         CastorPedestalWidth widthsp(detid);
00790         for (int icap=0; icap<4; icap++) {
00791           RawPedVals[icap]=fRawPedestals->getValues(detid)->getValue(icap);
00792           for (int icap2=icap; icap2<4; icap2++) {
00793             RawPedSigs[icap][icap2]=fRawPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00794             if(icap2!=icap)RawPedSigs[icap2][icap]=RawPedSigs[icap][icap2];
00795             widthsp.setSigma(icap2,icap,RawPedSigs[icap2][icap]);
00796           }
00797         }
00798         fValPedestalWidths->addValues(widthsp);
00799         CastorPedestal item(detid,RawPedVals[0],RawPedVals[1],RawPedVals[2],RawPedVals[3]);
00800         fValPedestals->addValues(item);
00801       }
00802     }
00803   }
00804   return erflag;
00805 }
00806