CMS 3D CMS Logo

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, DTOccupancyPointgetInitialPair ()
 
void sortClusters ()
 

Private Attributes

double maxMean
 
double maxRMS
 
std::vector< DTOccupancyClustertheClusters
 
std::map< double, std::pair< DTOccupancyPoint, DTOccupancyPoint > > theDistances
 
std::map< double, DTOccupancyPointtheDistancesFromTheCluster
 
std::set< DTOccupancyPointthePoints
 
std::set< DTLayerIdtheProblematicLayers
 

Detailed Description

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

Author
G. Cerminara - INFN Torino

Definition at line 20 of file DTOccupancyClusterBuilder.h.

Constructor & Destructor Documentation

DTOccupancyClusterBuilder::DTOccupancyClusterBuilder ( )

Constructor.

Definition at line 21 of file DTOccupancyClusterBuilder.cc.

DTOccupancyClusterBuilder::~DTOccupancyClusterBuilder ( )
virtual

Destructor.

Definition at line 25 of file DTOccupancyClusterBuilder.cc.

25 {}

Member Function Documentation

void DTOccupancyClusterBuilder::addPoint ( const DTOccupancyPoint point)

Add an occupancy point for a given layer.

Definition at line 29 of file DTOccupancyClusterBuilder.cc.

References EnergyCorrector::pt, theDistances, and thePoints.

29  {
30  // loop over points already stored
31  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
32  theDistances[(*pt).distance(point)] = make_pair(*pt, point);
33  }
34  thePoints.insert(point);
35 }
std::map< double, std::pair< DTOccupancyPoint, DTOccupancyPoint > > theDistances
std::set< DTOccupancyPoint > thePoints
void DTOccupancyClusterBuilder::buildClusters ( )

build the clusters

Definition at line 38 of file DTOccupancyClusterBuilder.cc.

References buildNewCluster(), LogTrace, DTOccupancyCluster::maxMean(), maxMean, DTOccupancyCluster::maxRMS(), maxRMS, EnergyCorrector::pt, sortClusters(), theClusters, and thePoints.

38  {
39  while(buildNewCluster()) {
40  if(thePoints.size() <= 1) break;
41  }
42 
43  // build single point clusters with the remaining points
44  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end();
45  ++pt) {
46  DTOccupancyCluster clusterCandidate(*pt);
47  theClusters.push_back(clusterCandidate);
48  // store the range for building the histograms later
49  if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean();
50  if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS();
51  }
52  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
53  << " # of valid clusters: " << theClusters.size() << endl;
54  sortClusters();
55 
56 }
std::set< DTOccupancyPoint > thePoints
#define LogTrace(id)
std::vector< DTOccupancyCluster > theClusters
bool DTOccupancyClusterBuilder::buildNewCluster ( )
private

Definition at line 108 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().

108  {
109  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
110  << "--------- New Cluster Candidate ----------------------" << endl;
111  pair<DTOccupancyPoint, DTOccupancyPoint> initialPair = getInitialPair();
112  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
113  << " Initial Pair: " << endl
114  << " point1: mean " << initialPair.first.mean()
115  << " rms " << initialPair.first.rms() << endl
116  << " point2: mean " << initialPair.second.mean()
117  << " rms " << initialPair.second.rms() << endl;
118  DTOccupancyCluster clusterCandidate(initialPair.first, initialPair.second);
119  if(clusterCandidate.isValid()) {
120  // remove already used pair
121  thePoints.erase(initialPair.first);
122  thePoints.erase(initialPair.second);
123  if(!thePoints.empty()) {
124  computeDistancesToCluster(clusterCandidate);
125  while(clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) {
126  thePoints.erase(theDistancesFromTheCluster.begin()->second);
127  if(thePoints.empty()) break;
128  computeDistancesToCluster(clusterCandidate);
129  }
130  }
131  } else {
132  return false;
133  }
134  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
135  << " # of layers: " << clusterCandidate.nPoints()
136  << " avrg. mean: " << clusterCandidate.averageMean() << " avrg. rms: " << clusterCandidate.averageRMS() << endl;
137  theClusters.push_back(clusterCandidate);
138  // store the range for building the histograms later
139  if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean();
140  if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS();
142  return true;
143 }
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 100 of file DTOccupancyClusterBuilder.cc.

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

Referenced by buildNewCluster().

100  {
102  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
104  }
105 }
std::map< double, DTOccupancyPoint > theDistancesFromTheCluster
std::set< DTOccupancyPoint > thePoints
double distance(const DTOccupancyPoint &point) const
void DTOccupancyClusterBuilder::computePointToPointDistances ( )
private

Definition at line 85 of file DTOccupancyClusterBuilder.cc.

References theDistances, and thePoints.

Referenced by buildNewCluster().

85  {
86  theDistances.clear();
87  for(set<DTOccupancyPoint>::const_iterator pt_i = thePoints.begin(); pt_i != thePoints.end();
88  ++pt_i) { // i loopo
89  for(set<DTOccupancyPoint>::const_iterator pt_j = thePoints.begin(); pt_j != thePoints.end();
90  ++pt_j) { // j loop
91  if(*pt_i != *pt_j) {
92  theDistances[pt_i->distance(*pt_j)] = make_pair(*pt_i, *pt_j);
93  }
94  }
95  }
96 }
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 59 of file DTOccupancyClusterBuilder.cc.

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

59  {
60  int nBinsX = 100;
61  int nBinsY = 100;
62  int colorMap[12] = {632, 600, 800, 400, 820, 416, 432, 880, 616, 860, 900, 920};
63 
64  TCanvas *canvas = new TCanvas(canvasName.c_str(),canvasName.c_str());
65  canvas->cd();
66  for(vector<DTOccupancyCluster>::const_iterator cluster = theClusters.begin();
67  cluster != theClusters.end(); ++cluster) {
68  stringstream stream;
69  stream << canvasName << "_" << cluster-theClusters.begin();
70  string histoName = stream.str();
71  TH2F *histo = (*cluster).getHisto(histoName, nBinsX, 0, maxMean+3*maxMean/100.,
72  nBinsY, 0, maxRMS+3*maxRMS/100., colorMap[cluster-theClusters.begin()]);
73  if(cluster == theClusters.begin())
74  histo->Draw("box");
75  else
76  histo->Draw("box,same");
77  }
78 }
std::vector< DTOccupancyCluster > theClusters
def canvas(sub, attr)
Definition: svgfig.py:482
DTOccupancyCluster DTOccupancyClusterBuilder::getBestCluster ( ) const

get the cluster correspondig to "normal" cell occupancy.

Definition at line 163 of file DTOccupancyClusterBuilder.cc.

References theClusters.

163  {
164  return theClusters.front();
165 }
std::vector< DTOccupancyCluster > theClusters
std::pair< DTOccupancyPoint, DTOccupancyPoint > DTOccupancyClusterBuilder::getInitialPair ( )
private

Definition at line 81 of file DTOccupancyClusterBuilder.cc.

References theDistances.

Referenced by buildNewCluster().

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

Definition at line 167 of file DTOccupancyClusterBuilder.cc.

References theProblematicLayers.

167  {
168  if(theProblematicLayers.find(layerId) != theProblematicLayers.end()) {
169  return true;
170  }
171  return false;
172 }
std::set< DTLayerId > theProblematicLayers
void DTOccupancyClusterBuilder::sortClusters ( )
private

Definition at line 147 of file DTOccupancyClusterBuilder.cc.

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

Referenced by buildClusters().

147  {
148  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl;
150  // we save the detid of the clusters which are not the best one
151  for(vector<DTOccupancyCluster>::const_iterator cluster = ++(theClusters.begin());
152  cluster != theClusters.end(); ++cluster) { // loop over clusters skipping the first
153  set<DTLayerId> clusterLayers = (*cluster).getLayerIDs();
154  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
155  << " # layers in the cluster: " << clusterLayers.size() << endl;
156  theProblematicLayers.insert(clusterLayers.begin(), clusterLayers.end());
157  }
158  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
159  << " # of problematic layers: " << theProblematicLayers.size() << endl;
160 }
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 62 of file DTOccupancyClusterBuilder.h.

Referenced by buildClusters(), buildNewCluster(), and drawClusters().

double DTOccupancyClusterBuilder::maxRMS
private

Definition at line 63 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 58 of file DTOccupancyClusterBuilder.h.

Referenced by buildNewCluster(), and computeDistancesToCluster().

std::set<DTOccupancyPoint> DTOccupancyClusterBuilder::thePoints
private
std::set<DTLayerId> DTOccupancyClusterBuilder::theProblematicLayers
private

Definition at line 60 of file DTOccupancyClusterBuilder.h.

Referenced by isProblematic(), and sortClusters().