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  thePoints.insert(point);
35 }
36 
37 
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 }
57 
58 
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 }
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();
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 }
97 
98 
99 
102  for(set<DTOccupancyPoint>::const_iterator pt = thePoints.begin(); pt != thePoints.end(); ++pt) {
104  }
105 }
106 
107 
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.size() != 0) {
124  computeDistancesToCluster(clusterCandidate);
125  while(clusterCandidate.addPoint(theDistancesFromTheCluster.begin()->second)) {
126  thePoints.erase(theDistancesFromTheCluster.begin()->second);
127  if(thePoints.size() ==0) 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 }
144 
145 
146 
148  LogTrace("DTDQM|DTMonitorClient|DTOccupancyTest|DTOccupancyClusterBuilder") << " sorting" << endl;
149  sort(theClusters.begin(), theClusters.end(), clusterIsLessThan);
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 }
161 
162 
164  return theClusters.front();
165 }
166 
168  if(theProblematicLayers.find(layerId) != theProblematicLayers.end()) {
169  return true;
170  }
171  return false;
172 }
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