test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 double const MaterialAccountingGroup::s_tolerance = 0.01; // 100um should be small enough that no elements from different layers/groups are so close
21 
23  m_name( name ),
24  m_elements(),
25  m_boundingbox(),
26  m_accounting(),
27  m_errors(),
28  m_tracks( 0 ),
29  m_counted( false ),
30  m_file( 0 )
31 {
32  // retrieve the elements from DDD
33  DDFilteredView fv( geometry );
35  filter.setCriteria(DDValue("TrackingMaterialGroup", name), DDCompOp::equals);
36  fv.addFilter(filter);
37  LogTrace("MaterialAccountingGroup") << "Elements within: " << name << std::endl;
38  while (fv.next()) {
39  // DD3Vector and DDTranslation are the same type as math::XYZVector
40  math::XYZVector position = fv.translation() / 10.; // mm -> cm
41  LogTrace("MaterialAccountingGroup") << "Adding element at(r,z): ("
42  << GlobalPoint(position.x(), position.y(), position.z()).perp()
43  << ", " << GlobalPoint(position.x(), position.y(), position.z()).z()
44  << ") cm" << std::endl;
45  LogTrace("MaterialAccountingGroup") << "Name of added element: "
46  << fv.logicalPart().toString() << std::endl;
47  m_elements.push_back( GlobalPoint(position.x(), position.y(), position.z()) );
48  }
49 
50  // grow the bounding box
51  for (unsigned int i = 0; i < m_elements.size(); ++i) {
53  }
55  LogTrace("MaterialAccountingGroup") << "Final BBox r_range: "
56  << m_boundingbox.range_r().first << ", " << m_boundingbox.range_r().second
57  << std::endl
58  << "Final BBox z_range: "
59  << m_boundingbox.range_z().first << ", " << m_boundingbox.range_z().second
60  << std::endl;
61 
62  // initialize the histograms
63  m_dedx_spectrum = new TH1F((m_name + "_dedx_spectrum").c_str(), "Energy loss spectrum", 1000, 0, 1);
64  m_radlen_spectrum = new TH1F((m_name + "_radlen_spectrum").c_str(), "Radiation lengths spectrum", 1000, 0, 1);
65  m_dedx_vs_eta = new TProfile((m_name + "_dedx_vs_eta").c_str(), "Energy loss vs. eta", 600, -3, 3);
66  m_dedx_vs_z = new TProfile((m_name + "_dedx_vs_z").c_str(), "Energy loss vs. Z", 6000, -300, 300);
67  m_dedx_vs_r = new TProfile((m_name + "_dedx_vs_r").c_str(), "Energy loss vs. R", 1200, 0, 120);
68  m_radlen_vs_eta = new TProfile((m_name + "_radlen_vs_eta").c_str(), "Radiation lengths vs. eta", 600, -3, 3);
69  m_radlen_vs_z = new TProfile((m_name + "_radlen_vs_z").c_str(), "Radiation lengths vs. Z", 6000, -300, 300);
70  m_radlen_vs_r = new TProfile((m_name + "_radlen_vs_r").c_str(), "Radiation lengths vs. R", 1200, 0, 120);
71  m_dedx_spectrum->SetDirectory( 0 );
72  m_radlen_spectrum->SetDirectory( 0 );
73  m_dedx_vs_eta->SetDirectory( 0 );
74  m_dedx_vs_z->SetDirectory( 0 );
75  m_dedx_vs_r->SetDirectory( 0 );
76  m_radlen_vs_eta->SetDirectory( 0 );
77  m_radlen_vs_z->SetDirectory( 0 );
78  m_radlen_vs_r->SetDirectory( 0 );
79 }
80 
82 {
83  delete m_dedx_spectrum;
84  delete m_dedx_vs_eta;
85  delete m_dedx_vs_z;
86  delete m_dedx_vs_r;
87  delete m_radlen_spectrum;
88  delete m_radlen_vs_eta;
89  delete m_radlen_vs_z;
90  delete m_radlen_vs_r;
91 }
92 
93 // TODO the inner check could be sped up in many ways
94 // (sorting the m_elements, partitioning the bounding box, ...)
95 // but is it worth?
96 // especially with the segmentation of the layers ?
98 {
99  const GlobalPoint & position = detector.position();
100  // first check to see if the point is inside the bounding box
101  LogTrace("MaterialAccountingGroup") << "Testing position: (x, y, z, r) = "
102  << position.x() << ", " << position.y()
103  << ", " << position.z() << ", " << position.perp()
104  << std::endl;
105  if (not m_boundingbox.inside(position.perp(), position.z())) {
106  LogTrace("MaterialAccountingGroup") << "r outside of: ("
107  << m_boundingbox.range_r().first << ", "
108  << m_boundingbox.range_r().second
109  << "), Z ouside of: ("
110  << m_boundingbox.range_z().first << ", "
111  << m_boundingbox.range_z().second << ")" << std::endl;
112  return false;
113  } else {
114  // now check if the point is actually close enough to any element
115  LogTrace("MaterialAccountingGroup") << "r within: ("
116  << m_boundingbox.range_r().first << ", "
117  << m_boundingbox.range_r().second
118  << "), Z within: ("
119  << m_boundingbox.range_z().first << ", "
120  << m_boundingbox.range_z().second << ")" << std::endl;
121  for (unsigned int i = 0; i < m_elements.size(); ++i) {
122  LogTrace("MaterialAccountingGroup") << "Closest testing agains(x, y, z, r): ("
123  << m_elements[i].x() << ", " << m_elements[i].y()
124  << ", " << m_elements[i].z() << ", " << m_elements[i].perp() << ") --> "
125  << (position - m_elements[i]).mag()
126  << " vs tolerance: " << s_tolerance
127  << std::endl;
128  if ((position - m_elements[i]).mag2() < (s_tolerance * s_tolerance))
129  return true;
130  }
131  return false;
132  }
133 }
134 
136 {
137  if (not inside(detector))
138  return false;
139 
140  // multiple hits in the same layer (overlaps, etc.) from a single track still count as one for averaging,
141  // since the energy deposits from the track have been already split between the different detectors
142  m_buffer += detector.material();
143  m_counted = true;
144 
145  return true;
146 }
147 
149  // add a detector
150  if (m_counted) {
153  ++m_tracks;
154 
155  GlobalPoint average( (m_buffer.in().x() + m_buffer.out().x()) / 2.,
156  (m_buffer.in().y() + m_buffer.out().y()) / 2.,
157  (m_buffer.in().z() + m_buffer.out().z()) / 2. );
158  m_dedx_spectrum->Fill( m_buffer.energyLoss() );
159  m_radlen_spectrum->Fill( m_buffer.radiationLengths() );
160  m_dedx_vs_eta->Fill( average.eta(), m_buffer.energyLoss(), 1. );
161  m_dedx_vs_z->Fill( average.z(), m_buffer.energyLoss(), 1. );
162  m_dedx_vs_r->Fill( average.perp(), m_buffer.energyLoss(), 1. );
163  m_radlen_vs_eta->Fill( average.eta(), m_buffer.radiationLengths(), 1. );
164  m_radlen_vs_z->Fill( average.z(), m_buffer.radiationLengths(), 1. );
165  m_radlen_vs_r->Fill( average.perp(), m_buffer.radiationLengths(), 1. );
166  }
167  m_counted = false;
169 }
170 
172 {
173  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
174  plot->SetFillColor(15); // grey
175  plot->SetLineColor(1); // black
176  plot->Draw("c e");
177  canvas.GetFrame()->SetFillColor(kWhite);
178  canvas.Draw();
179  canvas.SaveAs((name + ".png").c_str(), "");
180 
181  // store te plot into m_file
182  plot->SetDirectory(m_file);
183 }
184 
186 {
187  // Nota Bene:
188  // these "line" plots are not deleted explicitly since
189  // - deleting them before saving them to a TFile will not save them
190  // - deleting them after the TFile they're stored into results in a SEGV
191  // ROOT is probably "taking care" (read: messing things up) somehow...
192  TH1F * line = new TH1F((name + "_par").c_str(), "Parametrization", 1, plot->GetXaxis()->GetXmin(), plot->GetXaxis()->GetXmax());
193  line->SetBinContent(1, average);
194 
195  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
196  plot->SetFillColor(15); // grey
197  plot->SetLineColor(1); // black
198  plot->SetLineWidth(2);
199  plot->Draw("c e6");
200  line->SetLineColor(2); // red
201  line->SetLineWidth(2);
202  line->Draw("same");
203  canvas.GetFrame()->SetFillColor(kWhite);
204  canvas.Draw();
205  canvas.SaveAs((name + ".png").c_str(), "");
206 
207  // store te plots into m_file
208  plot->SetDirectory(m_file);
209  line->SetDirectory(m_file);
210 }
211 
214 {
215  std::stringstream out;
216  out << std::setw(48) << std::left << m_name << std::right << std::fixed;;
217  out << "BBox: " << std::setprecision(1) << std::setw(6) << m_boundingbox.range_z().first << " < Z < " << std::setprecision(1) << std::setw(6) << m_boundingbox.range_z().second;
218  out << ", " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().first << " < R < " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().second;
219  out << " Elements: " << std::setw(6) << m_elements.size();
220  return out.str();
221 }
222 
223 
225 {
226  m_file = new TFile((m_name + ".root").c_str(), "RECREATE");
227  savePlot(m_dedx_spectrum, m_name + "_dedx_spectrum");
228  savePlot(m_radlen_spectrum, m_name + "_radlen_spectrum");
229  savePlot(m_dedx_vs_eta, averageEnergyLoss(), m_name + "_dedx_vs_eta");
230  savePlot(m_dedx_vs_z, averageEnergyLoss(), m_name + "_dedx_vs_z");
231  savePlot(m_dedx_vs_r, averageEnergyLoss(), m_name + "_dedx_vs_r");
232  savePlot(m_radlen_vs_eta, averageRadiationLengths(), m_name + "_radlen_vs_eta");
235  m_file->Write();
236  m_file->Close();
237 
238  delete m_file;
239 }
int i
Definition: DBlmapReader.cc:9
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 &quot;already hit by this track&quot; flag
void addFilter(const DDFilter &, DDLogOp op=DDLogOp::AND)
MaterialAccountingStep m_errors
MaterialAccountingStep m_buffer
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
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
type of data representation of DDCompactView
Definition: DDCompactView.h:90
def canvas
Definition: svgfig.py:481
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
def plot
Definition: bigModule.py:19
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
std::string toString() const
Definition: DDBase.h:82
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void savePlots(void)
save the plots
const GlobalPoint & position() const
int average
Definition: PDRates.py:137
T eta() const
Definition: PV3DBase.h:76
ESHandle< TrackerGeometry > geometry
T perp() const
Magnitude of transverse component.
bool inside(double r, double z) const
static int position[264][3]
Definition: ReadPGInfo.cc:509
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
volatile std::atomic< bool > shutdown_flag false
T x() const
Definition: PV3DBase.h:62
bool inside(const MaterialAccountingDetector &detector) const
check if detector is inside any part of this layer
void setCriteria(const DDValue &nameVal, DDCompOp, DDLogOp l=DDLogOp::AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:253
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:33