#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 (std::string source) |
return the top folder | |
Private Attributes | |
std::map< uint32_t, std::map < std::string, MonitorElement * > > | chamberHistos |
DQMStore * | dbe |
bool | detailedPlots |
edm::InputTag | inputTagDCC |
edm::InputTag | inputTagDDU |
edm::InputTag | inputTagGMT |
edm::InputTag | inputTagMuons |
edm::InputTag | inputTagSEG |
int | maxBXDDU |
int | minBXDDU |
edm::ESHandle< DTGeometry > | muonGeom |
int | nevents |
int | nMinHitsPhi |
edm::ParameterSet | parameters |
float | phiAccRange |
bool | processDCC |
bool | processDDU |
std::vector< std::string > | processTags |
std::string | SegmArbitration |
DTTrigGeomUtils * | trigGeomUtils |
std::map< int, std::map < std::string, MonitorElement * > > | wheelHistos |
Definition at line 43 of file DTTriggerEfficiencyTask.h.
DTTriggerEfficiencyTask::DTTriggerEfficiencyTask | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 48 of file DTTriggerEfficiencyTask.cc.
References dbe, LogTrace, 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 132 of file DTTriggerEfficiencyTask.cc.
References reco::MuonSegmentMatch::BelongsToTrackByCleaning, reco::MuonSegmentMatch::BelongsToTrackByDR, reco::MuonSegmentMatch::BestInChamberByDR, chamberHistos, DTTrigGeomUtils::computeSCCoordinates(), detailedPlots, GeomDetEnumerators::DT, edm::Event::getByLabel(), hasRPCTriggers(), inputTagDCC, inputTagDDU, inputTagMuons, RPCpg::mu, patZpeak::muons, nevents, nMinHitsPhi, phiAccRange, processTags, DetId::rawId(), SegmArbitration, relativeConstraints::station, trigGeomUtils, wheelHistos, x, xdir, detailsBasic3DVector::y, and ydir.
{ nevents++; if (!hasRPCTriggers(e)) { return; } map<DTChamberId,const L1MuDTChambPhDigi*> phBestDCC; map<DTChamberId,const DTLocalTrigger*> phBestDDU; // 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(); DTChamberId chId(phwheel,phst,phsec); if( phcode < 7 && (phBestDCC.find(chId) == phBestDCC.end() || phcode>phBestDCC[chId]->code()) ) phBestDCC[chId] = &(*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; DTLocalTriggerCollection::const_iterator trigIt = range.first; DTLocalTriggerCollection::const_iterator trigEnd = range.second; for (; trigIt!= trigEnd;++trigIt){ int quality = trigIt->quality(); if(quality>-1 && quality<7 && (phBestDDU.find(id) == phBestDDU.end() || quality>phBestDDU[id]->quality()) ) phBestDDU[id] = &(*trigIt); } } //Getting Best Segments vector<const DTRecSegment4D*> best4DSegments; Handle<reco::MuonCollection> muons; e.getByLabel(inputTagMuons, muons); reco::MuonCollection::const_iterator mu; for( mu = muons->begin(); mu != muons->end(); ++mu ) { // Make sure that is standalone muon if( !((*mu).isStandAloneMuon()) ) {continue;} // Get the chambers compatible with the muon const vector<reco::MuonChamberMatch> matchedChambers = (*mu).matches(); vector<reco::MuonChamberMatch>::const_iterator chamber; for( chamber = matchedChambers.begin(); chamber != matchedChambers.end(); ++chamber ) { // look only in DTs if( chamber->detector() != MuonSubdetId::DT ) {continue;} // Get the matched segments in the chamber const vector<reco::MuonSegmentMatch> matchedSegments = (*chamber).segmentMatches; vector<reco::MuonSegmentMatch>::const_iterator segment; for( segment = matchedSegments.begin(); segment != matchedSegments.end(); ++segment ) { edm::Ref<DTRecSegment4DCollection> dtSegment = segment->dtSegmentRef; // Segment Arbitration if( SegmArbitration == "SegmentArbitration" && !((*segment).isMask(reco::MuonSegmentMatch::BestInChamberByDR)) ) {continue;} if( SegmArbitration == "SegmentAndTrackArbitration" && (!((*segment).isMask(reco::MuonSegmentMatch::BestInChamberByDR)) || !((*segment).isMask(reco::MuonSegmentMatch::BelongsToTrackByDR))) ) {continue;} if( SegmArbitration == "SegmentAndTrackArbitrationCleaned" && (!((*segment).isMask(reco::MuonSegmentMatch::BestInChamberByDR)) || !((*segment).isMask(reco::MuonSegmentMatch::BelongsToTrackByDR)) || !((*segment).isMask(reco::MuonSegmentMatch::BelongsToTrackByCleaning))) ) {continue;} if( (*dtSegment).hasPhi() ) { best4DSegments.push_back(&(*dtSegment)); } }// end loop on matched segments }// end loop on compatible chambers }// end loop on muons // 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)<phiAccRange && nHitsPhi>=nMinHitsPhi){ vector<string>::const_iterator tagIt = processTags.begin(); vector<string>::const_iterator tagEnd = processTags.end(); for (; tagIt!=tagEnd; ++tagIt) { int qual = (*tagIt) == "DCC" ? phBestDCC.find(dtChId) != phBestDCC.end() ? phBestDCC[dtChId]->code() : -1 : phBestDDU.find(dtChId) != phBestDDU.end() ? phBestDDU[dtChId]->quality() : -1; innerWhME.find((*tagIt) + "_TrigEffDenum")->second->Fill(scsector,station); if ( qual>=0 && qual<7 ) { innerWhME.find((*tagIt) + "_TrigEffNum")->second->Fill(scsector,station); if ( qual>=4 ) { innerWhME.find((*tagIt) + "_TrigEffCorrNum")->second->Fill(scsector,station); } } if (detailedPlots) { innerChME.find((*tagIt) + "_TrackPosvsAngle")->second->Fill(xdir,x); if ( qual>=0 && qual<7 ) { innerChME.find((*tagIt) + "_TrackPosvsAngleAnyQual")->second->Fill(xdir,x); if ( qual>=4 ) { innerChME.find((*tagIt) + "_TrackPosvsAngleCorr")->second->Fill(xdir,x); } } } } } } }
void DTTriggerEfficiencyTask::beginJob | ( | void | ) | [protected, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 65 of file DTTriggerEfficiencyTask.cc.
References detailedPlots, edm::ParameterSet::getUntrackedParameter(), inputTagDCC, inputTagDDU, inputTagGMT, inputTagMuons, inputTagSEG, LogTrace, maxBXDDU, minBXDDU, nMinHitsPhi, parameters, phiAccRange, processDCC, processDDU, processTags, SegmArbitration, and AlCaHLTBitMon_QueryRunRegistry::string.
{ LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: BeginJob" << endl; inputTagMuons = parameters.getUntrackedParameter<edm::InputTag>("inputTagMuons"); inputTagDCC = parameters.getUntrackedParameter<edm::InputTag>("inputTagDCC"); inputTagDDU = parameters.getUntrackedParameter<edm::InputTag>("inputTagDDU"); inputTagSEG = parameters.getUntrackedParameter<edm::InputTag>("inputTagSEG"); inputTagGMT = parameters.getUntrackedParameter<edm::InputTag>("inputTagGMT"); SegmArbitration = parameters.getUntrackedParameter<std::string>("SegmArbitration"); detailedPlots = parameters.getUntrackedParameter<bool>("detailedAnalysis"); processDCC = parameters.getUntrackedParameter<bool>("processDCC"); processDDU = parameters.getUntrackedParameter<bool>("processDDU"); minBXDDU = parameters.getUntrackedParameter<int>("minBXDDU"); maxBXDDU = parameters.getUntrackedParameter<int>("maxBXDDU"); nMinHitsPhi = parameters.getUntrackedParameter<int>("nMinHitsPhi"); phiAccRange = parameters.getUntrackedParameter<double>("phiAccRange"); if (processDCC) processTags.push_back("DCC"); if (processDDU) processTags.push_back("DDU"); }
void DTTriggerEfficiencyTask::beginLuminosityBlock | ( | const edm::LuminosityBlock & | lumiSeg, |
const edm::EventSetup & | context | ||
) | [protected, virtual] |
To reset the MEs.
Reimplemented from edm::EDAnalyzer.
Definition at line 118 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 92 of file DTTriggerEfficiencyTask.cc.
References bookChamberHistos(), bookWheelHistos(), detailedPlots, DTChamberId, edm::EventSetup::get(), LogTrace, muonGeom, nevents, processTags, and trigGeomUtils.
{ LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: BeginRun" << endl; nevents = 0; for (int wh=-2;wh<=2;++wh){ vector<string>::const_iterator tagIt = processTags.begin(); vector<string>::const_iterator tagEnd = processTags.end(); for (; tagIt!=tagEnd; ++tagIt) { bookWheelHistos(wh,(*tagIt),"Task"); if (detailedPlots) { for (int stat=1;stat<=4;++stat){ for (int sect=1;sect<=12;++sect){ bookChamberHistos(DTChamberId(wh,stat,sect),(*tagIt),"Segment"); } } } } } 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 201 of file DTRunConditionVar.cc.
References DQMStore::book1D(), DTRunConditionVar::chamberHistos, LogTrace, DetId::rawId(), DTChamberId::sector(), DQMStore::setCurrentFolder(), relativeConstraints::station, DTChamberId::station(), DTRunConditionVar::theDbe, and DTChamberId::wheel().
Referenced by beginRun().
{ 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 bookingFolder = "DT/02-Segments/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str(); string histoTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str(); theDbe->setCurrentFolder(bookingFolder); LogTrace ("DTDQM|DTMonitorModule|DTRunConditionVar") << "[DTRunConditionVar]: booking histos in " << bookingFolder << endl; string histoName = histoType + histoTag; string histoLabel = histoType; (chamberHistos[dtCh.rawId()])[histoType] = theDbe->book1D(histoName,histoLabel,nbins,min,max); }
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 python::rootplot::argparse::category, newFWLiteAna::fullName, LogTrace, MonitorElement::setAxisTitle(), and MonitorElement::setBinLabel().
Referenced by beginRun().
{ 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 125 of file DTTriggerEfficiencyTask.cc.
bool DTTriggerEfficiencyTask::hasRPCTriggers | ( | const edm::Event & | e | ) | [protected] |
checks for RPC Triggers
Definition at line 281 of file DTTriggerEfficiencyTask.cc.
References edm::Event::getByLabel(), and inputTagGMT.
Referenced by analyze().
{ edm::Handle<L1MuGMTReadoutCollection> gmtrc; e.getByLabel(inputTagGMT,gmtrc); std::vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords(); std::vector<L1MuGMTReadoutRecord>::const_iterator igmtrr = gmt_records.begin(); std::vector<L1MuGMTReadoutRecord>::const_iterator egmtrr = gmt_records.end(); for(; igmtrr!=egmtrr; igmtrr++) { std::vector<L1MuGMTExtendedCand> candsGMT = igmtrr->getGMTCands(); std::vector<L1MuGMTExtendedCand>::const_iterator candGMTIt = candsGMT.begin(); std::vector<L1MuGMTExtendedCand>::const_iterator candGMTEnd = candsGMT.end(); for(; candGMTIt!=candGMTEnd; ++candGMTIt){ if(!candGMTIt->empty()) { int quality = candGMTIt->quality(); if(candGMTIt->bx()==0 && (quality == 5 || quality == 7)){ return true; } } } } return false; }
std::string DTTriggerEfficiencyTask::topFolder | ( | std::string | source | ) | [inline, protected] |
return the top folder
Definition at line 71 of file DTTriggerEfficiencyTask.h.
{ return source=="DCC" ? "DT/03-LocalTrigger-DCC/" : "DT/04-LocalTrigger-DDU/"; }
std::map<uint32_t, std::map<std::string, MonitorElement*> > DTTriggerEfficiencyTask::chamberHistos [private] |
Definition at line 107 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().
DQMStore* DTTriggerEfficiencyTask::dbe [private] |
Definition at line 103 of file DTTriggerEfficiencyTask.h.
Referenced by DTTriggerEfficiencyTask().
bool DTTriggerEfficiencyTask::detailedPlots [private] |
Definition at line 88 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), beginJob(), and beginRun().
Definition at line 97 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 98 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 101 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob(), and hasRPCTriggers().
Definition at line 95 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 99 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob().
int DTTriggerEfficiencyTask::maxBXDDU [private] |
Definition at line 90 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob().
int DTTriggerEfficiencyTask::minBXDDU [private] |
Definition at line 90 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob().
Definition at line 105 of file DTTriggerEfficiencyTask.h.
Referenced by beginRun().
int DTTriggerEfficiencyTask::nevents [private] |
Definition at line 84 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), beginRun(), endJob(), and ~DTTriggerEfficiencyTask().
int DTTriggerEfficiencyTask::nMinHitsPhi [private] |
Definition at line 93 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 104 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob(), and DTTriggerEfficiencyTask().
float DTTriggerEfficiencyTask::phiAccRange [private] |
Definition at line 92 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
bool DTTriggerEfficiencyTask::processDCC [private] |
Definition at line 88 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob().
bool DTTriggerEfficiencyTask::processDDU [private] |
Definition at line 88 of file DTTriggerEfficiencyTask.h.
Referenced by beginJob().
std::vector<std::string> DTTriggerEfficiencyTask::processTags [private] |
Definition at line 89 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), beginJob(), and beginRun().
std::string DTTriggerEfficiencyTask::SegmArbitration [private] |
Definition at line 86 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginJob().
Definition at line 106 of file DTTriggerEfficiencyTask.h.
Referenced by analyze(), and beginRun().
std::map<int, std::map<std::string, MonitorElement*> > DTTriggerEfficiencyTask::wheelHistos [private] |
Definition at line 108 of file DTTriggerEfficiencyTask.h.
Referenced by analyze().