CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/DQM/HcalMonitorTasks/src/HcalNoiseMonitor.cc

Go to the documentation of this file.
00001 #include <cmath>
00002 #include <fstream>
00003 #include <algorithm>
00004 
00005 #include "DQM/HcalMonitorTasks/interface/HcalNoiseMonitor.h"
00006 
00007 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00008 #include "DataFormats/FEDRawData/interface/FEDTrailer.h"
00009 #include "DataFormats/FEDRawData/interface/FEDHeader.h"
00010 #include "EventFilter/HcalRawToDigi/interface/HcalDCCHeader.h"
00011 
00012 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00013 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00014 #include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
00015 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00016 
00017 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00018 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00019 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00020 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00021 
00022 #include "RecoMET/METAlgorithms/interface/HcalNoiseRBXArray.h"
00023 #include "CalibCalorimetry/HcalAlgos/interface/HcalPulseShapes.h"
00024 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00025 
00026 #include "FWCore/Common/interface/TriggerNames.h"
00027 
00028 HcalNoiseMonitor::HcalNoiseMonitor(const edm::ParameterSet& ps)
00029 {
00030    Online_                = ps.getUntrackedParameter<bool>("online",false);
00031    mergeRuns_             = ps.getUntrackedParameter<bool>("mergeRuns",false);
00032    enableCleanup_         = ps.getUntrackedParameter<bool>("enableCleanup",false);
00033    debug_                 = ps.getUntrackedParameter<int>("debug",0);
00034    prefixME_              = ps.getUntrackedParameter<std::string>("subSystemFolder","Hcal/");
00035    if(prefixME_.substr(prefixME_.size()-1,prefixME_.size())!="/")
00036       prefixME_.append("/");
00037    subdir_                = ps.getUntrackedParameter<std::string>("TaskFolder","NoiseMonitor_Hcal");
00038    if(subdir_.size()>0 && subdir_.substr(subdir_.size()-1,subdir_.size())!="/")
00039       subdir_.append("/");
00040    subdir_=prefixME_+subdir_;
00041    AllowedCalibTypes_     = ps.getUntrackedParameter<std::vector<int> > ("AllowedCalibTypes");
00042    skipOutOfOrderLS_      = ps.getUntrackedParameter<bool>("skipOutOfOrderLS",false);
00043    NLumiBlocks_           = ps.getUntrackedParameter<int>("NLumiBlocks",4000);
00044    makeDiagnostics_       = ps.getUntrackedParameter<bool>("makeDiagnostics",false);
00045 
00046    triggers_=ps.getUntrackedParameter<std::vector<std::string> >("nzsHLTnames");
00047       //["HLT_HcalPhiSym","HLT_HcalNoise_8E29]
00048    period_=ps.getUntrackedParameter<int>("NoiseeventPeriod",4096); //4096
00049    rawdataLabel_          = ps.getUntrackedParameter<edm::InputTag>("RawDataLabel");
00050    hltresultsLabel_       = ps.getUntrackedParameter<edm::InputTag>("HLTResultsLabel");
00051    hbheDigiLabel_         = ps.getUntrackedParameter<edm::InputTag>("hbheDigiLabel");
00052    hbheRechitLabel_       = ps.getUntrackedParameter<edm::InputTag>("hbheRechitLabel");
00053    noiseLabel_            = ps.getUntrackedParameter<edm::InputTag>("noiseLabel");
00054 
00055    mTrianglePeakTS        = 4;   // for now...
00056 
00057    mE2E10MinEnergy        = ps.getUntrackedParameter<double>("E2E10MinEnergy");
00058    mMinE2E10              = ps.getUntrackedParameter<double>("MinE2E10");
00059    mMaxE2E10              = ps.getUntrackedParameter<double>("MaxE2E10");
00060    mMaxHPDHitCount        = ps.getUntrackedParameter<int>("MaxHPDHitCount");
00061    mMaxHPDNoOtherHitCount = ps.getUntrackedParameter<int>("MaxHPDNoOtherHitCount");
00062    mMaxADCZeros           = ps.getUntrackedParameter<int>("MaxADCZeros");
00063    mTotalZeroMinEnergy    = ps.getUntrackedParameter<double>("TotalZeroMinEnergy");
00064 }
00065 
00066 HcalNoiseMonitor::~HcalNoiseMonitor() {}
00067 
00068 void HcalNoiseMonitor::reset()
00069 {
00070 }
00071 
00072 void HcalNoiseMonitor::cleanup()
00073 {
00074    if(dbe_)
00075    {
00076       dbe_->setCurrentFolder(subdir_);
00077       dbe_->removeContents();
00078    }
00079 }
00080 
00081 void HcalNoiseMonitor::beginRun(const edm::Run& run, const edm::EventSetup& c)
00082 {
00083    if(debug_ > 1)
00084       std::cout <<"HcalNoiseMonitor::beginRun"<< std::endl;
00085 
00086    HcalBaseDQMonitor::beginRun(run,c);
00087 
00088    if(tevt_ == 0)
00089       setup();
00090 
00091    if(mergeRuns_ == false)
00092       reset();
00093 
00094    return;
00095 }
00096 
00097 
00098 void HcalNoiseMonitor::setup()
00099 {
00100    HcalBaseDQMonitor::setup();
00101 
00102    if(debug_ > 1)
00103       std::cout << "<HcalNoiseMonitor::setup> Creating histograms" << std::endl;
00104 
00105    if(dbe_)
00106    {
00107       dbe_->setCurrentFolder(subdir_);
00108 
00109       // Fit-based
00110       dbe_->setCurrentFolder(subdir_ + "DoubleChi2/");
00111 
00112       hNominalChi2 = dbe_->book1D("Nominal_fit_chi2", "Nominal fit chi2, total charge > 20 fC", 100, 0, 200);
00113       hNominalChi2->setAxisTitle("Nominal fit #chi^{2}", 1);
00114 
00115       hLinearChi2 = dbe_->book1D("Linear_fit_chi2", "Linear fit chi2, total charge > 20 fC", 100, 0, 200);
00116       hLinearChi2->setAxisTitle("Linear fit #chi^{2}", 1);
00117       
00118       hLinearTestStatistics = dbe_->book1D("Lambda_linear", "#Lambda_{linear}, total charge > 20 fC", 100, -10, 10);
00119       hLinearTestStatistics->setAxisTitle("#Lambda_{linear}", 1);
00120 
00121       hRMS8OverMax = dbe_->book1D("RMS8_over_Max", "RMS8/max, total charge > 20 fC", 100, 0, 2);
00122       hRMS8OverMax->setAxisTitle("RMS8/max", 1);
00123 
00124       hRMS8OverMaxTestStatistics = dbe_->book1D("Lambda_RMS8_over_max", "#Lambda_{RMS8/Max}, total charge > 20 fC",
00125          100, -30, 10);
00126       hRMS8OverMaxTestStatistics->setAxisTitle("#Lambda_{RMS8/Max}", 1);
00127 
00128       hLambdaLinearVsTotalCharge = dbe_->book2D("Lambda_linear_vs_total_charge", "#Lambda_{Linear}",
00129          50, -5, 5, 25, 0, 500);
00130       hLambdaLinearVsTotalCharge->setAxisTitle("#Lambda_{linear}", 1);
00131       hLambdaLinearVsTotalCharge->setAxisTitle("Total charge", 2);
00132 
00133       hLambdaRMS8MaxVsTotalCharge = dbe_->book2D("Lambda_RMS8Max_vs_total_charge", "#Lambda_{RMS8/Max}",
00134          50, -15, 5, 25, 0, 500);
00135       hLambdaRMS8MaxVsTotalCharge->setAxisTitle("#Lambda_{RMS8/Max}", 1);
00136       hLambdaRMS8MaxVsTotalCharge->setAxisTitle("Total charge", 2);
00137 
00138       hTriangleLeftSlopeVsTS4 = dbe_->book2D("Triangle_fit_left_slope",
00139          "Triangle fit left distance vs. TS4", 50, 0, 10, 25, 0, 500);
00140       hTriangleLeftSlopeVsTS4->setAxisTitle("Left slope", 1);
00141       hTriangleLeftSlopeVsTS4->setAxisTitle("Peak time slice", 2);
00142 
00143       hTriangleRightSlopeVsTS4 = dbe_->book2D("Triangle_fit_right_slope",
00144          "Triangle fit right distance vs. peak time slice", 50, 0, 10, 25, 0, 500);
00145       hTriangleRightSlopeVsTS4->setAxisTitle("Left slope", 1);
00146       hTriangleRightSlopeVsTS4->setAxisTitle("Peak time slice", 2);
00147 
00148       SetupEtaPhiHists(hFailLinearEtaPhi, "Fail_linear_Eta_Phi_Map", "");
00149       SetupEtaPhiHists(hFailRMSMaxEtaPhi, "Fail_RMS8Max_Eta_Phi_Map", "");
00150       SetupEtaPhiHists(hFailTriangleEtaPhi, "Fail_triangle_Eta_Phi_Map", "");
00151 
00152       // High-level isolation filter
00153       dbe_->setCurrentFolder(subdir_ + "IsolationVariable/");
00154 
00155       SetupEtaPhiHists(hFailIsolationEtaPhi, "Fail_isolation_Eta_Phi_Map", "");
00156       
00157       // TS4 vs. TS5 variable
00158       dbe_->setCurrentFolder(subdir_ + "TS4TS5Variable/");
00159       
00160       hTS4TS5RelativeDifference = dbe_->book1D("TS4_TS5_relative_difference",
00161          "(TS4-TS5)/(TS4+TS5), total charge > 20 fC", 100, -1, 1);
00162       hTS4TS5RelativeDifference->setAxisTitle("(TS4 - TS5) / (TS4 + TS5)", 1);
00163 
00164       hTS4TS5RelativeDifferenceVsCharge = dbe_->book2D("TS4_TS5_relative_difference_charge",
00165          "(TS4-TS5)/(TS4+TS5) vs. Charge", 25, 0, 400, 75, -1, 1);
00166       hTS4TS5RelativeDifferenceVsCharge->setAxisTitle("Charge", 1);
00167       hTS4TS5RelativeDifferenceVsCharge->setAxisTitle("(TS4 - TS5) / (TS4 + TS5)", 2);
00168 
00169       // Noise summary object
00170       dbe_->setCurrentFolder(subdir_ + "NoiseMonitoring/");
00171    
00172       hMaxZeros = dbe_->book1D("Max_Zeros", "Max zeros", 15, -0.5, 14.5);
00173       
00174       hTotalZeros = dbe_->book1D("Total_Zeros", "Total zeros", 15, -0.5, 14.5);
00175       
00176       hE2OverE10Digi = dbe_->book1D("E2OverE10Digi", "E2/E10 of the highest digi in an HPD", 100, 0, 2);
00177       
00178       hE2OverE10Digi5 = dbe_->book1D("E2OverE10Digi5", "E2/E10 of the highest 5 digi in an HPD", 100, 0, 2);
00179       
00180       hE2OverE10RBX = dbe_->book1D("E2OverE10RBX", "E2/E10 of RBX", 100, 0, 2);
00181       
00182       hHPDHitCount = dbe_->book1D("HPDHitCount", "HPD hit count (1.5 GeV)", 19, -0.5, 18.5);
00183       
00184       hRBXHitCount = dbe_->book1D("RBXHitCount", "Number of hits in RBX", 74, -0.5, 73.5);
00185 
00186       hHcalNoiseCategory = dbe_->book1D("Hcal_noise_category", "Hcal noise category", 10, 0.5, 10.5);
00187       hHcalNoiseCategory->setBinLabel(1, "RBX noise", 1);
00188       hHcalNoiseCategory->setBinLabel(2, "RBX pedestal flatter", 1);
00189       hHcalNoiseCategory->setBinLabel(3, "RBX pedestal sharper", 1);
00190       hHcalNoiseCategory->setBinLabel(4, "RBX flash large hit count", 1);
00191       hHcalNoiseCategory->setBinLabel(5, "RBX flash small hit count", 1);
00192       hHcalNoiseCategory->setBinLabel(7, "HPD discharge", 1);
00193       hHcalNoiseCategory->setBinLabel(8, "HPD ion feedback", 1);
00194 
00195       hBadZeroRBX = dbe_->book1D("BadZeroRBX", "RBX with bad ADC zero counts", 72, 0.5, 72.5);
00196       hBadCountHPD = dbe_->book1D("BadCountHPD", "HPD with bad hit counts", 72 * 4, 0.5, 72 * 4 + 0.5);
00197       hBadNoOtherCountHPD = dbe_->book1D("BadNoOtherCountHPD", "HPD with bad \"no other\" hit counts", 72 * 4, 0.5, 72 * 4 + 0.5);
00198       hBadE2E10RBX = dbe_->book1D("BadE2E10RBX", "RBX with bad E2/E10 value", 72, 0.5, 72.5);
00199    }
00200 
00201    ReadHcalPulse();
00202 
00203    return;
00204 }
00205 
00206 
00207 void HcalNoiseMonitor::analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup)
00208 {
00209    edm::Handle<HBHEDigiCollection> hHBHEDigis;
00210    iEvent.getByLabel(edm::InputTag(hbheDigiLabel_),hHBHEDigis);
00211 
00212    edm::ESHandle<HcalDbService> hConditions;
00213    iSetup.get<HcalDbRecord>().get(hConditions);
00214 
00215    edm::Handle<HBHERecHitCollection> hRecHits;
00216    iEvent.getByLabel(edm::InputTag(hbheRechitLabel_), hRecHits);
00217 
00218    edm::Handle<reco::HcalNoiseRBXCollection> hRBXCollection;
00219    iEvent.getByLabel(edm::InputTag(noiseLabel_),hRBXCollection);
00220 
00221    HcalBaseDQMonitor::analyze(iEvent, iSetup);
00222 
00223    if(dbe_ == NULL)
00224    {
00225       if(debug_ > 0)
00226          std::cout << "HcalNoiseMonitor::processEvent DQMStore not instantiated!!!"<< std::endl;
00227       return;
00228    }
00229 
00230    // loop over digis
00231    for(HBHEDigiCollection::const_iterator iter = hHBHEDigis->begin(); iter != hHBHEDigis->end(); iter++)
00232    {
00233       HcalDetId id = iter->id();
00234       const HcalCalibrations &Calibrations = hConditions->getHcalCalibrations(id);
00235       const HcalQIECoder *ChannelCoder = hConditions->getHcalCoder(id);
00236       const HcalQIEShape *Shape = hConditions->getHcalShape();
00237       HcalCoderDb Coder(*ChannelCoder, *Shape);
00238       CaloSamples Tool;
00239       Coder.adc2fC(*iter, Tool);
00240 
00241       // int ieta = id.ieta();
00242       // int iphi = id.iphi();
00243       // int depth = id.depth();
00244 
00245       double Charge[10] = {0};
00246       for(int i = 0; i < iter->size(); i++)
00247          Charge[i] = Tool[i] - Calibrations.pedestal(iter->sample(i).capid());
00248 
00249       double TotalCharge = 0;
00250       for(int i = 0; i < 10; i++)
00251          TotalCharge = TotalCharge + Charge[i];
00252 
00253       if(TotalCharge > 20)
00254       {
00255          double NominalChi2 = 10000000;
00256          NominalChi2 = PerformNominalFit(Charge);
00257          
00258          double LinearChi2 = PerformLinearFit(Charge);
00259          double RMS8Max = CalculateRMS8Max(Charge);
00260          TriangleFitResult TriangleResult = PerformTriangleFit(Charge);
00261       
00262          double TS4LeftSlope = 100000;
00263          double TS4RightSlope = 100000;
00264          
00265          if(TriangleResult.LeftSlope > 1e-5)
00266             TS4LeftSlope = Charge[4] / fabs(TriangleResult.LeftSlope);
00267          if(TriangleResult.RightSlope < -1e-5)
00268             TS4RightSlope = Charge[4] / fabs(TriangleResult.RightSlope);
00269 
00270          if(TS4LeftSlope < -1000 || TS4LeftSlope > 1000)
00271             TS4LeftSlope = 1000;
00272          if(TS4RightSlope < -1000 || TS4RightSlope > 1000)
00273             TS4RightSlope = 1000;
00274 
00275          hNominalChi2->Fill(NominalChi2);
00276          hLinearChi2->Fill(LinearChi2);
00277          hLinearTestStatistics->Fill(log(LinearChi2) - log(NominalChi2));
00278          hRMS8OverMax->Fill(RMS8Max);
00279          hRMS8OverMaxTestStatistics->Fill(log(RMS8Max) - log(NominalChi2));
00280 
00281          hLambdaLinearVsTotalCharge->Fill(log(LinearChi2) - log(NominalChi2), TotalCharge);
00282          hLambdaRMS8MaxVsTotalCharge->Fill(log(RMS8Max) - log(NominalChi2), TotalCharge);
00283          hTriangleLeftSlopeVsTS4->Fill(TS4LeftSlope, Charge[4]);
00284          hTriangleRightSlopeVsTS4->Fill(TS4RightSlope, Charge[4]);
00285       }
00286 
00287       if(Charge[4] + Charge[5] > 1e-5)
00288       {
00289          hTS4TS5RelativeDifference->Fill((Charge[4] - Charge[5]) / (Charge[4] + Charge[5]));
00290          hTS4TS5RelativeDifferenceVsCharge->Fill(TotalCharge, (Charge[4] - Charge[5]) / (Charge[4] + Charge[5]));
00291       }
00292    }
00293 
00294    // loop over rechits - noise bits (fit-based, isolation)
00295    for(HBHERecHitCollection::const_iterator iter = hRecHits->begin(); iter != hRecHits->end(); iter++)
00296    {
00297       HcalDetId id = iter->id();
00298 
00299       int ieta = id.ieta();
00300       int iphi = id.iphi();
00301       int depth = id.depth();
00302 
00303       if(iter->flagField(HcalCaloFlagLabels::HBHEFlatNoise) == 1)
00304          hFailLinearEtaPhi.depth[depth-1]->Fill(ieta, iphi);
00305 
00306       if(iter->flagField(HcalCaloFlagLabels::HBHESpikeNoise) == 1)
00307          hFailRMSMaxEtaPhi.depth[depth-1]->Fill(ieta, iphi);
00308       
00309       if(iter->flagField(HcalCaloFlagLabels::HBHETriangleNoise) == 1)
00310          hFailTriangleEtaPhi.depth[depth-1]->Fill(ieta, iphi);
00311 
00312       if(iter->flagField(HcalCaloFlagLabels::HBHEIsolatedNoise) == 1)
00313          hFailIsolationEtaPhi.depth[depth-1]->Fill(ieta, iphi);
00314    }
00315 
00316    // Code analagous to Yifei's
00317    for(reco::HcalNoiseRBXCollection::const_iterator rbx = hRBXCollection->begin();
00318       rbx != hRBXCollection->end(); rbx++)
00319    {
00320       const reco::HcalNoiseRBX RBX = *rbx;
00321 
00322       int NumberRBXHits = RBX.numRecHits(1.5);
00323       double RBXEnergy = RBX.recHitEnergy(1.5);
00324       double RBXE2 = RBX.allChargeHighest2TS();
00325       double RBXE10 = RBX.allChargeTotal();
00326 
00327       std::vector<reco::HcalNoiseHPD> HPDs = RBX.HPDs();
00328       
00329       int RBXID = RBX.idnumber();
00330 
00331       if(RBXEnergy > mTotalZeroMinEnergy && RBX.totalZeros() >= mMaxADCZeros)
00332          hBadZeroRBX->Fill(RBXID);
00333       if(RBXEnergy > mE2E10MinEnergy && RBXE10 > 1e-5 && (RBXE2 / RBXE10 > mMaxE2E10 || RBXE2 / RBXE10 < mMinE2E10))
00334          hBadE2E10RBX->Fill(RBXID);
00335       for(std::vector<reco::HcalNoiseHPD>::const_iterator hpd = HPDs.begin(); hpd != HPDs.end(); hpd++)
00336       {
00337          reco::HcalNoiseHPD HPD = *hpd;
00338          int HPDHitCount = HPD.numRecHits(1.5);
00339          if(HPDHitCount >= mMaxHPDHitCount)
00340             hBadCountHPD->Fill(HPD.idnumber());
00341          if(HPDHitCount == NumberRBXHits && HPDHitCount >= mMaxHPDNoOtherHitCount)
00342             hBadNoOtherCountHPD->Fill(HPD.idnumber());
00343       }
00344       
00345       if(NumberRBXHits == 0 || RBXEnergy <= 10)
00346          continue;
00347 
00348       hRBXHitCount->Fill(NumberRBXHits);
00349 
00350       hMaxZeros->Fill(RBX.maxZeros());
00351       hTotalZeros->Fill(RBX.totalZeros());
00352    
00353       double HighestHPDEnergy = 0;
00354       int HighestHPDHits = 0;
00355 
00356       for(std::vector<reco::HcalNoiseHPD>::const_iterator hpd = HPDs.begin(); hpd != HPDs.end(); hpd++)
00357       {
00358          reco::HcalNoiseHPD HPD = *hpd;
00359 
00360          if(HPD.recHitEnergy(1.5) > HighestHPDEnergy)
00361          {
00362             HighestHPDEnergy = HPD.recHitEnergy(1.5);
00363             HighestHPDHits = HPD.numRecHits(1.5);
00364          }
00365 
00366          if(HPD.numRecHits(5) < 1)
00367             continue;
00368 
00369          if(HPD.bigChargeTotal() > 1e-5)
00370             hE2OverE10Digi->Fill(HPD.bigChargeHighest2TS() / HPD.bigChargeTotal());
00371          if(HPD.big5ChargeTotal() > 1e-5)
00372             hE2OverE10Digi5->Fill(HPD.big5ChargeHighest2TS() / HPD.big5ChargeTotal());
00373 
00374          hHPDHitCount->Fill(HPD.numRecHits(1.5));
00375       }
00376 
00377       int NoiseCategory = 0;
00378       bool IsRBXNoise = false;
00379       //      bool IsHPDNoise = false;
00380       //      bool IsHPDIonFeedback = false;
00381       //      bool IsHPDDischarge = false;
00382 
00383       if(RBXEnergy > 1e-5 && HighestHPDEnergy / RBXEnergy > 0.98)
00384       {
00385         //         IsHPDNoise = true;
00386 
00387          if(HighestHPDHits >= 9)
00388          {
00389            //            IsHPDDischarge = true;
00390             NoiseCategory = 7;
00391          }
00392          else
00393          {
00394            //            IsHPDIonFeedback = true;
00395             NoiseCategory = 8;
00396          }
00397       }
00398       else
00399       {
00400          IsRBXNoise = true;
00401          NoiseCategory = 1;
00402 
00403          if(RBXE10 > 1e-5)
00404          {
00405             if(RBXE2 / RBXE10 < 0.33)
00406                NoiseCategory = 2;
00407             else if(RBXE2 / RBXE10 < 0.8)
00408                NoiseCategory = 3;
00409             else if(RBXE2 / RBXE10 > 0.8 && NumberRBXHits > 10)
00410                NoiseCategory = 4;
00411             else if(RBXE2 / RBXE10 > 0.8 && NumberRBXHits < 10)  // [hic]
00412                NoiseCategory = 5;
00413          }
00414       }
00415 
00416       hHcalNoiseCategory->Fill(NoiseCategory);
00417 
00418       if(IsRBXNoise == true && RBXE10 > 1e-5)
00419          hE2OverE10RBX->Fill(RBXE2 / RBXE10);
00420    }
00421 
00422    return;
00423 }
00424 
00425 double HcalNoiseMonitor::PerformNominalFit(double Charge[10])
00426 {
00427    //
00428    // Performs a fit to the ideal pulse shape.  Returns best chi2
00429    //
00430    // A scan over different timing offset (for the ideal pulse) is carried out,
00431    //    and for each offset setting a one-parameter fit is performed
00432    //
00433 
00434    int DigiSize = 10;
00435 
00436    double MinimumChi2 = 100000;
00437 
00438    double F = 0;
00439 
00440    double SumF2 = 0;
00441    double SumTF = 0;
00442    double SumT2 = 0;
00443 
00444    for(int i = 0; i + 250 < (int)CumulativeIdealPulse.size(); i++)
00445    {
00446       if(CumulativeIdealPulse[i+250] - CumulativeIdealPulse[i] < 1e-5)
00447          continue;
00448 
00449       SumF2 = 0;
00450       SumTF = 0;
00451       SumT2 = 0;
00452 
00453       for(int j = 0; j < DigiSize; j++)
00454       {
00455          // get ideal pulse component for this time slice....
00456          F = CumulativeIdealPulse[i+j*25+25] - CumulativeIdealPulse[i+j*25];
00457 
00458          double Error2 = Charge[j];
00459          if(Error2 < 1)
00460             Error2 = 1;
00461 
00462          // ...and increment various summations
00463          SumF2 += F * F / Error2;
00464          SumTF += F * Charge[j] / Error2;
00465          SumT2 += Charge[j] * Charge[j] / Error2;
00466       }
00467 
00468       /* chi2= sum((Charge[j]-aF)^2/|Charge[j]|
00469          ( |Charge[j]| = assumed sigma^2 for Charge[j]; a bit wonky for Charge[j]<1 )
00470          chi2 = sum(|Charge[j]|) - 2*sum(aF*Charge[j]/|Charge[j]|) +sum( a^2*F^2/|Charge[j]|)
00471          chi2 minimimized when d(chi2)/da = 0:
00472          a = sum(F*Charge[j])/sum(F^2)
00473          ...
00474          chi2= sum(|Q[j]|) - sum(Q[j]/|Q[j]|*F)*sum(Q[j]/|Q[j]|*F)/sum(F^2/|Q[j]|), where Q = Charge
00475          chi2 = SumT2 - SumTF*SumTF/SumF2
00476       */
00477       
00478       double Chi2 = SumT2 - SumTF * SumTF / SumF2;
00479 
00480       if(Chi2 < MinimumChi2)
00481          MinimumChi2 = Chi2;
00482    }
00483 
00484    // safety protection in case of perfect fit - don't want the log(...) to explode
00485    if(MinimumChi2 < 1e-5)
00486       MinimumChi2 = 1e-5;
00487 
00488    return MinimumChi2;
00489 }
00490 
00491 double HcalNoiseMonitor::PerformDualNominalFit(double Charge[10])
00492 {
00493    //
00494    // Perform dual nominal fit and returns the chi2
00495    //
00496    // In this function we do a scan over possible "distance" (number of time slices between two components)
00497    //    and overall offset for the two components; first coarse, then finer
00498    // All the fitting is done in the DualNominalFitSingleTry function
00499    //                  
00500    
00501    double OverallMinimumChi2 = 1000000;
00502 
00503    int AvailableDistance[] = {-100, -75, -50, 50, 75, 100};
00504 
00505    // loop over possible pulse distances between two components
00506    for(int k = 0; k < 6; k++)
00507    {
00508       double SingleMinimumChi2 = 1000000;
00509       int MinOffset = 0;
00510 
00511       // scan coarsely through different offsets and find the minimum
00512       for(int i = 0; i + 250 < (int)CumulativeIdealPulse.size(); i += 10)
00513       {
00514          double Chi2 = DualNominalFitSingleTry(Charge, i, AvailableDistance[k]);
00515 
00516          if(Chi2 < SingleMinimumChi2)
00517          {
00518             SingleMinimumChi2 = Chi2;
00519             MinOffset = i;
00520          }
00521       }
00522 
00523       // around the minimum, scan finer for better a better minimum
00524       for(int i = MinOffset - 15; i + 250 < (int)CumulativeIdealPulse.size() && i < MinOffset + 15; i++)
00525       {
00526          double Chi2 = DualNominalFitSingleTry(Charge, i, AvailableDistance[k]);
00527          if(Chi2 < SingleMinimumChi2)
00528             SingleMinimumChi2 = Chi2;
00529       }
00530 
00531       // update overall minimum chi2
00532       if(SingleMinimumChi2 < OverallMinimumChi2)
00533          OverallMinimumChi2 = SingleMinimumChi2;
00534    }
00535 
00536    return OverallMinimumChi2;
00537 }
00538 
00539 double HcalNoiseMonitor::DualNominalFitSingleTry(double Charge[10], int Offset, int Distance)
00540 {
00541    //
00542    // Does a fit to dual signal pulse hypothesis given offset and distance of the two target pulses
00543    //
00544    // The only parameters to fit here are the two pulse heights of in-time and out-of-time components
00545    //    since offset is given
00546    // The calculation here is based from writing down the Chi2 formula and minimize against the two parameters,
00547    //    ie., Chi2 = Sum{((T[i] - a1 * F1[i] - a2 * F2[i]) / (Sigma[i]))^2}, where T[i] is the input pulse shape,
00548    //    and F1[i], F2[i] are the two ideal pulse components
00549    //
00550 
00551    int DigiSize = 10;
00552 
00553    if(Offset < 0 || Offset + 250 >= (int)CumulativeIdealPulse.size())
00554       return 1000000;
00555    if(CumulativeIdealPulse[Offset+250] - CumulativeIdealPulse[Offset] < 1e-5)
00556       return 1000000;
00557 
00558    static std::vector<double> F1;
00559    static std::vector<double> F2;
00560 
00561    F1.resize(DigiSize);
00562    F2.resize(DigiSize);
00563 
00564    double SumF1F1 = 0;
00565    double SumF1F2 = 0;
00566    double SumF2F2 = 0;
00567    double SumTF1 = 0;
00568    double SumTF2 = 0;
00569 
00570    double Error = 0;
00571 
00572    for(int j = 0; j < DigiSize; j++)
00573    {
00574       // this is the TS value for in-time component - no problem we can do a subtraction directly
00575       F1[j] = CumulativeIdealPulse[Offset+j*25+25] - CumulativeIdealPulse[Offset+j*25];
00576 
00577       // However for the out-of-time component the index might go out-of-bound.
00578       // Let's protect against this.
00579 
00580       int OffsetTemp = Offset + j * 25 + Distance;
00581       
00582       double C1 = 0;   // lower-indexed value in the cumulative pulse shape
00583       double C2 = 0;   // higher-indexed value in the cumulative pulse shape
00584       
00585       if(OffsetTemp + 25 < (int)CumulativeIdealPulse.size() && OffsetTemp + 25 >= 0)
00586          C1 = CumulativeIdealPulse[OffsetTemp+25];
00587       if(OffsetTemp + 25 >= (int)CumulativeIdealPulse.size())
00588          C1 = CumulativeIdealPulse[CumulativeIdealPulse.size()-1];
00589       if(OffsetTemp < (int)CumulativeIdealPulse.size() && OffsetTemp >= 0)
00590          C2 = CumulativeIdealPulse[OffsetTemp];
00591       if(OffsetTemp >= (int)CumulativeIdealPulse.size())
00592          C2 = CumulativeIdealPulse[CumulativeIdealPulse.size()-1];
00593       F2[j] = C1 - C2;
00594 
00595       Error = Charge[j];
00596       if(Error < 1)
00597          Error = 1;
00598 
00599       SumF1F1 += F1[j] * F1[j] / Error;
00600       SumF1F2 += F1[j] * F2[j] / Error; 
00601       SumF2F2 += F2[j] * F2[j] / Error;
00602       SumTF1  += F1[j] * Charge[j] / Error; 
00603       SumTF2  += F2[j] * Charge[j] / Error; 
00604    }
00605 
00606    double Height = (SumF1F2 * SumTF2 - SumF2F2 * SumTF1) / (SumF1F2 * SumF1F2 - SumF1F1 * SumF2F2);
00607    double Height2 = (SumF1F2 * SumTF1 - SumF1F1 * SumTF2) / (SumF1F2 * SumF1F2 - SumF1F1 * SumF2F2);
00608 
00609    double Chi2 = 0;
00610    for(int j = 0; j < DigiSize; j++)
00611    {
00612       double Error = Charge[j];
00613       if(Error < 1)
00614          Error = 1;
00615 
00616       double Residual = Height * F1[j] + Height2 * F2[j] - Charge[j];  
00617       Chi2 += Residual * Residual / Error;                             
00618    } 
00619 
00620    // Safety protection in case of zero
00621    if(Chi2 < 1e-5)
00622       Chi2 = 1e-5;
00623 
00624    return Chi2;
00625 }
00626 
00627 double HcalNoiseMonitor::PerformLinearFit(double Charge[10])
00628 {
00629    //
00630    // Performs a straight-line fit over all time slices, and returns the chi2 value
00631    //
00632    // The calculation here is based from writing down the formula for chi2 and minimize
00633    //    with respect to the parameters in the fit, ie., slope and intercept of the straight line
00634    // By doing two differentiation, we will get two equations, and the best parameters are determined by these
00635    //
00636 
00637    int DigiSize = 10;
00638 
00639    double SumTS2OverTi = 0;
00640    double SumTSOverTi = 0;
00641    double SumOverTi = 0;
00642    double SumTiTS = 0;
00643    double SumTi = 0;
00644 
00645    double Error2 = 0;
00646    for(int i = 0; i < DigiSize; i++)
00647    {
00648       Error2 = Charge[i];
00649       if(Charge[i] < 1)
00650          Error2 = 1;
00651 
00652       SumTS2OverTi += 1.* i * i / Error2;
00653       SumTSOverTi  += 1.* i / Error2;
00654       SumOverTi    += 1. / Error2;
00655       SumTiTS      += Charge[i] * i / Error2;
00656       SumTi        += Charge[i] / Error2;
00657    }
00658 
00659    double CM1 = SumTS2OverTi;   // Coefficient in front of slope in equation 1
00660    double CM2 = SumTSOverTi;   // Coefficient in front of slope in equation 2
00661    double CD1 = SumTSOverTi;   // Coefficient in front of intercept in equation 1
00662    double CD2 = SumOverTi;   // Coefficient in front of intercept in equation 2
00663    double C1 = SumTiTS;   // Constant coefficient in equation 1
00664    double C2 = SumTi;   // Constant coefficient in equation 2
00665 
00666    double Slope = (C1 * CD2 - C2 * CD1) / (CM1 * CD2 - CM2 * CD1);
00667    double Intercept = (C1 * CM2 - C2 * CM1) / (CD1 * CM2 - CD2 * CM1);
00668 
00669    // now that the best parameters are found, calculate chi2 from those
00670    double Chi2 = 0;
00671    for(int i = 0; i < DigiSize; i++)
00672    {
00673       double Deviation = Slope * i + Intercept - Charge[i];
00674       double Error2 = Charge[i];
00675       if(Charge[i] < 1)
00676          Error2 = 1;
00677       Chi2 += Deviation * Deviation / Error2;  
00678    }
00679 
00680    // safety protection in case of perfect fit
00681    if(Chi2 < 1e-5)
00682       Chi2 = 1e-5;
00683 
00684    return Chi2;
00685 }
00686 
00687 double HcalNoiseMonitor::CalculateRMS8Max(double Charge[10])
00688 {
00689    //
00690    // CalculateRMS8Max
00691    //
00692    // returns "RMS" divided by the largest charge in the time slices
00693    //    "RMS" is calculated using all but the two largest time slices.
00694    //    The "RMS" is not quite the actual RMS (see note below), but the value is only
00695    //    used for determining max values, and is not quoted as the actual RMS anywhere.
00696    //
00697 
00698    int DigiSize = 10;
00699 
00700    // Copy Charge vector again, since we are passing references around
00701    std::vector<double> TempCharge(Charge, Charge + 10);
00702 
00703    // Sort TempCharge vector from smallest to largest charge
00704    sort(TempCharge.begin(), TempCharge.end());
00705 
00706    double Total = 0;
00707    double Total2 = 0;
00708    for(int i = 0; i < DigiSize - 2; i++)
00709    {
00710       Total = Total + TempCharge[i];
00711       Total2 = Total2 + TempCharge[i] * TempCharge[i];
00712    }
00713 
00714    // This isn't quite the RMS (both Total2 and Total*Total need to be
00715    // divided by an extra (DigiSize-2) within the sqrt to get the RMS.)
00716    // We're only using this value for relative comparisons, though; we
00717    // aren't explicitly interpreting it as the RMS.  It might be nice
00718    // to either change the calculation or rename the variable in the future, though.
00719 
00720    double RMS = sqrt(Total2 - Total * Total / (DigiSize - 2));
00721 
00722    double RMS8Max = 99999;
00723    if(TempCharge[DigiSize-1] > 1e-5)
00724       RMS8Max = RMS / TempCharge[DigiSize-1];
00725    if(RMS8Max < 1e-5)   // protection against zero
00726       RMS8Max = 1e-5;
00727 
00728    return RMS8Max;
00729 }
00730 
00731 TriangleFitResult HcalNoiseMonitor::PerformTriangleFit(double Charge[10])
00732 {
00733    //
00734    // Perform a "triangle fit", and extract the slopes
00735    //
00736    // Left-hand side and right-hand side are not correlated to each other - do them separately
00737    //
00738 
00739    TriangleFitResult result;
00740    result.Chi2 = 0;
00741    result.LeftSlope = 0;
00742    result.RightSlope = 0;
00743 
00744    int DigiSize = 10;
00745 
00746    // right side, starting from TS4
00747    double MinimumRightChi2 = 1000000;
00748    double Numerator = 0;
00749    double Denominator = 0;
00750 
00751    for(int iTS = mTrianglePeakTS + 2; iTS <= DigiSize; iTS++)   // the place where first TS center in flat line
00752    {
00753       // fit a straight line to the triangle part
00754       Numerator = 0;
00755       Denominator = 0;
00756 
00757       for(int i = mTrianglePeakTS + 1; i < iTS; i++)
00758       {
00759          Numerator += (i - mTrianglePeakTS) * (Charge[i] - Charge[mTrianglePeakTS]);
00760          Denominator += (i - mTrianglePeakTS) * (i - mTrianglePeakTS);
00761       }
00762 
00763       double BestSlope = Numerator / Denominator;
00764       if(BestSlope > 0)
00765          BestSlope = 0;
00766 
00767       // check if the slope is reasonable
00768       if(iTS != DigiSize)
00769       {
00770          if(BestSlope > -1 * Charge[mTrianglePeakTS] / (iTS - mTrianglePeakTS))
00771             BestSlope = -1 * Charge[mTrianglePeakTS] / (iTS - mTrianglePeakTS);
00772          if(BestSlope < -1 * Charge[mTrianglePeakTS] / (iTS - 1 - mTrianglePeakTS))
00773             BestSlope = -1 * Charge[mTrianglePeakTS] / (iTS - 1 - mTrianglePeakTS);
00774       }
00775       else
00776       {
00777          if(BestSlope < -1 * Charge[mTrianglePeakTS] / (iTS - 1 - mTrianglePeakTS)) 
00778             BestSlope = -1 * Charge[mTrianglePeakTS] / (iTS - 1 - mTrianglePeakTS);
00779       }
00780 
00781       // calculate partial chi2
00782 
00783       // The shape I'm fitting is more like a tent than a triangle.
00784       // After the end of triangle edge assume a flat line
00785 
00786       double Chi2 = 0;
00787       for(int i = mTrianglePeakTS + 1; i < iTS; i++)
00788          Chi2 += (Charge[mTrianglePeakTS] - Charge[i] + (i - mTrianglePeakTS) * BestSlope)
00789             * (Charge[mTrianglePeakTS] - Charge[i] + (i - mTrianglePeakTS) * BestSlope);
00790       for(int i = iTS; i < DigiSize; i++)    // Assumes fit line = 0 for iTS > fit
00791          Chi2 += Charge[i] * Charge[i];
00792 
00793       if(Chi2 < MinimumRightChi2)
00794       {
00795          MinimumRightChi2 = Chi2;
00796          result.RightSlope = BestSlope;
00797       }
00798    }   // end of right-hand side loop
00799 
00800    // left side
00801    double MinimumLeftChi2 = 1000000;
00802 
00803    for(int iTS = 0; iTS < (int)mTrianglePeakTS; iTS++)   // the first time after linear fit ends
00804    {
00805       // fit a straight line to the triangle part
00806       Numerator = 0;
00807       Denominator = 0;
00808       for(int i = iTS; i < (int)mTrianglePeakTS; i++)
00809       {
00810          Numerator = Numerator + (i - mTrianglePeakTS) * (Charge[i] - Charge[mTrianglePeakTS]);
00811          Denominator = Denominator + (i - mTrianglePeakTS) * (i - mTrianglePeakTS);
00812       }
00813 
00814       double BestSlope = Numerator / Denominator;
00815       if(BestSlope < 0)
00816          BestSlope = 0;
00817 
00818       // check slope
00819       if(iTS != 0)
00820       {
00821          if(BestSlope > Charge[mTrianglePeakTS] / (mTrianglePeakTS - iTS))
00822             BestSlope = Charge[mTrianglePeakTS] / (mTrianglePeakTS - iTS);
00823          if(BestSlope < Charge[mTrianglePeakTS] / (mTrianglePeakTS + 1 - iTS))
00824             BestSlope = Charge[mTrianglePeakTS] / (mTrianglePeakTS + 1 - iTS);
00825       }
00826       else
00827       {
00828          if(BestSlope > Charge[mTrianglePeakTS] / (mTrianglePeakTS - iTS))
00829             BestSlope = Charge[mTrianglePeakTS] / (mTrianglePeakTS - iTS);
00830       }
00831 
00832       // calculate minimum chi2
00833       double Chi2 = 0;
00834       for(int i = 0; i < iTS; i++)
00835          Chi2 += Charge[i] * Charge[i];
00836       for(int i = iTS; i < (int)mTrianglePeakTS; i++)
00837          Chi2 += (Charge[mTrianglePeakTS] - Charge[i] + (i - mTrianglePeakTS) * BestSlope)
00838             * (Charge[mTrianglePeakTS] - Charge[i] + (i - mTrianglePeakTS) * BestSlope);
00839 
00840       if(MinimumLeftChi2 > Chi2)
00841       {
00842          MinimumLeftChi2 = Chi2;
00843          result.LeftSlope = BestSlope;
00844       }
00845    }   // end of left-hand side loop
00846 
00847    result.Chi2 = MinimumLeftChi2 + MinimumRightChi2;
00848 
00849    return result;
00850 }
00851 
00852 void HcalNoiseMonitor::ReadHcalPulse()
00853 {
00854    std::vector<double> PulseShape;
00855 
00856    HcalPulseShapes Shapes;
00857    HcalPulseShapes::Shape HPDShape = Shapes.hbShape();
00858 
00859    PulseShape.reserve(350);
00860    for(int i = 0; i < 200; i++)
00861       PulseShape.push_back(HPDShape.at(i));
00862    PulseShape.insert(PulseShape.begin(), 150, 0);   // Safety margin of a lot of zeros in the beginning
00863 
00864    CumulativeIdealPulse.reserve(350);
00865    CumulativeIdealPulse.clear();
00866    CumulativeIdealPulse.push_back(0);
00867    for(unsigned int i = 1; i < PulseShape.size(); i++)
00868       CumulativeIdealPulse.push_back(CumulativeIdealPulse[i-1] + PulseShape[i]);
00869 }
00870 
00871 DEFINE_FWK_MODULE(HcalNoiseMonitor);
00872