00001 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00002 #include "CondFormats/HcalObjects/interface/HcalQIECoder.h"
00003 #include "CondFormats/HcalObjects/interface/HcalPedestals.h"
00004 #include "CondFormats/HcalObjects/interface/HcalPedestalWidths.h"
00005
00006 #include "CalibCalorimetry/HcalAlgos/interface/HcalPedestalAnalysis.h"
00007 #include "TFile.h"
00008 #include <math.h>
00009
00010
00011
00012
00013
00014 using namespace std;
00015
00016 HcalPedestalAnalysis::HcalPedestalAnalysis(const edm::ParameterSet& ps)
00017 : fRefPedestals (0),
00018 fRefPedestalWidths (0),
00019 fRawPedestals (0),
00020 fRawPedestalWidths (0),
00021 fValPedestals (0),
00022 fValPedestalWidths (0)
00023 {
00024 evt=0;
00025 sample=0;
00026 m_file=0;
00027 m_AllPedsOK=0;
00028 for(int i=0; i<4; i++) m_stat[i]=0;
00029 for(int k=0;k<4;k++) state.push_back(true);
00030
00031
00032 m_outputFileMean = ps.getUntrackedParameter<string>("outputFileMeans", "");
00033 if ( m_outputFileMean.size() != 0 ) {
00034 cout << "Hcal pedestal means will be saved to " << m_outputFileMean.c_str() << endl;
00035 }
00036 m_outputFileWidth = ps.getUntrackedParameter<string>("outputFileWidths", "");
00037 if ( m_outputFileWidth.size() != 0 ) {
00038 cout << "Hcal pedestal widths will be saved to " << m_outputFileWidth.c_str() << endl;
00039 }
00040 m_outputFileROOT = ps.getUntrackedParameter<string>("outputFileHist", "");
00041 if ( m_outputFileROOT.size() != 0 ) {
00042 cout << "Hcal pedestal histograms will be saved to " << m_outputFileROOT.c_str() << endl;
00043 }
00044 m_nevtsample = ps.getUntrackedParameter<int>("nevtsample",0);
00045
00046 if(m_nevtsample==9999999) m_nevtsample=0;
00047 m_pedsinADC = ps.getUntrackedParameter<int>("pedsinADC",0);
00048 m_hiSaveflag = ps.getUntrackedParameter<int>("hiSaveflag",0);
00049 m_pedValflag = ps.getUntrackedParameter<int>("pedValflag",0);
00050 if(m_pedValflag<0) m_pedValflag=0;
00051 if (m_nevtsample>0 && m_pedValflag>0) {
00052 cout<<"WARNING - incompatible cfg options: nevtsample = "<<m_nevtsample<<", pedValflag = "<<m_pedValflag<<endl;
00053 cout<<"Setting pedValflag = 0"<<endl;
00054 m_pedValflag=0;
00055 }
00056 if(m_pedValflag>1) m_pedValflag=1;
00057 m_startTS = ps.getUntrackedParameter<int>("firstTS", 0);
00058 if(m_startTS<0) m_startTS=0;
00059 m_endTS = ps.getUntrackedParameter<int>("lastTS", 9);
00060
00061
00062
00063 hbHists.ALLPEDS = new TH1F("HBHE All Pedestals","HBHE All Peds",10,0,9);
00064 hbHists.PEDRMS= new TH1F("HBHE All Pedestal Widths","HBHE All Pedestal RMS",100,0,3);
00065 hbHists.PEDMEAN= new TH1F("HBHE All Pedestal Means","HBHE All Pedestal Means",100,0,9);
00066 hbHists.CHI2= new TH1F("HBHE Chi2/ndf for whole range Gauss fit","HBHE Chi2/ndf Gauss",200,0.,50.);
00067
00068 hoHists.ALLPEDS = new TH1F("HO All Pedestals","HO All Peds",10,0,9);
00069 hoHists.PEDRMS= new TH1F("HO All Pedestal Widths","HO All Pedestal RMS",100,0,3);
00070 hoHists.PEDMEAN= new TH1F("HO All Pedestal Means","HO All Pedestal Means",100,0,9);
00071 hoHists.CHI2= new TH1F("HO Chi2/ndf for whole range Gauss fit","HO Chi2/ndf Gauss",200,0.,50.);
00072
00073 hfHists.ALLPEDS = new TH1F("HF All Pedestals","HF All Peds",10,0,9);
00074 hfHists.PEDRMS= new TH1F("HF All Pedestal Widths","HF All Pedestal RMS",100,0,3);
00075 hfHists.PEDMEAN= new TH1F("HF All Pedestal Means","HF All Pedestal Means",100,0,9);
00076 hfHists.CHI2= new TH1F("HF Chi2/ndf for whole range Gauss fit","HF Chi2/ndf Gauss",200,0.,50.);
00077 }
00078
00079
00080 HcalPedestalAnalysis::~HcalPedestalAnalysis(){
00081
00082 for(_meot=hbHists.PEDTRENDS.begin(); _meot!=hbHists.PEDTRENDS.end(); _meot++){
00083 for(int i=0; i<16; i++) _meot->second[i].first->Delete();
00084 }
00085 for(_meot=hoHists.PEDTRENDS.begin(); _meot!=hoHists.PEDTRENDS.end(); _meot++){
00086 for(int i=0; i<16; i++) _meot->second[i].first->Delete();
00087 }
00088 for(_meot=hfHists.PEDTRENDS.begin(); _meot!=hfHists.PEDTRENDS.end(); _meot++){
00089 for(int i=0; i<16; i++) _meot->second[i].first->Delete();
00090 }
00091 hbHists.ALLPEDS->Delete();
00092 hbHists.PEDRMS->Delete();
00093 hbHists.PEDMEAN->Delete();
00094 hbHists.CHI2->Delete();
00095
00096 hoHists.ALLPEDS->Delete();
00097 hoHists.PEDRMS->Delete();
00098 hoHists.PEDMEAN->Delete();
00099 hoHists.CHI2->Delete();
00100
00101 hfHists.ALLPEDS->Delete();
00102 hfHists.PEDRMS->Delete();
00103 hfHists.PEDMEAN->Delete();
00104 hfHists.CHI2->Delete();
00105 }
00106
00107
00108 void HcalPedestalAnalysis::setup(const std::string& m_outputFileROOT) {
00109
00110 m_file=new TFile(m_outputFileROOT.c_str(),"RECREATE");
00111 m_file->mkdir("HB");
00112 m_file->cd();
00113 m_file->mkdir("HO");
00114 m_file->cd();
00115 m_file->mkdir("HF");
00116 m_file->cd();
00117 }
00118
00119
00120 void HcalPedestalAnalysis::processEvent(const HBHEDigiCollection& hbhe,
00121 const HODigiCollection& ho,
00122 const HFDigiCollection& hf,
00123 const HcalDbService& cond)
00124 {
00125 evt++;
00126 sample=1;
00127 evt_curr=evt;
00128 if(m_nevtsample>0) {
00129 sample = (evt-1)/m_nevtsample +1;
00130 evt_curr = evt%m_nevtsample;
00131 if(evt_curr==0)evt_curr=m_nevtsample;
00132 }
00133
00134 m_shape = cond.getHcalShape();
00135
00136
00137 try{
00138 if(!hbhe.size()) throw (int)hbhe.size();
00139 for (HBHEDigiCollection::const_iterator j=hbhe.begin(); j!=hbhe.end(); j++){
00140 const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00141 m_coder = cond.getHcalCoder(digi.id());
00142 for(int k=0; k<(int)state.size();k++) state[k]=true;
00143
00144
00145 for (int i=m_startTS; i<digi.size() && i<=m_endTS; i++) {
00146 for(int flag=0; flag<4; flag++){
00147 if(i+flag<digi.size() && i+flag<=m_endTS){
00148 per2CapsHists(flag,0,digi.id(),digi.sample(i),digi.sample(i+flag),hbHists.PEDTRENDS,cond);
00149 }
00150 }
00151 }
00152 if(m_startTS==0 && m_endTS>4){
00153 AllChanHists(digi.id(),digi.sample(0),digi.sample(1),digi.sample(2),digi.sample(3),digi.sample(4),digi.sample(5),hbHists.PEDTRENDS);
00154 }
00155 }
00156 }
00157 catch (int i ) {
00158
00159 }
00160
00161 try{
00162 if(!ho.size()) throw (int)ho.size();
00163 for (HODigiCollection::const_iterator j=ho.begin(); j!=ho.end(); j++){
00164 const HODataFrame digi = (const HODataFrame)(*j);
00165 m_coder = cond.getHcalCoder(digi.id());
00166 for (int i=m_startTS; i<digi.size() && i<=m_endTS; i++) {
00167 for(int flag=0; flag<4; flag++){
00168 if(i+flag<digi.size() && i+flag<=m_endTS){
00169 per2CapsHists(flag,1,digi.id(),digi.sample(i),digi.sample(i+flag),hoHists.PEDTRENDS,cond);
00170 }
00171 }
00172 }
00173 if(m_startTS==0 && m_endTS>4){
00174 AllChanHists(digi.id(),digi.sample(0),digi.sample(1),digi.sample(2),digi.sample(3),digi.sample(4),digi.sample(5),hoHists.PEDTRENDS);
00175 }
00176 }
00177 }
00178 catch (int i ) {
00179
00180 }
00181
00182 try{
00183 if(!hf.size()) throw (int)hf.size();
00184 for (HFDigiCollection::const_iterator j=hf.begin(); j!=hf.end(); j++){
00185 const HFDataFrame digi = (const HFDataFrame)(*j);
00186 m_coder = cond.getHcalCoder(digi.id());
00187 for (int i=m_startTS; i<digi.size() && i<=m_endTS; i++) {
00188 for(int flag=0; flag<4; flag++){
00189 if(i+flag<digi.size() && i+flag<=m_endTS){
00190 per2CapsHists(flag,2,digi.id(),digi.sample(i),digi.sample(i+flag),hfHists.PEDTRENDS,cond);
00191 }
00192 }
00193 }
00194 if(m_startTS==0 && m_endTS>4){
00195 AllChanHists(digi.id(),digi.sample(0),digi.sample(1),digi.sample(2),digi.sample(3),digi.sample(4),digi.sample(5),hfHists.PEDTRENDS);
00196 }
00197 }
00198 }
00199 catch (int i ) {
00200
00201 }
00202
00203 if(m_nevtsample>0) {
00204 if(evt%m_nevtsample==0) SampleAnalysis();
00205 }
00206 }
00207
00208
00209 void HcalPedestalAnalysis::per2CapsHists(int flag, int id, const HcalDetId detid, const HcalQIESample& qie1, const HcalQIESample& qie2, map<HcalDetId, map<int,PEDBUNCH> > &toolT, const HcalDbService& cond) {
00210
00211
00212
00213
00214 static const int bins=10;
00215 static const int bins2=100;
00216 float lo=-0.5; float hi=9.5;
00217 map<int,PEDBUNCH> _mei;
00218 static map<HcalDetId, map<int,float> > QieCalibMap;
00219 string type = "HBHE";
00220
00221 if(id==0){
00222 if(detid.ieta()<16) type = "HB";
00223 if(detid.ieta()>16) type = "HE";
00224 if(detid.ieta()==16){
00225 if(detid.depth()<3) type = "HB";
00226 if(detid.depth()==3) type = "HE";
00227 }
00228 }
00229 else if(id==1) type = "HO";
00230 else if(id==2) type = "HF";
00231
00232 _meot = toolT.find(detid);
00233
00234
00235 if (_meot==toolT.end()){
00236 map<int,PEDBUNCH> insert;
00237 map<int,float> qiecalib;
00238 char name[1024];
00239 for(int i=0; i<4; i++){
00240 lo=-0.5;
00241
00242
00243 if (m_pedsinADC) hi=9.5;
00244 else hi = 11.5;
00245 sprintf(name,"%s Pedestal, eta=%d phi=%d d=%d cap=%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth(),i);
00246 insert[i].first = new TH1F(name,name,bins,lo,hi);
00247 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);
00248 insert[4+i].first = new TH1F(name,name,bins2,0.,100.);
00249 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);
00250 insert[8+i].first = new TH1F(name,name,bins2,0.,100.);
00251 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);
00252 insert[12+i].first = new TH1F(name,name,bins2,0.,100.);
00253 }
00254 sprintf(name,"%s Signal in TS 4+5, eta=%d phi=%d d=%d",type.c_str(),detid.ieta(),detid.iphi(),detid.depth());
00255 insert[16].first = new TH1F(name,name,21,-0.5,20.5);
00256 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());
00257 insert[17].first = new TH1F(name,name,21,-10.5,10.5);
00258 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());
00259 insert[18].first = new TH1F(name,name,21,-10.5,10.5);
00260 toolT[detid] = insert;
00261 _meot = toolT.find(detid);
00262
00263 QieCalibMap[detid]=qiecalib;
00264 }
00265
00266 _mei = _meot->second;
00267
00268 const HcalQIECoder* coder = cond.getHcalCoder(detid);
00269 const HcalQIEShape* shape = cond.getHcalShape();
00270 float charge1 = coder->charge(*shape,qie1.adc(),qie1.capid());
00271 float charge2 = coder->charge(*shape,qie2.adc(),qie2.capid());
00272
00273
00274 if(flag==0){
00275 if(m_nevtsample>0) {
00276 if((evt-1)%m_nevtsample==0 && state[qie1.capid()]){
00277 state[qie1.capid()]=false;
00278 _mei[qie1.capid()].first->Reset();
00279 _mei[qie1.capid()+4].first->Reset();
00280 _mei[qie1.capid()+8].first->Reset();
00281 _mei[qie1.capid()+12].first->Reset();
00282 }
00283 }
00284 if (qie1.adc()<bins){
00285 if (m_pedsinADC) _mei[qie1.capid()].first->Fill(qie1.adc());
00286 else _mei[qie1.capid()].first->Fill(charge1);
00287 }
00288 else if(qie1.adc()>=bins){
00289 _mei[qie1.capid()].first->AddBinContent(bins+1,1);
00290 }
00291 }
00292
00293
00294 if(flag>0){
00295 map<int,float> qiecalib = QieCalibMap[detid];
00296
00297
00298 if (charge1*charge2<bins2){
00299 _mei[qie1.capid()+4*flag].first->Fill(charge1*charge2);
00300 }
00301 else{
00302 _mei[qie1.capid()+4*flag].first->Fill(bins2);
00303 }
00304 }
00305
00306 if(flag==0){
00307 if(id==0) hbHists.ALLPEDS->Fill(qie1.adc());
00308 else if(id==1) hoHists.ALLPEDS->Fill(qie1.adc());
00309 else if(id==2) hfHists.ALLPEDS->Fill(qie1.adc());
00310 }
00311 }
00312
00313
00314 void HcalPedestalAnalysis::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) {
00315
00316
00317
00318 _meot = toolT.find(detid);
00319 map<int,PEDBUNCH> _mei = _meot->second;
00320 _mei[16].first->Fill(qie4.adc()+qie5.adc()-1.);
00321 _mei[17].first->Fill(qie4.adc()+qie5.adc()-qie2.adc()-qie3.adc());
00322 _mei[18].first->Fill(qie4.adc()+qie5.adc()-(qie0.adc()+qie1.adc()+qie2.adc()+qie3.adc())/2.);
00323 }
00324
00325
00326 void HcalPedestalAnalysis::SampleAnalysis(){
00327
00328 char PedSampleNum[20];
00329
00330
00331 sprintf(PedSampleNum,"HB_Sample%d",sample);
00332 m_file->cd();
00333 m_file->mkdir(PedSampleNum);
00334 m_file->cd(PedSampleNum);
00335 GetPedConst(hbHists.PEDTRENDS,hbHists.PEDMEAN,hbHists.PEDRMS);
00336 sprintf(PedSampleNum,"HO_Sample%d",sample);
00337 m_file->cd();
00338 m_file->mkdir(PedSampleNum);
00339 m_file->cd(PedSampleNum);
00340 GetPedConst(hoHists.PEDTRENDS,hoHists.PEDMEAN,hoHists.PEDRMS);
00341 sprintf(PedSampleNum,"HF_Sample%d",sample);
00342 m_file->cd();
00343 m_file->mkdir(PedSampleNum);
00344 m_file->cd(PedSampleNum);
00345 GetPedConst(hfHists.PEDTRENDS,hfHists.PEDMEAN,hfHists.PEDRMS);
00346 }
00347
00348
00349 void HcalPedestalAnalysis::GetPedConst(map<HcalDetId, map<int,PEDBUNCH> > &toolT, TH1F* PedMeans, TH1F* PedWidths)
00350 {
00351
00352
00353 float cap[4]; float sig[4][4]; float dcap[4]; float dsig[4][4]; float chi2[4];
00354
00355 for(_meot=toolT.begin(); _meot!=toolT.end(); _meot++){
00356 HcalDetId detid = _meot->first;
00357
00358
00359 if(fitflag>0){
00360 for (int i=0; i<4; i++) {
00361 TF1 *fit = _meot->second[i].first->GetFunction("gaus");
00362 chi2[i]=0;
00363 if(fit->GetNDF()!=0) chi2[i]=fit->GetChisquare()/fit->GetNDF();
00364 cap[i]=fit->GetParameter(1);
00365 sig[i][i]=fit->GetParameter(2);
00366 dcap[i]=fit->GetParError(1);
00367 dsig[i][i]=fit->GetParError(2);
00368 }
00369 }
00370 else{
00371 for (int i=0; i<4; i++) {
00372 cap[i]=_meot->second[i].first->GetMean();
00373 sig[i][i]=_meot->second[i].first->GetRMS();
00374 m_stat[i]=0;
00375
00376 for(int j=m_startTS; j<m_endTS+1; j++){
00377 m_stat[i]+=_meot->second[i].first->GetBinContent(j+1);
00378 }
00379 dcap[i] = sig[i][i]/sqrt(m_stat[i]);
00380
00381 dsig[i][i] = sig[i][i]/sqrt(2.*m_stat[i]);
00382 chi2[i]=0.;
00383 }
00384 }
00385
00386 for (int i=0; i<4; i++) {
00387 if(m_hiSaveflag>0) {
00388 if (m_pedsinADC)
00389 _meot->second[i].first->GetXaxis()->SetTitle("ADC");
00390 else _meot->second[i].first->GetXaxis()->SetTitle("Charge, fC");
00391 _meot->second[i].first->GetYaxis()->SetTitle("CapID samplings");
00392 _meot->second[i].first->Write();
00393 }
00394 if(m_nevtsample>0) {
00395 _meot->second[i].second.first[0].push_back(cap[i]);
00396 _meot->second[i].second.first[1].push_back(dcap[i]);
00397 _meot->second[i].second.first[2].push_back(sig[i][i]);
00398 _meot->second[i].second.first[3].push_back(dsig[i][i]);
00399 _meot->second[i].second.first[4].push_back(chi2[i]);
00400 }
00401 PedMeans->Fill(cap[i]);
00402 PedWidths->Fill(sig[i][i]);
00403 }
00404
00405
00406 if(m_hiSaveflag==-100){
00407 for(int i=16; i<19; i++){
00408 if (m_pedsinADC)
00409 _meot->second[i].first->GetXaxis()->SetTitle("ADC");
00410 else _meot->second[i].first->GetXaxis()->SetTitle("Charge, fC");
00411 _meot->second[i].first->GetYaxis()->SetTitle("Events");
00412 _meot->second[i].first->Write();
00413 }
00414 }
00415
00416
00417 sig[0][0]=sig[0][0]*sig[0][0];
00418 sig[1][1]=sig[1][1]*sig[1][1];
00419 sig[2][2]=sig[2][2]*sig[2][2];
00420 sig[3][3]=sig[3][3]*sig[3][3];
00421
00422
00423
00424 sig[0][1]= _meot->second[4].first->GetMean()-cap[0]*cap[1];
00425 sig[0][2]= _meot->second[8].first->GetMean()-cap[0]*cap[2];
00426 sig[1][2]= _meot->second[5].first->GetMean()-cap[1]*cap[2];
00427 sig[1][3]= _meot->second[9].first->GetMean()-cap[1]*cap[3];
00428 sig[2][3]= _meot->second[6].first->GetMean()-cap[2]*cap[3];
00429 sig[0][3]= _meot->second[12].first->GetMean()-cap[0]*cap[3];
00430 sig[1][0]= _meot->second[13].first->GetMean()-cap[1]*cap[0];
00431 sig[2][0]= _meot->second[10].first->GetMean()-cap[2]*cap[0];
00432 sig[2][1]= _meot->second[14].first->GetMean()-cap[2]*cap[1];
00433 sig[3][1]= _meot->second[11].first->GetMean()-cap[3]*cap[1];
00434 sig[3][2]= _meot->second[15].first->GetMean()-cap[3]*cap[2];
00435 sig[3][0]= _meot->second[7].first->GetMean()-cap[3]*cap[0];
00436
00437
00438 for(int i=0; i<4; i++){
00439 if(m_nevtsample>0) {
00440 _meot->second[i].second.first[5].push_back(sig[i][(i+1)%4]);
00441 _meot->second[i].second.first[6].push_back(2*sig[i][i]*dsig[i][i]);
00442 _meot->second[i].second.first[7].push_back(sig[i][(i+2)%4]);
00443 _meot->second[i].second.first[8].push_back(2*sig[i][i]*dsig[i][i]);
00444 _meot->second[i].second.first[9].push_back(sig[i][(i+3)%4]);
00445 _meot->second[i].second.first[10].push_back(2*sig[i][i]*dsig[i][i]);
00446 }
00447
00448 if(m_hiSaveflag>10) {
00449 if (m_pedsinADC)
00450 _meot->second[i+4].first->GetXaxis()->SetTitle("ADC^2");
00451 else _meot->second[i+4].first->GetXaxis()->SetTitle("Charge^2, fC^2");
00452 _meot->second[i+4].first->GetYaxis()->SetTitle("2-CapID samplings");
00453 _meot->second[i+4].first->Write();
00454 if (m_pedsinADC)
00455 _meot->second[i+8].first->GetXaxis()->SetTitle("ADC^2");
00456 else _meot->second[i+8].first->GetXaxis()->SetTitle("Charge^2, fC^2");
00457 _meot->second[i+8].first->GetYaxis()->SetTitle("2-CapID samplings");
00458 _meot->second[i+8].first->Write();
00459 if (m_pedsinADC)
00460 _meot->second[i+12].first->GetXaxis()->SetTitle("ADC^2");
00461 else _meot->second[i+12].first->GetXaxis()->SetTitle("Charge^2, fC^2");
00462 _meot->second[i+12].first->GetYaxis()->SetTitle("2-CapID samplings");
00463 _meot->second[i+12].first->Write();
00464 }
00465 }
00466
00467
00468
00469 if (m_nevtsample<1) {
00470 sig[1][0]=sig[0][1];
00471 sig[2][0]=sig[0][2];
00472 sig[2][1]=sig[1][2];
00473 sig[3][1]=sig[1][3];
00474 sig[3][2]=sig[2][3];
00475 sig[0][3]=sig[3][0];
00476 if (fRawPedestals)
00477 {
00478 HcalPedestal item(detid,cap[0],cap[1],cap[2],cap[3]);
00479 fRawPedestals->addValues(item);
00480 }
00481 if (fRawPedestalWidths)
00482 {
00483 HcalPedestalWidth widthsp(detid);
00484 widthsp.setSigma(0,0,sig[0][0]);
00485 widthsp.setSigma(0,1,sig[0][1]);
00486 widthsp.setSigma(0,2,sig[0][2]);
00487 widthsp.setSigma(1,1,sig[1][1]);
00488 widthsp.setSigma(1,2,sig[1][2]);
00489 widthsp.setSigma(1,3,sig[1][3]);
00490 widthsp.setSigma(2,2,sig[2][2]);
00491 widthsp.setSigma(2,3,sig[2][3]);
00492 widthsp.setSigma(3,3,sig[3][3]);
00493 widthsp.setSigma(3,0,sig[0][3]);
00494 fRawPedestalWidths->addValues(widthsp);
00495 }
00496 }
00497 }
00498 }
00499
00500
00501 int HcalPedestalAnalysis::done(const HcalPedestals* fInputPedestals,
00502 const HcalPedestalWidths* fInputPedestalWidths,
00503 HcalPedestals* fOutputPedestals,
00504 HcalPedestalWidths* fOutputPedestalWidths)
00505 {
00506 int nstat[4];
00507
00508
00509
00510 fRefPedestals = fInputPedestals;
00511 fRefPedestalWidths = fInputPedestalWidths;
00512
00513
00514 if(m_pedValflag>0) {
00515 fValPedestals = fOutputPedestals;
00516 fValPedestalWidths = fOutputPedestalWidths;
00517 fRawPedestals = new HcalPedestals(m_pedsinADC);
00518 fRawPedestalWidths = new HcalPedestalWidths(m_pedsinADC);
00519 }
00520 else {
00521 fRawPedestals = fOutputPedestals;
00522 fRawPedestalWidths = fOutputPedestalWidths;
00523 fValPedestals = new HcalPedestals(m_pedsinADC);
00524 fValPedestalWidths = new HcalPedestalWidths(m_pedsinADC);
00525 }
00526
00527
00528 if(m_nevtsample<1) SampleAnalysis();
00529 if(m_nevtsample>0) {
00530 if(evt%m_nevtsample!=0) SampleAnalysis();
00531 }
00532
00533
00534 if(m_nevtsample>0){
00535 m_file->cd();
00536 m_file->cd("HB");
00537 Trendings(hbHists.PEDTRENDS,hbHists.CHI2,hbHists.CAPID_AVERAGE,hbHists.CAPID_CHI2);
00538 m_file->cd();
00539 m_file->cd("HO");
00540 Trendings(hoHists.PEDTRENDS,hoHists.CHI2,hoHists.CAPID_AVERAGE,hoHists.CAPID_CHI2);
00541 m_file->cd();
00542 m_file->cd("HF");
00543 Trendings(hfHists.PEDTRENDS,hfHists.CHI2,hfHists.CAPID_AVERAGE,hfHists.CAPID_CHI2);
00544 }
00545
00546 if (m_nevtsample<1) {
00547
00548
00549
00550
00551
00552 m_AllPedsOK=-1;
00553 if(m_pedValflag>0) {
00554 for (int i=0; i<4; i++) nstat[i]=(int)m_stat[i];
00555 int NPedErrors=HcalPedVal(nstat,fRefPedestals,fRefPedestalWidths,
00556 fRawPedestals,fRawPedestalWidths,
00557 fValPedestals,fValPedestalWidths);
00558 m_AllPedsOK=NPedErrors;
00559 }
00560
00561
00562
00563
00564
00565 }
00566
00567
00568
00569 m_file->cd();
00570 m_file->cd("HB");
00571 hbHists.ALLPEDS->Write();
00572 hbHists.PEDRMS->Write();
00573 hbHists.PEDMEAN->Write();
00574
00575 m_file->cd();
00576 m_file->cd("HO");
00577 hoHists.ALLPEDS->Write();
00578 hoHists.PEDRMS->Write();
00579 hoHists.PEDMEAN->Write();
00580
00581 m_file->cd();
00582 m_file->cd("HF");
00583 hfHists.ALLPEDS->Write();
00584 hfHists.PEDRMS->Write();
00585 hfHists.PEDMEAN->Write();
00586
00587 m_file->Close();
00588 cout << "Hcal histograms written to " << m_outputFileROOT.c_str() << endl;
00589 return (int)m_AllPedsOK;
00590 }
00591
00592
00593 void HcalPedestalAnalysis::Trendings(map<HcalDetId, map<int,PEDBUNCH> > &toolT, TH1F* Chi2, TH1F* CapidAverage, TH1F* CapidChi2){
00594
00595
00596
00597 map<int, std::vector<double> > AverageValues;
00598
00599 for(_meot=toolT.begin(); _meot!=toolT.end(); _meot++){
00600 for(int i=0; i<4; i++){
00601 char name[1024];
00602 HcalDetId detid = _meot->first;
00603 sprintf(name,"Pedestal trend, eta=%d phi=%d d=%d cap=%d",detid.ieta(),detid.iphi(),detid.depth(),i);
00604 int bins = _meot->second[i].second.first[0].size();
00605 float lo =0.5;
00606 float hi = (float)bins+0.5;
00607 _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00608 sprintf(name,"Width trend, eta=%d phi=%d d=%d cap=%d",detid.ieta(),detid.iphi(),detid.depth(),i);
00609 bins = _meot->second[i].second.first[2].size();
00610 hi = (float)bins+0.5;
00611 _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00612 sprintf(name,"Correlation trend, eta=%d phi=%d d=%d caps=%d*%d",detid.ieta(),detid.iphi(),detid.depth(),i,(i+1)%4);
00613 bins = _meot->second[i].second.first[5].size();
00614 hi = (float)bins+0.5;
00615 _meot->second[i].second.second.push_back(new TH1F(name,name,bins,lo,hi));
00616
00617
00618
00619
00620
00621
00622
00623
00624
00625 std::vector<double>::iterator sample_it;
00626
00627 int j=0;
00628 for(sample_it=_meot->second[i].second.first[0].begin();
00629 sample_it!=_meot->second[i].second.first[0].end();sample_it++){
00630 _meot->second[i].second.second[0]->SetBinContent(++j,*sample_it);
00631 }
00632 j=0;
00633 for(sample_it=_meot->second[i].second.first[1].begin();
00634 sample_it!=_meot->second[i].second.first[1].end();sample_it++){
00635 _meot->second[i].second.second[0]->SetBinError(++j,*sample_it);
00636 }
00637
00638 _meot->second[i].second.second[0]->Fit("pol0","Q");
00639 TF1 *fit = _meot->second[i].second.second[0]->GetFunction("pol0");
00640 AverageValues[0].push_back(fit->GetParameter(0));
00641 AverageValues[1].push_back(fit->GetParError(0));
00642 if(sample>1)
00643 AverageValues[2].push_back(fit->GetChisquare()/fit->GetNDF());
00644 else
00645 AverageValues[2].push_back(fit->GetChisquare());
00646 sprintf(name,"Sample (%d events)",m_nevtsample);
00647 _meot->second[i].second.second[0]->GetXaxis()->SetTitle(name);
00648 _meot->second[i].second.second[0]->GetYaxis()->SetTitle("Pedestal value");
00649 _meot->second[i].second.second[0]->Write();
00650
00651 j=0;
00652 for(sample_it=_meot->second[i].second.first[2].begin();
00653 sample_it!=_meot->second[i].second.first[2].end();sample_it++){
00654 _meot->second[i].second.second[1]->SetBinContent(++j,*sample_it);
00655 }
00656 j=0;
00657 for(sample_it=_meot->second[i].second.first[3].begin();
00658 sample_it!=_meot->second[i].second.first[3].end();sample_it++){
00659 _meot->second[i].second.second[1]->SetBinError(++j,*sample_it);
00660 }
00661 _meot->second[i].second.second[1]->GetXaxis()->SetTitle(name);
00662 _meot->second[i].second.second[1]->GetYaxis()->SetTitle("Pedestal width");
00663 _meot->second[i].second.second[1]->Write();
00664
00665 j=0;
00666 for(sample_it=_meot->second[i].second.first[5].begin();
00667 sample_it!=_meot->second[i].second.first[5].end();sample_it++){
00668 _meot->second[i].second.second[2]->SetBinContent(++j,*sample_it);
00669 }
00670 j=0;
00671 for(sample_it=_meot->second[i].second.first[6].begin();
00672 sample_it!=_meot->second[i].second.first[6].end();sample_it++){
00673 _meot->second[i].second.second[2]->SetBinError(++j,*sample_it);
00674 }
00675 _meot->second[i].second.second[2]->GetXaxis()->SetTitle(name);
00676 _meot->second[i].second.second[2]->GetYaxis()->SetTitle("Close correlation");
00677 _meot->second[i].second.second[2]->Write();
00678
00679
00680
00681
00682
00683
00684
00685
00686
00687
00688
00689
00690
00691
00692
00693
00694
00695
00696
00697
00698
00699
00700
00701
00702
00703
00704
00705 j=0;
00706 for(sample_it=_meot->second[i].second.first[4].begin();
00707 sample_it!=_meot->second[i].second.first[4].end();sample_it++){
00708 Chi2->Fill(*sample_it);
00709 }
00710 }
00711 }
00712 CapidAverage= new TH1F("Constant fit: Pedestal Values",
00713 "Constant fit: Pedestal Values",
00714 AverageValues[0].size(),0.,AverageValues[0].size());
00715 std::vector<double>::iterator sample_it;
00716 int j=0;
00717 for(sample_it=AverageValues[0].begin();
00718 sample_it!=AverageValues[0].end();sample_it++){
00719 CapidAverage->SetBinContent(++j,*sample_it);
00720 }
00721 j=0;
00722 for(sample_it=AverageValues[1].begin();
00723 sample_it!=AverageValues[1].end();sample_it++){
00724 CapidAverage->SetBinError(++j,*sample_it);
00725 }
00726 CapidChi2= new TH1F("Constant fit: Chi2/ndf",
00727 "Constant fit: Chi2/ndf",
00728 AverageValues[2].size(),0.,AverageValues[2].size());
00729 j=0;
00730 for(sample_it=AverageValues[2].begin();
00731 sample_it!=AverageValues[2].end();sample_it++){
00732 CapidChi2->SetBinContent(++j,*sample_it);
00733
00734 }
00735 Chi2->GetXaxis()->SetTitle("Chi2/ndf");
00736 Chi2->GetYaxis()->SetTitle("50 x [(16+2) x 4 x 4] `events`");
00737 Chi2->Write();
00738 CapidAverage->GetYaxis()->SetTitle("Pedestal value");
00739 CapidAverage->GetXaxis()->SetTitle("(16+2) x 4 x 4 `events`");
00740 CapidAverage->Write();
00741 CapidChi2->GetYaxis()->SetTitle("Chi2/ndf");
00742 CapidChi2->GetXaxis()->SetTitle("(16+2) x 4 x 4 `events`");
00743 CapidChi2->Write();
00744
00745 }
00746
00747
00748 int HcalPedestalAnalysis::HcalPedVal(int nstat[4], const HcalPedestals* fRefPedestals,
00749 const HcalPedestalWidths* fRefPedestalWidths,
00750 HcalPedestals* fRawPedestals,
00751 HcalPedestalWidths* fRawPedestalWidths,
00752 HcalPedestals* fValPedestals,
00753 HcalPedestalWidths* fValPedestalWidths)
00754 {
00755
00756
00757
00758 HcalDetId detid;
00759 float RefPedVals[4]; float RefPedSigs[4][4];
00760 float RawPedVals[4]; float RawPedSigs[4][4];
00761 map<HcalDetId,bool> isinRaw;
00762 map<HcalDetId,bool> isinRef;
00763 std::vector<DetId> RefChanns=fRefPedestals->getAllChannels();
00764 std::vector<DetId> RawChanns=fRawPedestals->getAllChannels();
00765 std::ofstream PedValLog;
00766 PedValLog.open("HcalPedVal.log");
00767
00768 if(nstat[0]+nstat[1]+nstat[2]+nstat[3]<2500) PedValLog<<"HcalPedVal: warning - low statistics"<<std::endl;
00769
00770 for (int i=0; i<(int)RawChanns.size(); i++){
00771 isinRef[HcalDetId(RawChanns[i])]=false;
00772 }
00773 for (int i=0; i<(int)RefChanns.size(); i++){
00774 detid=HcalDetId(RefChanns[i]);
00775 isinRaw[detid]=false;
00776 isinRef[detid]=true;
00777 }
00778 for (int i=0; i<(int)RawChanns.size(); i++){
00779 detid=HcalDetId(RawChanns[i]);
00780 isinRaw[detid]=true;
00781 if (isinRef[detid]==false) {
00782 PedValLog<<"HcalPedVal: channel "<<detid<<" not found in reference set"<<std::endl;
00783 std::cerr<<"HcalPedVal: channel "<<detid<<" not found in reference set"<<std::endl;
00784 }
00785 }
00786
00787
00788 int erflag=0;
00789 for (int i=0; i<(int)RefChanns.size(); i++){
00790 detid=HcalDetId(RefChanns[i]);
00791 for (int icap=0; icap<4; icap++) {
00792 RefPedVals[icap]=fRefPedestals->getValues(detid)->getValue(icap);
00793 for (int icap2=icap; icap2<4; icap2++) {
00794 RefPedSigs[icap][icap2]=fRefPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00795 if(icap2!=icap)RefPedSigs[icap2][icap]=RefPedSigs[icap][icap2];
00796 }
00797 }
00798
00799
00800 if(isinRaw[detid]) {
00801 for (int icap=0; icap<4; icap++) {
00802 RawPedVals[icap]=fRawPedestals->getValues(detid)->getValue(icap);
00803 for (int icap2=icap; icap2<4; icap2++) {
00804 RawPedSigs[icap][icap2]=fRawPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00805 if(icap2!=icap)RawPedSigs[icap2][icap]=RawPedSigs[icap][icap2];
00806 }
00807 }
00808
00809
00810 for (int icap=0; icap<4; icap++) {
00811 if(RawPedVals[icap]<1. || RawPedSigs[icap][icap]<0.01) isinRaw[detid]=false;
00812 for (int icap2=icap; icap2<4; icap2++){
00813 if(fabs(RawPedSigs[icap][icap2]/sqrt(RawPedSigs[icap][icap]*RawPedSigs[icap2][icap2]))>1.) isinRaw[detid]=false;
00814 }
00815 }
00816 }
00817
00818
00819 if(isinRaw[detid]) {
00820 for (int icap=0; icap<4; icap++) {
00821 int icap2=(icap+1)%4;
00822 float width=sqrt(RawPedSigs[icap][icap]);
00823 float erof1=width/sqrt((float)nstat[icap]);
00824 float erof2=sqrt(erof1*erof1+RawPedSigs[icap][icap]/(float)nstat[icap]);
00825 float erofwidth=width/sqrt(2.*nstat[icap]);
00826 float diffof1=RawPedVals[icap]-RefPedVals[icap];
00827 float diffof2=RawPedVals[icap]+RawPedVals[icap2]-RefPedVals[icap]-RefPedVals[icap2];
00828 float diffofw=width-sqrt(RefPedSigs[icap][icap]);
00829
00830
00831 int nTS=2;
00832 if(detid.subdet()==HcalForward) nTS=1;
00833 if(nTS==1 && fabs(diffof1)>0.5+erof1) {
00834 erflag+=1;
00835 PedValLog<<"HcalPedVal: drift in channel "<<detid<<" cap "<<icap<<": "<<RawPedVals[icap]<<" - "<<RefPedVals[icap]<<" = "<<diffof1<<std::endl;
00836 }
00837 if(nTS==2 && fabs(diffof2)>0.5+erof2) {
00838 erflag+=1;
00839 PedValLog<<"HcalPedVal: drift in channel "<<detid<<" caps "<<icap<<"+"<<icap2<<": "<<RawPedVals[icap]<<"+"<<RawPedVals[icap2]<<" - "<<RefPedVals[icap]<<"+"<<RefPedVals[icap2]<<" = "<<diffof2<<std::endl;
00840 }
00841 if(fabs(diffofw)>0.15*width+erofwidth) {
00842 erflag+=1;
00843 PedValLog<<"HcalPedVal: width changed in channel "<<detid<<" cap "<<icap<<": "<<width<<" - "<<sqrt(RefPedSigs[icap][icap])<<" = "<<diffofw<<std::endl;
00844 }
00845 }
00846 }
00847
00848
00849 else {
00850 PedValLog<<"HcalPedVal: no valid data from channel "<<detid<<std::endl;
00851 erflag+=100000;
00852 HcalPedestal item(detid,RefPedVals[0],RefPedVals[1],RefPedVals[2],RefPedVals[3]);
00853 fValPedestals->addValues(item);
00854 HcalPedestalWidth widthsp(detid);
00855 for (int icap=0; icap<4; icap++) {
00856 for (int icap2=icap; icap2<4; icap2++) widthsp.setSigma(icap2,icap,RefPedSigs[icap2][icap]);
00857 }
00858 fValPedestalWidths->addValues(widthsp);
00859 }
00860
00861
00862 }
00863
00864 if(erflag==0) PedValLog<<"HcalPedVal: all pedestals checked OK"<<std::endl;
00865
00866
00867
00868 if(erflag%100000 == 0) {
00869 for (int i=0; i<(int)RefChanns.size(); i++){
00870 detid=HcalDetId(RefChanns[i]);
00871 if (isinRaw[detid]) {
00872 HcalPedestalWidth widthsp(detid);
00873 for (int icap=0; icap<4; icap++) {
00874 RefPedVals[icap]=fRefPedestals->getValues(detid)->getValue(icap);
00875 for (int icap2=icap; icap2<4; icap2++) {
00876 RefPedSigs[icap][icap2]=fRefPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00877 if(icap2!=icap)RefPedSigs[icap2][icap]=RefPedSigs[icap][icap2];
00878 widthsp.setSigma(icap2,icap,RefPedSigs[icap2][icap]);
00879 }
00880 }
00881 fValPedestalWidths->addValues(widthsp);
00882 HcalPedestal item(detid,RefPedVals[0],RefPedVals[1],RefPedVals[2],RefPedVals[3]);
00883 fValPedestals->addValues(item);
00884 }
00885 }
00886 }
00887
00888
00889 else {
00890 for (int i=0; i<(int)RawChanns.size(); i++){
00891 detid=HcalDetId(RawChanns[i]);
00892 if (isinRaw[detid]) {
00893 HcalPedestalWidth widthsp(detid);
00894 for (int icap=0; icap<4; icap++) {
00895 RawPedVals[icap]=fRawPedestals->getValues(detid)->getValue(icap);
00896 for (int icap2=icap; icap2<4; icap2++) {
00897 RawPedSigs[icap][icap2]=fRawPedestalWidths->getValues(detid)->getSigma(icap,icap2);
00898 if(icap2!=icap)RawPedSigs[icap2][icap]=RawPedSigs[icap][icap2];
00899 widthsp.setSigma(icap2,icap,RawPedSigs[icap2][icap]);
00900 }
00901 }
00902 fValPedestalWidths->addValues(widthsp);
00903 HcalPedestal item(detid,RawPedVals[0],RawPedVals[1],RawPedVals[2],RawPedVals[3]);
00904 fValPedestals->addValues(item);
00905 }
00906 }
00907 }
00908 return erflag;
00909 }
00910