Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:32:07

0001 /*
0002  * \file EcalBarrelRecHitsValidation.cc
0003  *
0004  * \author C. Rovelli
0005  *
0006  */
0007 
0008 #include "DQMServices/Core/interface/DQMStore.h"
0009 #include <DataFormats/EcalDetId/interface/EBDetId.h>
0010 #include <Validation/EcalRecHits/interface/EcalBarrelRecHitsValidation.h>
0011 
0012 using namespace cms;
0013 using namespace edm;
0014 using namespace std;
0015 
0016 EcalBarrelRecHitsValidation::EcalBarrelRecHitsValidation(const ParameterSet &ps) : ecalPeds(esConsumes()) {
0017   // ----------------------
0018   EBdigiCollection_token_ = consumes<EBDigiCollection>(ps.getParameter<edm::InputTag>("EBdigiCollection"));
0019   EBuncalibrechitCollection_token_ =
0020       consumes<EBUncalibratedRecHitCollection>(ps.getParameter<edm::InputTag>("EBuncalibrechitCollection"));
0021 
0022   // ----------------------
0023   // verbosity switch
0024   verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
0025 
0026   // ----------------------
0027   meEBUncalibRecHitsOccupancy_ = nullptr;
0028   meEBUncalibRecHitsAmplitude_ = nullptr;
0029   meEBUncalibRecHitsPedestal_ = nullptr;
0030   meEBUncalibRecHitsJitter_ = nullptr;
0031   meEBUncalibRecHitsChi2_ = nullptr;
0032   meEBUncalibRecHitMaxSampleRatio_ = nullptr;
0033   meEBUncalibRecHitsOccupancyGt100adc_ = nullptr;
0034   meEBUncalibRecHitsAmplitudeGt100adc_ = nullptr;
0035   meEBUncalibRecHitsPedestalGt100adc_ = nullptr;
0036   meEBUncalibRecHitsJitterGt100adc_ = nullptr;
0037   meEBUncalibRecHitsChi2Gt100adc_ = nullptr;
0038   meEBUncalibRecHitMaxSampleRatioGt100adc_ = nullptr;
0039   meEBUncalibRecHitsAmpFullMap_ = nullptr;
0040   meEBUncalibRecHitsPedFullMap_ = nullptr;
0041   for (int i = 0; i < 36; i++) {
0042     meEBUncalibRecHitAmplMap_[i] = nullptr;
0043     meEBUncalibRecHitPedMap_[i] = nullptr;
0044   }
0045 }
0046 
0047 EcalBarrelRecHitsValidation::~EcalBarrelRecHitsValidation() {}
0048 
0049 void EcalBarrelRecHitsValidation::bookHistograms(DQMStore::IBooker &ibooker,
0050                                                  edm::Run const &,
0051                                                  edm::EventSetup const &) {
0052   Char_t histo[200];
0053 
0054   ibooker.setCurrentFolder("EcalRecHitsV/EcalBarrelRecHitsTask");
0055 
0056   sprintf(histo, "EB Occupancy");
0057   meEBUncalibRecHitsOccupancy_ = ibooker.book2D(histo, histo, 170, -85., 85., 360, 0., 360.);
0058 
0059   sprintf(histo, "EB Amplitude");
0060   meEBUncalibRecHitsAmplitude_ = ibooker.book1D(histo, histo, 201, -20., 4000.);
0061 
0062   sprintf(histo, "EB Pedestal");
0063   meEBUncalibRecHitsPedestal_ = ibooker.book1D(histo, histo, 50, 190., 210.);
0064 
0065   sprintf(histo, "EB Jitter");
0066   meEBUncalibRecHitsJitter_ = ibooker.book1D(histo, histo, 100, 0., 100.);
0067 
0068   sprintf(histo, "EB Chi2");
0069   meEBUncalibRecHitsChi2_ = ibooker.book1D(histo, histo, 100, 18000., 22000.);
0070 
0071   sprintf(histo, "EB RecHit Max Sample Ratio");
0072   meEBUncalibRecHitMaxSampleRatio_ = ibooker.book1D(histo, histo, 120, 0.90, 1.05);
0073 
0074   sprintf(histo, "EB Occupancy gt 100 adc counts");
0075   meEBUncalibRecHitsOccupancyGt100adc_ = ibooker.book2D(histo, histo, 170, -85., 85., 360, 0., 360.);
0076 
0077   sprintf(histo, "EB Amplitude gt 100 adc counts");
0078   meEBUncalibRecHitsAmplitudeGt100adc_ = ibooker.book1D(histo, histo, 200, 0., 4000.);
0079 
0080   sprintf(histo, "EB Pedestal gt 100 adc counts");
0081   meEBUncalibRecHitsPedestalGt100adc_ = ibooker.book1D(histo, histo, 50, 190., 210.);
0082 
0083   sprintf(histo, "EB Jitter gt 100 adc counts");
0084   meEBUncalibRecHitsJitterGt100adc_ = ibooker.book1D(histo, histo, 100, 0., 100.);
0085 
0086   sprintf(histo, "EB Chi2 gt 100 adc counts");
0087   meEBUncalibRecHitsChi2Gt100adc_ = ibooker.book1D(histo, histo, 100, 18000., 22000.);
0088 
0089   sprintf(histo, "EB RecHit Max Sample Ratio gt 100 adc counts");
0090   meEBUncalibRecHitMaxSampleRatioGt100adc_ = ibooker.book1D(histo, histo, 120, 0.90, 1.05);
0091 
0092   sprintf(histo, "EB Amplitude Full Map");
0093   meEBUncalibRecHitsAmpFullMap_ = ibooker.bookProfile2D(histo, histo, 170, -85., 85., 360, 0., 360., 200, 0., 4000.);
0094 
0095   sprintf(histo, "EB Pedestal Full Map");
0096   meEBUncalibRecHitsPedFullMap_ = ibooker.bookProfile2D(histo, histo, 170, -85., 85., 360, 0., 360., 50, 194., 201.);
0097 
0098   for (int i = 0; i < 36; i++) {
0099     sprintf(histo, "EB Amp SM%02d", i + 1);
0100     meEBUncalibRecHitAmplMap_[i] = ibooker.bookProfile2D(histo, histo, 85, 0., 85., 20, 0., 20., 200, 0., 4000.);
0101 
0102     sprintf(histo, "EB Ped SM%02d", i + 1);
0103     meEBUncalibRecHitPedMap_[i] = ibooker.bookProfile2D(histo, histo, 85, 0., 85., 20, 0., 20., 50, 194., 201.);
0104   }
0105 }
0106 
0107 void EcalBarrelRecHitsValidation::analyze(const Event &e, const EventSetup &c) {
0108   const EBUncalibratedRecHitCollection *EBUncalibRecHit = nullptr;
0109   Handle<EBUncalibratedRecHitCollection> EcalUncalibRecHitEB;
0110   e.getByToken(EBuncalibrechitCollection_token_, EcalUncalibRecHitEB);
0111   if (EcalUncalibRecHitEB.isValid()) {
0112     EBUncalibRecHit = EcalUncalibRecHitEB.product();
0113   } else {
0114     return;
0115   }
0116 
0117   bool skipDigis = false;
0118   const EBDigiCollection *EBDigi = nullptr;
0119   Handle<EBDigiCollection> EcalDigiEB;
0120   e.getByToken(EBdigiCollection_token_, EcalDigiEB);
0121   if (EcalDigiEB.isValid()) {
0122     EBDigi = EcalDigiEB.product();
0123   } else {
0124     skipDigis = true;
0125   }
0126 
0127   // ----------------------
0128   // loop over UncalibRecHits
0129   for (EcalUncalibratedRecHitCollection::const_iterator uncalibRecHit = EBUncalibRecHit->begin();
0130        uncalibRecHit != EBUncalibRecHit->end();
0131        ++uncalibRecHit) {
0132     EBDetId EBid = EBDetId(uncalibRecHit->id());
0133 
0134     // general checks
0135     if (meEBUncalibRecHitsOccupancy_)
0136       meEBUncalibRecHitsOccupancy_->Fill(EBid.ieta(), EBid.iphi());
0137     if (meEBUncalibRecHitsAmplitude_)
0138       meEBUncalibRecHitsAmplitude_->Fill(uncalibRecHit->amplitude());
0139     if (meEBUncalibRecHitsPedestal_)
0140       meEBUncalibRecHitsPedestal_->Fill(uncalibRecHit->pedestal());
0141     if (meEBUncalibRecHitsJitter_)
0142       meEBUncalibRecHitsJitter_->Fill(uncalibRecHit->jitter());
0143     if (meEBUncalibRecHitsChi2_)
0144       meEBUncalibRecHitsChi2_->Fill(uncalibRecHit->chi2());
0145     if (meEBUncalibRecHitsAmpFullMap_)
0146       meEBUncalibRecHitsAmpFullMap_->Fill(EBid.ieta(), EBid.iphi(), uncalibRecHit->amplitude());
0147     if (meEBUncalibRecHitsPedFullMap_)
0148       meEBUncalibRecHitsPedFullMap_->Fill(EBid.ieta(), EBid.iphi(), uncalibRecHit->pedestal());
0149 
0150     // general checks, with threshold at 3.5 GeV = 100 ADC counts
0151     if (uncalibRecHit->amplitude() > 100) {
0152       if (meEBUncalibRecHitsOccupancyGt100adc_)
0153         meEBUncalibRecHitsOccupancyGt100adc_->Fill(EBid.ieta(), EBid.iphi());
0154       if (meEBUncalibRecHitsAmplitudeGt100adc_)
0155         meEBUncalibRecHitsAmplitudeGt100adc_->Fill(uncalibRecHit->amplitude());
0156       if (meEBUncalibRecHitsPedestalGt100adc_)
0157         meEBUncalibRecHitsPedestalGt100adc_->Fill(uncalibRecHit->pedestal());
0158       if (meEBUncalibRecHitsJitterGt100adc_)
0159         meEBUncalibRecHitsJitterGt100adc_->Fill(uncalibRecHit->jitter());
0160       if (meEBUncalibRecHitsChi2Gt100adc_)
0161         meEBUncalibRecHitsChi2Gt100adc_->Fill(uncalibRecHit->chi2());
0162     }
0163 
0164     // supermodule maps
0165     int ic = EBid.ic();
0166     int ie = (ic - 1) / 20 + 1;
0167     int ip = (ic - 1) % 20 + 1;
0168     int ism = EBid.ism();
0169     float xie = ie - 0.5;
0170     float xip = ip - 0.5;
0171     if (meEBUncalibRecHitPedMap_[ism - 1])
0172       meEBUncalibRecHitPedMap_[ism - 1]->Fill(xie, xip, uncalibRecHit->pedestal());
0173     if (meEBUncalibRecHitAmplMap_[ism - 1])
0174       meEBUncalibRecHitAmplMap_[ism - 1]->Fill(xie, xip, uncalibRecHit->amplitude());
0175 
0176     if (!skipDigis) {
0177       // find the rechit corresponding digi and the max sample
0178       EBDigiCollection::const_iterator myDigi = EBDigi->find(EBid);
0179       // int sMax = -1; // UNUSED
0180       double eMax = 0.;
0181       if (myDigi != EBDigi->end()) {
0182         for (unsigned int sample = 0; sample < myDigi->size(); ++sample) {
0183           EcalMGPASample thisSample = (*myDigi)[sample];
0184           double analogSample = thisSample.adc();
0185           if (eMax < analogSample) {
0186             eMax = analogSample;
0187             // sMax = sample; // UNUSED
0188           }
0189         }
0190       } else
0191         continue;
0192 
0193       // ratio uncalibratedRecHit amplitude + ped / max energy digi
0194       const EcalPedestals *myped = &c.getData(ecalPeds);
0195       EcalPedestalsMap::const_iterator it = myped->getMap().find(EBid);
0196       if (it != myped->getMap().end()) {
0197         if (eMax > (*it).mean_x1 + 5 * (*it).rms_x1 && eMax != 0) {  // only real signal RecHit
0198 
0199           if (meEBUncalibRecHitMaxSampleRatio_)
0200             meEBUncalibRecHitMaxSampleRatio_->Fill((uncalibRecHit->amplitude() + uncalibRecHit->pedestal()) / eMax);
0201           if (meEBUncalibRecHitMaxSampleRatioGt100adc_ && (uncalibRecHit->amplitude() > 100))
0202             meEBUncalibRecHitMaxSampleRatioGt100adc_->Fill((uncalibRecHit->amplitude() + uncalibRecHit->pedestal()) /
0203                                                            eMax);
0204           LogDebug("EcalRecHitsTaskInfo")
0205               << "barrel, eMax = " << eMax << " Amplitude = " << uncalibRecHit->amplitude() + uncalibRecHit->pedestal();
0206         } else
0207           continue;
0208       } else
0209         continue;
0210     }
0211 
0212   }  // loop over the UncalibratedRecHitCollection
0213 }