CMS 3D CMS Logo

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...
 
const std::string & name (void) const
 get the layer name 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

 MaterialAccountingGroup (const MaterialAccountingGroup &layer)=delete
 stop default copy ctor More...
 
MaterialAccountingGroupoperator= (const MaterialAccountingGroup &layer)=delete
 stop default assignment operator More...
 
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 = 0.01
 

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 22 of file MaterialAccountingGroup.cc.

References ALCARECOTkAlBeamHalo_cff::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.

22  :
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( nullptr )
31 {
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()
42  << ") cm" << std::endl;
43  LogTrace("MaterialAccountingGroup") << "Name of added element: "
44  << fv.logicalPart().toString() << std::endl;
45  m_elements.push_back( GlobalPoint(position.x(), position.y(), position.z()) );
46  }
47 
48  // grow the bounding box
49  for (unsigned int i = 0; i < m_elements.size(); ++i) {
51  }
53  LogTrace("MaterialAccountingGroup") << "Final BBox r_range: "
54  << m_boundingbox.range_r().first << ", " << m_boundingbox.range_r().second
55  << std::endl
56  << "Final BBox z_range: "
57  << m_boundingbox.range_z().first << ", " << m_boundingbox.range_z().second
58  << std::endl;
59 
60  // initialize the histograms
61  m_dedx_spectrum = new TH1F((m_name + "_dedx_spectrum").c_str(), "Energy loss spectrum", 1000, 0, 1);
62  m_radlen_spectrum = new TH1F((m_name + "_radlen_spectrum").c_str(), "Radiation lengths spectrum", 1000, 0, 1);
63  m_dedx_vs_eta = new TProfile((m_name + "_dedx_vs_eta").c_str(), "Energy loss vs. eta", 600, -3, 3);
64  m_dedx_vs_z = new TProfile((m_name + "_dedx_vs_z").c_str(), "Energy loss vs. Z", 6000, -300, 300);
65  m_dedx_vs_r = new TProfile((m_name + "_dedx_vs_r").c_str(), "Energy loss vs. R", 1200, 0, 120);
66  m_radlen_vs_eta = new TProfile((m_name + "_radlen_vs_eta").c_str(), "Radiation lengths vs. eta", 600, -3, 3);
67  m_radlen_vs_z = new TProfile((m_name + "_radlen_vs_z").c_str(), "Radiation lengths vs. Z", 6000, -300, 300);
68  m_radlen_vs_r = new TProfile((m_name + "_radlen_vs_r").c_str(), "Radiation lengths vs. R", 1200, 0, 120);
69  m_dedx_spectrum->SetDirectory( nullptr );
70  m_radlen_spectrum->SetDirectory( nullptr );
71  m_dedx_vs_eta->SetDirectory( nullptr );
72  m_dedx_vs_z->SetDirectory( nullptr );
73  m_dedx_vs_r->SetDirectory( nullptr );
74  m_radlen_vs_eta->SetDirectory( nullptr );
75  m_radlen_vs_z->SetDirectory( nullptr );
76  m_radlen_vs_r->SetDirectory( nullptr );
77 }
MaterialAccountingStep m_errors
std::vector< GlobalPoint > m_elements
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::pair< double, double > range_r() const
std::pair< double, double > range_z() const
#define LogTrace(id)
const std::string & name(void) const
get the layer name
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
T perp() const
Magnitude of transverse component.
static int position[264][3]
Definition: ReadPGInfo.cc:509
MaterialAccountingStep m_accounting
MaterialAccountingGroup::~MaterialAccountingGroup ( void  )
MaterialAccountingGroup::MaterialAccountingGroup ( const MaterialAccountingGroup layer)
privatedelete

stop default copy ctor

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 133 of file MaterialAccountingGroup.cc.

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

134 {
135  if (not inside(detector))
136  return false;
137 
138  // multiple hits in the same layer (overlaps, etc.) from a single track still count as one for averaging,
139  // since the energy deposits from the track have been already split between the different detectors
140  m_buffer += detector.material();
141  m_counted = true;
142 
143  return true;
144 }
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 104 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 122 of file MaterialAccountingGroup.h.

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

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

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

return the average normalized layer thickness

Definition at line 110 of file MaterialAccountingGroup.h.

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

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

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

return the average normalized number of radiation lengths

Definition at line 116 of file MaterialAccountingGroup.h.

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

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

117  {
119  }
double radiationLengths(void) const
MaterialAccountingStep m_accounting
const std::vector<GlobalPoint>& MaterialAccountingGroup::elements ( void  ) const
inline

Definition at line 167 of file MaterialAccountingGroup.h.

References m_elements, name(), plotFactory::plot, savePlot(), and AlCaHLTBitMon_QueryRunRegistry::string.

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

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

Definition at line 146 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().

146  {
147  // add a detector
148  if (m_counted) {
151  ++m_tracks;
152 
153  GlobalPoint average( (m_buffer.in().x() + m_buffer.out().x()) / 2.,
154  (m_buffer.in().y() + m_buffer.out().y()) / 2.,
155  (m_buffer.in().z() + m_buffer.out().z()) / 2. );
156  m_dedx_spectrum->Fill( m_buffer.energyLoss() );
157  m_radlen_spectrum->Fill( m_buffer.radiationLengths() );
158  m_dedx_vs_eta->Fill( average.eta(), m_buffer.energyLoss(), 1. );
159  m_dedx_vs_z->Fill( average.z(), m_buffer.energyLoss(), 1. );
160  m_dedx_vs_r->Fill( average.perp(), m_buffer.energyLoss(), 1. );
161  m_radlen_vs_eta->Fill( average.eta(), m_buffer.radiationLengths(), 1. );
162  m_radlen_vs_z->Fill( average.z(), m_buffer.radiationLengths(), 1. );
163  m_radlen_vs_r->Fill( average.perp(), m_buffer.radiationLengths(), 1. );
164  }
165  m_counted = false;
166  m_buffer = MaterialAccountingStep();
167 }
MaterialAccountingStep m_errors
MaterialAccountingStep m_buffer
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 98 of file MaterialAccountingGroup.h.

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

98 {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 101 of file MaterialAccountingGroup.h.

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

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

get some infos

Definition at line 211 of file MaterialAccountingGroup.cc.

References alignBH_cfg::fixed, m_boundingbox, m_elements, m_name, MillePedeFileConverter_cfg::out, MaterialAccountingGroup::BoundingBox::range_r(), and MaterialAccountingGroup::BoundingBox::range_z().

Referenced by name().

212 {
213  std::stringstream out;
214  out << std::setw(48) << std::left << m_name << std::right << std::fixed;;
215  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;
216  out << ", " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().first << " < R < " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().second;
217  out << " Elements: " << std::setw(6) << m_elements.size();
218  return out.str();
219 }
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 95 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().

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

stop default assignment operator

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

Definition at line 169 of file MaterialAccountingGroup.cc.

References svgfig::canvas(), and m_file.

Referenced by elements(), and savePlots().

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

Definition at line 183 of file MaterialAccountingGroup.cc.

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

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

save the plots

Definition at line 222 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 name(), and TrackingMaterialAnalyser::saveLayerPlots().

223 {
224  m_file = new TFile((m_name + ".root").c_str(), "RECREATE");
225  savePlot(m_dedx_spectrum, m_name + "_dedx_spectrum");
226  savePlot(m_radlen_spectrum, m_name + "_radlen_spectrum");
227  savePlot(m_dedx_vs_eta, averageEnergyLoss(), m_name + "_dedx_vs_eta");
228  savePlot(m_dedx_vs_z, averageEnergyLoss(), m_name + "_dedx_vs_z");
229  savePlot(m_dedx_vs_r, averageEnergyLoss(), m_name + "_dedx_vs_r");
230  savePlot(m_radlen_vs_eta, averageRadiationLengths(), m_name + "_radlen_vs_eta");
233  m_file->Write();
234  m_file->Close();
235 
236  delete m_file;
237 }
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 140 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

141  {
143  }
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:18
double MaterialAccountingGroup::sigmaLength ( void  ) const
inline

return the sigma of the normalized layer thickness

Definition at line 128 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

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

return the sigma of the normalized number of radiation lengths

Definition at line 134 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

135  {
137  }
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:18
unsigned int MaterialAccountingGroup::tracks ( void  ) const
inline

return the number of tracks that hit this layer

Definition at line 146 of file MaterialAccountingGroup.h.

References m_tracks.

Referenced by TrackingMaterialAnalyser::saveParameters().

147  {
148  return m_tracks;
149  }

Member Data Documentation

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

Definition at line 182 of file MaterialAccountingGroup.h.

Referenced by addDetector(), and endOfTrack().

bool MaterialAccountingGroup::m_counted
private

Definition at line 181 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 176 of file MaterialAccountingGroup.h.

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

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

Definition at line 196 of file MaterialAccountingGroup.h.

Referenced by savePlot(), and savePlots().

std::string MaterialAccountingGroup::m_name
private

Definition at line 175 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 = 0.01
staticprivate

Definition at line 199 of file MaterialAccountingGroup.h.

Referenced by inside(), and MaterialAccountingGroup().