Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:34

0001 // Original Author:  Steven Won
0002 // Original Author:  Alan Campbell for castor
0003 //         Created:  Fri May  2 15:34:43 CEST 2008
0004 // Written to replace the combination of HcalPedestalAnalyzer and HcalPedestalAnalysis
0005 // This code runs 1000x faster and produces all outputs from a single run
0006 // (ADC, fC in .txt plus an .xml file)
0007 //
0008 
0009 #include "FWCore/Framework/interface/Frameworkfwd.h"
0010 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ServiceRegistry/interface/Service.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 
0018 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0019 #include "DataFormats/Common/interface/Handle.h"
0020 #include "DataFormats/HcalDetId/interface/HcalGenericDetId.h"
0021 #include "DataFormats/HcalDetId/interface/HcalElectronicsId.h"
0022 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0023 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0024 #include "DataFormats/DetId/interface/DetId.h"
0025 #include "CondFormats/CastorObjects/interface/CastorPedestals.h"
0026 #include "CondFormats/CastorObjects/interface/CastorPedestalWidths.h"
0027 #include "CondFormats/CastorObjects/interface/CastorQIECoder.h"
0028 #include "CondFormats/CastorObjects/interface/CastorQIEData.h"
0029 #include "CondFormats/CastorObjects/interface/CastorQIEShape.h"
0030 #include "CondFormats/CastorObjects/interface/CastorElectronicsMap.h"
0031 #include "CondFormats/CastorObjects/interface/AllObjects.h"
0032 
0033 #include "CalibFormats/CastorObjects/interface/CastorDbRecord.h"
0034 #include "CalibFormats/CastorObjects/interface/CastorDbService.h"
0035 #include "CalibFormats/CastorObjects/interface/CastorCalibrations.h"
0036 #include "CalibFormats/CastorObjects/interface/CastorCalibrationWidths.h"
0037 
0038 #include "CalibCalorimetry/CastorCalib/interface/CastorDbASCIIIO.h"
0039 #include "TBDataFormats/HcalTBObjects/interface/HcalTBTriggerData.h"
0040 
0041 #include "TProfile.h"
0042 #include "TH1.h"
0043 #include "TH2.h"
0044 #include "TCanvas.h"
0045 #include "TStyle.h"
0046 
0047 #include <cmath>
0048 #include <iostream>
0049 #include <map>
0050 #include <memory>
0051 #include <iomanip>
0052 #include <fstream>
0053 #include <vector>
0054 #include <string>
0055 
0056 namespace edm {
0057   class ParameterSet;
0058   class Event;
0059   class EventSetup;
0060 }  // namespace edm
0061 
0062 struct NewPedBunch {
0063   HcalCastorDetId detid;
0064   bool usedflag;
0065   float cap[4];
0066   float capfc[4];
0067   float sig[4][4];
0068   float sigfc[4][4];
0069   float prod[4][4];
0070   float prodfc[4][4];
0071   int num[4][4];
0072 };
0073 
0074 class CastorPedestalsAnalysis : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0075 public:
0076   //Constructor
0077   CastorPedestalsAnalysis(const edm::ParameterSet& ps);
0078   //Destructor
0079   ~CastorPedestalsAnalysis() override;
0080   //Analysis
0081   void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override;
0082 
0083 private:
0084   //Container for data, 1 per channel
0085   std::vector<NewPedBunch> Bunches;
0086   const edm::EDGetTokenT<CastorDigiCollection> castorDigiCollectionToken_;
0087   const edm::ESGetToken<CastorDbService, CastorDbRecord> tok_cond_;
0088   const edm::ESGetToken<CastorElectronicsMap, CastorElectronicsMapRcd> tok_map_;
0089   //Flag for saving histos
0090   const bool hiSaveFlag;
0091   const bool dumpXML;
0092   const bool verboseflag;
0093   const int firstTS;
0094   const int lastTS;
0095   bool firsttime;
0096   int runnum;
0097   std::string ROOTfilename;
0098   std::string pedsADCfilename;
0099   std::string pedsfCfilename;
0100   std::string widthsADCfilename;
0101   std::string widthsfCfilename;
0102   std::string XMLfilename;
0103   std::string XMLtag;
0104   std::string ZSfilename;
0105 
0106   TH1F* CASTORMeans;
0107   TH1F* CASTORWidths;
0108 
0109   // TH2F *dephist[4];
0110   TH2F* dephist;
0111 };
0112 
0113 CastorPedestalsAnalysis::CastorPedestalsAnalysis(const edm::ParameterSet& ps)
0114     : castorDigiCollectionToken_(
0115           consumes<CastorDigiCollection>(ps.getParameter<edm::InputTag>("castorDigiCollectionTag"))),
0116       tok_cond_(esConsumes<CastorDbService, CastorDbRecord>()),
0117       tok_map_(esConsumes<CastorElectronicsMap, CastorElectronicsMapRcd>()),
0118       hiSaveFlag(ps.getUntrackedParameter<bool>("hiSaveFlag", false)),
0119       dumpXML(ps.getUntrackedParameter<bool>("dumpXML", false)),
0120       verboseflag(ps.getUntrackedParameter<bool>("verbose", false)),
0121       firstTS(ps.getUntrackedParameter<int>("firstTS", 0)),
0122       lastTS(ps.getUntrackedParameter<int>("lastTS", 9)) {
0123   usesResource(TFileService::kSharedResource);
0124 
0125   firsttime = true;
0126 }
0127 
0128 CastorPedestalsAnalysis::~CastorPedestalsAnalysis() {
0129   CastorPedestals* rawPedsItem = new CastorPedestals(true);
0130   CastorPedestalWidths* rawWidthsItem = new CastorPedestalWidths(true);
0131   CastorPedestals* rawPedsItemfc = new CastorPedestals(false);
0132   CastorPedestalWidths* rawWidthsItemfc = new CastorPedestalWidths(false);
0133 
0134   //Calculate pedestal constants
0135   std::cout << "Calculating Pedestal constants...\n";
0136   std::vector<NewPedBunch>::iterator bunch_it;
0137   for (bunch_it = Bunches.begin(); bunch_it != Bunches.end(); ++bunch_it) {
0138     if (bunch_it->usedflag) {
0139       if (verboseflag)
0140         std::cout << "Analyzing channel sector= " << bunch_it->detid.sector()
0141                   << " module = " << bunch_it->detid.module() << std::endl;
0142       //pedestal constant is the mean
0143       bunch_it->cap[0] /= bunch_it->num[0][0];
0144       bunch_it->cap[1] /= bunch_it->num[1][1];
0145       bunch_it->cap[2] /= bunch_it->num[2][2];
0146       bunch_it->cap[3] /= bunch_it->num[3][3];
0147       bunch_it->capfc[0] /= bunch_it->num[0][0];
0148       bunch_it->capfc[1] /= bunch_it->num[1][1];
0149       bunch_it->capfc[2] /= bunch_it->num[2][2];
0150       bunch_it->capfc[3] /= bunch_it->num[3][3];
0151       //widths are the covariance matrix--assumed symmetric
0152       bunch_it->sig[0][0] = (bunch_it->prod[0][0] / bunch_it->num[0][0]) - (bunch_it->cap[0] * bunch_it->cap[0]);
0153       bunch_it->sig[0][1] = (bunch_it->prod[0][1] / bunch_it->num[0][1]) - (bunch_it->cap[0] * bunch_it->cap[1]);
0154       bunch_it->sig[0][2] = (bunch_it->prod[0][2] / bunch_it->num[0][2]) - (bunch_it->cap[0] * bunch_it->cap[2]);
0155       bunch_it->sig[0][3] = (bunch_it->prod[0][3] / bunch_it->num[0][3]) - (bunch_it->cap[0] * bunch_it->cap[3]);
0156       bunch_it->sig[1][0] = (bunch_it->prod[1][0] / bunch_it->num[1][0]) - (bunch_it->cap[1] * bunch_it->cap[0]);
0157       bunch_it->sig[1][1] = (bunch_it->prod[1][1] / bunch_it->num[1][1]) - (bunch_it->cap[1] * bunch_it->cap[1]);
0158       bunch_it->sig[1][2] = (bunch_it->prod[1][2] / bunch_it->num[1][2]) - (bunch_it->cap[1] * bunch_it->cap[2]);
0159       bunch_it->sig[1][3] = (bunch_it->prod[1][3] / bunch_it->num[1][3]) - (bunch_it->cap[1] * bunch_it->cap[3]);
0160       bunch_it->sig[2][0] = (bunch_it->prod[2][0] / bunch_it->num[2][0]) - (bunch_it->cap[2] * bunch_it->cap[0]);
0161       bunch_it->sig[2][1] = (bunch_it->prod[2][1] / bunch_it->num[2][1]) - (bunch_it->cap[2] * bunch_it->cap[1]);
0162       bunch_it->sig[2][2] = (bunch_it->prod[2][2] / bunch_it->num[2][2]) - (bunch_it->cap[2] * bunch_it->cap[2]);
0163       bunch_it->sig[2][3] = (bunch_it->prod[2][3] / bunch_it->num[2][3]) - (bunch_it->cap[2] * bunch_it->cap[3]);
0164       bunch_it->sig[3][0] = (bunch_it->prod[3][0] / bunch_it->num[3][0]) - (bunch_it->cap[3] * bunch_it->cap[0]);
0165       bunch_it->sig[3][1] = (bunch_it->prod[3][1] / bunch_it->num[3][1]) - (bunch_it->cap[3] * bunch_it->cap[1]);
0166       bunch_it->sig[3][2] = (bunch_it->prod[3][2] / bunch_it->num[3][2]) - (bunch_it->cap[3] * bunch_it->cap[2]);
0167       bunch_it->sig[3][3] = (bunch_it->prod[3][3] / bunch_it->num[3][3]) - (bunch_it->cap[3] * bunch_it->cap[3]);
0168 
0169       bunch_it->sigfc[0][0] =
0170           (bunch_it->prodfc[0][0] / bunch_it->num[0][0]) - (bunch_it->capfc[0] * bunch_it->capfc[0]);
0171       bunch_it->sigfc[0][1] =
0172           (bunch_it->prodfc[0][1] / bunch_it->num[0][1]) - (bunch_it->capfc[0] * bunch_it->capfc[1]);
0173       bunch_it->sigfc[0][2] =
0174           (bunch_it->prodfc[0][2] / bunch_it->num[0][2]) - (bunch_it->capfc[0] * bunch_it->capfc[2]);
0175       bunch_it->sigfc[0][3] =
0176           (bunch_it->prodfc[0][3] / bunch_it->num[0][3]) - (bunch_it->capfc[0] * bunch_it->capfc[3]);
0177       bunch_it->sigfc[1][0] =
0178           (bunch_it->prodfc[1][0] / bunch_it->num[1][0]) - (bunch_it->capfc[1] * bunch_it->capfc[0]);
0179       bunch_it->sigfc[1][1] =
0180           (bunch_it->prodfc[1][1] / bunch_it->num[1][1]) - (bunch_it->capfc[1] * bunch_it->capfc[1]);
0181       bunch_it->sigfc[1][2] =
0182           (bunch_it->prodfc[1][2] / bunch_it->num[1][2]) - (bunch_it->capfc[1] * bunch_it->capfc[2]);
0183       bunch_it->sigfc[1][3] =
0184           (bunch_it->prodfc[1][3] / bunch_it->num[1][3]) - (bunch_it->capfc[1] * bunch_it->capfc[3]);
0185       bunch_it->sigfc[2][0] =
0186           (bunch_it->prodfc[2][0] / bunch_it->num[2][0]) - (bunch_it->capfc[2] * bunch_it->capfc[0]);
0187       bunch_it->sigfc[2][1] =
0188           (bunch_it->prodfc[2][1] / bunch_it->num[2][1]) - (bunch_it->capfc[2] * bunch_it->capfc[1]);
0189       bunch_it->sigfc[2][2] =
0190           (bunch_it->prodfc[2][2] / bunch_it->num[2][2]) - (bunch_it->capfc[2] * bunch_it->capfc[2]);
0191       bunch_it->sigfc[2][3] =
0192           (bunch_it->prodfc[2][3] / bunch_it->num[2][3]) - (bunch_it->capfc[2] * bunch_it->capfc[3]);
0193       bunch_it->sigfc[3][0] =
0194           (bunch_it->prodfc[3][0] / bunch_it->num[3][0]) - (bunch_it->capfc[3] * bunch_it->capfc[0]);
0195       bunch_it->sigfc[3][1] =
0196           (bunch_it->prodfc[3][1] / bunch_it->num[3][1]) - (bunch_it->capfc[3] * bunch_it->capfc[1]);
0197       bunch_it->sigfc[3][2] =
0198           (bunch_it->prodfc[3][2] / bunch_it->num[3][2]) - (bunch_it->capfc[3] * bunch_it->capfc[2]);
0199       bunch_it->sigfc[3][3] =
0200           (bunch_it->prodfc[3][3] / bunch_it->num[3][3]) - (bunch_it->capfc[3] * bunch_it->capfc[3]);
0201 
0202       for (int i = 0; i != 3; i++) {
0203         CASTORMeans->Fill(bunch_it->cap[i]);
0204         CASTORWidths->Fill(bunch_it->sig[i][i]);
0205       }
0206 
0207       //if(bunch_it->detid.subdet() == 1){
0208 
0209       int fillphi = bunch_it->detid.sector();
0210       //if (bunch_it->detid.depth()==4) fillphi++;
0211 
0212       //    dephist[bunch_it->detid.module()-1]->Fill(bunch_it->detid.ieta(),fillphi,
0213       //             (bunch_it->cap[0]+bunch_it->cap[1]+bunch_it->cap[2]+bunch_it->cap[3])/4);
0214       dephist->Fill(bunch_it->detid.module(),
0215                     fillphi,
0216                     (bunch_it->cap[0] + bunch_it->cap[1] + bunch_it->cap[2] + bunch_it->cap[3]) / 4);
0217 
0218       const CastorPedestal item(bunch_it->detid,
0219                                 bunch_it->cap[0],
0220                                 bunch_it->cap[1],
0221                                 bunch_it->cap[2],
0222                                 bunch_it->cap[3],
0223                                 bunch_it->sig[0][0],
0224                                 bunch_it->sig[1][1],
0225                                 bunch_it->sig[2][2],
0226                                 bunch_it->sig[3][3]);
0227       rawPedsItem->addValues(item);
0228       CastorPedestalWidth widthsp(bunch_it->detid);
0229       widthsp.setSigma(0, 0, bunch_it->sig[0][0]);
0230       widthsp.setSigma(0, 1, bunch_it->sig[0][1]);
0231       widthsp.setSigma(0, 2, bunch_it->sig[0][2]);
0232       widthsp.setSigma(0, 3, bunch_it->sig[0][3]);
0233       widthsp.setSigma(1, 0, bunch_it->sig[1][0]);
0234       widthsp.setSigma(1, 1, bunch_it->sig[1][1]);
0235       widthsp.setSigma(1, 2, bunch_it->sig[1][2]);
0236       widthsp.setSigma(1, 3, bunch_it->sig[1][3]);
0237       widthsp.setSigma(2, 0, bunch_it->sig[2][0]);
0238       widthsp.setSigma(2, 1, bunch_it->sig[2][1]);
0239       widthsp.setSigma(2, 2, bunch_it->sig[2][2]);
0240       widthsp.setSigma(2, 3, bunch_it->sig[2][3]);
0241       widthsp.setSigma(3, 0, bunch_it->sig[3][0]);
0242       widthsp.setSigma(3, 1, bunch_it->sig[3][1]);
0243       widthsp.setSigma(3, 2, bunch_it->sig[3][2]);
0244       widthsp.setSigma(3, 3, bunch_it->sig[3][3]);
0245       rawWidthsItem->addValues(widthsp);
0246 
0247       const CastorPedestal itemfc(bunch_it->detid,
0248                                   bunch_it->capfc[0],
0249                                   bunch_it->capfc[1],
0250                                   bunch_it->capfc[2],
0251                                   bunch_it->capfc[3],
0252                                   bunch_it->sigfc[0][0],
0253                                   bunch_it->sigfc[1][1],
0254                                   bunch_it->sigfc[2][2],
0255                                   bunch_it->sigfc[3][3]);
0256       rawPedsItemfc->addValues(itemfc);
0257       CastorPedestalWidth widthspfc(bunch_it->detid);
0258       widthspfc.setSigma(0, 0, bunch_it->sigfc[0][0]);
0259       widthspfc.setSigma(0, 1, bunch_it->sigfc[0][1]);
0260       widthspfc.setSigma(0, 2, bunch_it->sigfc[0][2]);
0261       widthspfc.setSigma(0, 3, bunch_it->sigfc[0][3]);
0262       widthspfc.setSigma(1, 0, bunch_it->sigfc[1][0]);
0263       widthspfc.setSigma(1, 1, bunch_it->sigfc[1][1]);
0264       widthspfc.setSigma(1, 2, bunch_it->sigfc[1][2]);
0265       widthspfc.setSigma(1, 3, bunch_it->sigfc[1][3]);
0266       widthspfc.setSigma(2, 0, bunch_it->sigfc[2][0]);
0267       widthspfc.setSigma(2, 1, bunch_it->sigfc[2][1]);
0268       widthspfc.setSigma(2, 2, bunch_it->sigfc[2][2]);
0269       widthspfc.setSigma(2, 3, bunch_it->sigfc[2][3]);
0270       widthspfc.setSigma(3, 0, bunch_it->sigfc[3][0]);
0271       widthspfc.setSigma(3, 1, bunch_it->sigfc[3][1]);
0272       widthspfc.setSigma(3, 2, bunch_it->sigfc[3][2]);
0273       widthspfc.setSigma(3, 3, bunch_it->sigfc[3][3]);
0274       rawWidthsItemfc->addValues(widthspfc);
0275     }
0276   }
0277 
0278   // dump the resulting list of pedestals into a file
0279   std::ofstream outStream1(pedsADCfilename.c_str());
0280   CastorDbASCIIIO::dumpObject(outStream1, (*rawPedsItem));
0281   std::ofstream outStream2(widthsADCfilename.c_str());
0282   CastorDbASCIIIO::dumpObject(outStream2, (*rawWidthsItem));
0283 
0284   std::ofstream outStream3(pedsfCfilename.c_str());
0285   CastorDbASCIIIO::dumpObject(outStream3, (*rawPedsItemfc));
0286   std::ofstream outStream4(widthsfCfilename.c_str());
0287   CastorDbASCIIIO::dumpObject(outStream4, (*rawWidthsItemfc));
0288 
0289   if (dumpXML) {
0290     std::ofstream outStream5(XMLfilename.c_str());
0291     //  CastorCondXML::dumpObject (outStream5, runnum, runnum, runnum, XMLtag, 1, (*rawPedsItem), (*rawWidthsItem));
0292   }
0293 
0294   dephist->SetDrawOption("colz");
0295   dephist->GetXaxis()->SetTitle("module");
0296   dephist->GetYaxis()->SetTitle("sector");
0297 
0298   std::stringstream tempstringout;
0299   tempstringout << runnum;
0300   std::string name1 = tempstringout.str() + "_pedplots_1d.png";
0301   std::string name2 = tempstringout.str() + "_pedplots_2d.png";
0302 
0303   TStyle* theStyle = new TStyle("style", "null");
0304   theStyle->SetPalette(1, nullptr);
0305   theStyle->SetCanvasDefH(1200);  //Height of canvas
0306   theStyle->SetCanvasDefW(1600);  //Width of canvas
0307 
0308   gStyle = theStyle;
0309   /*
0310     TCanvas * c1 = new TCanvas("c1","graph",1);
0311     c1->Divide(2,2);
0312     c1->cd(1);
0313     CASTORMeans->Draw();
0314     c1->SaveAs(name1.c_str());   
0315 
0316     theStyle->SetOptStat("n");
0317     gStyle = theStyle;
0318 
0319     TCanvas * c2 = new TCanvas("c2","graph",1);
0320  //   c2->Divide(2,2);
0321     c2->cd(1);
0322     dephist->Draw();
0323     dephist->SetDrawOption("colz");
0324     //c2->cd(2);
0325     //dephist[1]->Draw();
0326     //dephist[1]->SetDrawOption("colz");
0327     //c2->cd(3);
0328     //dephist[2]->Draw();
0329     //dephist[2]->SetDrawOption("colz");
0330     //c2->cd(4);
0331     //dephist[3]->Draw();
0332     //dephist[3]->SetDrawOption("colz");
0333     c2->SaveAs(name2.c_str());
0334 */
0335   /*
0336   std::cout << "Writing ROOT file... ";
0337   theFile->Close();
0338   std::cout << "ROOT file closed.\n";
0339   */
0340 }
0341 
0342 // ------------ method called to for each event  ------------
0343 void CastorPedestalsAnalysis::analyze(const edm::Event& e, const edm::EventSetup& iSetup) {
0344   const edm::Handle<CastorDigiCollection>& castor = e.getHandle(castorDigiCollectionToken_);
0345 
0346   auto conditions = &iSetup.getData(tok_cond_);
0347   const CastorQIEShape* shape = conditions->getCastorShape();
0348 
0349   if (firsttime) {
0350     runnum = e.id().run();
0351     std::string runnum_string;
0352     std::stringstream tempstringout;
0353     tempstringout << runnum;
0354     runnum_string = tempstringout.str();
0355     ROOTfilename = runnum_string + "-peds_ADC.root";
0356     pedsADCfilename = runnum_string + "-peds_ADC.txt";
0357     pedsfCfilename = runnum_string + "-peds_fC.txt";
0358     widthsADCfilename = runnum_string + "-widths_ADC.txt";
0359     widthsfCfilename = runnum_string + "-widths_fC.txt";
0360     XMLfilename = runnum_string + "-peds_ADC_complete.xml";
0361     XMLtag = "Castor_pedestals_" + runnum_string;
0362 
0363     edm::Service<TFileService> fs;
0364     fs->mkdir("CASTOR");
0365     CASTORMeans = fs->make<TH1F>("All Ped Means CASTOR", "All Ped Means CASTOR", 100, 0, 9);
0366     CASTORWidths = fs->make<TH1F>("All Ped Widths CASTOR", "All Ped Widths CASTOR", 100, 0, 3);
0367 
0368     dephist = fs->make<TH2F>("Pedestals (ADC)", "All Castor", 14, 0., 14.5, 16, .5, 16.5);
0369 
0370     const CastorElectronicsMap* myRefEMap = &iSetup.getData(tok_map_);
0371     std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
0372     for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); ++it) {
0373       HcalGenericDetId mygenid(it->rawId());
0374       if (mygenid.isHcalCastorDetId()) {
0375         NewPedBunch a;
0376         HcalCastorDetId chanid(mygenid.rawId());
0377         a.detid = chanid;
0378         a.usedflag = false;
0379         std::string type = "CASTOR";
0380         for (int i = 0; i != 4; i++) {
0381           a.cap[i] = 0;
0382           a.capfc[i] = 0;
0383           for (int j = 0; j != 4; j++) {
0384             a.sig[i][j] = 0;
0385             a.sigfc[i][j] = 0;
0386             a.prod[i][j] = 0;
0387             a.prodfc[i][j] = 0;
0388             a.num[i][j] = 0;
0389           }
0390         }
0391         Bunches.push_back(a);
0392       }
0393     }
0394     firsttime = false;
0395   }
0396 
0397   std::vector<NewPedBunch>::iterator bunch_it;
0398 
0399   for (CastorDigiCollection::const_iterator j = castor->begin(); j != castor->end(); ++j) {
0400     const CastorDataFrame digi = (const CastorDataFrame)(*j);
0401     for (bunch_it = Bunches.begin(); bunch_it != Bunches.end(); ++bunch_it)
0402       if (bunch_it->detid.rawId() == digi.id().rawId())
0403         break;
0404     bunch_it->usedflag = true;
0405     for (int ts = firstTS; ts != lastTS + 1; ts++) {
0406       const CastorQIECoder* coder = conditions->getCastorCoder(digi.id().rawId());
0407       bunch_it->num[digi.sample(ts).capid()][digi.sample(ts).capid()] += 1;
0408       bunch_it->cap[digi.sample(ts).capid()] += digi.sample(ts).adc();
0409       double charge1 = coder->charge(*shape, digi.sample(ts).adc(), digi.sample(ts).capid());
0410       bunch_it->capfc[digi.sample(ts).capid()] += charge1;
0411       bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts).capid()] +=
0412           (digi.sample(ts).adc() * digi.sample(ts).adc());
0413       bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts).capid()] += charge1 * charge1;
0414       if ((ts + 1 < digi.size()) && (ts + 1 < lastTS)) {
0415         bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts + 1).capid()] +=
0416             digi.sample(ts).adc() * digi.sample(ts + 1).adc();
0417         double charge2 = coder->charge(*shape, digi.sample(ts + 1).adc(), digi.sample(ts + 1).capid());
0418         bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts + 1).capid()] += charge1 * charge2;
0419         bunch_it->num[digi.sample(ts).capid()][digi.sample(ts + 1).capid()] += 1;
0420       }
0421       if ((ts + 2 < digi.size()) && (ts + 2 < lastTS)) {
0422         bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts + 2).capid()] +=
0423             digi.sample(ts).adc() * digi.sample(ts + 2).adc();
0424         double charge2 = coder->charge(*shape, digi.sample(ts + 2).adc(), digi.sample(ts + 2).capid());
0425         bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts + 2).capid()] += charge1 * charge2;
0426         bunch_it->num[digi.sample(ts).capid()][digi.sample(ts + 2).capid()] += 1;
0427       }
0428       if ((ts + 3 < digi.size()) && (ts + 3 < lastTS)) {
0429         bunch_it->prod[digi.sample(ts).capid()][digi.sample(ts + 3).capid()] +=
0430             digi.sample(ts).adc() * digi.sample(ts + 3).adc();
0431         double charge2 = coder->charge(*shape, digi.sample(ts + 3).adc(), digi.sample(ts + 3).capid());
0432         bunch_it->prodfc[digi.sample(ts).capid()][digi.sample(ts + 3).capid()] += charge1 * charge2;
0433         bunch_it->num[digi.sample(ts).capid()][digi.sample(ts + 3).capid()] += 1;
0434       }
0435     }
0436   }
0437 
0438   //this is the last brace
0439 }
0440 
0441 //define this as a plug-in
0442 DEFINE_FWK_MODULE(CastorPedestalsAnalysis);