#include <DQM/DTMonitorClient/src/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 |
It's used by DTOccupancyTest.
Definition at line 22 of file DTOccupancyClusterBuilder.h.
DTOccupancyClusterBuilder::DTOccupancyClusterBuilder | ( | ) |
DTOccupancyClusterBuilder::~DTOccupancyClusterBuilder | ( | ) | [virtual] |
void DTOccupancyClusterBuilder::addPoint | ( | const DTOccupancyPoint & | point | ) |
Add an occupancy point for a given layer.
Definition at line 30 of file DTOccupancyClusterBuilder.cc.
References theDistances, and thePoints.
Referenced by DTOccupancyTest::runOccupancyTest().
00030 { 00031 // loop over points already stored 00032 for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) { 00033 theDistances[(*pt).distance(point)] = make_pair(*pt, point); 00034 } 00035 // cout << "[DTOccupancyClusterBuilder] Add point with mean: " << point.mean() 00036 // << " RMS: " << point.rms() << endl; 00037 thePoints.insert(point); 00038 }
void DTOccupancyClusterBuilder::buildClusters | ( | ) |
build the clusters
Definition at line 41 of file DTOccupancyClusterBuilder.cc.
References buildNewCluster(), lat::endl(), LogTrace, DTOccupancyCluster::maxMean(), maxMean, DTOccupancyCluster::maxRMS(), maxRMS, sortClusters(), theClusters, and thePoints.
Referenced by DTOccupancyTest::runOccupancyTest().
00041 { 00042 // cout << "[DTOccupancyClusterBuilder] buildClusters" << endl; 00043 while(buildNewCluster()) { 00044 // cout << "New cluster builded" << endl; 00045 // cout << "# of remaining points: " << thePoints.size() << endl; 00046 if(thePoints.size() <= 1) break; 00047 } 00048 00049 // build single point clusters with the remaining points 00050 for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); 00051 ++pt) { 00052 DTOccupancyCluster clusterCandidate(*pt); 00053 theClusters.push_back(clusterCandidate); 00054 // store the range for building the histograms later 00055 if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean(); 00056 if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS(); 00057 } 00058 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") 00059 << " # of valid clusters: " << theClusters.size() << endl; 00060 sortClusters(); 00061 00062 }
bool DTOccupancyClusterBuilder::buildNewCluster | ( | ) | [private] |
Definition at line 117 of file DTOccupancyClusterBuilder.cc.
References DTOccupancyCluster::addPoint(), DTOccupancyCluster::averageMean(), DTOccupancyCluster::averageRMS(), computeDistancesToCluster(), computePointToPointDistances(), lat::endl(), getInitialPair(), DTOccupancyCluster::isValid(), LogTrace, DTOccupancyCluster::maxMean(), maxMean, DTOccupancyCluster::maxRMS(), maxRMS, DTOccupancyCluster::nPoints(), theClusters, theDistancesFromTheCluster, and thePoints.
Referenced by buildClusters().
00117 { 00118 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") 00119 << "--------- New Cluster Candidate ----------------------" << endl; 00120 pair<DTOccupancyPoint, DTOccupancyPoint> initialPair = getInitialPair(); 00121 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") 00122 << " Initial Pair: " << endl 00123 << " point1: mean " << initialPair.first.mean() 00124 << " rms " << initialPair.first.rms() << endl 00125 << " point2: mean " << initialPair.second.mean() 00126 << " rms " << initialPair.second.rms() << endl; 00127 DTOccupancyCluster clusterCandidate(initialPair.first, initialPair.second); 00128 if(clusterCandidate.isValid()) { 00129 // cout << " cluster candidate is valid" << endl; 00130 // remove already used pair 00131 thePoints.erase(initialPair.first); 00132 thePoints.erase(initialPair.second); 00133 if(thePoints.size() != 0) { 00134 computeDistancesToCluster(clusterCandidate); 00135 while(clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) { 00136 thePoints.erase(theDistancesFromTheCluster.begin()->second); 00137 if(thePoints.size() ==0) break; 00138 computeDistancesToCluster(clusterCandidate); 00139 } 00140 } 00141 } else { 00142 return false; 00143 } 00144 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") 00145 << " # of layers: " << clusterCandidate.nPoints() 00146 << " avrg. mean: " << clusterCandidate.averageMean() << " avrg. rms: " << clusterCandidate.averageRMS() << endl; 00147 theClusters.push_back(clusterCandidate); 00148 // store the range for building the histograms later 00149 if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean(); 00150 if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS(); 00151 computePointToPointDistances(); 00152 return true; 00153 }
void DTOccupancyClusterBuilder::computeDistancesToCluster | ( | const DTOccupancyCluster & | cluster | ) | [private] |
Definition at line 109 of file DTOccupancyClusterBuilder.cc.
References DTOccupancyCluster::distance(), theDistancesFromTheCluster, and thePoints.
Referenced by buildNewCluster().
00109 { 00110 theDistancesFromTheCluster.clear(); 00111 for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) { 00112 theDistancesFromTheCluster[cluster.distance(*pt)] = *pt; 00113 } 00114 }
void DTOccupancyClusterBuilder::computePointToPointDistances | ( | ) | [private] |
Definition at line 94 of file DTOccupancyClusterBuilder.cc.
References theDistances, and thePoints.
Referenced by buildNewCluster().
00094 { 00095 theDistances.clear(); 00096 for(set<DTOccupancyPoint>::const_iterator pt_i = thePoints.begin(); pt_i != thePoints.end(); 00097 ++pt_i) { // i loopo 00098 for(set<DTOccupancyPoint>::const_iterator pt_j = thePoints.begin(); pt_j != thePoints.end(); 00099 ++pt_j) { // j loop 00100 if(*pt_i != *pt_j) { 00101 theDistances[pt_i->distance(*pt_j)] = make_pair(*pt_i, *pt_j); 00102 } 00103 } 00104 } 00105 }
void DTOccupancyClusterBuilder::drawClusters | ( | std::string | canvasName | ) |
draw a TH2F histograms showing the clusters
Definition at line 65 of file DTOccupancyClusterBuilder.cc.
References histo, maxMean, maxRMS, and theClusters.
00065 { 00066 int nBinsX = 100; 00067 int nBinsY = 100; 00068 int colorMap[12] = {632, 600, 800, 400, 820, 416, 432, 880, 616, 860, 900, 920}; 00069 00070 // cout << "Draw clusters: " << endl; 00071 // cout << " max mean: " << maxMean << " max rms: " << maxRMS << endl; 00072 00073 TCanvas *canvas = new TCanvas(canvasName.c_str(),canvasName.c_str()); 00074 canvas->cd(); 00075 for(vector<DTOccupancyCluster>::const_iterator cluster = theClusters.begin(); 00076 cluster != theClusters.end(); ++cluster) { 00077 stringstream stream; 00078 stream << canvasName << "_" << cluster-theClusters.begin(); 00079 string histoName = stream.str(); 00080 TH2F *histo = (*cluster).getHisto(histoName, nBinsX, 0, maxMean+3*maxMean/100., 00081 nBinsY, 0, maxRMS+3*maxRMS/100., colorMap[cluster-theClusters.begin()]); 00082 if(cluster == theClusters.begin()) 00083 histo->Draw("box"); 00084 else 00085 histo->Draw("box,same"); 00086 } 00087 }
DTOccupancyCluster DTOccupancyClusterBuilder::getBestCluster | ( | ) | const |
get the cluster correspondig to "normal" cell occupancy.
Definition at line 173 of file DTOccupancyClusterBuilder.cc.
References theClusters.
Referenced by DTOccupancyTest::runOccupancyTest().
00173 { 00174 return theClusters.front(); 00175 }
std::pair< DTOccupancyPoint, DTOccupancyPoint > DTOccupancyClusterBuilder::getInitialPair | ( | ) | [private] |
Definition at line 90 of file DTOccupancyClusterBuilder.cc.
References theDistances.
Referenced by buildNewCluster().
00090 { 00091 return theDistances.begin()->second; 00092 }
Definition at line 177 of file DTOccupancyClusterBuilder.cc.
References theProblematicLayers.
Referenced by DTOccupancyTest::runOccupancyTest().
00177 { 00178 if(theProblematicLayers.find(layerId) != theProblematicLayers.end()) { 00179 return true; 00180 } 00181 return false; 00182 }
void DTOccupancyClusterBuilder::sortClusters | ( | ) | [private] |
Definition at line 157 of file DTOccupancyClusterBuilder.cc.
References clusterIsLessThan(), lat::endl(), LogTrace, python::multivaluedict::sort(), theClusters, and theProblematicLayers.
Referenced by buildClusters().
00157 { 00158 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl; 00159 sort(theClusters.begin(), theClusters.end(), clusterIsLessThan); 00160 // we save the detid of the clusters which are not the best one 00161 for(vector<DTOccupancyCluster>::const_iterator cluster = ++(theClusters.begin()); 00162 cluster != theClusters.end(); ++cluster) { // loop over clusters skipping the first 00163 set<DTLayerId> clusterLayers = (*cluster).getLayerIDs(); 00164 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") 00165 << " # layers in the cluster: " << clusterLayers.size() << endl; 00166 theProblematicLayers.insert(clusterLayers.begin(), clusterLayers.end()); 00167 } 00168 LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") 00169 << " # of problematic layers: " << theProblematicLayers.size() << endl; 00170 }
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().