test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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...
 
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)
 stop default copy ctor More...
 
MaterialAccountingGroupoperator= (const MaterialAccountingGroup &layer)
 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 18 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 DDFilteredView::addFilter(), equals, alcazmumu_cfi::filter, MaterialAccountingGroup::BoundingBox::grow(), 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, DDSpecificsFilter::setCriteria(), 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( 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 }
int i
Definition: DBlmapReader.cc:9
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
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
MaterialAccountingGroup::~MaterialAccountingGroup ( void  )
MaterialAccountingGroup::MaterialAccountingGroup ( const MaterialAccountingGroup layer)
private

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

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

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 }
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 103 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 121 of file MaterialAccountingGroup.h.

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

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

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

return the average normalized layer thickness

Definition at line 109 of file MaterialAccountingGroup.h.

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

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

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

return the average normalized number of radiation lengths

Definition at line 115 of file MaterialAccountingGroup.h.

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

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

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

Definition at line 166 of file MaterialAccountingGroup.h.

References m_elements.

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

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

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

148  {
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;
168  m_buffer = MaterialAccountingStep();
169 }
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 97 of file MaterialAccountingGroup.h.

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

97 {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 100 of file MaterialAccountingGroup.h.

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

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

get some infos

Definition at line 213 of file MaterialAccountingGroup.cc.

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

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

References 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().

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 }
int i
Definition: DBlmapReader.cc:9
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)
private

stop default assignment operator

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

Definition at line 171 of file MaterialAccountingGroup.cc.

References svgfig::canvas(), and m_file.

Referenced by savePlots().

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 }
def canvas
Definition: svgfig.py:481
def plot
Definition: bigModule.py:19
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 185 of file MaterialAccountingGroup.cc.

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

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 }
def canvas
Definition: svgfig.py:481
def plot
Definition: bigModule.py:19
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 224 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().

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 }
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 139 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

140  {
142  }
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 127 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

128  {
130  }
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 133 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

134  {
136  }
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 145 of file MaterialAccountingGroup.h.

References m_tracks.

Referenced by TrackingMaterialAnalyser::saveParameters().

146  {
147  return m_tracks;
148  }

Member Data Documentation

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

Definition at line 181 of file MaterialAccountingGroup.h.

Referenced by addDetector(), and endOfTrack().

bool MaterialAccountingGroup::m_counted
private

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

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

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

Definition at line 195 of file MaterialAccountingGroup.h.

Referenced by savePlot(), and savePlots().

std::string MaterialAccountingGroup::m_name
private

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

Referenced by inside(), and MaterialAccountingGroup().