CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
MaterialAccountingGroup Class Reference

#include <MaterialAccountingGroup.h>

Classes

class  BoundingBox
 

Public Member Functions

bool addDetector (const MaterialAccountingDetector &detector)
 buffer material from a detector, if the detector is inside the DetLayer bounds More...
 
MaterialAccountingStep average (void) const
 return the average normalized material accounting informations More...
 
double averageEnergyLoss (void) const
 return the average normalized energy loss density factor for Bethe-Bloch More...
 
double averageLength (void) const
 return the average normalized layer thickness More...
 
double averageRadiationLengths (void) const
 return the average normalized number of radiation lengths More...
 
const std::vector< GlobalPoint > & elements (void) const
 
void endOfTrack (void)
 commit the buffer and reset the "already hit by this track" flag More...
 
std::pair< double, double > getBoundingR () const
 Return the bouding limit in R for the hosted Group. More...
 
std::pair< double, double > getBoundingZ () const
 Return the bouding limit in Z for the hosted Group. More...
 
std::string info (void) const
 get some infos More...
 
bool inside (const MaterialAccountingDetector &detector) const
 check if detector is inside any part of this layer More...
 
 MaterialAccountingGroup (const std::string &name, const DDCompactView &geometry)
 explicit constructors More...
 
 MaterialAccountingGroup (const MaterialAccountingGroup &layer)=delete
 stop default copy ctor More...
 
const std::string & name (void) const
 get the layer name More...
 
MaterialAccountingGroupoperator= (const MaterialAccountingGroup &layer)=delete
 stop default assignment operator More...
 
void savePlots (void)
 save the plots More...
 
double sigmaEnergyLoss (void) const
 return the sigma of the normalized energy loss density factor for Bethe-Bloch More...
 
double sigmaLength (void) const
 return the sigma of the normalized layer thickness More...
 
double sigmaRadiationLengths (void) const
 return the sigma of the normalized number of radiation lengths More...
 
unsigned int tracks (void) const
 return the number of tracks that hit this layer More...
 
 ~MaterialAccountingGroup (void)
 destructor More...
 

Private Member Functions

void savePlot (TH1F *plot, const std::string &name)
 
void savePlot (TProfile *plot, float average, const std::string &name)
 

Private Attributes

MaterialAccountingStep m_accounting
 
BoundingBox m_boundingbox
 
MaterialAccountingStep m_buffer
 
bool m_counted
 
TH1F * m_dedx_spectrum
 
TProfile * m_dedx_vs_eta
 
TProfile * m_dedx_vs_r
 
TProfile * m_dedx_vs_z
 
std::vector< GlobalPointm_elements
 
MaterialAccountingStep m_errors
 
TFile * m_file
 
std::string m_name
 
TH1F * m_radlen_spectrum
 
TProfile * m_radlen_vs_eta
 
TProfile * m_radlen_vs_r
 
TProfile * m_radlen_vs_z
 
unsigned int m_tracks
 

Static Private Attributes

static double const s_tolerance
 

Detailed Description

Definition at line 19 of file MaterialAccountingGroup.h.

Constructor & Destructor Documentation

MaterialAccountingGroup::MaterialAccountingGroup ( const std::string &  name,
const DDCompactView geometry 
)

explicit constructors

Definition at line 23 of file MaterialAccountingGroup.cc.

References alcazmumu_cfi::filter, MaterialAccountingGroup::BoundingBox::grow(), mps_fire::i, DDFilteredView::logicalPart(), LogTrace, m_boundingbox, m_dedx_spectrum, m_dedx_vs_eta, m_dedx_vs_r, m_dedx_vs_z, m_elements, m_name, m_radlen_spectrum, m_radlen_vs_eta, m_radlen_vs_r, m_radlen_vs_z, DDFilteredView::next(), perp(), position, MaterialAccountingGroup::BoundingBox::range_r(), MaterialAccountingGroup::BoundingBox::range_z(), s_tolerance, DDBase< N, C >::toString(), DDFilteredView::translation(), and z.

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  DDValue namevalue;
34  if (TString(name.c_str()).Contains("Tracker")) {
35  namevalue = DDValue("TrackingMaterialGroup", name);
36  } else if (TString(name.c_str()).Contains("HGCal")) {
37  namevalue = DDValue("Volume", name);
38  } else if (TString(name.c_str()).Contains("HFNose")) {
39  namevalue = DDValue("Volume", name);
40  } else {
41  LogTrace("MaterialAccountingGroup") << name << std::endl;
42  edm::LogError("MaterialAccountingGroup") << "Only Tracker , HGCal and HFNose are supported" << std::endl;
43  }
44 
46  DDFilteredView fv(geometry, filter);
47  LogTrace("MaterialAccountingGroup") << "Elements within: " << name << std::endl;
48  while (fv.next()) {
49  // DD3Vector and DDTranslation are the same type as math::XYZVector
50  math::XYZVector position = fv.translation() / 10.; // mm -> cm
51  LogTrace("MaterialAccountingGroup") << "Adding element at(r,z): ("
52  << GlobalPoint(position.x(), position.y(), position.z()).perp() << ", "
53  << GlobalPoint(position.x(), position.y(), position.z()).z() << ") cm"
54  << std::endl;
55  LogTrace("MaterialAccountingGroup") << "Name of added element: " << fv.logicalPart().toString() << std::endl;
56  m_elements.push_back(GlobalPoint(position.x(), position.y(), position.z()));
57  }
58 
59  // grow the bounding box
60  for (unsigned int i = 0; i < m_elements.size(); ++i) {
62  }
64  LogTrace("MaterialAccountingGroup") << "Final BBox r_range: " << m_boundingbox.range_r().first << ", "
65  << m_boundingbox.range_r().second << std::endl
66  << "Final BBox z_range: " << m_boundingbox.range_z().first << ", "
67  << m_boundingbox.range_z().second << std::endl;
68 
69  // initialize the histograms
70  m_dedx_spectrum = new TH1F((m_name + "_dedx_spectrum").c_str(), "Energy loss spectrum", 1000, 0, 1);
71  m_radlen_spectrum = new TH1F((m_name + "_radlen_spectrum").c_str(), "Radiation lengths spectrum", 1000, 0, 1);
72  m_dedx_vs_eta = new TProfile((m_name + "_dedx_vs_eta").c_str(), "Energy loss vs. eta", 600, -3, 3);
73  m_dedx_vs_z = new TProfile((m_name + "_dedx_vs_z").c_str(), "Energy loss vs. Z", 6000, -300, 300);
74  m_dedx_vs_r = new TProfile((m_name + "_dedx_vs_r").c_str(), "Energy loss vs. R", 1200, 0, 120);
75  m_radlen_vs_eta = new TProfile((m_name + "_radlen_vs_eta").c_str(), "Radiation lengths vs. eta", 600, -3, 3);
76  m_radlen_vs_z = new TProfile((m_name + "_radlen_vs_z").c_str(), "Radiation lengths vs. Z", 6000, -300, 300);
77  m_radlen_vs_r = new TProfile((m_name + "_radlen_vs_r").c_str(), "Radiation lengths vs. R", 1200, 0, 120);
78  m_dedx_spectrum->SetDirectory(nullptr);
79  m_radlen_spectrum->SetDirectory(nullptr);
80  m_dedx_vs_eta->SetDirectory(nullptr);
81  m_dedx_vs_z->SetDirectory(nullptr);
82  m_dedx_vs_r->SetDirectory(nullptr);
83  m_radlen_vs_eta->SetDirectory(nullptr);
84  m_radlen_vs_z->SetDirectory(nullptr);
85  m_radlen_vs_r->SetDirectory(nullptr);
86 }
MaterialAccountingStep m_errors
std::vector< GlobalPoint > m_elements
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::pair< double, double > range_r() const
Log< level::Error, false > LogError
#define LogTrace(id)
std::pair< double, double > range_z() const
const std::string & name(void) const
get the layer name
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
T perp() const
Magnitude of transverse component.
static int position[264][3]
Definition: ReadPGInfo.cc:289
MaterialAccountingStep m_accounting
MaterialAccountingGroup::MaterialAccountingGroup ( const MaterialAccountingGroup layer)
delete

stop default copy ctor

MaterialAccountingGroup::~MaterialAccountingGroup ( void  )

Member Function Documentation

bool MaterialAccountingGroup::addDetector ( const MaterialAccountingDetector detector)

buffer material from a detector, if the detector is inside the DetLayer bounds

Definition at line 132 of file MaterialAccountingGroup.cc.

References inside(), m_buffer, m_counted, and MaterialAccountingDetector::material().

132  {
133  if (not inside(detector))
134  return false;
135 
136  // multiple hits in the same layer (overlaps, etc.) from a single track still count as one for averaging,
137  // since the energy deposits from the track have been already split between the different detectors
138  m_buffer += detector.material();
139  m_counted = true;
140 
141  return true;
142 }
MaterialAccountingStep m_buffer
const MaterialAccountingStep & material() const
bool inside(const MaterialAccountingDetector &detector) const
check if detector is inside any part of this layer
MaterialAccountingStep MaterialAccountingGroup::average ( void  ) const
inline

return the average normalized material accounting informations

Definition at line 90 of file MaterialAccountingGroup.h.

References m_accounting, and m_tracks.

Referenced by endOfTrack().

double MaterialAccountingGroup::averageEnergyLoss ( void  ) const
inline

return the average normalized energy loss density factor for Bethe-Bloch

Definition at line 99 of file MaterialAccountingGroup.h.

References MaterialAccountingStep::energyLoss(), m_accounting, and m_tracks.

Referenced by TrackingMaterialAnalyser::saveParameters(), savePlots(), TrackingMaterialAnalyser::saveXml(), and sigmaEnergyLoss().

99 { return m_tracks ? m_accounting.energyLoss() / m_tracks : 0.; }
double energyLoss(void) const
MaterialAccountingStep m_accounting
double MaterialAccountingGroup::averageLength ( void  ) const
inline

return the average normalized layer thickness

Definition at line 93 of file MaterialAccountingGroup.h.

References MaterialAccountingStep::length(), m_accounting, and m_tracks.

Referenced by TrackingMaterialAnalyser::saveParameters(), and sigmaLength().

93 { return m_tracks ? m_accounting.length() / m_tracks : 0.; }
double length(void) const
MaterialAccountingStep m_accounting
double MaterialAccountingGroup::averageRadiationLengths ( void  ) const
inline

return the average normalized number of radiation lengths

Definition at line 96 of file MaterialAccountingGroup.h.

References m_accounting, m_tracks, and MaterialAccountingStep::radiationLengths().

Referenced by TrackingMaterialAnalyser::saveParameters(), savePlots(), TrackingMaterialAnalyser::saveXml(), and sigmaRadiationLengths().

96 { return m_tracks ? m_accounting.radiationLengths() / m_tracks : 0.; }
double radiationLengths(void) const
MaterialAccountingStep m_accounting
const std::vector<GlobalPoint>& MaterialAccountingGroup::elements ( void  ) const
inline

Definition at line 134 of file MaterialAccountingGroup.h.

References m_elements.

134 { return m_elements; }
std::vector< GlobalPoint > m_elements
void MaterialAccountingGroup::endOfTrack ( void  )

commit the buffer and reset the "already hit by this track" flag

Definition at line 144 of file MaterialAccountingGroup.cc.

References average(), PV3DBase< T, PVType, FrameType >::eta(), m_accounting, m_buffer, m_counted, m_dedx_spectrum, m_dedx_vs_eta, m_dedx_vs_r, m_dedx_vs_z, m_errors, m_radlen_spectrum, m_radlen_vs_eta, m_radlen_vs_r, m_radlen_vs_z, m_tracks, PV3DBase< T, PVType, FrameType >::perp(), and PV3DBase< T, PVType, FrameType >::z().

144  {
145  // add a detector
146  if (m_counted) {
149  ++m_tracks;
150 
151  GlobalPoint average((m_buffer.in().x() + m_buffer.out().x()) / 2.,
152  (m_buffer.in().y() + m_buffer.out().y()) / 2.,
153  (m_buffer.in().z() + m_buffer.out().z()) / 2.);
154  m_dedx_spectrum->Fill(m_buffer.energyLoss());
155  m_radlen_spectrum->Fill(m_buffer.radiationLengths());
156  m_dedx_vs_eta->Fill(average.eta(), m_buffer.energyLoss(), 1.);
157  m_dedx_vs_z->Fill(average.z(), m_buffer.energyLoss(), 1.);
158  m_dedx_vs_r->Fill(average.perp(), m_buffer.energyLoss(), 1.);
159  m_radlen_vs_eta->Fill(average.eta(), m_buffer.radiationLengths(), 1.);
160  m_radlen_vs_z->Fill(average.z(), m_buffer.radiationLengths(), 1.);
161  m_radlen_vs_r->Fill(average.perp(), m_buffer.radiationLengths(), 1.);
162  }
163  m_counted = false;
164  m_buffer = MaterialAccountingStep();
165 }
MaterialAccountingStep m_errors
MaterialAccountingStep m_buffer
double radiationLengths(void) const
double energyLoss(void) const
MaterialAccountingStep average(void) const
return the average normalized material accounting informations
MaterialAccountingStep m_accounting
std::pair<double, double> MaterialAccountingGroup::getBoundingR ( ) const
inline

Return the bouding limit in R for the hosted Group.

Definition at line 84 of file MaterialAccountingGroup.h.

References m_boundingbox, and MaterialAccountingGroup::BoundingBox::range_r().

84 { return m_boundingbox.range_r(); };
std::pair< double, double > range_r() const
std::pair<double, double> MaterialAccountingGroup::getBoundingZ ( ) const
inline

Return the bouding limit in Z for the hosted Group.

Definition at line 87 of file MaterialAccountingGroup.h.

References m_boundingbox, and MaterialAccountingGroup::BoundingBox::range_z().

87 { return m_boundingbox.range_z(); };
std::pair< double, double > range_z() const
std::string MaterialAccountingGroup::info ( void  ) const

get some infos

Definition at line 208 of file MaterialAccountingGroup.cc.

References m_boundingbox, m_elements, m_name, submitPVResolutionJobs::out, MaterialAccountingGroup::BoundingBox::range_r(), and MaterialAccountingGroup::BoundingBox::range_z().

208  {
209  std::stringstream out;
210  out << std::setw(48) << std::left << m_name << std::right << std::fixed;
211  ;
212  out << "BBox: " << std::setprecision(1) << std::setw(6) << m_boundingbox.range_z().first << " < Z < "
213  << std::setprecision(1) << std::setw(6) << m_boundingbox.range_z().second;
214  out << ", " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().first << " < R < "
215  << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().second;
216  out << " Elements: " << std::setw(6) << m_elements.size();
217  return out.str();
218 }
std::vector< GlobalPoint > m_elements
std::pair< double, double > range_r() const
std::pair< double, double > range_z() const
bool MaterialAccountingGroup::inside ( const MaterialAccountingDetector detector) const

check if detector is inside any part of this layer

Definition at line 103 of file MaterialAccountingGroup.cc.

References mps_fire::i, MaterialAccountingGroup::BoundingBox::inside(), LogTrace, m_boundingbox, m_elements, mag(), PV3DBase< T, PVType, FrameType >::perp(), MaterialAccountingDetector::position(), position, MaterialAccountingGroup::BoundingBox::range_r(), MaterialAccountingGroup::BoundingBox::range_z(), s_tolerance, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by addDetector().

103  {
104  const GlobalPoint& position = detector.position();
105  // first check to see if the point is inside the bounding box
106  LogTrace("MaterialAccountingGroup") << "Testing position: (x, y, z, r) = " << position.x() << ", " << position.y()
107  << ", " << position.z() << ", " << position.perp() << std::endl;
108  if (not m_boundingbox.inside(position.perp(), position.z())) {
109  LogTrace("MaterialAccountingGroup") << "r outside of: (" << m_boundingbox.range_r().first << ", "
110  << m_boundingbox.range_r().second << "), Z ouside of: ("
111  << m_boundingbox.range_z().first << ", " << m_boundingbox.range_z().second
112  << ")" << std::endl;
113  return false;
114  } else {
115  // now check if the point is actually close enough to any element
116  LogTrace("MaterialAccountingGroup") << "r within: (" << m_boundingbox.range_r().first << ", "
117  << m_boundingbox.range_r().second << "), Z within: ("
118  << m_boundingbox.range_z().first << ", " << m_boundingbox.range_z().second
119  << ")" << std::endl;
120  for (unsigned int i = 0; i < m_elements.size(); ++i) {
121  LogTrace("MaterialAccountingGroup")
122  << "Closest testing agains(x, y, z, r): (" << m_elements[i].x() << ", " << m_elements[i].y() << ", "
123  << m_elements[i].z() << ", " << m_elements[i].perp() << ") --> " << (position - m_elements[i]).mag()
124  << " vs tolerance: " << s_tolerance << std::endl;
125  if ((position - m_elements[i]).mag2() < (s_tolerance * s_tolerance))
126  return true;
127  }
128  return false;
129  }
130 }
T perp() const
Definition: PV3DBase.h:69
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< GlobalPoint > m_elements
T y() const
Definition: PV3DBase.h:60
std::pair< double, double > range_r() const
#define LogTrace(id)
T z() const
Definition: PV3DBase.h:61
std::pair< double, double > range_z() const
const GlobalPoint & position() const
bool inside(double r, double z) const
static int position[264][3]
Definition: ReadPGInfo.cc:289
T x() const
Definition: PV3DBase.h:59
const std::string& MaterialAccountingGroup::name ( void  ) const
inline
MaterialAccountingGroup& MaterialAccountingGroup::operator= ( const MaterialAccountingGroup layer)
delete

stop default assignment operator

void MaterialAccountingGroup::savePlot ( TH1F *  plot,
const std::string &  name 
)
private

Definition at line 167 of file MaterialAccountingGroup.cc.

References svgfig::canvas(), and m_file.

Referenced by savePlots().

167  {
168  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
169  plot->SetFillColor(15); // grey
170  plot->SetLineColor(1); // black
171  plot->Draw("c e");
172  canvas.GetFrame()->SetFillColor(kWhite);
173  canvas.Draw();
174  canvas.SaveAs((name + ".png").c_str(), "");
175 
176  // store te plot into m_file
177  plot->SetDirectory(m_file);
178 }
def canvas
Definition: svgfig.py:482
def plot
Definition: bigModule.py:18
const std::string & name(void) const
get the layer name
void MaterialAccountingGroup::savePlot ( TProfile *  plot,
float  average,
const std::string &  name 
)
private

Definition at line 180 of file MaterialAccountingGroup.cc.

References svgfig::canvas(), geometryCSVtoXML::line, and m_file.

180  {
181  // Nota Bene:
182  // these "line" plots are not deleted explicitly since
183  // - deleting them before saving them to a TFile will not save them
184  // - deleting them after the TFile they're stored into results in a SEGV
185  // ROOT is probably "taking care" (read: messing things up) somehow...
186  TH1F* line =
187  new TH1F((name + "_par").c_str(), "Parametrization", 1, plot->GetXaxis()->GetXmin(), plot->GetXaxis()->GetXmax());
188  line->SetBinContent(1, average);
189 
190  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
191  plot->SetFillColor(15); // grey
192  plot->SetLineColor(1); // black
193  plot->SetLineWidth(2);
194  plot->Draw("c e6");
195  line->SetLineColor(2); // red
196  line->SetLineWidth(2);
197  line->Draw("same");
198  canvas.GetFrame()->SetFillColor(kWhite);
199  canvas.Draw();
200  canvas.SaveAs((name + ".png").c_str(), "");
201 
202  // store te plots into m_file
203  plot->SetDirectory(m_file);
204  line->SetDirectory(m_file);
205 }
def canvas
Definition: svgfig.py:482
def plot
Definition: bigModule.py:18
MaterialAccountingStep average(void) const
return the average normalized material accounting informations
const std::string & name(void) const
get the layer name
void MaterialAccountingGroup::savePlots ( void  )

save the plots

Definition at line 220 of file MaterialAccountingGroup.cc.

References averageEnergyLoss(), averageRadiationLengths(), m_dedx_spectrum, m_dedx_vs_eta, m_dedx_vs_r, m_dedx_vs_z, m_file, m_name, m_radlen_spectrum, m_radlen_vs_eta, m_radlen_vs_r, m_radlen_vs_z, and savePlot().

Referenced by TrackingMaterialAnalyser::saveLayerPlots().

220  {
221  m_file = new TFile((m_name + ".root").c_str(), "RECREATE");
222  savePlot(m_dedx_spectrum, m_name + "_dedx_spectrum");
223  savePlot(m_radlen_spectrum, m_name + "_radlen_spectrum");
224  savePlot(m_dedx_vs_eta, averageEnergyLoss(), m_name + "_dedx_vs_eta");
225  savePlot(m_dedx_vs_z, averageEnergyLoss(), m_name + "_dedx_vs_z");
226  savePlot(m_dedx_vs_r, averageEnergyLoss(), m_name + "_dedx_vs_r");
227  savePlot(m_radlen_vs_eta, averageRadiationLengths(), m_name + "_radlen_vs_eta");
230  m_file->Write();
231  m_file->Close();
232 
233  delete m_file;
234 }
double averageRadiationLengths(void) const
return the average normalized number of radiation lengths
double averageEnergyLoss(void) const
return the average normalized energy loss density factor for Bethe-Bloch
void savePlot(TH1F *plot, const std::string &name)
double MaterialAccountingGroup::sigmaEnergyLoss ( void  ) const
inline

return the sigma of the normalized energy loss density factor for Bethe-Bloch

Definition at line 114 of file MaterialAccountingGroup.h.

References averageEnergyLoss(), MaterialAccountingStep::energyLoss(), m_errors, m_tracks, and mathSSE::sqrt().

Referenced by TrackingMaterialAnalyser::saveParameters().

114  {
116  }
MaterialAccountingStep m_errors
double energyLoss(void) const
double averageEnergyLoss(void) const
return the average normalized energy loss density factor for Bethe-Bloch
T sqrt(T t)
Definition: SSEVec.h:19
double MaterialAccountingGroup::sigmaLength ( void  ) const
inline

return the sigma of the normalized layer thickness

Definition at line 102 of file MaterialAccountingGroup.h.

References averageLength(), MaterialAccountingStep::length(), m_errors, m_tracks, and mathSSE::sqrt().

Referenced by TrackingMaterialAnalyser::saveParameters().

102  {
104  }
MaterialAccountingStep m_errors
double length(void) const
double averageLength(void) const
return the average normalized layer thickness
T sqrt(T t)
Definition: SSEVec.h:19
double MaterialAccountingGroup::sigmaRadiationLengths ( void  ) const
inline

return the sigma of the normalized number of radiation lengths

Definition at line 107 of file MaterialAccountingGroup.h.

References averageRadiationLengths(), m_errors, m_tracks, MaterialAccountingStep::radiationLengths(), and mathSSE::sqrt().

Referenced by TrackingMaterialAnalyser::saveParameters().

107  {
110  : 0.;
111  }
double averageRadiationLengths(void) const
return the average normalized number of radiation lengths
MaterialAccountingStep m_errors
double radiationLengths(void) const
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int MaterialAccountingGroup::tracks ( void  ) const
inline

return the number of tracks that hit this layer

Definition at line 119 of file MaterialAccountingGroup.h.

References m_tracks.

Referenced by TrackingMaterialAnalyser::saveParameters().

119 { return m_tracks; }

Member Data Documentation

MaterialAccountingStep MaterialAccountingGroup::m_accounting
private
BoundingBox MaterialAccountingGroup::m_boundingbox
private
MaterialAccountingStep MaterialAccountingGroup::m_buffer
private

Definition at line 147 of file MaterialAccountingGroup.h.

Referenced by addDetector(), and endOfTrack().

bool MaterialAccountingGroup::m_counted
private

Definition at line 146 of file MaterialAccountingGroup.h.

Referenced by addDetector(), and endOfTrack().

TH1F* MaterialAccountingGroup::m_dedx_spectrum
private
TProfile* MaterialAccountingGroup::m_dedx_vs_eta
private
TProfile* MaterialAccountingGroup::m_dedx_vs_r
private
TProfile* MaterialAccountingGroup::m_dedx_vs_z
private
std::vector<GlobalPoint> MaterialAccountingGroup::m_elements
private

Definition at line 141 of file MaterialAccountingGroup.h.

Referenced by elements(), info(), inside(), and MaterialAccountingGroup().

MaterialAccountingStep MaterialAccountingGroup::m_errors
private
TFile* MaterialAccountingGroup::m_file
mutableprivate

Definition at line 161 of file MaterialAccountingGroup.h.

Referenced by savePlot(), and savePlots().

std::string MaterialAccountingGroup::m_name
private

Definition at line 140 of file MaterialAccountingGroup.h.

Referenced by info(), MaterialAccountingGroup(), name(), and savePlots().

TH1F* MaterialAccountingGroup::m_radlen_spectrum
private
TProfile* MaterialAccountingGroup::m_radlen_vs_eta
private
TProfile* MaterialAccountingGroup::m_radlen_vs_r
private
TProfile* MaterialAccountingGroup::m_radlen_vs_z
private
unsigned int MaterialAccountingGroup::m_tracks
private
double const MaterialAccountingGroup::s_tolerance
staticprivate
Initial value:
=
0.01

Definition at line 164 of file MaterialAccountingGroup.h.

Referenced by inside(), and MaterialAccountingGroup().