CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
TrackingMaterialPlotter Class Reference

#include <TrackingMaterialPlotter.h>

Public Types

typedef std::pair< double, double > Range
 

Public Member Functions

void draw (void)
 
void normalize (void)
 
void plotSegmentInLayer (const MaterialAccountingStep &step, int layer)
 
void plotSegmentUnassigned (const MaterialAccountingStep &step)
 
 TrackingMaterialPlotter (float maxZ, float maxR, float resolution)
 

Private Member Functions

void fill_color ()
 
unsigned int fill_gradient (const TColor &first, const TColor &last, unsigned int steps=100, unsigned int index=0)
 
unsigned int fill_gradient (unsigned int first, unsigned int last, unsigned int steps=100, unsigned int index=0)
 

Private Attributes

std::vector< int > m_color
 
std::vector< int > m_gradient
 
XHistogram m_tracker
 

Detailed Description

Definition at line 16 of file TrackingMaterialPlotter.h.

Member Typedef Documentation

typedef std::pair<double, double> TrackingMaterialPlotter::Range

Definition at line 19 of file TrackingMaterialPlotter.h.

Constructor & Destructor Documentation

TrackingMaterialPlotter::TrackingMaterialPlotter ( float  maxZ,
float  maxR,
float  resolution 
)

Definition at line 103 of file TrackingMaterialPlotter.cc.

References fill_color(), fill_gradient(), createfilelist::int, m_color, m_tracker, SiStripPI::max, and pvSelector_cfi::maxZ.

104 {
105  const float rzMinZ = -maxZ;
106  const float rzMaxZ = maxZ;
107  const float rzMinR = 0.;
108  const float rzMaxR = maxR;
109  const int rzBinsZ = (int) (2. * maxZ * resolution);
110  const int rzBinsR = (int) ( maxR * resolution);
111 
112  std::vector<double> max;
113  max.push_back( 0.08 );
114  max.push_back( 0.00016 );
115  m_tracker = XHistogram( 2, rzBinsZ, rzBinsR, std::make_pair(rzMinZ, rzMaxZ), std::make_pair(rzMinR, rzMaxR), m_color.size(), max);
116 
117  TColor::InitializeColors();
118  fill_color();
119  fill_gradient( kWhite, kBlack, 100); // 100-steps gradient from white to black
120 }
unsigned int fill_gradient(const TColor &first, const TColor &last, unsigned int steps=100, unsigned int index=0)

Member Function Documentation

void TrackingMaterialPlotter::draw ( void  )

Definition at line 143 of file TrackingMaterialPlotter.cc.

References svgfig::canvas(), XHistogram::colormap(), XHistogram::get(), cmsRelvalreport::green(), mps_fire::i, m_color, m_gradient, m_tracker, cmsRelvalreport::red(), and Scenarios_cff::scale.

Referenced by normalize().

144 {
145  const double scale = 10.;
146  TCanvas* canvas;
147 
148  XHistogram::Histogram* radlen = m_tracker.get(0);
149  canvas = new TCanvas("radlen_rz", "RadiationLengths - RZ view", (int) (600 * scale * 1.25), (int) (120 * scale * 1.50));
150  gStyle->SetOptStat(0);
151  gStyle->SetPalette( m_gradient.size(), & m_gradient.front() );
152  gStyle->SetNumberContours( m_gradient.size() );
153  canvas->GetFrame()->SetFillColor(kWhite);
154  radlen->Draw("colz");
155  radlen->Draw("same axis y+");
156  radlen->SaveAs("radlen.root");
157  canvas->SaveAs("radlen.png");
158  // Replicate RainBow palette, with White in the first white_slots
159  // positions
160  int white_slots = 1;
161  int MyPalette[100];
162  double stops[9] = { 0.0000, 0.1250, 0.2500, 0.3750, 0.5000, 0.6250, 0.7500, 0.8750, 1.0000};
163  double red[9] = { 0./255., 5./255., 15./255., 35./255., 102./255., 196./255., 208./255., 199./255., 110./255.};
164  double green[9] = { 0./255., 48./255., 124./255., 192./255., 206./255., 226./255., 97./255., 16./255., 0./255.};
165  double blue[9] = { 99./255., 142./255., 198./255., 201./255., 90./255., 22./255., 13./255., 8./255., 2./255.};
166  int palette_index = TColor::CreateGradientColorTable(9, stops, red, green, blue, 100 - white_slots);
167  for (int i=0; i<white_slots; i++) MyPalette[i] = kWhite;
168  for (int i=0; i<100-white_slots; i++) MyPalette[i + white_slots] = palette_index + i;
169  canvas->Clear();
170  gStyle->SetNumberContours( 100 );
171  gStyle->SetPalette(100, MyPalette); // ROOT Rainbow color palette
172  radlen->Draw("colz");
173  radlen->Draw("same axis y+");
174  canvas->SaveAs("radlenColor.png");
175 
176  delete canvas;
177 
179  canvas = new TCanvas("dedx_rz", "-dE/dx term - RZ view", (int) (600 * scale * 1.25), (int) (120 * scale * 1.50));
180  canvas->GetFrame()->SetFillColor(kWhite);
181  gStyle->SetOptStat(0);
182  gStyle->SetPalette( m_gradient.size(), & m_gradient.front() );
183  gStyle->SetNumberContours( m_gradient.size() );
184  dedx->Draw("colz");
185  dedx->Draw("same axis y+");
186  dedx->SaveAs("dedx.root");
187  canvas->SaveAs("dedx.png");
188  canvas->Clear();
189  gStyle->SetNumberContours( 100 );
190  gStyle->SetPalette(100, MyPalette); // ROOT Rainbow color palette
191  dedx->Draw("colz");
192  dedx->Draw("same axis y+");
193  canvas->SaveAs("dedxColor.png");
194  delete canvas;
195 
197  canvas = new TCanvas("layer_rz", "Layers - RZ view", (int) (600 * scale * 1.25), (int) (120 * scale * 1.50));
198  canvas->GetFrame()->SetFillColor(kWhite);
199  gStyle->SetOptStat(0);
200  gStyle->SetPalette( m_color.size(), & m_color.front() );
201  gStyle->SetNumberContours( m_color.size() );
202  colormap->SetMinimum( 1 );
203  colormap->SetMaximum( m_color.size() );
204  colormap->Draw("col");
205  colormap->Draw("same axis y+");
206  colormap->SaveAs("layers.root");
207  canvas->SaveAs("layers.png");
208  delete canvas;
209 }
ColorMap * colormap(void) const
access the colormap
Definition: XHistogram.h:104
TH2F Histogram
Definition: XHistogram.h:17
Histogram * get(size_t h=0) const
access one of the histograms
Definition: XHistogram.h:89
def green(string)
TH2I ColorMap
Definition: XHistogram.h:16
def canvas(sub, attr)
Definition: svgfig.py:481
void TrackingMaterialPlotter::fill_color ( void  )
private

Definition at line 14 of file TrackingMaterialPlotter.cc.

References m_color.

Referenced by TrackingMaterialPlotter().

15 {
16  m_color.push_back(kBlack); // unassigned
17  m_color.push_back(kAzure); // PixelBarrel
18  m_color.push_back(kAzure + 1) ; //
19  m_color.push_back(kAzure + 1) ; //
20  m_color.push_back(kAzure + 3) ; //
21  m_color.push_back(kAzure + 3) ; //
22  m_color.push_back(kGreen); // TIB
23  m_color.push_back(kGreen); //
24  m_color.push_back(kGreen + 2); //
25  m_color.push_back(kGreen + 2); //
26  m_color.push_back(kGreen - 3); //
27  m_color.push_back(kGreen - 3); //
28  m_color.push_back(kGreen - 1); //
29  m_color.push_back(kGreen - 1); //
30  m_color.push_back(kRed); // TOB
31  m_color.push_back(kRed); //
32  m_color.push_back(kRed); //
33  m_color.push_back(kRed + 3); //
34  m_color.push_back(kRed + 3); //
35  m_color.push_back(kRed + 3); //
36  m_color.push_back(kRed - 3); //
37  m_color.push_back(kRed - 3); //
38  m_color.push_back(kRed - 3); //
39  m_color.push_back(kOrange + 9); //
40  m_color.push_back(kOrange + 9); //
41  m_color.push_back(kOrange + 9); //
42  m_color.push_back(kOrange + 7); //
43  m_color.push_back(kOrange + 7); //
44  m_color.push_back(kOrange + 7); //
45  m_color.push_back(kOrange + 5); //
46  m_color.push_back(kOrange + 5); //
47  m_color.push_back(kOrange + 5); //
48  m_color.push_back(kOrange + 8); // PixelEndcap Z-
49  m_color.push_back(kOrange + 10); //
50  m_color.push_back(kOrange - 3); //
51  m_color.push_back(kOrange - 1); // PixelEndcap Z+
52  m_color.push_back(kOrange - 8); //
53  m_color.push_back(kYellow); // TID Z-
54  m_color.push_back(kYellow); //
55  m_color.push_back(kYellow + 2); //
56  m_color.push_back(kYellow + 2); //
57  m_color.push_back(kYellow + 2); //
58  m_color.push_back(kYellow + 3); //
59  m_color.push_back(kMagenta); //
60  m_color.push_back(kMagenta); //
61  m_color.push_back(kMagenta); //
62  m_color.push_back(kMagenta); //
63  m_color.push_back(kMagenta); //
64  m_color.push_back(kMagenta + 1); //
65  m_color.push_back(kMagenta + 2); //
66  m_color.push_back(kMagenta + 3); //
67  m_color.push_back(kMagenta + 4); //
68  m_color.push_back(kMagenta + 5); //
69  m_color.push_back(kMagenta + 6); //
70  m_color.push_back(kMagenta + 7); //
71  m_color.push_back(kMagenta + 8); //
72 }
unsigned int TrackingMaterialPlotter::fill_gradient ( const TColor &  first,
const TColor &  last,
unsigned int  steps = 100,
unsigned int  index = 0 
)
private

Definition at line 75 of file TrackingMaterialPlotter.cc.

References hitfit::delta_r(), diffTwoXMLs::g1, diffTwoXMLs::g2, mps_fire::i, m_gradient, diffTwoXMLs::r1, diffTwoXMLs::r2, and customisers::steps.

Referenced by fill_gradient(), and TrackingMaterialPlotter().

76 {
77  if (index == 0) {
78  // if no index was given, find the highest used one and start from that plus one
79  index = ((TObjArray*) gROOT->GetListOfColors())->GetLast() + 1;
80  }
81 
82  float r1, g1, b1, r2, g2, b2;
83  first.GetRGB(r1, g1, b1);
84  last.GetRGB(r2, g2, b2);
85  float delta_r = (r2 - r1) / (steps - 1);
86  float delta_g = (g2 - g1) / (steps - 1);
87  float delta_b = (b2 - b1) / (steps - 1);
88 
89  m_gradient.resize(steps);
90  for (unsigned int i = 0; i < steps; ++i) {
91  new TColor(static_cast<Int_t>(index + i), r1 + delta_r * i, g1 + delta_g * i, b1 + delta_b * i);
92  m_gradient[i] = index + i;
93  }
94 
95  return index;
96 }
double delta_r(const Fourvec &a, const Fourvec &b)
Find the distance between two four-vectors in the two-dimensional space .
Definition: fourvec.cc:249
unsigned int TrackingMaterialPlotter::fill_gradient ( unsigned int  first,
unsigned int  last,
unsigned int  steps = 100,
unsigned int  index = 0 
)
private

Definition at line 98 of file TrackingMaterialPlotter.cc.

References fill_gradient(), and customisers::steps.

99 {
100  return fill_gradient(* (TColor *) gROOT->GetListOfColors()->At(first), * (TColor *) gROOT->GetListOfColors()->At(last), steps, index);
101 }
unsigned int fill_gradient(const TColor &first, const TColor &last, unsigned int steps=100, unsigned int index=0)
void TrackingMaterialPlotter::normalize ( void  )
inline

Definition at line 25 of file TrackingMaterialPlotter.h.

References draw(), m_tracker, and XHistogram::normalize().

25  {
27  }
void normalize(void)
normalize the histograms
Definition: XHistogram.cc:119
void TrackingMaterialPlotter::plotSegmentInLayer ( const MaterialAccountingStep step,
int  layer 
)

Definition at line 132 of file TrackingMaterialPlotter.cc.

References MaterialAccountingStep::energyLoss(), XHistogram::fill(), MaterialAccountingStep::in(), MaterialAccountingStep::length(), m_tracker, MaterialAccountingStep::out(), PV3DBase< T, PVType, FrameType >::perp(), MaterialAccountingStep::radiationLengths(), w, and PV3DBase< T, PVType, FrameType >::z().

Referenced by TrackingMaterialAnalyser::split().

133 {
134  std::vector<double> w(2);
135  w[0] = step.radiationLengths();
136  w[1] = step.energyLoss();
137  m_tracker.fill( std::make_pair(step.in().z(), step.out().z()),
138  std::make_pair(step.in().perp(), step.out().perp()),
139  w, step.length(), layer+1 ); // layer is 1-based, but plot uses: 0 is empty, 1 is unassigned
140 }
T perp() const
Definition: PV3DBase.h:72
const double w
Definition: UKUtility.cc:23
double length(void) const
double radiationLengths(void) const
double energyLoss(void) const
T z() const
Definition: PV3DBase.h:64
const GlobalPoint & out(void) const
void fill(double x, double y, const std::vector< double > &weight, double norm)
fill one point
Definition: XHistogram.cc:71
const GlobalPoint & in(void) const
void TrackingMaterialPlotter::plotSegmentUnassigned ( const MaterialAccountingStep step)

Definition at line 122 of file TrackingMaterialPlotter.cc.

References MaterialAccountingStep::energyLoss(), XHistogram::fill(), MaterialAccountingStep::in(), MaterialAccountingStep::length(), m_tracker, MaterialAccountingStep::out(), PV3DBase< T, PVType, FrameType >::perp(), MaterialAccountingStep::radiationLengths(), w, and PV3DBase< T, PVType, FrameType >::z().

Referenced by TrackingMaterialAnalyser::split().

123 {
124  std::vector<double> w(2);
125  w[0] = step.radiationLengths();
126  w[1] = step.energyLoss();
127  m_tracker.fill( std::make_pair(step.in().z(), step.out().z()),
128  std::make_pair(step.in().perp(), step.out().perp()),
129  w, step.length(), 1 ); // 0 is empty, 1 is unassigned
130 }
T perp() const
Definition: PV3DBase.h:72
const double w
Definition: UKUtility.cc:23
double length(void) const
double radiationLengths(void) const
double energyLoss(void) const
T z() const
Definition: PV3DBase.h:64
const GlobalPoint & out(void) const
void fill(double x, double y, const std::vector< double > &weight, double norm)
fill one point
Definition: XHistogram.cc:71
const GlobalPoint & in(void) const

Member Data Documentation

std::vector<int> TrackingMaterialPlotter::m_color
private

Definition at line 34 of file TrackingMaterialPlotter.h.

Referenced by draw(), fill_color(), and TrackingMaterialPlotter().

std::vector<int> TrackingMaterialPlotter::m_gradient
private

Definition at line 35 of file TrackingMaterialPlotter.h.

Referenced by draw(), and fill_gradient().

XHistogram TrackingMaterialPlotter::m_tracker
private