Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:54:18

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  */
0005 
0006 #include "DTOccupancyEfficiency.h"
0007 
0008 // Framework
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/ServiceRegistry/interface/Service.h"
0014 
0015 #include "DQMServices/Core/interface/DQMStore.h"
0016 
0017 #include <iterator>
0018 
0019 using namespace edm;
0020 using namespace std;
0021 
0022 DTOccupancyEfficiency::DTOccupancyEfficiency(const ParameterSet& pset) {
0023   debug = pset.getUntrackedParameter<bool>("debug", false);
0024   // label for dtDigis
0025   dtDigiToken_ = consumes<DTDigiCollection>(edm::InputTag(pset.getUntrackedParameter<string>("digiLabel")));
0026   // the name of the 4D rec hits collection
0027   recHits4DToken_ =
0028       consumes<DTRecSegment4DCollection>(edm::InputTag(pset.getUntrackedParameter<string>("recHits4DLabel")));
0029   // the name of the rechits collection
0030   recHitToken_ = consumes<DTRecHitCollection>(edm::InputTag(pset.getUntrackedParameter<string>("recHitLabel")));
0031 
0032   parameters = pset;
0033 }
0034 
0035 DTOccupancyEfficiency::~DTOccupancyEfficiency() {}
0036 
0037 void DTOccupancyEfficiency::bookHistograms(DQMStore::IBooker& ibooker,
0038                                            edm::Run const& iRun,
0039                                            edm::EventSetup const& context) {
0040   ibooker.setCurrentFolder("DT/04-OccupancyEfficiency/digisPerRing");
0041   // Digis per ring
0042   for (int station = 1; station < 5; station++) {
0043     string station_s = to_string(station);
0044     for (int wheel = -2; wheel < 3; wheel++) {
0045       string wheel_s = to_string(wheel);
0046       if (wheel > 0)
0047         wheel_s = "+" + wheel_s;
0048       string histoName = "digisPerMB" + station_s + "W" + wheel_s;
0049       string histoTitle = "Number of digis in MB" + station_s + "YB" + wheel_s;
0050       (digisPerRing[station])[wheel] = ibooker.book1D(histoName, histoTitle, 100, 0, 150);
0051     }
0052   }
0053 
0054   ibooker.setCurrentFolder("DT/04-OccupancyEfficiency/timeBoxesPerRing");
0055   // TimeBoxes per ring
0056   for (int station = 1; station < 5; station++) {
0057     string station_s = to_string(station);
0058     for (int wheel = -2; wheel < 3; wheel++) {
0059       string wheel_s = to_string(wheel);
0060       if (wheel > 0)
0061         wheel_s = "+" + wheel_s;
0062       string histoName = "timeBoxesPerMB" + station_s + "W" + wheel_s;
0063       string histoTitle = "Number of TDC counts in MB" + station_s + "YB" + wheel_s;
0064       (timeBoxesPerRing[station])[wheel] = ibooker.book1D(histoName, histoTitle, 400, 0, 1600);
0065     }
0066   }
0067 
0068   ibooker.setCurrentFolder("DT/04-OccupancyEfficiency");
0069 
0070   // TimeBoxes
0071   timeBoxesPerEvent = ibooker.book1D("timeBoxesPerEvent", "TDC counts per event", 400, 0, 1600);
0072 
0073   // Digis
0074   digisPerEvent = ibooker.book1D("digisPerEvent", "Number of digis per event", 100, 0, 900);
0075 
0076   // RecHits
0077   recHitsPerEvent = ibooker.book1D("recHitsPerEvent", "Number of RecHits per event", 100, 0, 250);
0078 
0079   // 4D segments
0080   segments4DPerEvent = ibooker.book1D("segments4DPerEvent", "Number of 4D Segments per event", 50, 0, 50);
0081 
0082   recHitsPer4DSegment = ibooker.book1D("recHitsPer4DSegment", "Number of RecHits per segment", 16, 0.5, 16.5);
0083 
0084   // T0 from segements
0085   t0From4DPhiSegment = ibooker.book1D("t0From4DPhiSegment", "T0 from 4D Phi segments", 100, -150, 150);
0086   t0From4DZSegment = ibooker.book1D("t0From4DZSegment", "T0 from 4D Z segments", 100, -150, 150);
0087 }
0088 
0089 void DTOccupancyEfficiency::analyze(const edm::Event& event, const edm::EventSetup& setup) {
0090   if (debug)
0091     cout << "[DTOccupancyEfficiency] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event() << endl;
0092 
0093   // Digi collection
0094   edm::Handle<DTDigiCollection> dtdigis;
0095   event.getByToken(dtDigiToken_, dtdigis);
0096 
0097   int numberOfDigis = 0;
0098   std::map<int, std::map<int, int>> numberOfDigisPerRing;
0099 
0100   DTDigiCollection::DigiRangeIterator dtLayerId_It;
0101   for (dtLayerId_It = dtdigis->begin(); dtLayerId_It != dtdigis->end(); ++dtLayerId_It) {  // Loop over layers
0102     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
0103          digiIt != ((*dtLayerId_It).second).second;
0104          ++digiIt) {
0105       numberOfDigis++;
0106       int wheelId = ((*dtLayerId_It).first).wheel();
0107       int stationId = ((*dtLayerId_It).first).station();
0108       (numberOfDigisPerRing[stationId])[wheelId] += 1;
0109 
0110       timeBoxesPerEvent->Fill((*digiIt).countsTDC());
0111       (timeBoxesPerRing[stationId])[wheelId]->Fill((*digiIt).countsTDC());
0112     }
0113   }
0114 
0115   // Total number of Digis per Event
0116   digisPerEvent->Fill(numberOfDigis);
0117 
0118   // Digis per Ring in every wheel
0119   for (int station = 1; station < 5; station++)
0120     for (int wheel = -2; wheel < 3; wheel++)
0121       (digisPerRing[station])[wheel]->Fill((numberOfDigisPerRing[station])[wheel]);
0122 
0123   // Get the 4D segment collection from the event
0124   edm::Handle<DTRecSegment4DCollection> all4DSegments;
0125   event.getByToken(recHits4DToken_, all4DSegments);
0126 
0127   segments4DPerEvent->Fill(all4DSegments->size());
0128 
0129   // Get the rechit collection from the event
0130   Handle<DTRecHitCollection> dtRecHits;
0131   event.getByToken(recHitToken_, dtRecHits);
0132 
0133   recHitsPerEvent->Fill(dtRecHits->size());
0134 
0135   // Number of RecHits per segment and t0 from Segment
0136   // Loop on all segments
0137   for (DTRecSegment4DCollection::const_iterator segment = all4DSegments->begin(); segment != all4DSegments->end();
0138        ++segment) {
0139     unsigned int nHits = (segment->hasPhi() ? (segment->phiSegment())->recHits().size() : 0);
0140     nHits += (segment->hasZed() ? (segment->zSegment())->recHits().size() : 0);
0141     recHitsPer4DSegment->Fill(nHits);
0142 
0143     if (segment->hasPhi()) {
0144       double segmentPhiT0 = segment->phiSegment()->t0();
0145       if (segment->phiSegment()->ist0Valid())
0146         t0From4DPhiSegment->Fill(segmentPhiT0);
0147     }
0148     if (segment->hasZed()) {
0149       double segmentZT0 = segment->zSegment()->t0();
0150       if (segment->zSegment()->ist0Valid())
0151         t0From4DZSegment->Fill(segmentZT0);
0152     }
0153   }
0154 }
0155 
0156 // Local Variables:
0157 // show-trailing-whitespace: t
0158 // truncate-lines: t
0159 // End: