00001
00002
00003
00004
00005
00006
00007 #include <memory>
00008 #include "CalibCalorimetry/HcalStandardModules/interface/HcalPedestalWidthsValidation.h"
00009
00010 HcalPedestalWidthsValidation::HcalPedestalWidthsValidation(const edm::ParameterSet& ps)
00011 {
00012 firstTS = ps.getUntrackedParameter<int>("firstTS", 0);
00013 lastTS = ps.getUntrackedParameter<int>("lastTS", 9);
00014 firsttime = true;
00015
00016 rawPedsItem = new HcalPedestals();
00017 rawWidthsItem = new HcalPedestalWidths();
00018 rawPedsItemfc = new HcalPedestals();
00019 rawWidthsItemfc = new HcalPedestalWidths();
00020 }
00021
00022
00023 HcalPedestalWidthsValidation::~HcalPedestalWidthsValidation()
00024 {
00025
00026 std::cout << "Calculating Pedestal constants...\n";
00027 std::vector<NewPedBunchVal>::iterator bunch_it;
00028 for(bunch_it=BunchVales.begin(); bunch_it != BunchVales.end(); bunch_it++)
00029 {
00030 if(bunch_it->usedflag){
00031
00032 bunch_it->cap[0] /= bunch_it->num[0][0];
00033 bunch_it->cap[1] /= bunch_it->num[1][1];
00034 bunch_it->cap[2] /= bunch_it->num[2][2];
00035 bunch_it->cap[3] /= bunch_it->num[3][3];
00036 bunch_it->capfc[0] /= bunch_it->num[0][0];
00037 bunch_it->capfc[1] /= bunch_it->num[1][1];
00038 bunch_it->capfc[2] /= bunch_it->num[2][2];
00039 bunch_it->capfc[3] /= bunch_it->num[3][3];
00040
00041 bunch_it->sig[0][0] = (bunch_it->prod[0][0]/bunch_it->num[0][0])-(bunch_it->cap[0])*(bunch_it->cap[0]);
00042 bunch_it->sig[1][1] = (bunch_it->prod[1][1]/bunch_it->num[1][1])-(bunch_it->cap[1])*(bunch_it->cap[1]);
00043 bunch_it->sig[2][2] = (bunch_it->prod[2][2]/bunch_it->num[2][2])-(bunch_it->cap[2])*(bunch_it->cap[2]);
00044 bunch_it->sig[3][3] = (bunch_it->prod[3][3]/bunch_it->num[3][3])-(bunch_it->cap[3])*(bunch_it->cap[3]);
00045 bunch_it->sig[0][1] = (bunch_it->prod[0][1])/(bunch_it->num[0][1])-(bunch_it->cap[0]*bunch_it->cap[1]);
00046 bunch_it->sig[0][2] = (bunch_it->prod[0][2])/(bunch_it->num[0][2])-(bunch_it->cap[0]*bunch_it->cap[2]);
00047 bunch_it->sig[0][3] = (bunch_it->prod[0][3])/(bunch_it->num[0][3])-(bunch_it->cap[0]*bunch_it->cap[3]);
00048 bunch_it->sig[1][2] = (bunch_it->prod[1][2])/(bunch_it->num[1][2])-(bunch_it->cap[1]*bunch_it->cap[2]);
00049 bunch_it->sig[1][3] = (bunch_it->prod[1][3])/(bunch_it->num[1][3])-(bunch_it->cap[1]*bunch_it->cap[3]);
00050 bunch_it->sig[2][3] = (bunch_it->prod[2][3])/(bunch_it->num[2][3])-(bunch_it->cap[2]*bunch_it->cap[3]);
00051
00052 bunch_it->sig[1][0] = (bunch_it->prod[1][0]/bunch_it->num[1][0])-(bunch_it->cap[1])*(bunch_it->cap[0]);
00053 bunch_it->sig[2][0] = (bunch_it->prod[2][0]/bunch_it->num[2][0])-(bunch_it->cap[2])*(bunch_it->cap[0]);
00054 bunch_it->sig[2][1] = (bunch_it->prod[2][1]/bunch_it->num[2][1])-(bunch_it->cap[2])*(bunch_it->cap[1]);
00055 bunch_it->sig[3][0] = (bunch_it->prod[3][0]/bunch_it->num[3][0])-(bunch_it->cap[3])*(bunch_it->cap[0]);
00056 bunch_it->sig[3][1] = (bunch_it->prod[3][1]/bunch_it->num[3][1])-(bunch_it->cap[3])*(bunch_it->cap[1]);
00057 bunch_it->sig[3][2] = (bunch_it->prod[3][2]/bunch_it->num[3][2])-(bunch_it->cap[3])*(bunch_it->cap[2]);
00058
00059
00060
00061 bunch_it->sigfc[0][0] = (bunch_it->prodfc[0][0]/bunch_it->num[0][0])-(bunch_it->capfc[0])*(bunch_it->capfc[0]);
00062 bunch_it->sigfc[1][1] = (bunch_it->prodfc[1][1]/bunch_it->num[1][1])-(bunch_it->capfc[1])*(bunch_it->capfc[1]);
00063 bunch_it->sigfc[2][2] = (bunch_it->prodfc[2][2]/bunch_it->num[2][2])-(bunch_it->capfc[2])*(bunch_it->capfc[2]);
00064 bunch_it->sigfc[3][3] = (bunch_it->prodfc[3][3]/bunch_it->num[3][3])-(bunch_it->capfc[3])*(bunch_it->capfc[3]);
00065 bunch_it->sigfc[0][1] = (bunch_it->prodfc[0][1]/(bunch_it->num[0][1]))-(bunch_it->capfc[0]*bunch_it->capfc[1]);
00066 bunch_it->sigfc[0][2] = (bunch_it->prodfc[0][2]/(bunch_it->num[0][2]))-(bunch_it->capfc[0]*bunch_it->capfc[2]);
00067 bunch_it->sigfc[0][3] = (bunch_it->prodfc[0][3]/(bunch_it->num[0][3]))-(bunch_it->capfc[0]*bunch_it->capfc[3]);
00068 bunch_it->sigfc[1][2] = (bunch_it->prodfc[1][2]/(bunch_it->num[1][2]))-(bunch_it->capfc[1]*bunch_it->capfc[2]);
00069 bunch_it->sigfc[1][3] = (bunch_it->prodfc[1][3]/(bunch_it->num[1][3]))-(bunch_it->capfc[1]*bunch_it->capfc[3]);
00070 bunch_it->sigfc[2][3] = (bunch_it->prodfc[2][3]/(bunch_it->num[2][3]))-(bunch_it->capfc[2]*bunch_it->capfc[3]);
00071
00072 bunch_it->sigfc[1][0] = (bunch_it->prodfc[1][0]/(bunch_it->num[1][0]))-(bunch_it->capfc[1]*bunch_it->capfc[0]);
00073 bunch_it->sigfc[2][0] = (bunch_it->prodfc[2][0]/(bunch_it->num[2][0]))-(bunch_it->capfc[2]*bunch_it->capfc[0]);
00074 bunch_it->sigfc[2][1] = (bunch_it->prodfc[2][1]/(bunch_it->num[2][1]))-(bunch_it->capfc[2]*bunch_it->capfc[1]);
00075 bunch_it->sigfc[3][0] = (bunch_it->prodfc[3][0]/(bunch_it->num[3][0]))-(bunch_it->capfc[3]*bunch_it->capfc[0]);
00076 bunch_it->sigfc[3][1] = (bunch_it->prodfc[3][1]/(bunch_it->num[3][1]))-(bunch_it->capfc[3]*bunch_it->capfc[1]);
00077 bunch_it->sigfc[3][2] = (bunch_it->prodfc[3][2]/(bunch_it->num[3][2]))-(bunch_it->capfc[3]*bunch_it->capfc[2]);
00078
00079 for(int i = 0; i != 4; i++){
00080 for(int j = 0; j != 4; j++){
00081 bunch_it->covarhistADC->Fill(i,j,bunch_it->sig[i][j]);
00082 bunch_it->covarhistfC->Fill(i,j,bunch_it->sigfc[i][j]);
00083 }
00084 }
00085
00086 if(bunch_it->detid.subdet() == 1){
00087 theFile->cd("HB");
00088 bunch_it->covarhistADC->Write();
00089 bunch_it->covarhistfC->Write();
00090 for(int i = 0; i != 4; i++){
00091 bunch_it->hist[i]->Write();
00092 HBMeans->Fill(bunch_it->cap[i]);
00093 HBWidths->Fill(sqrt(bunch_it->sig[i][i]));
00094 }
00095 }
00096 if(bunch_it->detid.subdet() == 2){
00097 theFile->cd("HE");
00098 bunch_it->covarhistADC->Write();
00099 bunch_it->covarhistfC->Write();
00100 for(int i = 0; i != 4; i++){
00101 bunch_it->hist[i]->Write();
00102 HEMeans->Fill(bunch_it->cap[i]);
00103 HEWidths->Fill(sqrt(bunch_it->sig[i][i]));
00104 }
00105 }
00106 if(bunch_it->detid.subdet() == 3){
00107 theFile->cd("HO");
00108 bunch_it->covarhistADC->Write();
00109 bunch_it->covarhistfC->Write();
00110 for(int i = 0; i != 4; i++){
00111 bunch_it->hist[i]->Write();
00112 HOMeans->Fill(bunch_it->cap[i]);
00113 HOWidths->Fill(sqrt(bunch_it->sig[i][i]));
00114 }
00115 }
00116 if(bunch_it->detid.subdet() == 4){
00117 theFile->cd("HF");
00118 bunch_it->covarhistADC->Write();
00119 bunch_it->covarhistfC->Write();
00120 for(int i = 0; i != 4; i++){
00121 bunch_it->hist[i]->Write();
00122 HFMeans->Fill(bunch_it->cap[i]);
00123 HFWidths->Fill(sqrt(bunch_it->sig[i][i]));
00124 }
00125 }
00126
00127 theFile->cd();
00128
00129 const HcalPedestal item(bunch_it->detid, bunch_it->cap[0], bunch_it->cap[1], bunch_it->cap[2], bunch_it->cap[3]);
00130 rawPedsItem->addValues(item);
00131 HcalPedestalWidth widthsp(bunch_it->detid);
00132 widthsp.setSigma(0,0,bunch_it->sig[0][0]);
00133 widthsp.setSigma(0,1,bunch_it->sig[0][1]);
00134 widthsp.setSigma(0,2,bunch_it->sig[0][2]);
00135 widthsp.setSigma(0,3,bunch_it->sig[0][3]);
00136 widthsp.setSigma(1,1,bunch_it->sig[1][1]);
00137 widthsp.setSigma(1,2,bunch_it->sig[1][2]);
00138 widthsp.setSigma(1,3,bunch_it->sig[1][3]);
00139 widthsp.setSigma(2,2,bunch_it->sig[2][2]);
00140 widthsp.setSigma(2,3,bunch_it->sig[2][3]);
00141 widthsp.setSigma(3,3,bunch_it->sig[3][3]);
00142 rawWidthsItem->addValues(widthsp);
00143
00144 const HcalPedestal itemfc(bunch_it->detid, bunch_it->capfc[0], bunch_it->capfc[1],
00145 bunch_it->capfc[2], bunch_it->capfc[3]);
00146 rawPedsItemfc->addValues(itemfc);
00147 HcalPedestalWidth widthspfc(bunch_it->detid);
00148 widthspfc.setSigma(0,0,bunch_it->sigfc[0][0]);
00149 widthspfc.setSigma(0,1,bunch_it->sigfc[0][1]);
00150 widthspfc.setSigma(0,2,bunch_it->sigfc[0][2]);
00151 widthspfc.setSigma(0,3,bunch_it->sigfc[0][3]);
00152 widthspfc.setSigma(1,1,bunch_it->sigfc[1][1]);
00153 widthspfc.setSigma(1,2,bunch_it->sigfc[1][2]);
00154 widthspfc.setSigma(1,3,bunch_it->sigfc[1][3]);
00155 widthspfc.setSigma(2,2,bunch_it->sigfc[2][2]);
00156 widthspfc.setSigma(2,3,bunch_it->sigfc[2][3]);
00157 widthspfc.setSigma(3,3,bunch_it->sigfc[3][3]);
00158 rawWidthsItemfc->addValues(widthspfc);
00159 }
00160 }
00161
00162
00163 std::ofstream outStream1(pedsADCfilename.c_str());
00164 HcalDbASCIIIO::dumpObject (outStream1, (*rawPedsItem) );
00165 std::ofstream outStream2(widthsADCfilename.c_str());
00166 HcalDbASCIIIO::dumpObject (outStream2, (*rawWidthsItem) );
00167
00168 std::ofstream outStream3(pedsfCfilename.c_str());
00169 HcalDbASCIIIO::dumpObject (outStream3, (*rawPedsItemfc) );
00170 std::ofstream outStream4(widthsfCfilename.c_str());
00171 HcalDbASCIIIO::dumpObject (outStream4, (*rawWidthsItemfc) );
00172
00173 theFile->cd();
00174 theFile->cd("HB");
00175 HBMeans->Write();
00176 HBWidths->Write();
00177 theFile->cd();
00178 theFile->cd("HF");
00179 HFMeans->Write();
00180 HFWidths->Write();
00181 theFile->cd();
00182 theFile->cd("HE");
00183 HEMeans->Write();
00184 HEWidths->Write();
00185 theFile->cd();
00186 theFile->cd("HO");
00187 HOMeans->Write();
00188 HOWidths->Write();
00189
00190 std::cout << "Writing ROOT file... ";
00191 theFile->Close();
00192 std::cout << "ROOT file closed.\n";
00193
00194 delete rawPedsItem;
00195 delete rawWidthsItem;
00196 delete rawPedsItemfc;
00197 delete rawWidthsItemfc;
00198 }
00199
00200
00201 void
00202 HcalPedestalWidthsValidation::analyze(const edm::Event& e, const edm::EventSetup& iSetup)
00203 {
00204
00205 edm::Handle<HBHEDigiCollection> hbhe; e.getByType(hbhe);
00206 edm::Handle<HODigiCollection> ho; e.getByType(ho);
00207 edm::Handle<HFDigiCollection> hf; e.getByType(hf);
00208 edm::ESHandle<HcalDbService> conditions;
00209 iSetup.get<HcalDbRecord>().get(conditions);
00210
00211 const HcalQIEShape* shape = conditions->getHcalShape();
00212
00213 if(firsttime)
00214 {
00215 runnum = e.id().run();
00216 std::string runnum_string;
00217 std::stringstream tempstringout;
00218 tempstringout << runnum;
00219 runnum_string = tempstringout.str();
00220 ROOTfilename = runnum_string + "-val_peds_ADC.root";
00221 pedsADCfilename = runnum_string + "-val_peds_ADC.txt";
00222 pedsfCfilename = runnum_string + "-val_peds_fC.txt";
00223 widthsADCfilename = runnum_string + "-val_widths_ADC.txt";
00224 widthsfCfilename = runnum_string + "-val_widths_fC.txt";
00225
00226 theFile = new TFile(ROOTfilename.c_str(), "RECREATE");
00227 theFile->cd();
00228
00229 theFile->mkdir("HB");
00230 theFile->mkdir("HE");
00231 theFile->mkdir("HF");
00232 theFile->mkdir("HO");
00233 theFile->cd();
00234
00235 HBMeans = new TH1F("All Ped Means HB","All Ped Means HB", 100, 0, 9);
00236 HBWidths = new TH1F("All Ped Widths HB","All Ped Widths HB", 100, 0, 3);
00237 HEMeans = new TH1F("All Ped Means HE","All Ped Means HE", 100, 0, 9);
00238 HEWidths = new TH1F("All Ped Widths HE","All Ped Widths HE", 100, 0, 3);
00239 HFMeans = new TH1F("All Ped Means HF","All Ped Means HF", 100, 0, 9);
00240 HFWidths = new TH1F("All Ped Widths HF","All Ped Widths HF", 100, 0, 3);
00241 HOMeans = new TH1F("All Ped Means HO","All Ped Means HO", 100, 0, 9);
00242 HOWidths = new TH1F("All Ped Widths HO","All Ped Widths HO", 100, 0, 3);
00243
00244 edm::ESHandle<HcalElectronicsMap> refEMap;
00245 iSetup.get<HcalElectronicsMapRcd>().get(refEMap);
00246 const HcalElectronicsMap* myRefEMap = refEMap.product();
00247 std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
00248 for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); it++)
00249 {
00250 HcalGenericDetId mygenid(it->rawId());
00251 if(mygenid.isHcalDetId())
00252 {
00253 NewPedBunchVal a;
00254 HcalDetId chanid(mygenid.rawId());
00255 a.detid = chanid;
00256 a.usedflag = false;
00257 std::ostringstream s1;
00258 s1 << mygenid;
00259 std::string histname = s1.str();
00260 histname += " Covariance Matrix";
00261 a.covarhistADC = new TH2F(histname.c_str(), histname.c_str(), 4, -.5, 3.5, 4, -.5, 3.5);
00262 a.covarhistADC->SetOption("TEXT");
00263 a.covarhistADC->SetStats(0);
00264 histname += " fC";
00265 a.covarhistfC = new TH2F(histname.c_str(), histname.c_str(), 4, -.5, 3.5, 4, -.5, 3.5);
00266 a.covarhistfC->SetOption("TEXT");
00267 a.covarhistfC->SetStats(0);
00268 for(int i = 0; i != 4; i++)
00269 {
00270 a.cap[i] = 0;
00271 a.capfc[i] = 0;
00272 std::ostringstream s3;
00273 s3 << mygenid;
00274 std::string histname = s3.str() + " Cap ";
00275 std::ostringstream s4;
00276 s4 << i;
00277 histname += s4.str();
00278 std::ostringstream s5;
00279 s5 << std::hex << (mygenid.rawId()) << std::dec;
00280 std::string histnamedetid = s5.str();
00281 a.hist[i] = new TH1F(histname.c_str(), histnamedetid.c_str(), 16, -.5, 15.5);
00282 for(int j = 0; j != 4; j++)
00283 {
00284 a.sig[i][j] = 0;
00285 a.sigfc[i][j] = 0;
00286 a.prod[i][j] = 0;
00287 a.prodfc[i][j] = 0;
00288 a.num[i][j] = 0;
00289 }
00290 }
00291 BunchVales.push_back(a);
00292
00293 }
00294 }
00295 firsttime = false;
00296 }
00297
00298 std::vector<NewPedBunchVal>::iterator bunch_it;
00299
00300 for(HBHEDigiCollection::const_iterator j = hbhe->begin(); j != hbhe->end(); j++)
00301 {
00302 const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00303 for(bunch_it = BunchVales.begin(); bunch_it != BunchVales.end(); bunch_it++)
00304 if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00305 bunch_it->usedflag = true;
00306 for(int ts = firstTS; ts != lastTS+1; ts++)
00307 {
00308 if(digi.sample(ts).adc() > 15) continue;
00309 bunch_it->hist[digi.sample(ts).capid()]->Fill(digi.sample(ts).adc());
00310 const HcalQIECoder* coder = conditions->getHcalCoder(digi.id().rawId());
00311 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
00312 bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
00313 double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
00314 bunch_it->capfc[digi.sample(ts).capid()] += charge1;
00315 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] += (digi.sample(ts).adc() * digi.sample(ts).adc());
00316 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
00317 if((ts+1 < digi.size()) && (ts+1 < lastTS)){
00318 if(digi.sample(ts+1).adc() > 15) continue;
00319 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += digi.sample(ts).adc()*digi.sample(ts+1).adc();
00320 double charge2 = coder->charge(*shape, digi.sample(ts+1).adc(), digi.sample(ts+1).capid());
00321 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += charge1*charge2;
00322 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += 1;
00323 }
00324 if((ts+2 < digi.size()) && (ts+2 < lastTS)){
00325 if(digi.sample(ts+2).adc() > 15) continue;
00326 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += digi.sample(ts).adc()*digi.sample(ts+2).adc();
00327 double charge2 = coder->charge(*shape, digi.sample(ts+2).adc(), digi.sample(ts+2).capid());
00328 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += charge1*charge2;
00329 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += 1;
00330 }
00331 if((ts+3 < digi.size()) && (ts+3 < lastTS)){
00332 if(digi.sample(ts+3).adc() > 15) continue;
00333 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += digi.sample(ts).adc()*digi.sample(ts+3).adc();
00334 double charge2 = coder->charge(*shape, digi.sample(ts+3).adc(), digi.sample(ts+3).capid());
00335 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += charge1*charge2;
00336 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += 1;
00337 }
00338 }
00339 }
00340
00341 for(HODigiCollection::const_iterator j = ho->begin(); j != ho->end(); j++)
00342 {
00343 const HODataFrame digi = (const HODataFrame)(*j);
00344 for(bunch_it = BunchVales.begin(); bunch_it != BunchVales.end(); bunch_it++)
00345 if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00346 bunch_it->usedflag = true;
00347 for(int ts = firstTS; ts <= lastTS; ts++)
00348 {
00349 if(digi.sample(ts).adc() > 15) continue;
00350 bunch_it->hist[digi.sample(ts).capid()]->Fill(digi.sample(ts).adc());
00351 const HcalQIECoder* coder = conditions->getHcalCoder(digi.id().rawId());
00352 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
00353 bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
00354 double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
00355 bunch_it->capfc[digi.sample(ts).capid()] += charge1;
00356 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] += (digi.sample(ts).adc() * digi.sample(ts).adc());
00357 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
00358 if((ts+1 < digi.size()) && (ts+1 < lastTS)){
00359 if(digi.sample(ts+1).adc() > 15) continue;
00360 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += digi.sample(ts).adc()*digi.sample(ts+1).adc();
00361 double charge2 = coder->charge(*shape, digi.sample(ts+1).adc(), digi.sample(ts+1).capid());
00362 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += charge1*charge2;
00363 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += 1;
00364 }
00365 if((ts+2 < digi.size()) && (ts+2 < lastTS)){
00366 if(digi.sample(ts+2).adc() > 15) continue;
00367 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += digi.sample(ts).adc()*digi.sample(ts+2).adc();
00368 double charge2 = coder->charge(*shape, digi.sample(ts+2).adc(), digi.sample(ts+2).capid());
00369 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += charge1*charge2;
00370 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += 1;
00371 }
00372 if((ts+3 < digi.size()) && (ts+3 < lastTS)){
00373 if(digi.sample(ts+3).adc() > 15) continue;
00374 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += digi.sample(ts).adc()*digi.sample(ts+3).adc();
00375 double charge2 = coder->charge(*shape, digi.sample(ts+3).adc(), digi.sample(ts+3).capid());
00376 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += charge1*charge2;
00377 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += 1;
00378 }
00379 }
00380 }
00381
00382 for(HFDigiCollection::const_iterator j = hf->begin(); j != hf->end(); j++)
00383 {
00384 const HFDataFrame digi = (const HFDataFrame)(*j);
00385 for(bunch_it = BunchVales.begin(); bunch_it != BunchVales.end(); bunch_it++)
00386 if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00387 bunch_it->usedflag = true;
00388 for(int ts = firstTS; ts <= lastTS; ts++)
00389 {
00390 if(digi.sample(ts).adc() > 15) continue;
00391 bunch_it->hist[digi.sample(ts).capid()]->Fill(digi.sample(ts).adc());
00392 const HcalQIECoder* coder = conditions->getHcalCoder(digi.id().rawId());
00393 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
00394 bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
00395 double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
00396 bunch_it->capfc[digi.sample(ts).capid()] += charge1;
00397 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] += (digi.sample(ts).adc() * digi.sample(ts).adc());
00398 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
00399 if((ts+1 < digi.size()) && (ts+1 < lastTS)){
00400 if(digi.sample(ts+1).adc() > 15) continue;
00401 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += digi.sample(ts).adc()*digi.sample(ts+1).adc();
00402 double charge2 = coder->charge(*shape, digi.sample(ts+1).adc(), digi.sample(ts+1).capid());
00403 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += charge1*charge2;
00404 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+1).capid()] += 1;
00405 }
00406 if((ts+2 < digi.size()) && (ts+2 < lastTS)){
00407 if(digi.sample(ts+2).adc() > 15) continue;
00408 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += digi.sample(ts).adc()*digi.sample(ts+2).adc();
00409 double charge2 = coder->charge(*shape, digi.sample(ts+2).adc(), digi.sample(ts+2).capid());
00410 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += charge1*charge2;
00411 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+2).capid()] += 1;
00412 }
00413 if((ts+3 < digi.size()) && (ts+3 < lastTS)){
00414 if(digi.sample(ts+3).adc() > 15) continue;
00415 bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += digi.sample(ts).adc()*digi.sample(ts+3).adc();
00416 double charge2 = coder->charge(*shape, digi.sample(ts+3).adc(), digi.sample(ts+3).capid());
00417 bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += charge1*charge2;
00418 bunch_it->num[digi.sample(ts).capid()][digi.sample(ts+3).capid()] += 1;
00419 }
00420 }
00421 }
00422
00423
00424 }
00425
00426
00427 DEFINE_FWK_MODULE(HcalPedestalWidthsValidation);