CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  maxRMS(-1.) {
23 }
24 
26 
27 
28 
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  // cout << "[DTOccupancyClusterBuilder] Add point with mean: " << point.mean()
35  // << " RMS: " << point.rms() << endl;
36  thePoints.insert(point);
37 }
38 
39 
41  // cout << "[DTOccupancyClusterBuilder] buildClusters" << endl;
42  while(buildNewCluster()) {
43  // cout << "New cluster builded" << endl;
44  // cout << "# of remaining points: " << thePoints.size() << endl;
45  if(thePoints.size() <= 1) break;
46  }
47 
48  // build single point clusters with the remaining points
49  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end();
50  ++pt) {
51  DTOccupancyCluster clusterCandidate(*pt);
52  theClusters.push_back(clusterCandidate);
53  // store the range for building the histograms later
54  if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean();
55  if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS();
56  }
57  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
58  << " # of valid clusters: " << theClusters.size() << endl;
59  sortClusters();
60 
61 }
62 
63 
65  int nBinsX = 100;
66  int nBinsY = 100;
67  int colorMap[12] = {632, 600, 800, 400, 820, 416, 432, 880, 616, 860, 900, 920};
68 
69  // cout << "Draw clusters: " << endl;
70  // cout << " max mean: " << maxMean << " max rms: " << maxRMS << endl;
71 
72  TCanvas *canvas = new TCanvas(canvasName.c_str(),canvasName.c_str());
73  canvas->cd();
74  for(vector<DTOccupancyCluster>::const_iterator cluster = theClusters.begin();
75  cluster != theClusters.end(); ++cluster) {
76  stringstream stream;
77  stream << canvasName << "_" << cluster-theClusters.begin();
78  string histoName = stream.str();
79  TH2F *histo = (*cluster).getHisto(histoName, nBinsX, 0, maxMean+3*maxMean/100.,
80  nBinsY, 0, maxRMS+3*maxRMS/100., colorMap[cluster-theClusters.begin()]);
81  if(cluster == theClusters.begin())
82  histo->Draw("box");
83  else
84  histo->Draw("box,same");
85  }
86 }
87 
88 
89 std::pair<DTOccupancyPoint, DTOccupancyPoint> DTOccupancyClusterBuilder::getInitialPair() {
90  return theDistances.begin()->second;
91 }
92 
94  theDistances.clear();
95  for(set<DTOccupancyPoint>::const_iterator pt_i = thePoints.begin(); pt_i != thePoints.end();
96  ++pt_i) { // i loopo
97  for(set<DTOccupancyPoint>::const_iterator pt_j = thePoints.begin(); pt_j != thePoints.end();
98  ++pt_j) { // j loop
99  if(*pt_i != *pt_j) {
100  theDistances[pt_i->distance(*pt_j)] = make_pair(*pt_i, *pt_j);
101  }
102  }
103  }
104 }
105 
106 
107 
110  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
112  }
113 }
114 
115 
117  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
118  << "--------- New Cluster Candidate ----------------------" << endl;
119  pair<DTOccupancyPoint, DTOccupancyPoint> initialPair = getInitialPair();
120  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
121  << " Initial Pair: " << endl
122  << " point1: mean " << initialPair.first.mean()
123  << " rms " << initialPair.first.rms() << endl
124  << " point2: mean " << initialPair.second.mean()
125  << " rms " << initialPair.second.rms() << endl;
126  DTOccupancyCluster clusterCandidate(initialPair.first, initialPair.second);
127  if(clusterCandidate.isValid()) {
128  // cout << " cluster candidate is valid" << endl;
129  // remove already used pair
130  thePoints.erase(initialPair.first);
131  thePoints.erase(initialPair.second);
132  if(thePoints.size() != 0) {
133  computeDistancesToCluster(clusterCandidate);
134  while(clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) {
135  thePoints.erase(theDistancesFromTheCluster.begin()->second);
136  if(thePoints.size() ==0) break;
137  computeDistancesToCluster(clusterCandidate);
138  }
139  }
140  } else {
141  return false;
142  }
143  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
144  << " # of layers: " << clusterCandidate.nPoints()
145  << " avrg. mean: " << clusterCandidate.averageMean() << " avrg. rms: " << clusterCandidate.averageRMS() << endl;
146  theClusters.push_back(clusterCandidate);
147  // store the range for building the histograms later
148  if(clusterCandidate.maxMean() > maxMean) maxMean = clusterCandidate.maxMean();
149  if(clusterCandidate.maxRMS() > maxRMS) maxRMS = clusterCandidate.maxRMS();
151  return true;
152 }
153 
154 
155 
157  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl;
159  // we save the detid of the clusters which are not the best one
160  for(vector<DTOccupancyCluster>::const_iterator cluster = ++(theClusters.begin());
161  cluster != theClusters.end(); ++cluster) { // loop over clusters skipping the first
162  set<DTLayerId> clusterLayers = (*cluster).getLayerIDs();
163  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
164  << " # layers in the cluster: " << clusterLayers.size() << endl;
165  theProblematicLayers.insert(clusterLayers.begin(), clusterLayers.end());
166  }
167  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder")
168  << " # of problematic layers: " << theProblematicLayers.size() << endl;
169 }
170 
171 
173  return theClusters.front();
174 }
175 
177  if(theProblematicLayers.find(layerId) != theProblematicLayers.end()) {
178  return true;
179  }
180  return false;
181 }
bool isValid() const
Check if the cluster candidate satisfies the quality requirements.
DTOccupancyCluster getBestCluster() const
get the cluster correspondig to &quot;normal&quot; 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
def canvas
Definition: svgfig.py:481
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
std::vector< DTOccupancyCluster > theClusters
bool isProblematic(DTLayerId layerId) const
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