CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTOccupancyClusterBuilder Class Reference

#include <DTOccupancyClusterBuilder.h>

Public Member Functions

void addPoint (const DTOccupancyPoint &point)
 Add an occupancy point for a given layer. More...
 
void buildClusters ()
 build the clusters More...
 
void drawClusters (std::string canvasName)
 draw a TH2F histograms showing the clusters More...
 
 DTOccupancyClusterBuilder ()
 Constructor. More...
 
DTOccupancyCluster getBestCluster () const
 get the cluster correspondig to "normal" cell occupancy. More...
 
bool isProblematic (DTLayerId layerId) const
 
virtual ~DTOccupancyClusterBuilder ()
 Destructor. More...
 

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< DTOccupancyClustertheClusters
 
std::map< double, std::pair
< DTOccupancyPoint,
DTOccupancyPoint > > 
theDistances
 
std::map< double,
DTOccupancyPoint
theDistancesFromTheCluster
 
std::set< DTOccupancyPointthePoints
 
std::set< DTLayerIdtheProblematicLayers
 

Detailed Description

Build clusters of layer occupancies (DTOccupancyCluster) to spot problematic layers. It's used by DTOccupancyTest.

Date:
2008/10/16 09:33:39
Revision:
1.3
Author
G. Cerminara - INFN Torino

Definition at line 22 of file DTOccupancyClusterBuilder.h.

Constructor & Destructor Documentation

DTOccupancyClusterBuilder::DTOccupancyClusterBuilder ( )

Constructor.

Definition at line 23 of file DTOccupancyClusterBuilder.cc.

DTOccupancyClusterBuilder::~DTOccupancyClusterBuilder ( )
virtual

Destructor.

Definition at line 27 of file DTOccupancyClusterBuilder.cc.

27 {}

Member Function Documentation

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().

31  {
32  // loop over points already stored
33  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
34  theDistances[(*pt).distance(point)] = make_pair(*pt, point);
35  }
36  // cout << "[DTOccupancyClusterBuilder] Add point with mean: " << point.mean()
37  // << " RMS: " << point.rms() << endl;
38  thePoints.insert(point);
39 }
std::map< double, std::pair< DTOccupancyPoint, DTOccupancyPoint > > theDistances
std::set< DTOccupancyPoint > thePoints
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().

42  {
43  // cout << "[DTOccupancyClusterBuilder] buildClusters" << endl;
44  while(buildNewCluster()) {
45  // cout << "New cluster builded" << endl;
46  // cout << "# of remaining points: " << thePoints.size() << endl;
47  if(thePoints.size() <= 1) break;
48  }
49 
50  // build single point clusters with the remaining points
51  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end();
52  ++pt) {
53  DTOccupancyCluster clusterCandidate(*pt);
54  theClusters.push_back(clusterCandidate);
55  // store the range for building the histograms later
56  if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean();
57  if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS();
58  }
59  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
60  << " # of valid clusters: " << theClusters.size() << endl;
61  sortClusters();
62 
63 }
std::set< DTOccupancyPoint > thePoints
#define LogTrace(id)
std::vector< DTOccupancyCluster > theClusters
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().

118  {
119  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
120  << "--------- New Cluster Candidate ----------------------" << endl;
121  pair<DTOccupancyPoint, DTOccupancyPoint> initialPair = getInitialPair();
122  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
123  << " Initial Pair: " << endl
124  << " point1: mean " << initialPair.first.mean()
125  << " rms " << initialPair.first.rms() << endl
126  << " point2: mean " << initialPair.second.mean()
127  << " rms " << initialPair.second.rms() << endl;
128  DTOccupancyCluster clusterCandidate(initialPair.first, initialPair.second);
129  if(clusterCandidate.isValid()) {
130  // cout << " cluster candidate is valid" << endl;
131  // remove already used pair
132  thePoints.erase(initialPair.first);
133  thePoints.erase(initialPair.second);
134  if(thePoints.size() != 0) {
135  computeDistancesToCluster(clusterCandidate);
136  while(clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) {
137  thePoints.erase(theDistancesFromTheCluster.begin()->second);
138  if(thePoints.size() ==0) break;
139  computeDistancesToCluster(clusterCandidate);
140  }
141  }
142  } else {
143  return false;
144  }
145  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
146  << " # of layers: " << clusterCandidate.nPoints()
147  << " avrg. mean: " << clusterCandidate.averageMean() << " avrg. rms: " << clusterCandidate.averageRMS() << endl;
148  theClusters.push_back(clusterCandidate);
149  // store the range for building the histograms later
150  if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean();
151  if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS();
153  return true;
154 }
std::map< double, DTOccupancyPoint > theDistancesFromTheCluster
std::set< DTOccupancyPoint > thePoints
#define LogTrace(id)
void computeDistancesToCluster(const DTOccupancyCluster &cluster)
std::pair< DTOccupancyPoint, DTOccupancyPoint > getInitialPair()
std::vector< DTOccupancyCluster > theClusters
void DTOccupancyClusterBuilder::computeDistancesToCluster ( const DTOccupancyCluster cluster)
private

Definition at line 110 of file DTOccupancyClusterBuilder.cc.

References DTOccupancyCluster::distance(), theDistancesFromTheCluster, and thePoints.

Referenced by buildNewCluster().

110  {
112  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
113  theDistancesFromTheCluster[cluster.distance(*pt)] = *pt;
114  }
115 }
std::map< double, DTOccupancyPoint > theDistancesFromTheCluster
std::set< DTOccupancyPoint > thePoints
double distance(const DTOccupancyPoint &point) const
void DTOccupancyClusterBuilder::computePointToPointDistances ( )
private

Definition at line 95 of file DTOccupancyClusterBuilder.cc.

References theDistances, and thePoints.

Referenced by buildNewCluster().

95  {
96  theDistances.clear();
97  for(set<DTOccupancyPoint>::const_iterator pt_i = thePoints.begin(); pt_i != thePoints.end();
98  ++pt_i) { // i loopo
99  for(set<DTOccupancyPoint>::const_iterator pt_j = thePoints.begin(); pt_j != thePoints.end();
100  ++pt_j) { // j loop
101  if(*pt_i != *pt_j) {
102  theDistances[pt_i->distance(*pt_j)] = make_pair(*pt_i, *pt_j);
103  }
104  }
105  }
106 }
std::map< double, std::pair< DTOccupancyPoint, DTOccupancyPoint > > theDistances
std::set< DTOccupancyPoint > thePoints
void DTOccupancyClusterBuilder::drawClusters ( std::string  canvasName)

draw a TH2F histograms showing the clusters

Definition at line 66 of file DTOccupancyClusterBuilder.cc.

References svgfig::canvas(), timingPdfMaker::histo, maxMean, maxRMS, and theClusters.

66  {
67  int nBinsX = 100;
68  int nBinsY = 100;
69  int colorMap[12] = {632, 600, 800, 400, 820, 416, 432, 880, 616, 860, 900, 920};
70 
71  // cout << "Draw clusters: " << endl;
72  // cout << " max mean: " << maxMean << " max rms: " << maxRMS << endl;
73 
74  TCanvas *canvas = new TCanvas(canvasName.c_str(),canvasName.c_str());
75  canvas->cd();
76  for(vector<DTOccupancyCluster>::const_iterator cluster = theClusters.begin();
77  cluster != theClusters.end(); ++cluster) {
78  stringstream stream;
79  stream << canvasName << "_" << cluster-theClusters.begin();
80  string histoName = stream.str();
81  TH2F *histo = (*cluster).getHisto(histoName, nBinsX, 0, maxMean+3*maxMean/100.,
82  nBinsY, 0, maxRMS+3*maxRMS/100., colorMap[cluster-theClusters.begin()]);
83  if(cluster == theClusters.begin())
84  histo->Draw("box");
85  else
86  histo->Draw("box,same");
87  }
88 }
def canvas
Definition: svgfig.py:481
std::vector< DTOccupancyCluster > theClusters
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().

174  {
175  return theClusters.front();
176 }
std::vector< DTOccupancyCluster > theClusters
std::pair< DTOccupancyPoint, DTOccupancyPoint > DTOccupancyClusterBuilder::getInitialPair ( )
private

Definition at line 91 of file DTOccupancyClusterBuilder.cc.

References theDistances.

Referenced by buildNewCluster().

91  {
92  return theDistances.begin()->second;
93 }
std::map< double, std::pair< DTOccupancyPoint, DTOccupancyPoint > > theDistances
bool DTOccupancyClusterBuilder::isProblematic ( DTLayerId  layerId) const

Definition at line 178 of file DTOccupancyClusterBuilder.cc.

References theProblematicLayers.

Referenced by DTOccupancyTest::runOccupancyTest().

178  {
179  if(theProblematicLayers.find(layerId) != theProblematicLayers.end()) {
180  return true;
181  }
182  return false;
183 }
std::set< DTLayerId > theProblematicLayers
void DTOccupancyClusterBuilder::sortClusters ( )
private

Definition at line 158 of file DTOccupancyClusterBuilder.cc.

References clusterIsLessThan(), LogTrace, python.multivaluedict::sort(), theClusters, and theProblematicLayers.

Referenced by buildClusters().

158  {
159  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl;
161  // we save the detid of the clusters which are not the best one
162  for(vector<DTOccupancyCluster>::const_iterator cluster = ++(theClusters.begin());
163  cluster != theClusters.end(); ++cluster) { // loop over clusters skipping the first
164  set<DTLayerId> clusterLayers = (*cluster).getLayerIDs();
165  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
166  << " # layers in the cluster: " << clusterLayers.size() << endl;
167  theProblematicLayers.insert(clusterLayers.begin(), clusterLayers.end());
168  }
169  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
170  << " # of problematic layers: " << theProblematicLayers.size() << endl;
171 }
std::set< DTLayerId > theProblematicLayers
#define LogTrace(id)
std::vector< DTOccupancyCluster > theClusters
bool clusterIsLessThan(const DTOccupancyCluster &clusterOne, const DTOccupancyCluster &clusterTwo)
for DTOccupancyCluster sorting

Member Data Documentation

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
std::map<double, std::pair<DTOccupancyPoint, DTOccupancyPoint> > DTOccupancyClusterBuilder::theDistances
private
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
std::set<DTLayerId> DTOccupancyClusterBuilder::theProblematicLayers
private

Definition at line 62 of file DTOccupancyClusterBuilder.h.

Referenced by isProblematic(), and sortClusters().