#include <DTOccupancyTest.h>
Public Member Functions | |
DTOccupancyTest (const edm::ParameterSet &ps) | |
Constructor. | |
virtual | ~DTOccupancyTest () |
Destructor. | |
Protected Member Functions | |
void | analyze (const edm::Event &event, const edm::EventSetup &context) |
Analyze. | |
void | beginJob () |
BeginJob. | |
void | beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) |
void | beginRun (edm::Run const &run, edm::EventSetup const &context) |
BeginRun. | |
void | endJob () |
Endjob. | |
void | endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context) |
DQM Client Diagnostic. | |
Private Member Functions | |
void | bookHistos (const int wheelId, std::string folder, std::string histoTag) |
book the summary histograms | |
std::string | getMEName (std::string histoTag, const DTChamberId &chId) |
Get the ME name. | |
int | runOccupancyTest (TH2F *histo, const DTChamberId &chId, float &chamberPercentage) |
std::string | topFolder () const |
Private Attributes | |
DQMStore * | dbe |
MonitorElement * | glbSummaryHisto |
int | lsCounter |
std::set< DTLayerId > | monitoredLayers |
edm::ESHandle< DTGeometry > | muonGeom |
std::string | nameMonitoredHisto |
int | nevents |
int | nMinEvts |
TNtuple * | ntuple |
TFile * | rootFile |
bool | runOnAllHitsOccupancies |
bool | runOnInTimeOccupancies |
bool | runOnNoiseOccupancies |
MonitorElement * | summaryHisto |
bool | tpMode |
std::map< int, MonitorElement * > | wheelHistos |
bool | writeRootFile |
* DQM Test Client
Definition at line 35 of file DTOccupancyTest.h.
DTOccupancyTest::DTOccupancyTest | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 33 of file DTOccupancyTest.cc.
References edm::ParameterSet::getUntrackedParameter(), and cppFunctionSkipper::operator.
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest]: Constructor"; // Get the DQM service dbe = Service<DQMStore>().operator->(); lsCounter = 0; writeRootFile = ps.getUntrackedParameter<bool>("writeRootFile", false); if(writeRootFile) { rootFile = new TFile("DTOccupancyTest.root","RECREATE"); ntuple = new TNtuple("OccupancyNtuple", "OccupancyNtuple", "ls:wh:st:se:lay1MeanCell:lay1RMS:lay2MeanCell:lay2RMS:lay3MeanCell:lay3RMS:lay4MeanCell:lay4RMS:lay5MeanCell:lay5RMS:lay6MeanCell:lay6RMS:lay7MeanCell:lay7RMS:lay8MeanCell:lay8RMS:lay9MeanCell:lay9RMS:lay10MeanCell:lay10RMS:lay11MeanCell:lay11RMS:lay12MeanCell:lay12RMS"); } // switch on the mode for running on test pulses (different top folder) tpMode = ps.getUntrackedParameter<bool>("testPulseMode", false); runOnAllHitsOccupancies = ps.getUntrackedParameter<bool>("runOnAllHitsOccupancies", true); runOnNoiseOccupancies = ps.getUntrackedParameter<bool>("runOnNoiseOccupancies", false); runOnInTimeOccupancies = ps.getUntrackedParameter<bool>("runOnInTimeOccupancies", false); nMinEvts = ps.getUntrackedParameter<int>("nEventsCert", 5000); }
DTOccupancyTest::~DTOccupancyTest | ( | ) | [virtual] |
Destructor.
Definition at line 60 of file DTOccupancyTest.cc.
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << " destructor called" << endl; }
void DTOccupancyTest::analyze | ( | const edm::Event & | event, |
const edm::EventSetup & | context | ||
) | [protected, virtual] |
Analyze.
Implements edm::EDAnalyzer.
Definition at line 132 of file DTOccupancyTest.cc.
References nevents.
{ nevents++; // if(nevents%1000) // LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest]: "<<nevents<<" events"; }
void DTOccupancyTest::beginJob | ( | void | ) | [protected, virtual] |
BeginJob.
Reimplemented from edm::EDAnalyzer.
Definition at line 69 of file DTOccupancyTest.cc.
References bookHistos(), nevents, and indexGen::title.
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest]: BeginJob"; // Event counter nevents = 0; // Book the summary histos // - one summary per wheel for(int wh = -2; wh <= 2; ++wh) { // loop over wheels bookHistos(wh, string("Occupancies"), "OccupancySummary"); } dbe->setCurrentFolder(topFolder()); string title = "Occupancy Summary"; if(tpMode) { title = "Test Pulse Occupancy Summary"; } // - global summary with alarms summaryHisto = dbe->book2D("OccupancySummary",title.c_str(),12,1,13,5,-2,3); summaryHisto->setAxisTitle("sector",1); summaryHisto->setAxisTitle("wheel",2); // - global summary with percentages glbSummaryHisto = dbe->book2D("OccupancyGlbSummary",title.c_str(),12,1,13,5,-2,3); glbSummaryHisto->setAxisTitle("sector",1); glbSummaryHisto->setAxisTitle("wheel",2); // assign the name of the input histogram if(runOnAllHitsOccupancies) { nameMonitoredHisto = "OccupancyAllHits_perCh"; } else if(runOnNoiseOccupancies) { nameMonitoredHisto = "OccupancyNoise_perCh"; } else if(runOnInTimeOccupancies) { nameMonitoredHisto = "OccupancyInTimeHits_perCh"; } else { // default is AllHits histo nameMonitoredHisto = "OccupancyAllHits_perCh"; } }
void DTOccupancyTest::beginLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, |
edm::EventSetup const & | context | ||
) | [protected, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 125 of file DTOccupancyTest.cc.
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") <<"[DTOccupancyTest]: Begin of LS transition"; }
void DTOccupancyTest::beginRun | ( | edm::Run const & | run, |
edm::EventSetup const & | context | ||
) | [protected, virtual] |
BeginRun.
Reimplemented from edm::EDAnalyzer.
Definition at line 113 of file DTOccupancyTest.cc.
References edm::EventSetup::get().
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest]: BeginRun"; // Get the geometry context.get<MuonGeometryRecord>().get(muonGeom); }
void DTOccupancyTest::bookHistos | ( | const int | wheelId, |
std::string | folder, | ||
std::string | histoTag | ||
) | [private] |
book the summary histograms
void DTOccupancyTest::endJob | ( | void | ) | [protected, virtual] |
Endjob.
Reimplemented from edm::EDAnalyzer.
Definition at line 222 of file DTOccupancyTest.cc.
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest] endjob called!"; if(writeRootFile) { rootFile->cd(); ntuple->Write(); rootFile->Close(); } }
void DTOccupancyTest::endLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, |
edm::EventSetup const & | context | ||
) | [protected, virtual] |
DQM Client Diagnostic.
Reimplemented from edm::EDAnalyzer.
Definition at line 141 of file DTOccupancyTest.cc.
References chambers, MonitorElement::getFloatValue(), MonitorElement::getTH2F(), interpolateCardsSimple::histo, query::result, DTChamberId::sector(), DTChamberId::station(), and DTChamberId::wheel().
{ LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") <<"[DTOccupancyTest]: End of LS transition, performing the DQM client operation"; lsCounter++; // Reset the global summary summaryHisto->Reset(); glbSummaryHisto->Reset(); // Get all the DT chambers vector<DTChamber*> chambers = muonGeom->chambers(); for(vector<DTChamber*>::const_iterator chamber = chambers.begin(); chamber != chambers.end(); ++chamber) { // Loop over all chambers DTChamberId chId = (*chamber)->id(); MonitorElement * chamberOccupancyHisto = dbe->get(getMEName(nameMonitoredHisto, chId)); // Run the tests on the plot for the various granularities if(chamberOccupancyHisto != 0) { // Get the 2D histo TH2F* histo = chamberOccupancyHisto->getTH2F(); float chamberPercentage = 1.; int result = runOccupancyTest(histo, chId, chamberPercentage); int sector = chId.sector(); if(sector == 13) { sector = 4; float resultSect4 = wheelHistos[chId.wheel()]->getBinContent(sector, chId.station()); if(resultSect4 > result) { result = (int)resultSect4; } } else if(sector == 14) { sector = 10; float resultSect10 = wheelHistos[chId.wheel()]->getBinContent(sector, chId.station()); if(resultSect10 > result) { result = (int)resultSect10; } } // the 2 MB4 of Sect 4 and 10 count as half a chamber if((sector == 4 || sector == 10) && chId.station() == 4) chamberPercentage = chamberPercentage/2.; wheelHistos[chId.wheel()]->setBinContent(sector, chId.station(),result); if(result > summaryHisto->getBinContent(sector, chId.wheel()+3)) { summaryHisto->setBinContent(sector, chId.wheel()+3, result); } glbSummaryHisto->Fill(sector, chId.wheel(), chamberPercentage*1./4.); } else { LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest] ME: " << getMEName(nameMonitoredHisto, chId) << " not found!" << endl; } } string nEvtsName = "DT/EventInfo/Counters/nProcessedEventsDigi"; MonitorElement * meProcEvts = dbe->get(nEvtsName); if (meProcEvts) { int nProcEvts = meProcEvts->getFloatValue(); glbSummaryHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts); summaryHisto->setEntries(nProcEvts < nMinEvts ? 10. : nProcEvts); } else { glbSummaryHisto->setEntries(nMinEvts +1); summaryHisto->setEntries(nMinEvts + 1); LogVerbatim ("DTDQM|DTMonitorClient|DTOccupancyTest") << "[DTOccupancyTest] ME: " << nEvtsName << " not found!" << endl; } // Fill the global summary // Check for entire sectors off and report them on the global summary //FIXME: TODO if(writeRootFile) ntuple->AutoSave("SaveSelf"); }
std::string DTOccupancyTest::getMEName | ( | std::string | histoTag, |
const DTChamberId & | chId | ||
) | [private] |
Get the ME name.
int DTOccupancyTest::runOccupancyTest | ( | TH2F * | histo, |
const DTChamberId & | chId, | ||
float & | chamberPercentage | ||
) | [private] |
Definition at line 299 of file DTOccupancyTest.cc.
References DTOccupancyClusterBuilder::addPoint(), DTOccupancyCluster::averageMean(), DTOccupancyCluster::averageRMS(), DTOccupancyClusterBuilder::buildClusters(), DTOccupancyClusterBuilder::getBestCluster(), getHLTprescales::index, DTOccupancyClusterBuilder::isProblematic(), LogTrace, point, plotscripts::rms(), DTChamberId::sector(), mathSSE::sqrt(), DTChamberId::station(), makeHLTPrescaleTable::values, and DTChamberId::wheel().
{ int nBinsX = histo->GetNbinsX(); // Reset the error flags bool failSL = false; bool failLayer = false; bool failCells = false; // Check that the chamber has digis if(histo->Integral() == 0) { chamberPercentage = 0; return 4; } LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << "--- Occupancy test for chamber: " << chId << endl; // set the # of SLs int nSL = 3; if(chId.station() == 4) nSL = 2; // float values[28]; if(writeRootFile) { values[0] = lsCounter; values[1] = chId.wheel(); values[2] = chId.station(); values[3] = chId.sector(); } // Compute the average occupancy per layer and its RMS // we also look of the layer with the smallest RMS in order to find a reference value // for the cell occupancy double totalChamberOccupp = 0; double squaredLayerOccupSum = 0; map<DTLayerId, pair<double, double> > averageCellOccupAndRMS; map<DTLayerId, double> layerOccupancyMap; int index = 3; for(int slay = 1; slay <= 3; ++slay) { // loop over SLs // Skip layer 2 on MB4 if(chId.station() == 4 && slay == 2) { if(writeRootFile) { values[12] = -1; values[13] = -1; values[14] = -1; values[15] = -1; values[16] = -1; values[17] = -1; values[18] = -1; values[19] = -1; } index = 19; continue; } // check the SL occupancy int binYlow = ((slay-1)*4)+1; int binYhigh = binYlow+3; double slInteg = histo->Integral(1,nBinsX,binYlow,binYhigh); if(slInteg == 0) { chamberPercentage = 1.-1./(float)nSL; return 3; } for(int lay = 1; lay <= 4; ++lay) { // loop over layers DTLayerId layID(chId,slay,lay); int binY = binYlow+(lay-1); double layerInteg = histo->Integral(1,nBinsX,binY,binY); squaredLayerOccupSum += layerInteg*layerInteg; totalChamberOccupp+= layerInteg; layerOccupancyMap[layID] = layerInteg; // We look for the distribution of hits within the layer int nWires = muonGeom->layer(layID)->specificTopology().channels(); int firstWire = muonGeom->layer(layID)->specificTopology().firstChannel(); double layerSquaredSum = 0; // reset the alert bit in the plot (used by render plugins) histo->SetBinContent(nBinsX+1,binY,0.); for(int cell = firstWire; cell != (nWires+firstWire); ++cell) { // loop over cells double cellOccup = histo->GetBinContent(cell,binY); layerSquaredSum+=cellOccup*cellOccup; } // compute the average cell occpuancy and RMS double averageCellOccup = layerInteg/nWires; double averageSquaredCellOccup = layerSquaredSum/nWires; double rmsCellOccup = sqrt(averageSquaredCellOccup - averageCellOccup*averageCellOccup); averageCellOccupAndRMS[layID] = make_pair(averageCellOccup, rmsCellOccup); LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " " << layID << " average cell occ.: " << averageCellOccup << " RMS: " << rmsCellOccup << endl; if(writeRootFile) { index++; values[index] = averageCellOccup; index++; values[index] = rmsCellOccup; } } } if(writeRootFile) ntuple->Fill(values); // double averageLayerOcc = totalChamberOccupp/(nSL*4); // double averageSquaredLayeroccup = squaredLayerOccupSum/(nSL*4); // double layerOccupRMS = sqrt(averageSquaredLayeroccup - averageLayerOcc*averageLayerOcc); double minCellRMS = 99999999; double referenceCellOccup = -1; DTOccupancyClusterBuilder builder; // find the cell reference value for(map<DTLayerId, pair<double, double> >::const_iterator layAndValues = averageCellOccupAndRMS.begin(); layAndValues != averageCellOccupAndRMS.end(); layAndValues++) { DTLayerId lid = (*layAndValues).first; double rms = (*layAndValues).second.second; double lOcc = layerOccupancyMap[lid]; // FIXME: useless double avCellOcc = (*layAndValues).second.first; LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " " << lid << " tot. occ: " << lOcc << " average cell occ: " << avCellOcc << " RMS: " << rms << endl; if(avCellOcc != 0) { DTOccupancyPoint point(avCellOcc, rms, lid); builder.addPoint(point); } else { if(monitoredLayers.find(lid) == monitoredLayers.end()) monitoredLayers.insert(lid); } } builder.buildClusters(); referenceCellOccup = builder.getBestCluster().averageMean(); minCellRMS = builder.getBestCluster().averageRMS(); // set<DTLayerId> bestLayers getLayerIDs() double safeFactor = 3.; // if(minCellRMS > referenceCellOccup) safeFactor = 5; LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " Reference cell occup.: " << referenceCellOccup << " RMS: " << minCellRMS << endl; // Set a warning for particularly high RMS: noise can "mask" dead channels // bool rmsWarning = false; // if(layerOccupRMS > averageLayerOcc) { // cout << " Warning RMS is too big: monitoring all layers" << endl; // rmsWarning = true; // } int nFailingSLs = 0; // Check the layer occupancy for(int slay = 1; slay <= 3; ++slay) { // loop over SLs // Skip layer 2 on MB4 if(chId.station() == 4 && slay == 2) continue; int binYlow = ((slay-1)*4)+1; // int binYhigh = binYlow+3; int nFailingLayers = 0; for(int lay = 1; lay <= 4; ++lay) { // loop over layers DTLayerId layID(chId,slay,lay); int nWires = muonGeom->layer(layID)->specificTopology().channels(); int firstWire = muonGeom->layer(layID)->specificTopology().firstChannel(); int binY = binYlow+(lay-1); // compute the integral of the layer occupancy double layerInteg = histo->Integral(1,nBinsX,binY,binY); LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " layer: " << layID << " integral: " << layerInteg << endl; // Check if in the list of layers which are monitored bool alreadyMonitored = false; if(monitoredLayers.find(layID) != monitoredLayers.end()) alreadyMonitored = true; if(layerInteg == 0) { // layer is dead (no need to go further LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " fail layer: no entries" << endl; // Add it to the list of of monitored layers if(!alreadyMonitored) monitoredLayers.insert(layID); nFailingLayers++; failLayer = true; histo->SetBinContent(nBinsX+1,binY,-1.); // go to next layer continue; } // double avCellOccInLayer = averageCellOccupAndRMS[layID].first; // double cellOccupRMS = averageCellOccupAndRMS[layID].second; // if(monitoredLayers.find(layID) != monitoredLayers.end() || // layerInteg == 0 || // layerInteg < (averageLayerOcc - 3*layerOccupRMS) || // cellOccupRMS > avCellOccInLayer || // avCellOccInLayer < referenceCellOccup/3.) { // check the layer if(alreadyMonitored || builder.isProblematic(layID)) { // check the layer // Add it to the list of of monitored layers if(monitoredLayers.find(layID) == monitoredLayers.end()) monitoredLayers.insert(layID); // if(layerInteg != 0) { // check # of dead cells int totalDeadCells = 0; int nDeadCellsInARow = 1; int nDeadCellsInARowMax = 0; int nCellsZeroCount = 0; bool previousIsDead = false; int interDeadCells = 0; for(int cell = firstWire; cell != (nWires+firstWire); ++cell) { // loop over cells double cellOccup = histo->GetBinContent(cell,binY); LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " cell occup: " << cellOccup; if(cellOccup == 0 || cellOccup < (referenceCellOccup-safeFactor*sqrt(referenceCellOccup))) { if(cellOccup == 0) nCellsZeroCount++; totalDeadCells++; if(previousIsDead) nDeadCellsInARow++; previousIsDead = true; interDeadCells = 0; LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " below reference" << endl; } else { previousIsDead = false; interDeadCells++; // 3 cells not dead between a group of dead cells don't break the count if(interDeadCells > 3) { if(nDeadCellsInARow > nDeadCellsInARowMax) nDeadCellsInARowMax = nDeadCellsInARow; nDeadCellsInARow = 1; } } } if(nDeadCellsInARow > nDeadCellsInARowMax) nDeadCellsInARowMax = nDeadCellsInARow; LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " # wires: " << nWires << " # cells 0 count: " << nCellsZeroCount << " # dead cells in a row: " << nDeadCellsInARowMax << " total # of dead cells: " << totalDeadCells; // Count dead cells // if(TMath::Erf(referenceCellOccup/sqrt(referenceCellOccup)) > 2./3. && // nDeadCellsInARowMax > nWires/3. // && nDeadCellsInARowMax < 2*nWires/3.) { // cout << " -> fail cells!" << endl; // failCells = true; // histo->SetBinContent(nBinsX+1,binY,-1.); // } else if((TMath::Erfc(referenceCellOccup/sqrt(referenceCellOccup)) < 10./(double)nWires && nDeadCellsInARowMax>= 10.) || (TMath::Erfc(referenceCellOccup/sqrt(referenceCellOccup)) < 0.5 && totalDeadCells > nWires/2.)) { LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " -> fail layer!" << endl; nFailingLayers++; failLayer = true; histo->SetBinContent(nBinsX+1,binY,-1.); } else if(referenceCellOccup > 10 && nCellsZeroCount > nWires/3. && (double)nCellsZeroCount/(double)nWires > 2.*TMath::Erfc(referenceCellOccup/sqrt(referenceCellOccup))) { LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " -> would fail cells!" << endl; LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " # of cells with 0 count: " << nCellsZeroCount << " # wires: " << nWires << " erfc: " << TMath::Erfc(referenceCellOccup/sqrt(referenceCellOccup)) << endl; // failCells = true; // histo->SetBinContent(nBinsX+1,binY,-1.); } // } else { // all layer is dead // LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest") << " fail layer: no entries" << endl; // nFailingLayers++; // failLayer = true; // histo->SetBinContent(nBinsX+1,binY,-1.); // } } } // Check if the whole layer is off if( nFailingLayers == 4) { nFailingSLs++; failSL = true; } } // All the chamber is off if(nFailingSLs == nSL) { chamberPercentage = 0; return 4; } else { chamberPercentage = 1.-(float)nFailingSLs/(float)nSL; } // FIXME add check on cells if(failSL) return 3; if(failLayer) return 2; if(failCells) return 1; return 0; }
string DTOccupancyTest::topFolder | ( | ) | const [private] |
Definition at line 610 of file DTOccupancyTest.cc.
{ if(tpMode) return string("DT/10-TestPulses/"); return string("DT/01-Digi/"); }
DQMStore* DTOccupancyTest::dbe [private] |
Definition at line 80 of file DTOccupancyTest.h.
MonitorElement* DTOccupancyTest::glbSummaryHisto [private] |
Definition at line 87 of file DTOccupancyTest.h.
int DTOccupancyTest::lsCounter [private] |
Definition at line 91 of file DTOccupancyTest.h.
std::set<DTLayerId> DTOccupancyTest::monitoredLayers [private] |
Definition at line 89 of file DTOccupancyTest.h.
edm::ESHandle<DTGeometry> DTOccupancyTest::muonGeom [private] |
Definition at line 82 of file DTOccupancyTest.h.
std::string DTOccupancyTest::nameMonitoredHisto [private] |
Definition at line 102 of file DTOccupancyTest.h.
int DTOccupancyTest::nevents [private] |
Definition at line 78 of file DTOccupancyTest.h.
int DTOccupancyTest::nMinEvts [private] |
Definition at line 92 of file DTOccupancyTest.h.
TNtuple* DTOccupancyTest::ntuple [private] |
Definition at line 96 of file DTOccupancyTest.h.
TFile* DTOccupancyTest::rootFile [private] |
Definition at line 95 of file DTOccupancyTest.h.
bool DTOccupancyTest::runOnAllHitsOccupancies [private] |
Definition at line 99 of file DTOccupancyTest.h.
bool DTOccupancyTest::runOnInTimeOccupancies [private] |
Definition at line 101 of file DTOccupancyTest.h.
bool DTOccupancyTest::runOnNoiseOccupancies [private] |
Definition at line 100 of file DTOccupancyTest.h.
MonitorElement* DTOccupancyTest::summaryHisto [private] |
Definition at line 86 of file DTOccupancyTest.h.
bool DTOccupancyTest::tpMode [private] |
Definition at line 97 of file DTOccupancyTest.h.
std::map< int, MonitorElement* > DTOccupancyTest::wheelHistos [private] |
Definition at line 85 of file DTOccupancyTest.h.
bool DTOccupancyTest::writeRootFile [private] |
Definition at line 94 of file DTOccupancyTest.h.