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...
 
void endOfTrack (void)
 commit the buffer and reset the "already hit by this track" flag 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 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(), DDSpecificsFilter::equals, align_tpl::filter, MaterialAccountingGroup::BoundingBox::grow(), i, 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, s_tolerance, DDSpecificsFilter::setCriteria(), 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), DDSpecificsFilter::equals);
36  fv.addFilter(filter);
37  while (fv.next()) {
38  // DD3Vector and DDTranslation are the same type as math::XYZVector
39  math::XYZVector position = fv.translation() / 10.; // mm -> cm
40  m_elements.push_back( GlobalPoint(position.x(), position.y(), position.z()) );
41  }
42 
43  // grow the bounding box
44  for (unsigned int i = 0; i < m_elements.size(); ++i) {
46  }
48 
49  // initialize the histograms
50  m_dedx_spectrum = new TH1F((m_name + "_dedx_spectrum").c_str(), "Energy loss spectrum", 1000, 0, 1);
51  m_radlen_spectrum = new TH1F((m_name + "_radlen_spectrum").c_str(), "Radiation lengths spectrum", 1000, 0, 1);
52  m_dedx_vs_eta = new TProfile((m_name + "_dedx_vs_eta").c_str(), "Energy loss vs. eta", 600, -3, 3);
53  m_dedx_vs_z = new TProfile((m_name + "_dedx_vs_z").c_str(), "Energy loss vs. Z", 6000, -300, 300);
54  m_dedx_vs_r = new TProfile((m_name + "_dedx_vs_r").c_str(), "Energy loss vs. R", 1200, 0, 120);
55  m_radlen_vs_eta = new TProfile((m_name + "_radlen_vs_eta").c_str(), "Radiation lengths vs. eta", 600, -3, 3);
56  m_radlen_vs_z = new TProfile((m_name + "_radlen_vs_z").c_str(), "Radiation lengths vs. Z", 6000, -300, 300);
57  m_radlen_vs_r = new TProfile((m_name + "_radlen_vs_r").c_str(), "Radiation lengths vs. R", 1200, 0, 120);
58  m_dedx_spectrum->SetDirectory( 0 );
59  m_radlen_spectrum->SetDirectory( 0 );
60  m_dedx_vs_eta->SetDirectory( 0 );
61  m_dedx_vs_z->SetDirectory( 0 );
62  m_dedx_vs_r->SetDirectory( 0 );
63  m_radlen_vs_eta->SetDirectory( 0 );
64  m_radlen_vs_z->SetDirectory( 0 );
65  m_radlen_vs_r->SetDirectory( 0 );
66 }
int i
Definition: DBlmapReader.cc:9
MaterialAccountingStep m_errors
T perp() const
Magnitude of transverse component.
std::vector< GlobalPoint > m_elements
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
Definition: DDAxes.h:10
const std::string & name(void) const
get the layer name
tuple filter
USE THIS FOR SKIMMED TRACKS process.p = cms.Path(process.hltLevel1GTSeed*process.skimming*process.offlineBeamSpot*process.TrackRefitter2) OTHERWISE USE THIS.
Definition: align_tpl.py:86
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
MaterialAccountingStep m_accounting
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37
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 99 of file MaterialAccountingGroup.cc.

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

100 {
101  if (not inside(detector))
102  return false;
103 
104  // multiple hits in the same layer (overlaps, etc.) from a single track still count as one for averaging,
105  // since the energy deposits from the track have been already split between the different detectors
106  m_buffer += detector.material();
107  m_counted = true;
108 
109  return true;
110 }
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 97 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 115 of file MaterialAccountingGroup.h.

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

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

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

return the average normalized layer thickness

Definition at line 103 of file MaterialAccountingGroup.h.

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

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

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

return the average normalized number of radiation lengths

Definition at line 109 of file MaterialAccountingGroup.h.

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

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

110  {
112  }
double radiationLengths(void) const
MaterialAccountingStep m_accounting
void MaterialAccountingGroup::endOfTrack ( void  )

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

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

112  {
113  // add a detector
114  if (m_counted) {
117  ++m_tracks;
118 
119  GlobalPoint average( (m_buffer.in().x() + m_buffer.out().x()) / 2.,
120  (m_buffer.in().y() + m_buffer.out().y()) / 2.,
121  (m_buffer.in().z() + m_buffer.out().z()) / 2. );
122  m_dedx_spectrum->Fill( m_buffer.energyLoss() );
123  m_radlen_spectrum->Fill( m_buffer.radiationLengths() );
124  m_dedx_vs_eta->Fill( average.eta(), m_buffer.energyLoss(), 1. );
125  m_dedx_vs_z->Fill( average.z(), m_buffer.energyLoss(), 1. );
126  m_dedx_vs_r->Fill( average.perp(), m_buffer.energyLoss(), 1. );
127  m_radlen_vs_eta->Fill( average.eta(), m_buffer.radiationLengths(), 1. );
128  m_radlen_vs_z->Fill( average.z(), m_buffer.radiationLengths(), 1. );
129  m_radlen_vs_r->Fill( average.perp(), m_buffer.radiationLengths(), 1. );
130  }
131  m_counted = false;
132  m_buffer = MaterialAccountingStep();
133 }
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::string MaterialAccountingGroup::info ( void  ) const

get some infos

Definition at line 177 of file MaterialAccountingGroup.cc.

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

178 {
179  std::stringstream out;
180  out << std::setw(48) << std::left << m_name << std::right << std::fixed;;
181  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;
182  out << ", " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().first << " < R < " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().second;
183  out << " Elements: " << std::setw(6) << m_elements.size();
184  return out.str();
185 }
std::vector< GlobalPoint > m_elements
std::pair< double, double > range_r() const
std::pair< double, double > range_z() const
tuple out
Definition: dbtoconf.py:99
bool MaterialAccountingGroup::inside ( const MaterialAccountingDetector detector) const

check if detector is inside any part of this layer

Definition at line 84 of file MaterialAccountingGroup.cc.

References i, MaterialAccountingGroup::BoundingBox::inside(), m_boundingbox, m_elements, PV3DBase< T, PVType, FrameType >::perp(), MaterialAccountingDetector::position(), position, s_tolerance, and PV3DBase< T, PVType, FrameType >::z().

Referenced by addDetector().

85 {
86  const GlobalPoint & position = detector.position();
87  // first check to see if the point is inside the bounding box
88  if (not m_boundingbox.inside(position.perp(), position.z())) {
89  return false;
90  } else {
91  // now check if the point is actually close enough to any element
92  for (unsigned int i = 0; i < m_elements.size(); ++i)
93  if ((position - m_elements[i]).mag2() < (s_tolerance * s_tolerance))
94  return true;
95  return false;
96  }
97 }
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:66
std::vector< GlobalPoint > m_elements
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
T z() const
Definition: PV3DBase.h:58
const GlobalPoint & position() const
bool inside(double r, double z) const
const std::string& MaterialAccountingGroup::name ( void  ) const
inline

get the layer name

Definition at line 145 of file MaterialAccountingGroup.h.

References m_name.

Referenced by BeautifulSoup.Tag::_invert(), TrackingMaterialAnalyser::saveParameters(), and TrackingMaterialAnalyser::saveXml().

146  {
147  return m_name;
148  }
MaterialAccountingGroup& MaterialAccountingGroup::operator= ( const MaterialAccountingGroup layer)
private

stop default assignment operator

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

Definition at line 135 of file MaterialAccountingGroup.cc.

References svgfig::canvas(), and m_file.

Referenced by savePlots().

136 {
137  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
138  plot->SetFillColor(15); // grey
139  plot->SetLineColor(1); // black
140  plot->Draw("c e");
141  canvas.GetFrame()->SetFillColor(kWhite);
142  canvas.Draw();
143  canvas.SaveAs((name + ".png").c_str(), "");
144 
145  // store te plot into m_file
146  plot->SetDirectory(m_file);
147 }
def canvas
Definition: svgfig.py:481
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 149 of file MaterialAccountingGroup.cc.

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

150 {
151  // Nota Bene:
152  // these "line" plots are not deleted explicitly since
153  // - deleting them before saving them to a TFile will not save them
154  // - deleting them after the TFile they're stored into results in a SEGV
155  // ROOT is probably "taking care" (read: messing things up) somehow...
156  TH1F * line = new TH1F((name + "_par").c_str(), "Parametrization", 1, plot->GetXaxis()->GetXmin(), plot->GetXaxis()->GetXmax());
157  line->SetBinContent(1, average);
158 
159  TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024);
160  plot->SetFillColor(15); // grey
161  plot->SetLineColor(1); // black
162  plot->SetLineWidth(2);
163  plot->Draw("c e6");
164  line->SetLineColor(2); // red
165  line->SetLineWidth(2);
166  line->Draw("same");
167  canvas.GetFrame()->SetFillColor(kWhite);
168  canvas.Draw();
169  canvas.SaveAs((name + ".png").c_str(), "");
170 
171  // store te plots into m_file
172  plot->SetDirectory(m_file);
173  line->SetDirectory(m_file);
174 }
def canvas
Definition: svgfig.py:481
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 188 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().

189 {
190  m_file = new TFile((m_name + ".root").c_str(), "RECREATE");
191  savePlot(m_dedx_spectrum, m_name + "_dedx_spectrum");
192  savePlot(m_radlen_spectrum, m_name + "_radlen_spectrum");
193  savePlot(m_dedx_vs_eta, averageEnergyLoss(), m_name + "_dedx_vs_eta");
194  savePlot(m_dedx_vs_z, averageEnergyLoss(), m_name + "_dedx_vs_z");
195  savePlot(m_dedx_vs_r, averageEnergyLoss(), m_name + "_dedx_vs_r");
196  savePlot(m_radlen_vs_eta, averageRadiationLengths(), m_name + "_radlen_vs_eta");
199  m_file->Write();
200  m_file->Close();
201 
202  delete m_file;
203 }
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 133 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

134  {
136  }
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:28
double MaterialAccountingGroup::sigmaLength ( void  ) const
inline

return the sigma of the normalized layer thickness

Definition at line 121 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

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

return the sigma of the normalized number of radiation lengths

Definition at line 127 of file MaterialAccountingGroup.h.

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

Referenced by TrackingMaterialAnalyser::saveParameters().

128  {
130  }
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:28
unsigned int MaterialAccountingGroup::tracks ( void  ) const
inline

return the number of tracks that hit this layer

Definition at line 139 of file MaterialAccountingGroup.h.

References m_tracks.

Referenced by TrackingMaterialAnalyser::saveParameters().

140  {
141  return m_tracks;
142  }

Member Data Documentation

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

Definition at line 162 of file MaterialAccountingGroup.h.

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

MaterialAccountingStep MaterialAccountingGroup::m_buffer
private

Definition at line 167 of file MaterialAccountingGroup.h.

Referenced by addDetector(), and endOfTrack().

bool MaterialAccountingGroup::m_counted
private

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

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

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

Definition at line 181 of file MaterialAccountingGroup.h.

Referenced by savePlot(), and savePlots().

std::string MaterialAccountingGroup::m_name
private

Definition at line 160 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 MaterialAccountingGroup::s_tolerance = 0.01
staticprivate

Definition at line 184 of file MaterialAccountingGroup.h.

Referenced by inside(), and MaterialAccountingGroup().