CMS 3D CMS Logo

DTSegmentsTask.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2008/10/21 12:07:47 $
00006  *  $Revision: 1.3 $
00007  *  \author G. Mila - INFN Torino
00008  */
00009 
00010 #include "DQMOffline/Muon/interface/DTSegmentsTask.h"
00011 
00012 // Framework
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 
00019 #include "DQMServices/Core/interface/DQMStore.h"
00020 #include "DQMServices/Core/interface/MonitorElement.h"
00021 
00022 //Geometry
00023 #include "DataFormats/GeometryVector/interface/Pi.h"
00024 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
00025 
00026 //RecHit
00027 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00028 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00029 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00030 
00031 #include <iterator>
00032 
00033 using namespace edm;
00034 using namespace std;
00035 
00036 DTSegmentsTask::DTSegmentsTask(const edm::ParameterSet& pset) {
00037 
00038   debug = pset.getUntrackedParameter<bool>("debug","false");
00039 
00040   // Get the DQM needed services
00041   theDbe = edm::Service<DQMStore>().operator->();
00042   theDbe->setVerbose(1);
00043 
00044   parameters = pset;
00045 
00046 }
00047 
00048 
00049 DTSegmentsTask::~DTSegmentsTask(){
00050 }
00051 
00052 
00053 void DTSegmentsTask::beginJob(const edm::EventSetup& context){
00054  
00055  // the name of the 4D rec hits collection
00056   theRecHits4DLabel = parameters.getParameter<string>("recHits4DLabel");
00057 
00058   theDbe->setCurrentFolder("Muons/DTSegmentsMonitor");
00059 
00060   // histos for phi segments
00061   phiHistos.push_back(theDbe->book2D("phiSegments_numHitsVsWheel", "phiSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
00062   phiHistos[0]->setBinLabel(1,"W-2",1);
00063   phiHistos[0]->setBinLabel(2,"W-1",1);
00064   phiHistos[0]->setBinLabel(3,"W0",1);
00065   phiHistos[0]->setBinLabel(4,"W1",1);
00066   phiHistos[0]->setBinLabel(5,"W2",1);
00067   phiHistos.push_back(theDbe->book2D("phiSegments_numHitsVsSector", "phiSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
00068   phiHistos[1]->setBinLabel(1,"Sec1",1);
00069   phiHistos[1]->setBinLabel(2,"Sec2",1);
00070   phiHistos[1]->setBinLabel(3,"Sec3",1);
00071   phiHistos[1]->setBinLabel(4,"Sec4",1);
00072   phiHistos[1]->setBinLabel(5,"Sec5",1);
00073   phiHistos[1]->setBinLabel(6,"Sec6",1);
00074   phiHistos[1]->setBinLabel(7,"Sec7",1);
00075   phiHistos[1]->setBinLabel(8,"Sec8",1);
00076   phiHistos[1]->setBinLabel(9,"Sec9",1);
00077   phiHistos[1]->setBinLabel(10,"Sec10",1);
00078   phiHistos[1]->setBinLabel(11,"Sec11",1);
00079   phiHistos[1]->setBinLabel(12,"Sec12",1);
00080   phiHistos[1]->setBinLabel(13,"Sec13",1);
00081   phiHistos[1]->setBinLabel(14,"Sec14",1);
00082   phiHistos.push_back(theDbe->book2D("phiSegments_numHitsVsStation", "phiSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
00083   phiHistos[2]->setBinLabel(1,"St1",1);
00084   phiHistos[2]->setBinLabel(2,"St2",1);
00085   phiHistos[2]->setBinLabel(3,"St3",1);
00086   phiHistos[2]->setBinLabel(4,"St4",1);
00087 
00088   // histos for theta segments
00089   thetaHistos.push_back(theDbe->book2D("thetaSegments_numHitsVsWheel", "thetaSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
00090   thetaHistos[0]->setBinLabel(1,"W-2",1);
00091   thetaHistos[0]->setBinLabel(2,"W-1",1);
00092   thetaHistos[0]->setBinLabel(3,"W0",1);
00093   thetaHistos[0]->setBinLabel(4,"W1",1);
00094   thetaHistos[0]->setBinLabel(5,"W2",1);
00095   thetaHistos.push_back(theDbe->book2D("thetaSegments_numHitsVsSector", "thetaSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
00096   thetaHistos[1]->setBinLabel(1,"Sec1",1);
00097   thetaHistos[1]->setBinLabel(2,"Sec2",1);
00098   thetaHistos[1]->setBinLabel(3,"Sec3",1);
00099   thetaHistos[1]->setBinLabel(4,"Sec4",1);
00100   thetaHistos[1]->setBinLabel(5,"Sec5",1);
00101   thetaHistos[1]->setBinLabel(6,"Sec6",1);
00102   thetaHistos[1]->setBinLabel(7,"Sec7",1);
00103   thetaHistos[1]->setBinLabel(8,"Sec8",1);
00104   thetaHistos[1]->setBinLabel(9,"Sec9",1);
00105   thetaHistos[1]->setBinLabel(10,"Sec10",1);
00106   thetaHistos[1]->setBinLabel(11,"Sec11",1);
00107   thetaHistos[1]->setBinLabel(12,"Sec12",1);
00108   thetaHistos[1]->setBinLabel(13,"Sec13",1);
00109   thetaHistos[1]->setBinLabel(14,"Sec14",1);
00110   thetaHistos.push_back(theDbe->book2D("thetaSegments_numHitsVsStation", "thetaSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
00111   thetaHistos[2]->setBinLabel(1,"St1",1);
00112   thetaHistos[2]->setBinLabel(2,"St2",1);
00113   thetaHistos[2]->setBinLabel(3,"St3",1);
00114   thetaHistos[2]->setBinLabel(4,"St4",1);
00115 
00116 }
00117 
00118 
00119 void DTSegmentsTask::endJob(){
00120   bool outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
00121   std::string outputFileName = parameters.getParameter<std::string>("OutputFileName");
00122   if(outputMEsInRootFile){
00123     theDbe->save(outputFileName);
00124   }
00125 
00126   theDbe->rmdir("DT/DTSegmentsTask");
00127 }
00128   
00129 void DTSegmentsTask::analyze(const edm::Event& event, const edm::EventSetup& setup) {
00130 
00131 //  if(!(event.id().event()%1000) && debug)
00132 //    {
00133 //      cout << "[DTSegmentsTask] Analyze #Run: " << event.id().run()
00134 //         << " #Event: " << event.id().event() << endl;
00135 //    }
00136 
00137 
00138   // Get the map of noisy channels
00139   bool checkNoisyChannels = parameters.getUntrackedParameter<bool>("checkNoisyChannels","false");
00140   ESHandle<DTStatusFlag> statusMap;
00141   if(checkNoisyChannels) {
00142     setup.get<DTStatusFlagRcd>().get(statusMap);
00143   } 
00144 
00145   // Get the 4D segment collection from the event
00146   edm::Handle<DTRecSegment4DCollection> all4DSegments;
00147   event.getByLabel(theRecHits4DLabel, all4DSegments);
00148 
00149   // Loop over all chambers containing a segment
00150   DTRecSegment4DCollection::id_iterator chamberId;
00151   for (chamberId = all4DSegments->id_begin();
00152        chamberId != all4DSegments->id_end();
00153        ++chamberId){
00154     // Get the range for the corresponding ChamerId
00155     DTRecSegment4DCollection::range  range = all4DSegments->get(*chamberId);
00156     int nsegm = distance(range.first, range.second);
00157     if(debug)
00158       cout << "   Chamber: " << *chamberId << " has " << nsegm
00159            << " 4D segments" << endl;
00160    
00161     
00162      // Loop over the rechits of this ChamerId
00163     for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00164          segment4D!=range.second;
00165            ++segment4D){
00166 
00167       //FOR NOISY CHANNELS////////////////////////////////
00168      bool segmNoisy = false;
00169      if((*segment4D).hasPhi()){
00170        const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
00171        vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
00172        map<DTSuperLayerId,vector<DTRecHit1D> > hitsBySLMap; 
00173        for(vector<DTRecHit1D>::const_iterator hit = phiHits.begin();
00174            hit != phiHits.end(); ++hit) {
00175          DTWireId wireId = (*hit).wireId();
00176          
00177          // Check for noisy channels to skip them
00178          if(checkNoisyChannels) {
00179            bool isNoisy = false;
00180            bool isFEMasked = false;
00181            bool isTDCMasked = false;
00182            bool isTrigMask = false;
00183            bool isDead = false;
00184            bool isNohv = false;
00185            statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00186            if(isNoisy) {
00187              if(debug)
00188                cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
00189              segmNoisy = true;
00190           }      
00191          }
00192        }
00193      }
00194 
00195      if((*segment4D).hasZed()) {
00196        const DTSLRecSegment2D* zSeg = (*segment4D).zSegment();  // zSeg lives in the SL RF
00197        // Check for noisy channels to skip them
00198        vector<DTRecHit1D> zHits = zSeg->specificRecHits();
00199        for(vector<DTRecHit1D>::const_iterator hit = zHits.begin();
00200            hit != zHits.end(); ++hit) {
00201          DTWireId wireId = (*hit).wireId();
00202          if(checkNoisyChannels) {
00203            bool isNoisy = false;
00204             bool isFEMasked = false;
00205             bool isTDCMasked = false;
00206             bool isTrigMask = false;
00207             bool isDead = false;
00208             bool isNohv = false;
00209             //cout<<"wire id "<<wireId<<endl;
00210             statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00211             if(isNoisy) {
00212               if(debug)
00213                 cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
00214               segmNoisy = true;
00215             }      
00216          }
00217        }
00218      } 
00219      
00220      if (segmNoisy) {
00221        if(debug)
00222          cout<<"skipping the segment: it contains noisy cells"<<endl;
00223        continue;
00224      }
00225      //END FOR NOISY CHANNELS////////////////////////////////
00226 
00227      // Fill the histos
00228      int nHits=0;
00229      if((*segment4D).hasPhi()){
00230        nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
00231        if(debug)
00232          cout<<"Phi segment with number of hits: "<<nHits<<endl;
00233        phiHistos[0]->Fill((*chamberId).wheel(), nHits);
00234        phiHistos[1]->Fill((*chamberId).sector(), nHits);
00235        phiHistos[2]->Fill((*chamberId).station(), nHits);
00236      }
00237      if((*segment4D).hasZed()) {
00238        nHits = (((*segment4D).zSegment())->specificRecHits()).size();
00239        if(debug)
00240          cout<<"Zed segment with number of hits: "<<nHits<<endl;
00241        thetaHistos[0]->Fill((*chamberId).wheel(), nHits);
00242        thetaHistos[1]->Fill((*chamberId).sector(), nHits);
00243        thetaHistos[2]->Fill((*chamberId).station(), nHits);
00244      }
00245 
00246     } //loop over segments
00247   } // loop over chambers
00248 
00249 }
00250 

Generated on Tue Jun 9 17:33:49 2009 for CMSSW by  doxygen 1.5.4