#include <DTLocalTriggerSynchTest.h>
Public Member Functions | |
DTLocalTriggerSynchTest (const edm::ParameterSet &ps) | |
Constructor. | |
virtual | ~DTLocalTriggerSynchTest () |
Destructor. | |
Protected Member Functions | |
void | beginJob () |
Begin Job. | |
void | beginRun (const edm::Run &run, const edm::EventSetup &c) |
begin Run | |
void | bookChambHistos (DTChamberId chambId, std::string htype, std::string subfolder="") |
Book the new MEs (for each chamber) | |
void | endJob () |
End Job. | |
float | getFloatFromME (DTChamberId chId, std::string meType) |
Get float MEs. | |
void | makeRatioME (TH1F *numerator, TH1F *denominator, MonitorElement *result) |
Compute efficiency plots. | |
void | runClientDiagnostic () |
DQM Client Diagnostic. | |
Private Attributes | |
double | bxTime |
std::map< uint32_t, std::map < std::string, MonitorElement * > > | chambME |
std::string | denHistoTag |
int | minEntries |
int | nBXHigh |
int | nBXLow |
std::string | numHistoTag |
bool | rangeInBX |
std::string | ratioHistoTag |
DTTPGParameters | wPhaseMap |
bool | writeDB |
* DQM Test Client
Definition at line 21 of file DTLocalTriggerSynchTest.h.
DTLocalTriggerSynchTest::DTLocalTriggerSynchTest | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 44 of file DTLocalTriggerSynchTest.cc.
{ setConfig(ps,"DTLocalTriggerSynch"); baseFolderDCC = "DT/90-LocalTriggerSynch/"; baseFolderDDU = "DT/90-LocalTriggerSynch/"; }
DTLocalTriggerSynchTest::~DTLocalTriggerSynchTest | ( | ) | [virtual] |
void DTLocalTriggerSynchTest::beginJob | ( | void | ) | [protected, virtual] |
Begin Job.
Reimplemented from DTLocalTriggerBaseTest.
Definition at line 58 of file DTLocalTriggerSynchTest.cc.
References Parameters::parameters.
{ numHistoTag = parameters.getParameter<string>("numHistoTag"); denHistoTag = parameters.getParameter<string>("denHistoTag"); ratioHistoTag = parameters.getParameter<string>("ratioHistoTag"); bxTime = parameters.getParameter<double>("bxTimeInterval"); rangeInBX = parameters.getParameter<bool>("rangeWithinBX"); nBXLow = parameters.getParameter<int>("nBXLow"); nBXHigh = parameters.getParameter<int>("nBXHigh"); minEntries = parameters.getParameter<int>("minEntries"); }
void DTLocalTriggerSynchTest::beginRun | ( | const edm::Run & | run, |
const edm::EventSetup & | c | ||
) | [protected, virtual] |
begin Run
Reimplemented from DTLocalTriggerBaseTest.
Definition at line 71 of file DTLocalTriggerSynchTest.cc.
References DTLocalTriggerBaseTest::beginRun(), argparse::category, edm::EventSetup::get(), and Parameters::parameters.
{ DTLocalTriggerBaseTest::beginRun(run,c); vector<string>::const_iterator iTr = trigSources.begin(); vector<string>::const_iterator trEnd = trigSources.end(); vector<string>::const_iterator iHw = hwSources.begin(); vector<string>::const_iterator hwEnd = hwSources.end(); //Booking if(parameters.getUntrackedParameter<bool>("staticBooking", true)){ for (; iTr != trEnd; ++iTr){ trigSource = (*iTr); for (; iHw != hwEnd; ++iHw){ hwSource = (*iHw); std::vector<DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin(); std::vector<DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end(); for (; chambIt!=chambEnd; ++chambIt) { DTChamberId chId = ((*chambIt)->id()); bookChambHistos(chId,ratioHistoTag); } } } } LogVerbatim(category()) << "[" << testName << "Test]: beginRun" << endl; if (parameters.getParameter<bool>("fineParamDiff")) { ESHandle<DTTPGParameters> wPhaseHandle; c.get<DTTPGParametersRcd>().get(wPhaseHandle); wPhaseMap = (*wPhaseHandle); } }
void DTLocalTriggerSynchTest::bookChambHistos | ( | DTChamberId | chambId, |
std::string | htype, | ||
std::string | subfolder = "" |
||
) | [protected] |
Book the new MEs (for each chamber)
void DTLocalTriggerSynchTest::endJob | ( | void | ) | [protected, virtual] |
End Job.
Reimplemented from DTLocalTriggerBaseTest.
Definition at line 155 of file DTLocalTriggerSynchTest.cc.
References DTTPGParameters::begin(), argparse::category, DTTPGParameters::end(), DTLocalTriggerBaseTest::endJob(), DTTimeUnits::ns, Parameters::parameters, DTTPGParameters::set(), and DTCalibDBUtils::writeToDB().
{ DTLocalTriggerBaseTest::endJob(); if ( parameters.getParameter<bool>("writeDB")) { LogVerbatim(category()) << "[" << testName << "Test]: writeDB flag set to true. Producing peak position database." << endl; DTTPGParameters* delayMap = new DTTPGParameters(); hwSource = parameters.getParameter<bool>("dbFromDCC") ? "DCC" : "DDU"; std::vector<DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin(); std::vector<DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end(); for (; chambIt!=chambEnd; ++chambIt) { DTChamberId chId = (*chambIt)->id(); float fineDelay = 0; int coarseDelay = static_cast<int>((getFloatFromME(chId,"tTrig_SL1") + getFloatFromME(chId,"tTrig_SL3"))*0.5/bxTime); bool fineDiff = parameters.getParameter<bool>("fineParamDiff"); bool coarseDiff = parameters.getParameter<bool>("coarseParamDiff"); TH1F *ratioH = getHisto<TH1F>(dbe->get(getMEName(ratioHistoTag,"", chId))); if (ratioH->GetEntries()>minEntries) { TF1 *fitF=ratioH->GetFunction("pol8"); if (fitF) { fineDelay=fitF->GetMaximumX(0,bxTime); } } else { LogInfo(category()) << "[" << testName << "Test]: Ratio histogram for chamber " << chId << " is empty. Worst Phase value set to 0." << endl; } if (fineDiff || coarseDiff) { float wFine; int wCoarse; wPhaseMap.get(chId,wCoarse,wFine,DTTimeUnits::ns); if (fineDiff) { fineDelay = wFine - fineDelay; } if (coarseDiff) { coarseDelay = wCoarse - coarseDelay; } } delayMap->set(chId,coarseDelay,fineDelay,DTTimeUnits::ns); } std::vector< std::pair<DTTPGParametersId,DTTPGParametersData> >::const_iterator dbIt = delayMap->begin(); std::vector< std::pair<DTTPGParametersId,DTTPGParametersData> >::const_iterator dbEnd = delayMap->end(); for (; dbIt!=dbEnd; ++dbIt) { LogVerbatim(category()) << "[" << testName << "Test]: DB entry for Wh " << (*dbIt).first.wheelId << " Sec " << (*dbIt).first.sectorId << " St " << (*dbIt).first.stationId << " has coarse " << (*dbIt).second.nClock << " and phase " << (*dbIt).second.tPhase << std::endl; } string delayRecord = "DTTPGParametersRcd"; DTCalibDBUtils::writeToDB(delayRecord,delayMap); } }
float DTLocalTriggerSynchTest::getFloatFromME | ( | DTChamberId | chId, |
std::string | meType | ||
) | [protected] |
Get float MEs.
Definition at line 223 of file DTLocalTriggerSynchTest.cc.
References argparse::category, MonitorElement::getFloatValue(), DTChamberId::sector(), relativeConstraints::station, DTChamberId::station(), and DTChamberId::wheel().
{ stringstream wheel; wheel << chId.wheel(); stringstream station; station << chId.station(); stringstream sector; sector << chId.sector(); string folderName = topFolder(hwSource=="DCC") + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() + "/" ; string histoname = sourceFolder + folderName + meType + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str(); MonitorElement* me = dbe->get(histoname); if (me) { return me->getFloatValue(); } else { LogProblem(category()) << "[" << testName << "Test]: " << histoname << " is not a valid ME. 0 returned" << std::endl; } return 0; }
void DTLocalTriggerSynchTest::makeRatioME | ( | TH1F * | numerator, |
TH1F * | denominator, | ||
MonitorElement * | result | ||
) | [protected] |
Compute efficiency plots.
Definition at line 216 of file DTLocalTriggerSynchTest.cc.
References postValidation_cfi::efficiency, and MonitorElement::getTH1F().
{ TH1F* efficiency = result->getTH1F(); efficiency->Divide(numerator,denominator,1,1,""); }
void DTLocalTriggerSynchTest::runClientDiagnostic | ( | ) | [protected, virtual] |
DQM Client Diagnostic.
Implements DTLocalTriggerBaseTest.
Definition at line 108 of file DTLocalTriggerSynchTest.cc.
References argparse::category, FitTarget::Fit, newFWLiteAna::fullName, MonitorElement::getName(), and DetId::rawId().
{ // Loop over Trig & Hw sources for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr){ trigSource = (*iTr); for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){ hwSource = (*iHw); std::vector<DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin(); std::vector<DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end(); for (; chambIt!=chambEnd; ++chambIt) { DTChamberId chId = (*chambIt)->id(); uint32_t indexCh = chId.rawId(); // Perform peak finding TH1F *numH = getHisto<TH1F>(dbe->get(getMEName(numHistoTag,"", chId))); TH1F *denH = getHisto<TH1F>(dbe->get(getMEName(denHistoTag,"", chId))); if (numH && denH && numH->GetEntries()>minEntries && denH->GetEntries()>minEntries) { std::map<std::string,MonitorElement*> innerME = chambME[indexCh]; MonitorElement* ratioH = innerME.find(fullName(ratioHistoTag))->second; makeRatioME(numH,denH,ratioH); try { getHisto<TH1F>(ratioH)->Fit("pol8","CQO"); } catch (...) { edm::LogPrint(category()) << "[" << testName << "Test]: Error fitting " << ratioH->getName() << " returned 0" << endl; } } else { if (!numH || !denH) { LogPrint(category()) << "[" << testName << "Test]: At least one of the required Histograms was not found for chamber " << chId << ". Peaks not computed" << endl; } else { LogPrint(category()) << "[" << testName << "Test]: Number of plots entries for " << chId << " is less than minEntries=" << minEntries <<". Peaks not computed" << endl; } } } } } }
double DTLocalTriggerSynchTest::bxTime [private] |
Definition at line 60 of file DTLocalTriggerSynchTest.h.
std::map<uint32_t,std::map<std::string,MonitorElement*> > DTLocalTriggerSynchTest::chambME [private] |
Definition at line 56 of file DTLocalTriggerSynchTest.h.
std::string DTLocalTriggerSynchTest::denHistoTag [private] |
Definition at line 58 of file DTLocalTriggerSynchTest.h.
int DTLocalTriggerSynchTest::minEntries [private] |
Definition at line 64 of file DTLocalTriggerSynchTest.h.
int DTLocalTriggerSynchTest::nBXHigh [private] |
Definition at line 63 of file DTLocalTriggerSynchTest.h.
int DTLocalTriggerSynchTest::nBXLow [private] |
Definition at line 62 of file DTLocalTriggerSynchTest.h.
std::string DTLocalTriggerSynchTest::numHistoTag [private] |
Definition at line 57 of file DTLocalTriggerSynchTest.h.
bool DTLocalTriggerSynchTest::rangeInBX [private] |
Definition at line 61 of file DTLocalTriggerSynchTest.h.
std::string DTLocalTriggerSynchTest::ratioHistoTag [private] |
Definition at line 59 of file DTLocalTriggerSynchTest.h.
Definition at line 66 of file DTLocalTriggerSynchTest.h.
bool DTLocalTriggerSynchTest::writeDB [private] |
Definition at line 65 of file DTLocalTriggerSynchTest.h.