#include <DTOccupancyClusterBuilder.h>
Public Member Functions | |
void | addPoint (const DTOccupancyPoint &point) |
Add an occupancy point for a given layer. | |
void | buildClusters () |
build the clusters | |
void | drawClusters (std::string canvasName) |
draw a TH2F histograms showing the clusters | |
DTOccupancyClusterBuilder () | |
Constructor. | |
DTOccupancyCluster | getBestCluster () const |
get the cluster correspondig to "normal" cell occupancy. | |
bool | isProblematic (DTLayerId layerId) const |
virtual | ~DTOccupancyClusterBuilder () |
Destructor. | |
Private Member Functions | |
bool | buildNewCluster () |
void | computeDistancesToCluster (const DTOccupancyCluster &cluster) |
void | computePointToPointDistances () |
std::pair< DTOccupancyPoint, DTOccupancyPoint > | getInitialPair () |
void | sortClusters () |
Private Attributes | |
double | maxMean |
double | maxRMS |
std::vector< DTOccupancyCluster > | theClusters |
std::map< double, std::pair < DTOccupancyPoint, DTOccupancyPoint > > | theDistances |
std::map< double, DTOccupancyPoint > | theDistancesFromTheCluster |
std::set< DTOccupancyPoint > | thePoints |
std::set< DTLayerId > | theProblematicLayers |
Build clusters of layer occupancies (DTOccupancyCluster) to spot problematic layers. It's used by DTOccupancyTest.
Definition at line 22 of file DTOccupancyClusterBuilder.h.
DTOccupancyClusterBuilder::DTOccupancyClusterBuilder | ( | ) |
Constructor.
Definition at line 23 of file DTOccupancyClusterBuilder.cc.
DTOccupancyClusterBuilder::~DTOccupancyClusterBuilder | ( | ) | [virtual] |
void DTOccupancyClusterBuilder::addPoint | ( | const DTOccupancyPoint & | point | ) |
Add an occupancy point for a given layer.
Definition at line 31 of file DTOccupancyClusterBuilder.cc.
References theDistances, and thePoints.
Referenced by DTOccupancyTest::runOccupancyTest().
{ // loop over points already stored for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) { theDistances[(*pt).distance(point)] = make_pair(*pt, point); } // cout << "[DTOccupancyClusterBuilder] Add point with mean: " << point.mean() // << " RMS: " << point.rms() << endl; thePoints.insert(point); }
void DTOccupancyClusterBuilder::buildClusters | ( | ) |
build the clusters
Definition at line 42 of file DTOccupancyClusterBuilder.cc.
References buildNewCluster(), LogTrace, DTOccupancyCluster::maxMean(), maxMean, DTOccupancyCluster::maxRMS(), maxRMS, sortClusters(), theClusters, and thePoints.
Referenced by DTOccupancyTest::runOccupancyTest().
{ // cout << "[DTOccupancyClusterBuilder] buildClusters" << endl; while(buildNewCluster()) { // cout << "New cluster builded" << endl; // cout << "# of remaining points: " << thePoints.size() << endl; if(thePoints.size() <= 1) break; } // build single point clusters with the remaining points for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) { DTOccupancyCluster clusterCandidate(*pt); theClusters.push_back(clusterCandidate); // store the range for building the histograms later if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean(); if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS(); } LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " # of valid clusters: " << theClusters.size() << endl; sortClusters(); }
bool DTOccupancyClusterBuilder::buildNewCluster | ( | ) | [private] |
Definition at line 118 of file DTOccupancyClusterBuilder.cc.
References DTOccupancyCluster::addPoint(), DTOccupancyCluster::averageMean(), DTOccupancyCluster::averageRMS(), computeDistancesToCluster(), computePointToPointDistances(), getInitialPair(), DTOccupancyCluster::isValid(), LogTrace, DTOccupancyCluster::maxMean(), maxMean, DTOccupancyCluster::maxRMS(), maxRMS, DTOccupancyCluster::nPoints(), theClusters, theDistancesFromTheCluster, and thePoints.
Referenced by buildClusters().
{ LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << "--------- New Cluster Candidate ----------------------" << endl; pair<DTOccupancyPoint, DTOccupancyPoint> initialPair = getInitialPair(); LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " Initial Pair: " << endl << " point1: mean " << initialPair.first.mean() << " rms " << initialPair.first.rms() << endl << " point2: mean " << initialPair.second.mean() << " rms " << initialPair.second.rms() << endl; DTOccupancyCluster clusterCandidate(initialPair.first, initialPair.second); if(clusterCandidate.isValid()) { // cout << " cluster candidate is valid" << endl; // remove already used pair thePoints.erase(initialPair.first); thePoints.erase(initialPair.second); if(thePoints.size() != 0) { computeDistancesToCluster(clusterCandidate); while(clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) { thePoints.erase(theDistancesFromTheCluster.begin()->second); if(thePoints.size() ==0) break; computeDistancesToCluster(clusterCandidate); } } } else { return false; } LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " # of layers: " << clusterCandidate.nPoints() << " avrg. mean: " << clusterCandidate.averageMean() << " avrg. rms: " << clusterCandidate.averageRMS() << endl; theClusters.push_back(clusterCandidate); // store the range for building the histograms later if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean(); if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS(); computePointToPointDistances(); return true; }
void DTOccupancyClusterBuilder::computeDistancesToCluster | ( | const DTOccupancyCluster & | cluster | ) | [private] |
Definition at line 110 of file DTOccupancyClusterBuilder.cc.
References DTOccupancyCluster::distance(), theDistancesFromTheCluster, and thePoints.
Referenced by buildNewCluster().
{ theDistancesFromTheCluster.clear(); for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) { theDistancesFromTheCluster[cluster.distance(*pt)] = *pt; } }
void DTOccupancyClusterBuilder::computePointToPointDistances | ( | ) | [private] |
Definition at line 95 of file DTOccupancyClusterBuilder.cc.
References theDistances, and thePoints.
Referenced by buildNewCluster().
{ theDistances.clear(); for(set<DTOccupancyPoint>::const_iterator pt_i = thePoints.begin(); pt_i != thePoints.end(); ++pt_i) { // i loopo for(set<DTOccupancyPoint>::const_iterator pt_j = thePoints.begin(); pt_j != thePoints.end(); ++pt_j) { // j loop if(*pt_i != *pt_j) { theDistances[pt_i->distance(*pt_j)] = make_pair(*pt_i, *pt_j); } } } }
void DTOccupancyClusterBuilder::drawClusters | ( | std::string | canvasName | ) |
draw a TH2F histograms showing the clusters
Definition at line 66 of file DTOccupancyClusterBuilder.cc.
References MultipleCompare::canvas, trackerHits::histo, maxMean, maxRMS, and theClusters.
{ int nBinsX = 100; int nBinsY = 100; int colorMap[12] = {632, 600, 800, 400, 820, 416, 432, 880, 616, 860, 900, 920}; // cout << "Draw clusters: " << endl; // cout << " max mean: " << maxMean << " max rms: " << maxRMS << endl; TCanvas *canvas = new TCanvas(canvasName.c_str(),canvasName.c_str()); canvas->cd(); for(vector<DTOccupancyCluster>::const_iterator cluster = theClusters.begin(); cluster != theClusters.end(); ++cluster) { stringstream stream; stream << canvasName << "_" << cluster-theClusters.begin(); string histoName = stream.str(); TH2F *histo = (*cluster).getHisto(histoName, nBinsX, 0, maxMean+3*maxMean/100., nBinsY, 0, maxRMS+3*maxRMS/100., colorMap[cluster-theClusters.begin()]); if(cluster == theClusters.begin()) histo->Draw("box"); else histo->Draw("box,same"); } }
DTOccupancyCluster DTOccupancyClusterBuilder::getBestCluster | ( | ) | const |
get the cluster correspondig to "normal" cell occupancy.
Definition at line 174 of file DTOccupancyClusterBuilder.cc.
References theClusters.
Referenced by DTOccupancyTest::runOccupancyTest().
{ return theClusters.front(); }
std::pair< DTOccupancyPoint, DTOccupancyPoint > DTOccupancyClusterBuilder::getInitialPair | ( | ) | [private] |
Definition at line 91 of file DTOccupancyClusterBuilder.cc.
References theDistances.
Referenced by buildNewCluster().
{ return theDistances.begin()->second; }
bool DTOccupancyClusterBuilder::isProblematic | ( | DTLayerId | layerId | ) | const |
Definition at line 178 of file DTOccupancyClusterBuilder.cc.
References theProblematicLayers.
Referenced by DTOccupancyTest::runOccupancyTest().
{ if(theProblematicLayers.find(layerId) != theProblematicLayers.end()) { return true; } return false; }
void DTOccupancyClusterBuilder::sortClusters | ( | ) | [private] |
Definition at line 158 of file DTOccupancyClusterBuilder.cc.
References clusterIsLessThan(), LogTrace, python::multivaluedict::sort(), theClusters, and theProblematicLayers.
Referenced by buildClusters().
{ LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl; sort(theClusters.begin(), theClusters.end(), clusterIsLessThan); // we save the detid of the clusters which are not the best one for(vector<DTOccupancyCluster>::const_iterator cluster = ++(theClusters.begin()); cluster != theClusters.end(); ++cluster) { // loop over clusters skipping the first set<DTLayerId> clusterLayers = (*cluster).getLayerIDs(); LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " # layers in the cluster: " << clusterLayers.size() << endl; theProblematicLayers.insert(clusterLayers.begin(), clusterLayers.end()); } LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " # of problematic layers: " << theProblematicLayers.size() << endl; }
double DTOccupancyClusterBuilder::maxMean [private] |
Definition at line 64 of file DTOccupancyClusterBuilder.h.
Referenced by buildClusters(), buildNewCluster(), and drawClusters().
double DTOccupancyClusterBuilder::maxRMS [private] |
Definition at line 65 of file DTOccupancyClusterBuilder.h.
Referenced by buildClusters(), buildNewCluster(), and drawClusters().
std::vector<DTOccupancyCluster> DTOccupancyClusterBuilder::theClusters [private] |
Definition at line 61 of file DTOccupancyClusterBuilder.h.
Referenced by buildClusters(), buildNewCluster(), drawClusters(), getBestCluster(), and sortClusters().
std::map<double, std::pair<DTOccupancyPoint, DTOccupancyPoint> > DTOccupancyClusterBuilder::theDistances [private] |
Definition at line 59 of file DTOccupancyClusterBuilder.h.
Referenced by addPoint(), computePointToPointDistances(), and getInitialPair().
std::map<double, DTOccupancyPoint> DTOccupancyClusterBuilder::theDistancesFromTheCluster [private] |
Definition at line 60 of file DTOccupancyClusterBuilder.h.
Referenced by buildNewCluster(), and computeDistancesToCluster().
std::set<DTOccupancyPoint> DTOccupancyClusterBuilder::thePoints [private] |
Definition at line 58 of file DTOccupancyClusterBuilder.h.
Referenced by addPoint(), buildClusters(), buildNewCluster(), computeDistancesToCluster(), and computePointToPointDistances().
std::set<DTLayerId> DTOccupancyClusterBuilder::theProblematicLayers [private] |
Definition at line 62 of file DTOccupancyClusterBuilder.h.
Referenced by isProblematic(), and sortClusters().