Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:50

0001 /*
0002  *  See header file for a description of this class.
0003  *
0004  *  \author G. Mila - INFN Torino
0005  */
0006 
0007 #include "DQMOffline/Muon/interface/DTSegmentsTask.h"
0008 
0009 // Framework
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015 
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 
0018 //Geometry
0019 #include "DataFormats/GeometryVector/interface/Pi.h"
0020 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
0021 
0022 #include <iterator>
0023 
0024 using namespace edm;
0025 using namespace std;
0026 
0027 DTSegmentsTask::DTSegmentsTask(const edm::ParameterSet& pset)
0028     : statusMapToken_(esConsumes<DTStatusFlag, DTStatusFlagRcd>()) {
0029   debug = pset.getUntrackedParameter<bool>("debug", false);
0030   parameters = pset;
0031 
0032   // should be init from pset, but it was commented out...
0033   // better false than undefined
0034   checkNoisyChannels = false;
0035 
0036   // the name of the 4D rec hits collection
0037   theRecHits4DLabel_ = consumes<DTRecSegment4DCollection>(parameters.getParameter<std::string>("recHits4DLabel"));
0038 }
0039 
0040 DTSegmentsTask::~DTSegmentsTask() {}
0041 
0042 void DTSegmentsTask::bookHistograms(DQMStore::IBooker& ibooker,
0043                                     edm::Run const& /*iRun*/,
0044                                     edm::EventSetup const& /* iSetup */) {
0045   ibooker.cd();
0046   ibooker.setCurrentFolder("Muons/DTSegmentsMonitor");
0047 
0048   // histos for phi segments
0049   phiHistos.push_back(
0050       ibooker.book2D("phiSegments_numHitsVsWheel", "phiSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
0051   phiHistos[0]->setBinLabel(1, "W-2", 1);
0052   phiHistos[0]->setBinLabel(2, "W-1", 1);
0053   phiHistos[0]->setBinLabel(3, "W0", 1);
0054   phiHistos[0]->setBinLabel(4, "W1", 1);
0055   phiHistos[0]->setBinLabel(5, "W2", 1);
0056   phiHistos.push_back(
0057       ibooker.book2D("phiSegments_numHitsVsSector", "phiSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
0058   phiHistos[1]->setBinLabel(1, "Sec1", 1);
0059   phiHistos[1]->setBinLabel(2, "Sec2", 1);
0060   phiHistos[1]->setBinLabel(3, "Sec3", 1);
0061   phiHistos[1]->setBinLabel(4, "Sec4", 1);
0062   phiHistos[1]->setBinLabel(5, "Sec5", 1);
0063   phiHistos[1]->setBinLabel(6, "Sec6", 1);
0064   phiHistos[1]->setBinLabel(7, "Sec7", 1);
0065   phiHistos[1]->setBinLabel(8, "Sec8", 1);
0066   phiHistos[1]->setBinLabel(9, "Sec9", 1);
0067   phiHistos[1]->setBinLabel(10, "Sec10", 1);
0068   phiHistos[1]->setBinLabel(11, "Sec11", 1);
0069   phiHistos[1]->setBinLabel(12, "Sec12", 1);
0070   phiHistos[1]->setBinLabel(13, "Sec13", 1);
0071   phiHistos[1]->setBinLabel(14, "Sec14", 1);
0072   phiHistos.push_back(
0073       ibooker.book2D("phiSegments_numHitsVsStation", "phiSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
0074   phiHistos[2]->setBinLabel(1, "St1", 1);
0075   phiHistos[2]->setBinLabel(2, "St2", 1);
0076   phiHistos[2]->setBinLabel(3, "St3", 1);
0077   phiHistos[2]->setBinLabel(4, "St4", 1);
0078 
0079   // histos for theta segments
0080   thetaHistos.push_back(
0081       ibooker.book2D("thetaSegments_numHitsVsWheel", "thetaSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
0082   thetaHistos[0]->setBinLabel(1, "W-2", 1);
0083   thetaHistos[0]->setBinLabel(2, "W-1", 1);
0084   thetaHistos[0]->setBinLabel(3, "W0", 1);
0085   thetaHistos[0]->setBinLabel(4, "W1", 1);
0086   thetaHistos[0]->setBinLabel(5, "W2", 1);
0087   thetaHistos.push_back(
0088       ibooker.book2D("thetaSegments_numHitsVsSector", "thetaSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
0089   thetaHistos[1]->setBinLabel(1, "Sec1", 1);
0090   thetaHistos[1]->setBinLabel(2, "Sec2", 1);
0091   thetaHistos[1]->setBinLabel(3, "Sec3", 1);
0092   thetaHistos[1]->setBinLabel(4, "Sec4", 1);
0093   thetaHistos[1]->setBinLabel(5, "Sec5", 1);
0094   thetaHistos[1]->setBinLabel(6, "Sec6", 1);
0095   thetaHistos[1]->setBinLabel(7, "Sec7", 1);
0096   thetaHistos[1]->setBinLabel(8, "Sec8", 1);
0097   thetaHistos[1]->setBinLabel(9, "Sec9", 1);
0098   thetaHistos[1]->setBinLabel(10, "Sec10", 1);
0099   thetaHistos[1]->setBinLabel(11, "Sec11", 1);
0100   thetaHistos[1]->setBinLabel(12, "Sec12", 1);
0101   thetaHistos[1]->setBinLabel(13, "Sec13", 1);
0102   thetaHistos[1]->setBinLabel(14, "Sec14", 1);
0103   thetaHistos.push_back(
0104       ibooker.book2D("thetaSegments_numHitsVsStation", "thetaSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
0105   thetaHistos[2]->setBinLabel(1, "St1", 1);
0106   thetaHistos[2]->setBinLabel(2, "St2", 1);
0107   thetaHistos[2]->setBinLabel(3, "St3", 1);
0108   thetaHistos[2]->setBinLabel(4, "St4", 1);
0109 }
0110 
0111 void DTSegmentsTask::analyze(const edm::Event& event, const edm::EventSetup& setup) {
0112   // Get the map of noisy channels
0113   //  bool checkNoisyChannels = parameters.getUntrackedParameter<bool>("checkNoisyChannels",false);
0114   ESHandle<DTStatusFlag> statusMap;
0115   if (checkNoisyChannels) {
0116     statusMap = setup.getHandle(statusMapToken_);
0117   }
0118 
0119   // Get the 4D segment collection from the event
0120   edm::Handle<DTRecSegment4DCollection> all4DSegments;
0121   event.getByToken(theRecHits4DLabel_, all4DSegments);
0122 
0123   // Loop over all chambers containing a segment
0124   DTRecSegment4DCollection::id_iterator chamberId;
0125   for (chamberId = all4DSegments->id_begin(); chamberId != all4DSegments->id_end(); ++chamberId) {
0126     // Get the range for the corresponding ChamerId
0127     DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
0128     int nsegm = distance(range.first, range.second);
0129     if (debug)
0130       cout << "   Chamber: " << *chamberId << " has " << nsegm << " 4D segments" << endl;
0131 
0132     // Loop over the rechits of this ChamerId
0133     for (DTRecSegment4DCollection::const_iterator segment4D = range.first; segment4D != range.second; ++segment4D) {
0134       //FOR NOISY CHANNELS////////////////////////////////
0135       bool segmNoisy = false;
0136       if ((*segment4D).hasPhi()) {
0137         const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
0138         vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
0139         map<DTSuperLayerId, vector<DTRecHit1D> > hitsBySLMap;
0140         for (vector<DTRecHit1D>::const_iterator hit = phiHits.begin(); hit != phiHits.end(); ++hit) {
0141           DTWireId wireId = (*hit).wireId();
0142 
0143           // Check for noisy channels to skip them
0144           if (checkNoisyChannels) {
0145             bool isNoisy = false;
0146             bool isFEMasked = false;
0147             bool isTDCMasked = false;
0148             bool isTrigMask = false;
0149             bool isDead = false;
0150             bool isNohv = false;
0151             statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
0152             if (isNoisy) {
0153               if (debug)
0154                 cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
0155               segmNoisy = true;
0156             }
0157           }
0158         }
0159       }
0160 
0161       if ((*segment4D).hasZed()) {
0162         const DTSLRecSegment2D* zSeg = (*segment4D).zSegment();  // zSeg lives in the SL RF
0163         // Check for noisy channels to skip them
0164         vector<DTRecHit1D> zHits = zSeg->specificRecHits();
0165         for (vector<DTRecHit1D>::const_iterator hit = zHits.begin(); hit != zHits.end(); ++hit) {
0166           DTWireId wireId = (*hit).wireId();
0167           if (checkNoisyChannels) {
0168             bool isNoisy = false;
0169             bool isFEMasked = false;
0170             bool isTDCMasked = false;
0171             bool isTrigMask = false;
0172             bool isDead = false;
0173             bool isNohv = false;
0174             //cout<<"wire id "<<wireId<<endl;
0175             statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
0176             if (isNoisy) {
0177               if (debug)
0178                 cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
0179               segmNoisy = true;
0180             }
0181           }
0182         }
0183       }
0184 
0185       if (segmNoisy) {
0186         if (debug)
0187           cout << "skipping the segment: it contains noisy cells" << endl;
0188         continue;
0189       }
0190       //END FOR NOISY CHANNELS////////////////////////////////
0191 
0192       // Fill the histos
0193       int nHits = 0;
0194       if ((*segment4D).hasPhi()) {
0195         nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
0196         if (debug)
0197           cout << "Phi segment with number of hits: " << nHits << endl;
0198         phiHistos[0]->Fill((*chamberId).wheel(), nHits);
0199         phiHistos[1]->Fill((*chamberId).sector(), nHits);
0200         phiHistos[2]->Fill((*chamberId).station(), nHits);
0201       }
0202       if ((*segment4D).hasZed()) {
0203         nHits = (((*segment4D).zSegment())->specificRecHits()).size();
0204         if (debug)
0205           cout << "Zed segment with number of hits: " << nHits << endl;
0206         thetaHistos[0]->Fill((*chamberId).wheel(), nHits);
0207         thetaHistos[1]->Fill((*chamberId).sector(), nHits);
0208         thetaHistos[2]->Fill((*chamberId).station(), nHits);
0209       }
0210 
0211     }  //loop over segments
0212   }  // loop over chambers
0213 }