CMS 3D CMS Logo

DTOccupancyClusterBuilder.cc
Go to the documentation of this file.
1 
2 /*
3  * See header file for a description of this class.
4  *
5  * \author G. Cerminara - INFN Torino
6  */
7 
10 
11 #include "TCanvas.h"
12 #include "TH2F.h"
13 
14 #include <algorithm>
15 #include <sstream>
16 #include <iostream>
17 
18 using namespace std;
19 using namespace edm;
20 
22 
24 
26  // loop over points already stored
27  for (set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
28  theDistances[(*pt).distance(point)] = make_pair(*pt, point);
29  }
30  thePoints.insert(point);
31 }
32 
34  while (buildNewCluster()) {
35  if (thePoints.size() <= 1)
36  break;
37  }
38 
39  // build single point clusters with the remaining points
40  for (set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
41  DTOccupancyCluster clusterCandidate(*pt);
42  theClusters.push_back(clusterCandidate);
43  // store the range for building the histograms later
44  if (clusterCandidate.maxMean() > maxMean)
45  maxMean = clusterCandidate.maxMean();
46  if (clusterCandidate.maxRMS() > maxRMS)
47  maxRMS = clusterCandidate.maxRMS();
48  }
49  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
50  << " # of valid clusters: " << theClusters.size() << endl;
51  sortClusters();
52 }
53 
55  int nBinsX = 100;
56  int nBinsY = 100;
57  int colorMap[12] = {632, 600, 800, 400, 820, 416, 432, 880, 616, 860, 900, 920};
58 
59  TCanvas* canvas = new TCanvas(canvasName.c_str(), canvasName.c_str());
60  canvas->cd();
61  for (vector<DTOccupancyCluster>::const_iterator cluster = theClusters.begin(); cluster != theClusters.end();
62  ++cluster) {
63  stringstream stream;
64  stream << canvasName << "_" << cluster - theClusters.begin();
65  string histoName = stream.str();
66  TH2F* histo = (*cluster).getHisto(histoName,
67  nBinsX,
68  0,
69  maxMean + 3 * maxMean / 100.,
70  nBinsY,
71  0,
72  maxRMS + 3 * maxRMS / 100.,
73  colorMap[cluster - theClusters.begin()]);
74  if (cluster == theClusters.begin())
75  histo->Draw("box");
76  else
77  histo->Draw("box,same");
78  }
79 }
80 
81 std::pair<DTOccupancyPoint, DTOccupancyPoint> DTOccupancyClusterBuilder::getInitialPair() {
82  return theDistances.begin()->second;
83 }
84 
86  theDistances.clear();
87  for (set<DTOccupancyPoint>::const_iterator pt_i = thePoints.begin(); pt_i != thePoints.end(); ++pt_i) { // i loopo
88  for (set<DTOccupancyPoint>::const_iterator pt_j = thePoints.begin(); pt_j != thePoints.end(); ++pt_j) { // j loop
89  if (*pt_i != *pt_j) {
90  theDistances[pt_i->distance(*pt_j)] = make_pair(*pt_i, *pt_j);
91  }
92  }
93  }
94 }
95 
98  for (set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
100  }
101 }
102 
104  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
105  << "--------- New Cluster Candidate ----------------------" << endl;
106  pair<DTOccupancyPoint, DTOccupancyPoint> initialPair = getInitialPair();
107  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
108  << " Initial Pair: " << endl
109  << " point1: mean " << initialPair.first.mean() << " rms " << initialPair.first.rms() << endl
110  << " point2: mean " << initialPair.second.mean() << " rms " << initialPair.second.rms() << endl;
111  DTOccupancyCluster clusterCandidate(initialPair.first, initialPair.second);
112  if (clusterCandidate.isValid()) {
113  // remove already used pair
114  thePoints.erase(initialPair.first);
115  thePoints.erase(initialPair.second);
116  if (!thePoints.empty()) {
117  computeDistancesToCluster(clusterCandidate);
118  while (clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) {
119  thePoints.erase(theDistancesFromTheCluster.begin()->second);
120  if (thePoints.empty())
121  break;
122  computeDistancesToCluster(clusterCandidate);
123  }
124  }
125  } else {
126  return false;
127  }
128  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
129  << " # of layers: " << clusterCandidate.nPoints() << " avrg. mean: " << clusterCandidate.averageMean()
130  << " avrg. rms: " << clusterCandidate.averageRMS() << endl;
131  theClusters.push_back(clusterCandidate);
132  // store the range for building the histograms later
133  if (clusterCandidate.maxMean() > maxMean)
134  maxMean = clusterCandidate.maxMean();
135  if (clusterCandidate.maxRMS() > maxRMS)
136  maxRMS = clusterCandidate.maxRMS();
138  return true;
139 }
140 
142  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl;
143  sort(theClusters.begin(), theClusters.end(), clusterIsLessThan);
144  // we save the detid of the clusters which are not the best one
145  for (vector<DTOccupancyCluster>::const_iterator cluster = ++(theClusters.begin()); cluster != theClusters.end();
146  ++cluster) { // loop over clusters skipping the first
147  set<DTLayerId> clusterLayers = (*cluster).getLayerIDs();
148  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
149  << " # layers in the cluster: " << clusterLayers.size() << endl;
150  theProblematicLayers.insert(clusterLayers.begin(), clusterLayers.end());
151  }
152  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
153  << " # of problematic layers: " << theProblematicLayers.size() << endl;
154 }
155 
157 
159  if (theProblematicLayers.find(layerId) != theProblematicLayers.end()) {
160  return true;
161  }
162  return false;
163 }
bool isValid() const
Check if the cluster candidate satisfies the quality requirements.
DTOccupancyCluster getBestCluster() const
get the cluster correspondig to "normal" cell occupancy.
std::set< DTLayerId > theProblematicLayers
std::map< double, std::pair< DTOccupancyPoint, DTOccupancyPoint > > theDistances
std::map< double, DTOccupancyPoint > theDistancesFromTheCluster
double maxMean() const
max average cell occupancy of the layers in the cluster
double averageRMS() const
average RMS of the cell occpuancy distributions of the layers in the cluster
double maxRMS() const
max RMS of the cell occpuancy distributions of the layers in the cluster
std::set< DTOccupancyPoint > thePoints
bool addPoint(const DTOccupancyPoint &anotherPoint)
void drawClusters(std::string canvasName)
draw a TH2F histograms showing the clusters
#define LogTrace(id)
virtual ~DTOccupancyClusterBuilder()
Destructor.
double averageMean() const
average cell occupancy of the layers in the cluster
void computeDistancesToCluster(const DTOccupancyCluster &cluster)
void addPoint(const DTOccupancyPoint &point)
Add an occupancy point for a given layer.
void buildClusters()
build the clusters
std::pair< DTOccupancyPoint, DTOccupancyPoint > getInitialPair()
int nPoints() const
of layers belonging to the cluster
HLT enums.
std::vector< DTOccupancyCluster > theClusters
bool isProblematic(DTLayerId layerId) const
def canvas(sub, attr)
Definition: svgfig.py:482
bool clusterIsLessThan(const DTOccupancyCluster &clusterOne, const DTOccupancyCluster &clusterTwo)
for DTOccupancyCluster sorting
double distance(const DTOccupancyPoint &point) const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5