CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/DQM/EcalBarrelMonitorClient/src/OccupancyClient.cc

Go to the documentation of this file.
00001 #include "../interface/OccupancyClient.h"
00002 
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00005 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00006 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00007 #include "DataFormats/EcalDetId/interface/EcalTrigTowerDetId.h"
00008 
00009 #include "DQM/EcalBarrelMonitorTasks/interface/OccupancyTask.h"
00010 
00011 #include "DQM/EcalCommon/interface/EcalDQMCommonUtils.h"
00012 
00013 namespace ecaldqm {
00014 
00015   OccupancyClient::OccupancyClient(const edm::ParameterSet& _params, const edm::ParameterSet& _paths) :
00016     DQWorkerClient(_params, _paths, "OccupancyClient"),
00017     geometry_(0),
00018     minHits_(0),
00019     deviationThreshold_(0.)
00020   {
00021     edm::ParameterSet const& taskParams(_params.getUntrackedParameterSet(name_));
00022     minHits_ = taskParams.getUntrackedParameter<int>("minHits");
00023     deviationThreshold_ = taskParams.getUntrackedParameter<double>("deviationThreshold");
00024 
00025     edm::ParameterSet const& sources(_params.getUntrackedParameterSet("sources"));
00026     source_(sDigi, "OccupancyTask", OccupancyTask::kDigi, sources);
00027     source_(sRecHitThr, "OccupancyTask", OccupancyTask::kRecHitThr, sources);
00028     source_(sTPDigiThr, "OccupancyTask", OccupancyTask::kTPDigiThr, sources);
00029   }
00030 
00031   void
00032   OccupancyClient::beginRun(const edm::Run &, const edm::EventSetup &_es)
00033   {
00034     edm::ESHandle<CaloGeometry> geomHndl;
00035     _es.get<CaloGeometryRecord>().get(geomHndl);
00036     geometry_ = geomHndl.product();
00037     if(!geometry_)
00038       throw cms::Exception("EventSetup") << "CaloGeometry invalid";
00039   }
00040 
00041   void
00042   OccupancyClient::bookMEs()
00043   {
00044     DQWorker::bookMEs();
00045 
00046     MEs_[kQualitySummary]->resetAll(-1.);
00047   }
00048 
00049   void
00050   OccupancyClient::producePlots()
00051   {
00052     using namespace std;
00053 
00054     MEs_[kHotDigi]->reset();
00055     MEs_[kHotRecHitThr]->reset();
00056     MEs_[kHotTPDigiThr]->reset();
00057 
00058     uint32_t mask(1 << EcalDQMStatusHelper::PEDESTAL_ONLINE_HIGH_GAIN_RMS_ERROR);
00059 
00060     vector<double> digiPhiRingMean(28, 0.);
00061     vector<double> rechitPhiRingMean(28, 0.);
00062     vector<int> numCrystals(28, 0); // this is static, but is easier to count now
00063 
00064     for(unsigned dccid(1); dccid <= 54; dccid++){
00065       for(unsigned tower(1); tower <= getNSuperCrystals(dccid); tower++){
00066         vector<DetId> ids(getElectronicsMap()->dccTowerConstituents(dccid, tower));
00067 
00068         if(ids.size() == 0) continue;
00069 
00070         for(vector<DetId>::iterator idItr(ids.begin()); idItr != ids.end(); ++idItr){
00071           float entries(sources_[sDigi]->getBinContent(*idItr));
00072           float rhentries(sources_[sRecHitThr]->getBinContent(*idItr));
00073 
00074           int ieta(getTrigTowerMap()->towerOf(*idItr).ietaAbs());
00075           digiPhiRingMean.at(ieta - 1) += entries;
00076           rechitPhiRingMean.at(ieta - 1) += rhentries;
00077 
00078           numCrystals.at(ieta - 1) += 1;
00079         }
00080       }
00081     }
00082 
00083     for(int ie(0); ie < 28; ie++){
00084       digiPhiRingMean[ie] /= numCrystals[ie];
00085       rechitPhiRingMean[ie] /= numCrystals[ie];
00086     }
00087 
00088     // second round to find hot towers
00089     for(unsigned dccid(1); dccid <= 54; dccid++){
00090       for(unsigned tower(1); tower <= getNSuperCrystals(dccid); tower++){
00091         vector<DetId> ids(getElectronicsMap()->dccTowerConstituents(dccid, tower));
00092 
00093         if(ids.size() == 0) continue;
00094 
00095         float quality(1.);
00096         for(vector<DetId>::iterator idItr(ids.begin()); idItr != ids.end(); ++idItr){
00097           float entries(sources_[sDigi]->getBinContent(*idItr));
00098           float rhentries(sources_[sRecHitThr]->getBinContent(*idItr));
00099 
00100           int ieta(getTrigTowerMap()->towerOf(*idItr).ietaAbs());
00101 
00102           if(entries > minHits_ && entries > digiPhiRingMean.at(ieta - 1) * deviationThreshold_){
00103             MEs_[kHotDigi]->fill(*idItr);
00104             quality = 0.;
00105           }
00106           if(rhentries > minHits_ && rhentries > rechitPhiRingMean.at(ieta - 1) * deviationThreshold_){
00107             MEs_[kHotRecHitThr]->fill(*idItr);
00108             quality = 0.;
00109           }
00110         }
00111         if(dccid <= 9 || dccid >= 46){
00112           vector<EcalScDetId> scs(getElectronicsMap()->getEcalScDetId(dccid, tower));
00113           for(vector<EcalScDetId>::iterator scItr(scs.begin()); scItr != scs.end(); ++scItr)
00114             fillQuality_(kQualitySummary, *scItr, mask, quality);
00115         }
00116         else
00117           fillQuality_(kQualitySummary, ids[0], mask, quality);
00118       }
00119     }
00120 
00121     vector<double> tpdigiPhiRingMean(28, 0.);
00122 
00123     for(unsigned iTT(0); iTT < EcalTrigTowerDetId::kSizeForDenseIndexing; iTT++){
00124       EcalTrigTowerDetId ttid(EcalTrigTowerDetId::detIdFromDenseIndex(iTT));
00125       float entries(sources_[sTPDigiThr]->getBinContent(ttid));
00126 
00127       tpdigiPhiRingMean.at(ttid.ietaAbs() - 1) += entries;
00128     }
00129 
00130     for(int ie(0); ie < 28; ie++){
00131       float denom(-1.);
00132       if(ie < 27) denom = 72.;
00133       else denom = 36.;
00134       tpdigiPhiRingMean[ie] /= denom;
00135     }
00136 
00137     for(unsigned iTT(0); iTT < EcalTrigTowerDetId::kSizeForDenseIndexing; iTT++){
00138       EcalTrigTowerDetId ttid(EcalTrigTowerDetId::detIdFromDenseIndex(iTT));
00139       float entries(sources_[sTPDigiThr]->getBinContent(ttid));
00140 
00141       if(entries > minHits_ && entries > tpdigiPhiRingMean.at(ttid.ietaAbs() - 1) * deviationThreshold_){
00142         MEs_[kHotTPDigiThr]->fill(ttid);
00143         vector<DetId> ids(getTrigTowerMap()->constituentsOf(ttid));
00144         for(vector<DetId>::iterator idItr(ids.begin()); idItr != ids.end(); ++idItr){
00145           if(MEs_[kQualitySummary]->getBinContent(*idItr) > 0.)
00146             fillQuality_(kQualitySummary, *idItr, mask, 0.);
00147         }
00148       }   
00149     }
00150 
00151   }
00152 
00153   /*static*/
00154   void
00155   OccupancyClient::setMEData(std::vector<MEData>& _data)
00156   {
00157     _data[kHotDigi] = MEData("HotDigi", BinService::kChannel, BinService::kCrystal, MonitorElement::DQM_KIND_TH1F);
00158     _data[kHotRecHitThr] = MEData("HotRecHitThr", BinService::kChannel, BinService::kCrystal, MonitorElement::DQM_KIND_TH1F);
00159     _data[kHotTPDigiThr] = MEData("HotTPDigiThr", BinService::kChannel, BinService::kTriggerTower, MonitorElement::DQM_KIND_TH1F);
00160     _data[kQualitySummary] = MEData("QualitySummary", BinService::kEcal2P, BinService::kSuperCrystal, MonitorElement::DQM_KIND_TH2F);
00161   }
00162 
00163   DEFINE_ECALDQM_WORKER(OccupancyClient);
00164 }
00165