#include <DTTriggerEfficiencyTask.h>
Public Member Functions | |
DTTriggerEfficiencyTask (const edm::ParameterSet &ps) | |
Constructor. | |
virtual | ~DTTriggerEfficiencyTask () |
Destructor. | |
Protected Member Functions | |
void | analyze (const edm::Event &e, const edm::EventSetup &c) |
Analyze. | |
void | beginJob () |
void | beginLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &context) |
To reset the MEs. | |
void | beginRun (const edm::Run &run, const edm::EventSetup &context) |
BeginRun. | |
void | bookChamberHistos (const DTChamberId &dtCh, std::string histoTag, std::string folder="") |
Book chamber granularity histograms. | |
void | bookWheelHistos (int wheel, std::string histoTag, std::string folder="") |
Book wheel granularity histograms. | |
void | endJob (void) |
EndJob. | |
bool | hasRPCTriggers (const edm::Event &e) |
checks for RPC Triggers | |
std::string | topFolder (bool source) |
return the top folder 0=DDU 1=DCC | |
Private Attributes | |
std::map< uint32_t, std::map < std::string, MonitorElement * > > | chamberHistos |
DQMStore * | dbe |
bool | detailedPlots |
int | maxBXDDU |
int | minBXDDU |
edm::ESHandle< DTGeometry > | muonGeom |
int | nevents |
edm::ParameterSet | parameters |
const L1MuDTChambPhDigi * | phBestDCC [6][5][13] |
const DTLocalTrigger * | phBestDDU [6][5][13] |
int | phCodeBestDCC [6][5][13] |
int | phCodeBestDDU [6][5][13] |
bool | processDCC |
bool | processDDU |
DTTrigGeomUtils * | trigGeomUtils |
std::map< int, std::map < std::string, MonitorElement * > > | wheelHistos |
Definition at line 45 of file DTTriggerEfficiencyTask.h.
DTTriggerEfficiencyTask::DTTriggerEfficiencyTask | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 46 of file DTTriggerEfficiencyTask.cc.
References dbe, LogTrace, cmsCodeRules::cppFunctionSkipper::operator, and parameters.
: trigGeomUtils(0) { LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: Constructor" << endl; parameters = ps; dbe = edm::Service<DQMStore>().operator->(); }
DTTriggerEfficiencyTask::~DTTriggerEfficiencyTask | ( | ) | [virtual] |
void DTTriggerEfficiencyTask::analyze | ( | const edm::Event & | e, |
const edm::EventSetup & | c | ||
) | [protected, virtual] |
Analyze.
Implements edm::EDAnalyzer.
Definition at line 135 of file DTTriggerEfficiencyTask.cc.
References DTLocalTrigger::bx(), chamberHistos, DTRecSegment4D::chamberId(), DTTrigGeomUtils::computeSCCoordinates(), detailedPlots, edm::Event::getByLabel(), edm::ParameterSet::getUntrackedParameter(), hasRPCTriggers(), i, j, gen::k, maxBXDDU, minBXDDU, nevents, parameters, phBestDCC, phBestDDU, phCodeBestDCC, phCodeBestDDU, processDCC, processDDU, DetId::rawId(), DTChamberId::sector(), relativeConstraints::station, trigGeomUtils, DTChamberId::wheel(), wheelHistos, x, xdir, detailsBasic3DVector::y, and ydir.
{ nevents++; if (!hasRPCTriggers(e)) { return; } InputTag inputTagDCC = parameters.getUntrackedParameter<edm::InputTag>("inputTagDCC"); InputTag inputTagDDU = parameters.getUntrackedParameter<edm::InputTag>("inputTagDDU"); InputTag inputTagSEG = parameters.getUntrackedParameter<edm::InputTag>("inputTagSEG"); for (int i=0;i<5;++i){ for (int j=0;j<6;++j){ for (int k=0;k<13;++k){ phCodeBestDCC[j][i][k] = -1; phCodeBestDDU[j][i][k] = -1; } } } // Getting best DCC Stuff edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh; e.getByLabel(inputTagDCC,l1DTTPGPh); vector<L1MuDTChambPhDigi>* phTrigs = l1DTTPGPh->getContainer(); vector<L1MuDTChambPhDigi>::const_iterator iph = phTrigs->begin(); vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end(); for(; iph !=iphe ; ++iph) { int phwheel = iph->whNum(); int phsec = iph->scNum() + 1; // DTTF numbering [0:11] -> DT numbering [1:12] int phst = iph->stNum(); int phcode = iph->code(); if(phcode>phCodeBestDCC[phwheel+3][phst][phsec] && phcode<7) { phCodeBestDCC[phwheel+3][phst][phsec]=phcode; phBestDCC[phwheel+3][phst][phsec] = &(*iph); } } //Getting Best DDU Stuff Handle<DTLocalTriggerCollection> trigsDDU; e.getByLabel(inputTagDDU,trigsDDU); DTLocalTriggerCollection::DigiRangeIterator detUnitIt; for (detUnitIt=trigsDDU->begin();detUnitIt!=trigsDDU->end();++detUnitIt){ const DTChamberId& id = (*detUnitIt).first; const DTLocalTriggerCollection::Range& range = (*detUnitIt).second; int wh = id.wheel(); int sec = id.sector(); int st = id.station(); for (DTLocalTriggerCollection::const_iterator trigIt = range.first; trigIt!=range.second;++trigIt){ int quality = trigIt->quality(); if(quality>-1 && quality<7 && quality>phCodeBestDDU[wh+3][st][sec]) { phCodeBestDDU[wh+3][st][sec]=quality; phBestDDU[wh+3][st][sec] = &(*trigIt); } } } //Getting Best Segments vector<const DTRecSegment4D*> best4DSegments; Handle<DTRecSegment4DCollection> segments4D; e.getByLabel(inputTagSEG, segments4D); DTRecSegment4DCollection::const_iterator track; DTRecSegment4DCollection::id_iterator chamberId; for (chamberId = segments4D->id_begin(); chamberId != segments4D->id_end(); ++chamberId){ DTRecSegment4DCollection::range range = segments4D->get(*chamberId); const DTRecSegment4D* tmpBest=0; int tmpHit = 0; int hit = 0; for ( track = range.first; track != range.second; ++track){ if( (*track).hasPhi() ) { hit = (*track).phiSegment()->degreesOfFreedom()+2; if ( hit>tmpHit ){ tmpBest = &(*track); tmpHit = hit; int sec = (*track).chamberId().sector(); if (sec==13){ sec=4; } else if (sec==14){ sec=10; } } } } if (tmpBest) { best4DSegments.push_back(tmpBest); } } // Plot filling vector<const DTRecSegment4D*>::const_iterator btrack; for ( btrack = best4DSegments.begin(); btrack != best4DSegments.end(); ++btrack ){ int wheel = (*btrack)->chamberId().wheel(); int station = (*btrack)->chamberId().station(); int scsector = 0; float x, xdir, y, ydir; trigGeomUtils->computeSCCoordinates((*btrack),scsector,x,xdir,y,ydir); int nHitsPhi = (*btrack)->phiSegment()->degreesOfFreedom()+2; DTChamberId dtChId(wheel,station,scsector); uint32_t indexCh = dtChId.rawId(); map<string, MonitorElement*> &innerChME = chamberHistos[indexCh]; map<string, MonitorElement*> &innerWhME = wheelHistos[wheel]; if (fabs(xdir)<30. && nHitsPhi>=7){ if (processDCC) { int qual = phCodeBestDCC[wheel+3][station][scsector]; innerWhME.find("DCC_TrigEffDenum")->second->Fill(scsector,station); if ( qual>=0 && qual<7 ) { innerWhME.find("DCC_TrigEffNum")->second->Fill(scsector,station); if ( qual>=4 ) { innerWhME.find("DCC_TrigEffCorrNum")->second->Fill(scsector,station); // if (qual==7 ) { // innerWhME.find("DCC_TrackPosvsAngleHH")->second->Fill(scsector,wheel); // } } } if (detailedPlots) { innerChME.find("DCC_TrackPosvsAngle")->second->Fill(xdir,x); if ( qual>=0 && qual<7 ) { innerChME.find("DCC_TrackPosvsAngleAnyQual")->second->Fill(xdir,x); if ( qual>=4 ) { innerChME.find("DCC_TrackPosvsAngleCorr")->second->Fill(xdir,x); // if (qual==7 ) { // innerChME.find("DCC_TrackPosvsAngleHH")->second->Fill(xdir,x); // } } } } } if (processDDU) { int qual = phCodeBestDDU[wheel+3][station][scsector]; innerWhME.find("DDU_TrigEffDenum")->second->Fill(scsector,station); bool qualOK = qual>=0 && qual<7; int bx = qualOK ? phBestDDU[wheel+3][station][scsector]->bx() : -10; if ( qualOK && bx>=minBXDDU && bx<=maxBXDDU ) { innerWhME.find("DDU_TrigEffNum")->second->Fill(scsector,station); if ( qual>=4 ) { innerWhME.find("DDU_TrigEffCorrNum")->second->Fill(scsector,station); // if (qual==7 ) { // innerWhME.find("DDU_TrackPosvsAngleHH")->second->Fill(scsector,wheel); // } } } if (detailedPlots) { innerChME.find("DDU_TrackPosvsAngle")->second->Fill(xdir,x); if ( qualOK && bx>+minBXDDU && bx<maxBXDDU ) { innerChME.find("DDU_TrackPosvsAngleAnyQual")->second->Fill(xdir,x); if ( qual>=4 ) { innerChME.find("DDU_TrackPosvsAngleCorr")->second->Fill(xdir,x); // if (qual==7 ) { // innerChME.find("DDU_TrackPosvsAngleHH")->second->Fill(xdir,x); // } } } } } } } }
void DTTriggerEfficiencyTask::beginJob | ( | void | ) | [protected, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 63 of file DTTriggerEfficiencyTask.cc.
References bookChamberHistos(), bookWheelHistos(), detailedPlots, edm::ParameterSet::getUntrackedParameter(), LogTrace, maxBXDDU, minBXDDU, nevents, parameters, processDCC, and processDDU.
{ LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: BeginJob" << endl; detailedPlots = parameters.getUntrackedParameter<bool>("detailedAnalysis",false); processDCC = parameters.getUntrackedParameter<bool>("processDCC",true); processDDU = parameters.getUntrackedParameter<bool>("processDDU",true); minBXDDU = parameters.getUntrackedParameter<int>("minBXDDU",0); maxBXDDU = parameters.getUntrackedParameter<int>("maxBXDDU",20); nevents = 0; for (int wh=-2;wh<=2;++wh){ if (processDCC) { bookWheelHistos(wh,"DCC_TrigEffDenum"); bookWheelHistos(wh,"DCC_TrigEffNum"); bookWheelHistos(wh,"DCC_TrigEffCorrNum"); if (detailedPlots) { for (int stat=1;stat<=4;++stat){ for (int sect=1;sect<=12;++sect){ DTChamberId dtChId(wh,stat,sect); bookChamberHistos(dtChId,"DCC_TrackPosvsAngleAnyQual","Segment"); bookChamberHistos(dtChId,"DCC_TrackPosvsAngleCorr","Segment"); bookChamberHistos(dtChId,"DCC_TrackPosvsAngle","Segment"); } } } } if (processDDU) { bookWheelHistos(wh,"DDU_TrigEffDenum"); bookWheelHistos(wh,"DDU_TrigEffNum"); bookWheelHistos(wh,"DDU_TrigEffCorrNum"); if (detailedPlots) { for (int stat=1;stat<=4;++stat){ for (int sect=1;sect<=12;++sect){ DTChamberId dtChId(wh,stat,sect); bookChamberHistos(dtChId,"DDU_TrackPosvsAngleAnyQual","Segment"); bookChamberHistos(dtChId,"DDU_TrackPosvsAngleCorr","Segment"); bookChamberHistos(dtChId,"DDU_TrackPosvsAngle","Segment"); } } } } } // end of wheel loop }
void DTTriggerEfficiencyTask::beginLuminosityBlock | ( | const edm::LuminosityBlock & | lumiSeg, |
const edm::EventSetup & | context | ||
) | [protected, virtual] |
To reset the MEs.
Reimplemented from edm::EDAnalyzer.
Definition at line 119 of file DTTriggerEfficiencyTask.cc.
References LogTrace.
{ LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") <<"[DTTriggerEfficiencyTask]: Begin of LS transition"<<endl; }
void DTTriggerEfficiencyTask::beginRun | ( | const edm::Run & | run, |
const edm::EventSetup & | context | ||
) | [protected, virtual] |
BeginRun.
Reimplemented from edm::EDAnalyzer.
Definition at line 110 of file DTTriggerEfficiencyTask.cc.
References edm::EventSetup::get(), LogTrace, muonGeom, and trigGeomUtils.
{ LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: BeginRun" << endl; context.get<MuonGeometryRecord>().get(muonGeom); trigGeomUtils = new DTTrigGeomUtils(muonGeom); }
void DTTriggerEfficiencyTask::bookChamberHistos | ( | const DTChamberId & | dtCh, |
std::string | histoTag, | ||
std::string | folder = "" |
||
) | [protected] |
Book chamber granularity histograms.
Definition at line 336 of file DTTriggerEfficiencyTask.cc.
References DQMStore::book2D(), chamberHistos, dbe, LogTrace, max(), min, RecoTauCommonJetSelections_cfi::nbins, DTTrigGeomUtils::phiRange(), DetId::rawId(), DTChamberId::sector(), DQMStore::setCurrentFolder(), relativeConstraints::station, DTChamberId::station(), topFolder(), trigGeomUtils, and DTChamberId::wheel().
Referenced by beginJob().
{ int wh = dtCh.wheel(); int sc = dtCh.sector(); int st = dtCh.station(); stringstream wheel; wheel << wh; stringstream station; station << st; stringstream sector; sector << sc; string histoType = histoTag.substr(4,histoTag.find("_",4)-4); string hwFolder = topFolder(histoTag.substr(0,3)=="DCC"); string bookingFolder = hwFolder + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() + "/" + folder; string histoName = histoTag + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str(); dbe->setCurrentFolder(bookingFolder); LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: booking " << bookingFolder << "/" << histoName << endl; if( histoType.find("TrackPosvsAngle") == 0 ){ float min, max; int nbins; trigGeomUtils->phiRange(dtCh,min,max,nbins,20); string histoLabel = "Position vs Angle (phi)"; if (histoType.find("Corr") != string::npos) histoLabel += " for correlated triggers"; else if (histoType.find("AnyQual") != string::npos) histoLabel += " for any qual triggers"; (chamberHistos[dtCh.rawId()])[histoTag] = dbe->book2D(histoName,histoLabel,12,-30.,30.,nbins,min,max); return ; } }
void DTTriggerEfficiencyTask::bookWheelHistos | ( | int | wheel, |
std::string | histoTag, | ||
std::string | folder = "" |
||
) | [protected] |
Book wheel granularity histograms.
Definition at line 282 of file DTLocalTriggerBaseTest.cc.
References category(), newFWLiteAna::fullName, LogTrace, MonitorElement::setAxisTitle(), and MonitorElement::setBinLabel().
Referenced by beginJob().
{ stringstream wh; wh << wheel; string basedir; bool isDCC = hwSource=="DCC" ; if (hTag.find("Summary") != string::npos ) { basedir = topFolder(isDCC); //Book summary histo outside wheel directories } else { basedir = topFolder(isDCC) + "Wheel" + wh.str() + "/" ; } if (folder != "") { basedir += folder +"/" ; } dbe->setCurrentFolder(basedir); string fullTag = fullName(hTag); string hname = fullTag+ "_W" + wh.str(); LogTrace(category()) << "[" << testName << "Test]: booking "<< basedir << hname; if (hTag.find("Phi")!= string::npos || hTag.find("Summary") != string::npos ){ MonitorElement* me = dbe->book2D(hname.c_str(),hname.c_str(),12,1,13,4,1,5); // setLabelPh(me); me->setBinLabel(1,"MB1",2); me->setBinLabel(2,"MB2",2); me->setBinLabel(3,"MB3",2); me->setBinLabel(4,"MB4",2); me->setAxisTitle("Sector",1); whME[wheel][fullTag] = me; return; } if (hTag.find("Theta") != string::npos){ MonitorElement* me =dbe->book2D(hname.c_str(),hname.c_str(),12,1,13,3,1,4); // setLabelTh(me); me->setBinLabel(1,"MB1",2); me->setBinLabel(2,"MB2",2); me->setBinLabel(3,"MB3",2); me->setAxisTitle("Sector",1); whME[wheel][fullTag] = me; return; } }
void DTTriggerEfficiencyTask::endJob | ( | void | ) | [protected, virtual] |
EndJob.
Reimplemented from edm::EDAnalyzer.
Definition at line 126 of file DTTriggerEfficiencyTask.cc.
bool DTTriggerEfficiencyTask::hasRPCTriggers | ( | const edm::Event & | e | ) | [protected] |
checks for RPC Triggers
Definition at line 311 of file DTTriggerEfficiencyTask.cc.
References edm::Event::getByLabel(), edm::ParameterSet::getUntrackedParameter(), and parameters.
Referenced by analyze().
{ InputTag inputTagGMT = parameters.getUntrackedParameter<edm::InputTag>("inputTagGMT"); edm::Handle<L1MuGMTReadoutCollection> gmtrc; e.getByLabel(inputTagGMT,gmtrc); std::vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords(); std::vector<L1MuGMTReadoutRecord>::const_iterator igmtrr; for(igmtrr=gmt_records.begin(); igmtrr!=gmt_records.end(); igmtrr++) { if( (*igmtrr).getBxInEvent()==0 ) { std::vector<L1MuRegionalCand>::const_iterator iter1; std::vector<L1MuRegionalCand> rmc = (*igmtrr).getBrlRPCCands(); for(iter1=rmc.begin(); iter1!=rmc.end(); iter1++) { if ( !(*iter1).empty() ) { return true; } } } } return false; }
std::string DTTriggerEfficiencyTask::topFolder | ( | bool | source | ) | [inline, protected] |
return the top folder 0=DDU 1=DCC
Definition at line 73 of file DTTriggerEfficiencyTask.h.
Referenced by bookChamberHistos().
{ return source ? "DT/03-LocalTrigger-DCC/" : "DT/04-LocalTrigger-DDU/"; }
std::map<uint32_t, std::map<std::string, MonitorElement*> > DTTriggerEfficiencyTask::chamberHistos [private] |
Definition at line 100 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and bookChamberHistos().
DQMStore* DTTriggerEfficiencyTask::dbe [private] |
Definition at line 96 of file DTTriggerEfficiencyTask.h.
Referenced by bookChamberHistos(), and DTTriggerEfficiencyTask().
bool DTTriggerEfficiencyTask::detailedPlots [private] |
Definition at line 88 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
int DTTriggerEfficiencyTask::maxBXDDU [private] |
Definition at line 89 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
int DTTriggerEfficiencyTask::minBXDDU [private] |
Definition at line 89 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 98 of file DTTriggerEfficiencyTask.h.
Referenced by beginRun().
int DTTriggerEfficiencyTask::nevents [private] |
Definition at line 86 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), beginJob(), endJob(), and ~DTTriggerEfficiencyTask().
Definition at line 97 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), beginJob(), DTTriggerEfficiencyTask(), and hasRPCTriggers().
const L1MuDTChambPhDigi* DTTriggerEfficiencyTask::phBestDCC[6][5][13] [private] |
Definition at line 93 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().
const DTLocalTrigger* DTTriggerEfficiencyTask::phBestDDU[6][5][13] [private] |
Definition at line 94 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().
int DTTriggerEfficiencyTask::phCodeBestDCC[6][5][13] [private] |
Definition at line 91 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().
int DTTriggerEfficiencyTask::phCodeBestDDU[6][5][13] [private] |
Definition at line 92 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().
bool DTTriggerEfficiencyTask::processDCC [private] |
Definition at line 88 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
bool DTTriggerEfficiencyTask::processDDU [private] |
Definition at line 88 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 99 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), beginRun(), and bookChamberHistos().
std::map<int, std::map<std::string, MonitorElement*> > DTTriggerEfficiencyTask::wheelHistos [private] |
Definition at line 101 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().