CMS 3D CMS Logo

MaterialAccountingGroup.cc
Go to the documentation of this file.
1 #include <sstream>
2 #include <iomanip>
3 #include <string>
4 #include <stdexcept>
5 
6 #include <TFile.h>
7 #include <TH1F.h>
8 #include <TProfile.h>
9 #include <TCanvas.h>
10 #include <TFrame.h>
11 
19 
21  0.01; // 100um should be small enough that no elements from different layers/groups are so close
22 
24  : m_name(name),
25  m_elements(),
26  m_boundingbox(),
27  m_accounting(),
28  m_errors(),
29  m_tracks(0),
30  m_counted(false),
31  m_file(nullptr) {
32  // retrieve the elements from DDD
33  DDSpecificsMatchesValueFilter filter{DDValue("TrackingMaterialGroup", name)};
34  DDFilteredView fv(geometry, filter);
35  LogTrace("MaterialAccountingGroup") << "Elements within: " << name << std::endl;
36  while (fv.next()) {
37  // DD3Vector and DDTranslation are the same type as math::XYZVector
38  math::XYZVector position = fv.translation() / 10.; // mm -> cm
39  LogTrace("MaterialAccountingGroup") << "Adding element at(r,z): ("
40  << GlobalPoint(position.x(), position.y(), position.z()).perp() << ", "
41  << GlobalPoint(position.x(), position.y(), position.z()).z() << ") cm"
42  << std::endl;
43  LogTrace("MaterialAccountingGroup") << "Name of added element: " << fv.logicalPart().toString() << std::endl;
44  m_elements.push_back(GlobalPoint(position.x(), position.y(), position.z()));
45  }
46 
47  // grow the bounding box
48  for (unsigned int i = 0; i < m_elements.size(); ++i) {
50  }
52  LogTrace("MaterialAccountingGroup") << "Final BBox r_range: " << m_boundingbox.range_r().first << ", "
53  << m_boundingbox.range_r().second << std::endl
54  << "Final BBox z_range: " << m_boundingbox.range_z().first << ", "
55  << m_boundingbox.range_z().second << std::endl;
56 
57  // initialize the histograms
58  m_dedx_spectrum = new TH1F((m_name + "_dedx_spectrum").c_str(), "Energy loss spectrum", 1000, 0, 1);
59  m_radlen_spectrum = new TH1F((m_name + "_radlen_spectrum").c_str(), "Radiation lengths spectrum", 1000, 0, 1);
60  m_dedx_vs_eta = new TProfile((m_name + "_dedx_vs_eta").c_str(), "Energy loss vs. eta", 600, -3, 3);
61  m_dedx_vs_z = new TProfile((m_name + "_dedx_vs_z").c_str(), "Energy loss vs. Z", 6000, -300, 300);
62  m_dedx_vs_r = new TProfile((m_name + "_dedx_vs_r").c_str(), "Energy loss vs. R", 1200, 0, 120);
63  m_radlen_vs_eta = new TProfile((m_name + "_radlen_vs_eta").c_str(), "Radiation lengths vs. eta", 600, -3, 3);
64  m_radlen_vs_z = new TProfile((m_name + "_radlen_vs_z").c_str(), "Radiation lengths vs. Z", 6000, -300, 300);
65  m_radlen_vs_r = new TProfile((m_name + "_radlen_vs_r").c_str(), "Radiation lengths vs. R", 1200, 0, 120);
66  m_dedx_spectrum->SetDirectory(nullptr);
67  m_radlen_spectrum->SetDirectory(nullptr);
68  m_dedx_vs_eta->SetDirectory(nullptr);
69  m_dedx_vs_z->SetDirectory(nullptr);
70  m_dedx_vs_r->SetDirectory(nullptr);
71  m_radlen_vs_eta->SetDirectory(nullptr);
72  m_radlen_vs_z->SetDirectory(nullptr);
73  m_radlen_vs_r->SetDirectory(nullptr);
74 }
75 
77  delete m_dedx_spectrum;
78  delete m_dedx_vs_eta;
79  delete m_dedx_vs_z;
80  delete m_dedx_vs_r;
81  delete m_radlen_spectrum;
82  delete m_radlen_vs_eta;
83  delete m_radlen_vs_z;
84  delete m_radlen_vs_r;
85 }
86 
87 // TODO the inner check could be sped up in many ways
88 // (sorting the m_elements, partitioning the bounding box, ...)
89 // but is it worth?
90 // especially with the segmentation of the layers ?
92  const GlobalPoint& position = detector.position();
93  // first check to see if the point is inside the bounding box
94  LogTrace("MaterialAccountingGroup") << "Testing position: (x, y, z, r) = " << position.x() << ", " << position.y()
95  << ", " << position.z() << ", " << position.perp() << std::endl;
96  if (not m_boundingbox.inside(position.perp(), position.z())) {
97  LogTrace("MaterialAccountingGroup") << "r outside of: (" << m_boundingbox.range_r().first << ", "
98  << m_boundingbox.range_r().second << "), Z ouside of: ("
99  << m_boundingbox.range_z().first << ", " << m_boundingbox.range_z().second
100  << ")" << std::endl;
101  return false;
102  } else {
103  // now check if the point is actually close enough to any element
104  LogTrace("MaterialAccountingGroup") << "r within: (" << m_boundingbox.range_r().first << ", "
105  << m_boundingbox.range_r().second << "), Z within: ("
106  << m_boundingbox.range_z().first << ", " << m_boundingbox.range_z().second
107  << ")" << std::endl;
108  for (unsigned int i = 0; i < m_elements.size(); ++i) {
109  LogTrace("MaterialAccountingGroup")
110  << "Closest testing agains(x, y, z, r): (" << m_elements[i].x() << ", " << m_elements[i].y() << ", "
111  << m_elements[i].z() << ", " << m_elements[i].perp() << ") --> " << (position - m_elements[i]).mag()
112  << " vs tolerance: " << s_tolerance << std::endl;
113  if ((position - m_elements[i]).mag2() < (s_tolerance * s_tolerance))
114  return true;
115  }
116  return false;
117  }
118 }
119 
121  if (not inside(detector))
122  return false;
123 
124  // multiple hits in the same layer (overlaps, etc.) from a single track still count as one for averaging,
125  // since the energy deposits from the track have been already split between the different detectors
126  m_buffer += detector.material();
127  m_counted = true;
128 
129  return true;
130 }
131 
133  // add a detector
134  if (m_counted) {
137  ++m_tracks;
138 
139  GlobalPoint average((m_buffer.in().x() + m_buffer.out().x()) / 2.,
140  (m_buffer.in().y() + m_buffer.out().y()) / 2.,
141  (m_buffer.in().z() + m_buffer.out().z()) / 2.);
142  m_dedx_spectrum->Fill(m_buffer.energyLoss());
143  m_radlen_spectrum->Fill(m_buffer.radiationLengths());
144  m_dedx_vs_eta->Fill(average.eta(), m_buffer.energyLoss(), 1.);
145  m_dedx_vs_z->Fill(average.z(), m_buffer.energyLoss(), 1.);
146  m_dedx_vs_r->Fill(average.perp(), m_buffer.energyLoss(), 1.);
147  m_radlen_vs_eta->Fill(average.eta(), m_buffer.radiationLengths(), 1.);
148  m_radlen_vs_z->Fill(average.z(), m_buffer.radiationLengths(), 1.);
149  m_radlen_vs_r->Fill(average.perp(), m_buffer.radiationLengths(), 1.);
150  }
151  m_counted = false;
153 }
154 
156  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
157  plot->SetFillColor(15); // grey
158  plot->SetLineColor(1); // black
159  plot->Draw("c e");
160  canvas.GetFrame()->SetFillColor(kWhite);
161  canvas.Draw();
162  canvas.SaveAs((name + ".png").c_str(), "");
163 
164  // store te plot into m_file
165  plot->SetDirectory(m_file);
166 }
167 
169  // Nota Bene:
170  // these "line" plots are not deleted explicitly since
171  // - deleting them before saving them to a TFile will not save them
172  // - deleting them after the TFile they're stored into results in a SEGV
173  // ROOT is probably "taking care" (read: messing things up) somehow...
174  TH1F* line =
175  new TH1F((name + "_par").c_str(), "Parametrization", 1, plot->GetXaxis()->GetXmin(), plot->GetXaxis()->GetXmax());
176  line->SetBinContent(1, average);
177 
178  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
179  plot->SetFillColor(15); // grey
180  plot->SetLineColor(1); // black
181  plot->SetLineWidth(2);
182  plot->Draw("c e6");
183  line->SetLineColor(2); // red
184  line->SetLineWidth(2);
185  line->Draw("same");
186  canvas.GetFrame()->SetFillColor(kWhite);
187  canvas.Draw();
188  canvas.SaveAs((name + ".png").c_str(), "");
189 
190  // store te plots into m_file
191  plot->SetDirectory(m_file);
192  line->SetDirectory(m_file);
193 }
194 
197  std::stringstream out;
198  out << std::setw(48) << std::left << m_name << std::right << std::fixed;
199  ;
200  out << "BBox: " << std::setprecision(1) << std::setw(6) << m_boundingbox.range_z().first << " < Z < "
201  << std::setprecision(1) << std::setw(6) << m_boundingbox.range_z().second;
202  out << ", " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().first << " < R < "
203  << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().second;
204  out << " Elements: " << std::setw(6) << m_elements.size();
205  return out.str();
206 }
207 
209  m_file = new TFile((m_name + ".root").c_str(), "RECREATE");
210  savePlot(m_dedx_spectrum, m_name + "_dedx_spectrum");
211  savePlot(m_radlen_spectrum, m_name + "_radlen_spectrum");
212  savePlot(m_dedx_vs_eta, averageEnergyLoss(), m_name + "_dedx_vs_eta");
213  savePlot(m_dedx_vs_z, averageEnergyLoss(), m_name + "_dedx_vs_z");
214  savePlot(m_dedx_vs_r, averageEnergyLoss(), m_name + "_dedx_vs_r");
215  savePlot(m_radlen_vs_eta, averageRadiationLengths(), m_name + "_radlen_vs_eta");
218  m_file->Write();
219  m_file->Close();
220 
221  delete m_file;
222 }
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
MaterialAccountingGroup(const std::string &name, const DDCompactView &geometry)
explicit constructors
T perp() const
Definition: PV3DBase.h:72
double averageRadiationLengths(void) const
return the average normalized number of radiation lengths
void endOfTrack(void)
commit the buffer and reset the "already hit by this track" flag
MaterialAccountingStep m_errors
MaterialAccountingStep m_buffer
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
#define nullptr
std::vector< GlobalPoint > m_elements
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
std::pair< double, double > range_r() const
std::string info(void) const
get some infos
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
double averageEnergyLoss(void) const
return the average normalized energy loss density factor for Bethe-Bloch
bool next()
set current node to the next node in the filtered tree
void savePlot(TH1F *plot, const std::string &name)
T z() const
Definition: PV3DBase.h:64
MaterialAccountingStep average(void) const
return the average normalized material accounting informations
std::pair< double, double > range_z() const
#define LogTrace(id)
const MaterialAccountingStep & material() const
const std::string & name(void) const
get the layer name
std::string toString() const
Definition: DDBase.h:78
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void savePlots(void)
save the plots
const GlobalPoint & position() const
T eta() const
Definition: PV3DBase.h:76
T perp() const
Magnitude of transverse component.
bool inside(double r, double z) const
static int position[264][3]
Definition: ReadPGInfo.cc:509
def canvas(sub, attr)
Definition: svgfig.py:482
bool addDetector(const MaterialAccountingDetector &detector)
buffer material from a detector, if the detector is inside the DetLayer bounds
~MaterialAccountingGroup(void)
destructor
const DDTranslation & translation() const
The absolute translation of the current node.
MaterialAccountingStep m_accounting
T x() const
Definition: PV3DBase.h:62
bool inside(const MaterialAccountingDetector &detector) const
check if detector is inside any part of this layer