Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:16

0001 /*
0002  * \file DQMSourceEleCalib.cc
0003  *
0004  * \author Andrea Gozzelino - Universita  e INFN Torino
0005  * \author Stefano Argiro
0006  *
0007  *
0008  *
0009  * Description: Monitoring of Phi Symmetry Calibration Stream
0010  */
0011 
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/MakerMacros.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 
0020 // DQM include files
0021 
0022 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0023 #include "DQMServices/Core/interface/DQMStore.h"
0024 
0025 // work on collections
0026 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0027 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0028 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
0029 #include "DataFormats/DetId/interface/DetId.h"
0030 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0031 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0032 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0033 
0034 class DQMSourceEleCalib : public DQMEDAnalyzer {
0035 public:
0036   DQMSourceEleCalib(const edm::ParameterSet &);
0037   ~DQMSourceEleCalib() override;
0038 
0039 protected:
0040   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0041 
0042   void analyze(const edm::Event &e, const edm::EventSetup &c) override;
0043 
0044 private:
0045   //! find the MOX
0046   DetId findMaxHit(const std::vector<std::pair<DetId, float>> &,
0047                    const EcalRecHitCollection *,
0048                    const EcalRecHitCollection *);
0049   //! fills local occupancy graphs
0050   void fillAroundBarrel(const EcalRecHitCollection *, int, int);
0051   void fillAroundEndcap(const EcalRecHitCollection *, int, int);
0052 
0053   int eventCounter_;
0054 
0055   //! Number of recHits per electron
0056   MonitorElement *recHitsPerElectron_;
0057   //! Number of electrons
0058   MonitorElement *ElectronsNumber_;
0059   //! ESCoP
0060   MonitorElement *ESCoP_;
0061   //! Occupancy
0062   MonitorElement *OccupancyEB_;
0063   MonitorElement *OccupancyEEP_;
0064   MonitorElement *OccupancyEEM_;
0065   MonitorElement *LocalOccupancyEB_;
0066   MonitorElement *LocalOccupancyEE_;
0067 
0068   //! recHits over associated recHits
0069   MonitorElement *HitsVsAssociatedHits_;
0070 
0071   /// object to monitor
0072   edm::EDGetTokenT<EcalRecHitCollection> productMonitoredEB_;
0073 
0074   /// object to monitor
0075   edm::EDGetTokenT<EcalRecHitCollection> productMonitoredEE_;
0076   //! electrons to monitor
0077   edm::EDGetTokenT<reco::GsfElectronCollection> productMonitoredElectrons_;
0078 
0079   /// Monitor every prescaleFactor_ events
0080   unsigned int prescaleFactor_;
0081 
0082   /// DQM folder name
0083   std::string folderName_;
0084 
0085   /// Write to file
0086   bool saveToFile_;
0087 
0088   /// Output file name if required
0089   std::string fileName_;
0090 };
0091 
0092 // ******************************************
0093 // constructors
0094 // *****************************************
0095 
0096 DQMSourceEleCalib::DQMSourceEleCalib(const edm::ParameterSet &ps) : eventCounter_(0) {
0097   folderName_ = ps.getUntrackedParameter<std::string>("FolderName", "ALCAStreamEcalSingleEle");
0098   productMonitoredEB_ = consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("AlCaStreamEBTag"));
0099   productMonitoredEE_ = consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("AlCaStreamEETag"));
0100 
0101   saveToFile_ = ps.getUntrackedParameter<bool>("SaveToFile", false);
0102   fileName_ = ps.getUntrackedParameter<std::string>("FileName", "MonitorAlCaEcalSingleEle.root");
0103   productMonitoredElectrons_ =
0104       consumes<reco::GsfElectronCollection>(ps.getParameter<edm::InputTag>("electronCollection"));
0105   prescaleFactor_ = ps.getUntrackedParameter<unsigned int>("prescaleFactor", 1);
0106 }
0107 
0108 DQMSourceEleCalib::~DQMSourceEleCalib() {}
0109 
0110 //--------------------------------------------------------
0111 void DQMSourceEleCalib::bookHistograms(DQMStore::IBooker &ibooker,
0112                                        edm::Run const &irun,
0113                                        edm::EventSetup const &isetup) {
0114   // create and cd into new folder
0115   ibooker.setCurrentFolder(folderName_);
0116 
0117   recHitsPerElectron_ = ibooker.book1D("recHitsPerElectron_", "recHitPerElectron", 200, 0, 200);
0118   ElectronsNumber_ = ibooker.book1D("ElectronsNumber_", "electrons in the event", 40, 0, 40);
0119   ESCoP_ = ibooker.book1D("ESCoP", "ESCoP", 50, 0, 5);
0120 
0121   OccupancyEB_ = ibooker.book2D("OccupancyEB_", "OccupancyEB", 360, 1, 361, 171, -85, 86);
0122   OccupancyEEP_ = ibooker.book2D("OccupancyEEP_", "Occupancy EE Plus", 100, 1, 101, 100, 1, 101);
0123   OccupancyEEM_ = ibooker.book2D("OccupancyEEM_", "Occupancy EE Minus", 100, 1, 101, 100, 1, 101);
0124   HitsVsAssociatedHits_ = ibooker.book1D("HitsVsAssociatedHits_", "HitsVsAssociatedHits", 100, 0, 5);
0125   LocalOccupancyEB_ = ibooker.book2D("LocalOccupancyEB_", "Local occupancy Barrel", 9, -4, 5, 9, -4, 5);
0126   LocalOccupancyEE_ = ibooker.book2D("LocalOccupancyEE_", "Local occupancy Endcap", 9, -4, 5, 9, -4, 5);
0127 }
0128 
0129 //--------------------------------------------------------
0130 
0131 //-------------------------------------------------------------
0132 
0133 void DQMSourceEleCalib::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0134   //  if (eventCounter_% prescaleFactor_ ) return; //FIXME
0135   eventCounter_++;
0136   int numberOfHits = 0;
0137   int numberOfElectrons = 0;
0138   int numberOfAssociatedHits = 0;
0139   // reads the recHits
0140   edm::Handle<EcalRecHitCollection> rhEB;
0141   edm::Handle<EcalRecHitCollection> rhEE;
0142 
0143   iEvent.getByToken(productMonitoredEB_, rhEB);
0144   iEvent.getByToken(productMonitoredEE_, rhEE);
0145 
0146   EcalRecHitCollection::const_iterator itb;
0147 
0148   // reads the electrons
0149   edm::Handle<reco::GsfElectronCollection> pElectrons;
0150   iEvent.getByToken(productMonitoredElectrons_, pElectrons);
0151 
0152   if (pElectrons.isValid()) {
0153     ElectronsNumber_->Fill(pElectrons->size() + 0.1);
0154     numberOfElectrons = pElectrons->size();
0155     for (reco::GsfElectronCollection::const_iterator eleIt = pElectrons->begin(); eleIt != pElectrons->end(); ++eleIt) {
0156       ESCoP_->Fill(eleIt->eSuperClusterOverP());
0157       numberOfAssociatedHits += eleIt->superCluster()->size();
0158       DetId Max = findMaxHit(eleIt->superCluster()->hitsAndFractions(), rhEB.product(), rhEE.product());
0159       if (!Max.det())
0160         continue;
0161       if (Max.subdetId() == EcalBarrel) {
0162         EBDetId EBMax(Max);
0163         fillAroundBarrel(rhEB.product(), EBMax.ieta(), EBMax.iphi());
0164       }
0165       if (Max.subdetId() == EcalEndcap) {
0166         EEDetId EEMax(Max);
0167         fillAroundEndcap(rhEE.product(), EEMax.ix(), EEMax.iy());
0168       }
0169     }
0170   }  // is valid electron
0171 
0172   // fill EB histos
0173   if (rhEB.isValid()) {
0174     numberOfHits += rhEB->size();
0175     for (itb = rhEB->begin(); itb != rhEB->end(); ++itb) {
0176       EBDetId id(itb->id());
0177       OccupancyEB_->Fill(id.iphi(), id.ieta());
0178     }  // Eb rechits
0179   }    // is Valid
0180   if (rhEE.isValid()) {
0181     numberOfHits += rhEE->size();
0182     for (itb = rhEE->begin(); itb != rhEE->end(); ++itb) {
0183       EEDetId id(itb->id());
0184       if (id.zside() > 0) {
0185         OccupancyEEP_->Fill(id.ix(), id.iy());
0186       }  // zside>0
0187       else if (id.zside() < 0) {
0188         OccupancyEEM_->Fill(id.ix(), id.iy());
0189       }  // zside<0
0190 
0191     }  // EE reChit
0192   }    // is Valid
0193   if (numberOfElectrons)
0194     recHitsPerElectron_->Fill((double)numberOfHits / ((double)numberOfElectrons));
0195   if (numberOfHits)
0196     HitsVsAssociatedHits_->Fill((double)numberOfAssociatedHits / ((double)numberOfHits));
0197 }  // end of the analyzer
0198 
0199 //--------------------------------------------------------
0200 
0201 //------------------------------------------------
0202 
0203 DetId DQMSourceEleCalib::findMaxHit(const std::vector<std::pair<DetId, float>> &v1,
0204                                     const EcalRecHitCollection *EBhits,
0205                                     const EcalRecHitCollection *EEhits) {
0206   double currEnergy = 0.;
0207   DetId maxHit;
0208   for (std::vector<std::pair<DetId, float>>::const_iterator idsIt = v1.begin(); idsIt != v1.end(); ++idsIt) {
0209     if (idsIt->first.subdetId() == EcalBarrel) {
0210       EcalRecHitCollection::const_iterator itrechit;
0211       itrechit = EBhits->find((*idsIt).first);
0212       if (itrechit == EBhits->end()) {
0213         edm::LogInfo("reading") << "[findMaxHit] rechit not found! ";
0214         continue;
0215       }
0216       // FIXME: wnat to use the fraction i.e. .second??
0217       if (itrechit->energy() > currEnergy) {
0218         currEnergy = itrechit->energy();
0219         maxHit = (*idsIt).first;
0220       }
0221     } else {
0222       EcalRecHitCollection::const_iterator itrechit;
0223       itrechit = EEhits->find((*idsIt).first);
0224       if (itrechit == EEhits->end()) {
0225         edm::LogInfo("reading") << "[findMaxHit] rechit not found! ";
0226         continue;
0227       }
0228 
0229       // FIXME: wnat to use the fraction i.e. .second??
0230       if (itrechit->energy() > currEnergy) {
0231         currEnergy = itrechit->energy();
0232         maxHit = (*idsIt).first;
0233       }
0234     }
0235   }
0236   return maxHit;
0237 }
0238 
0239 void DQMSourceEleCalib::fillAroundBarrel(const EcalRecHitCollection *recHits, int eta, int phi) {
0240   for (EcalRecHitCollection::const_iterator elem = recHits->begin(); elem != recHits->end(); ++elem) {
0241     EBDetId elementId = elem->id();
0242     LocalOccupancyEB_->Fill(elementId.ieta() - eta, elementId.iphi() - phi, elem->energy());
0243   }
0244   return;
0245 }
0246 
0247 // ----------------------------------------------------------------
0248 
0249 void DQMSourceEleCalib::fillAroundEndcap(const EcalRecHitCollection *recHits, int ics, int ips) {
0250   for (EcalRecHitCollection::const_iterator elem = recHits->begin(); elem != recHits->end(); ++elem) {
0251     EEDetId elementId = elem->id();
0252     LocalOccupancyEE_->Fill(elementId.ix() - ics, elementId.iy() - ips, elem->energy());
0253   }
0254   return;
0255 }
0256 
0257 DEFINE_FWK_MODULE(DQMSourceEleCalib);