#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 | |
MaterialAccountingStep | average (void) const |
return the average normalized material accounting informations | |
double | averageEnergyLoss (void) const |
return the average normalized energy loss density factor for Bethe-Bloch | |
double | averageLength (void) const |
return the average normalized layer thickness | |
double | averageRadiationLengths (void) const |
return the average normalized number of radiation lengths | |
void | endOfTrack (void) |
commit the buffer and reset the "already hit by this track" flag | |
std::string | info (void) const |
get some infos | |
bool | inside (const MaterialAccountingDetector &detector) const |
check if detector is inside any part of this layer | |
MaterialAccountingGroup (const std::string &name, const DDCompactView &geometry) | |
explicit constructors | |
const std::string & | name (void) const |
get the layer name | |
void | savePlots (void) |
save the plots | |
double | sigmaEnergyLoss (void) const |
return the sigma of the normalized energy loss density factor for Bethe-Bloch | |
double | sigmaLength (void) const |
return the sigma of the normalized layer thickness | |
double | sigmaRadiationLengths (void) const |
return the sigma of the normalized number of radiation lengths | |
unsigned int | tracks (void) const |
return the number of tracks that hit this layer | |
~MaterialAccountingGroup (void) | |
destructor | |
Private Member Functions | |
MaterialAccountingGroup (const MaterialAccountingGroup &layer) | |
stop default copy ctor | |
MaterialAccountingGroup & | operator= (const MaterialAccountingGroup &layer) |
stop default assignment operator | |
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< GlobalPoint > | m_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 |
Definition at line 18 of file MaterialAccountingGroup.h.
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.
: m_name( name ), m_elements(), m_boundingbox(), m_accounting(), m_errors(), m_tracks( 0 ), m_counted( false ), m_file( 0 ) { // retrieve the elements from DDD DDFilteredView fv( geometry ); DDSpecificsFilter filter; filter.setCriteria(DDValue("TrackingMaterialGroup", name), DDSpecificsFilter::equals); fv.addFilter(filter); while (fv.next()) { // DD3Vector and DDTranslation are the same type as math::XYZVector math::XYZVector position = fv.translation() / 10.; // mm -> cm m_elements.push_back( GlobalPoint(position.x(), position.y(), position.z()) ); } // grow the bounding box for (unsigned int i = 0; i < m_elements.size(); ++i) { m_boundingbox.grow(m_elements[i].perp(), m_elements[i].z()); } m_boundingbox.grow(s_tolerance); // initialize the histograms m_dedx_spectrum = new TH1F((m_name + "_dedx_spectrum").c_str(), "Energy loss spectrum", 1000, 0, 1); m_radlen_spectrum = new TH1F((m_name + "_radlen_spectrum").c_str(), "Radiation lengths spectrum", 1000, 0, 1); m_dedx_vs_eta = new TProfile((m_name + "_dedx_vs_eta").c_str(), "Energy loss vs. eta", 600, -3, 3); m_dedx_vs_z = new TProfile((m_name + "_dedx_vs_z").c_str(), "Energy loss vs. Z", 6000, -300, 300); m_dedx_vs_r = new TProfile((m_name + "_dedx_vs_r").c_str(), "Energy loss vs. R", 1200, 0, 120); m_radlen_vs_eta = new TProfile((m_name + "_radlen_vs_eta").c_str(), "Radiation lengths vs. eta", 600, -3, 3); m_radlen_vs_z = new TProfile((m_name + "_radlen_vs_z").c_str(), "Radiation lengths vs. Z", 6000, -300, 300); m_radlen_vs_r = new TProfile((m_name + "_radlen_vs_r").c_str(), "Radiation lengths vs. R", 1200, 0, 120); m_dedx_spectrum->SetDirectory( 0 ); m_radlen_spectrum->SetDirectory( 0 ); m_dedx_vs_eta->SetDirectory( 0 ); m_dedx_vs_z->SetDirectory( 0 ); m_dedx_vs_r->SetDirectory( 0 ); m_radlen_vs_eta->SetDirectory( 0 ); m_radlen_vs_z->SetDirectory( 0 ); m_radlen_vs_r->SetDirectory( 0 ); }
MaterialAccountingGroup::~MaterialAccountingGroup | ( | void | ) |
destructor
Definition at line 68 of file MaterialAccountingGroup.cc.
References m_dedx_spectrum, m_dedx_vs_eta, m_dedx_vs_r, m_dedx_vs_z, m_radlen_spectrum, m_radlen_vs_eta, m_radlen_vs_r, and m_radlen_vs_z.
{ delete m_dedx_spectrum; delete m_dedx_vs_eta; delete m_dedx_vs_z; delete m_dedx_vs_r; delete m_radlen_spectrum; delete m_radlen_vs_eta; delete m_radlen_vs_z; delete m_radlen_vs_r; }
MaterialAccountingGroup::MaterialAccountingGroup | ( | const MaterialAccountingGroup & | layer | ) | [private] |
stop default copy ctor
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().
{ if (not inside(detector)) return false; // multiple hits in the same layer (overlaps, etc.) from a single track still count as one for averaging, // since the energy deposits from the track have been already split between the different detectors m_buffer += detector.material(); m_counted = true; return true; }
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().
{ return m_tracks ? m_accounting / m_tracks : MaterialAccountingStep(); }
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().
{ return m_tracks ? m_accounting.energyLoss() / m_tracks : 0.; }
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().
{ return m_tracks ? m_accounting.length() / m_tracks : 0.; }
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().
{ return m_tracks ? m_accounting.radiationLengths() / m_tracks : 0.; }
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().
{ // add a detector if (m_counted) { m_accounting += m_buffer; m_errors += m_buffer * m_buffer; ++m_tracks; GlobalPoint average( (m_buffer.in().x() + m_buffer.out().x()) / 2., (m_buffer.in().y() + m_buffer.out().y()) / 2., (m_buffer.in().z() + m_buffer.out().z()) / 2. ); m_dedx_spectrum->Fill( m_buffer.energyLoss() ); m_radlen_spectrum->Fill( m_buffer.radiationLengths() ); m_dedx_vs_eta->Fill( average.eta(), m_buffer.energyLoss(), 1. ); m_dedx_vs_z->Fill( average.z(), m_buffer.energyLoss(), 1. ); m_dedx_vs_r->Fill( average.perp(), m_buffer.energyLoss(), 1. ); m_radlen_vs_eta->Fill( average.eta(), m_buffer.radiationLengths(), 1. ); m_radlen_vs_z->Fill( average.z(), m_buffer.radiationLengths(), 1. ); m_radlen_vs_r->Fill( average.perp(), m_buffer.radiationLengths(), 1. ); } m_counted = false; m_buffer = MaterialAccountingStep(); }
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().
{ std::stringstream out; out << std::setw(48) << std::left << m_name << std::right << std::fixed;; 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; out << ", " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().first << " < R < " << std::setprecision(1) << std::setw(5) << m_boundingbox.range_r().second; out << " Elements: " << std::setw(6) << m_elements.size(); return out.str(); }
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(), position, MaterialAccountingDetector::position(), s_tolerance, and PV3DBase< T, PVType, FrameType >::z().
Referenced by addDetector().
{ const GlobalPoint & position = detector.position(); // first check to see if the point is inside the bounding box if (not m_boundingbox.inside(position.perp(), position.z())) { return false; } else { // now check if the point is actually close enough to any element for (unsigned int i = 0; i < m_elements.size(); ++i) if ((position - m_elements[i]).mag2() < (s_tolerance * s_tolerance)) return true; return false; } }
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 TrackingMaterialAnalyser::saveParameters(), and TrackingMaterialAnalyser::saveXml().
{ return m_name; }
MaterialAccountingGroup& MaterialAccountingGroup::operator= | ( | const MaterialAccountingGroup & | layer | ) | [private] |
stop default assignment operator
void MaterialAccountingGroup::savePlot | ( | TProfile * | plot, |
float | average, | ||
const std::string & | name | ||
) | [private] |
Definition at line 149 of file MaterialAccountingGroup.cc.
References MultipleCompare::canvas, geometryCSVtoXML::line, and m_file.
{ // Nota Bene: // these "line" plots are not deleted explicitly since // - deleting them before saving them to a TFile will not save them // - deleting them after the TFile they're stored into results in a SEGV // ROOT is probably "taking care" (read: messing things up) somehow... TH1F * line = new TH1F((name + "_par").c_str(), "Parametrization", 1, plot->GetXaxis()->GetXmin(), plot->GetXaxis()->GetXmax()); line->SetBinContent(1, average); TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024); plot->SetFillColor(15); // grey plot->SetLineColor(1); // black plot->SetLineWidth(2); plot->Draw("c e6"); line->SetLineColor(2); // red line->SetLineWidth(2); line->Draw("same"); canvas.GetFrame()->SetFillColor(kWhite); canvas.Draw(); canvas.SaveAs((name + ".png").c_str(), ""); // store te plots into m_file plot->SetDirectory(m_file); line->SetDirectory(m_file); }
void MaterialAccountingGroup::savePlot | ( | TH1F * | plot, |
const std::string & | name | ||
) | [private] |
Definition at line 135 of file MaterialAccountingGroup.cc.
References MultipleCompare::canvas, and m_file.
Referenced by savePlots().
{ TCanvas canvas(name.c_str(), plot->GetTitle(), 1280, 1024); plot->SetFillColor(15); // grey plot->SetLineColor(1); // black plot->Draw("c e"); canvas.GetFrame()->SetFillColor(kWhite); canvas.Draw(); canvas.SaveAs((name + ".png").c_str(), ""); // store te plot into m_file plot->SetDirectory(m_file); }
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().
{ m_file = new TFile((m_name + ".root").c_str(), "RECREATE"); savePlot(m_dedx_spectrum, m_name + "_dedx_spectrum"); savePlot(m_radlen_spectrum, m_name + "_radlen_spectrum"); savePlot(m_dedx_vs_eta, averageEnergyLoss(), m_name + "_dedx_vs_eta"); savePlot(m_dedx_vs_z, averageEnergyLoss(), m_name + "_dedx_vs_z"); savePlot(m_dedx_vs_r, averageEnergyLoss(), m_name + "_dedx_vs_r"); savePlot(m_radlen_vs_eta, averageRadiationLengths(), m_name + "_radlen_vs_eta"); savePlot(m_radlen_vs_z, averageRadiationLengths(), m_name + "_radlen_vs_z"); savePlot(m_radlen_vs_r, averageRadiationLengths(), m_name + "_radlen_vs_r"); m_file->Write(); m_file->Close(); delete m_file; }
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().
{ return m_tracks ? std::sqrt(m_errors.energyLoss() / m_tracks - averageEnergyLoss()*averageEnergyLoss()) : 0.; }
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().
{ return m_tracks ? std::sqrt(m_errors.length() / m_tracks - averageLength()*averageLength()) : 0.; }
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().
{ return m_tracks ? std::sqrt(m_errors.radiationLengths() / m_tracks - averageRadiationLengths()*averageRadiationLengths()) : 0.; }
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().
{ return m_tracks; }
Definition at line 163 of file MaterialAccountingGroup.h.
Referenced by average(), averageEnergyLoss(), averageLength(), averageRadiationLengths(), and endOfTrack().
Definition at line 162 of file MaterialAccountingGroup.h.
Referenced by info(), inside(), and MaterialAccountingGroup().
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] |
Definition at line 170 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
TProfile* MaterialAccountingGroup::m_dedx_vs_eta [private] |
Definition at line 173 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
TProfile* MaterialAccountingGroup::m_dedx_vs_r [private] |
Definition at line 175 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
TProfile* MaterialAccountingGroup::m_dedx_vs_z [private] |
Definition at line 174 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
std::vector<GlobalPoint> MaterialAccountingGroup::m_elements [private] |
Definition at line 161 of file MaterialAccountingGroup.h.
Referenced by info(), inside(), and MaterialAccountingGroup().
Definition at line 164 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), sigmaEnergyLoss(), sigmaLength(), and sigmaRadiationLengths().
TFile* MaterialAccountingGroup::m_file [mutable, private] |
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] |
Definition at line 171 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
TProfile* MaterialAccountingGroup::m_radlen_vs_eta [private] |
Definition at line 176 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
TProfile* MaterialAccountingGroup::m_radlen_vs_r [private] |
Definition at line 178 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
TProfile* MaterialAccountingGroup::m_radlen_vs_z [private] |
Definition at line 177 of file MaterialAccountingGroup.h.
Referenced by endOfTrack(), MaterialAccountingGroup(), savePlots(), and ~MaterialAccountingGroup().
unsigned int MaterialAccountingGroup::m_tracks [private] |
Definition at line 165 of file MaterialAccountingGroup.h.
Referenced by average(), averageEnergyLoss(), averageLength(), averageRadiationLengths(), endOfTrack(), sigmaEnergyLoss(), sigmaLength(), sigmaRadiationLengths(), and tracks().
double MaterialAccountingGroup::s_tolerance = 0.01 [static, private] |
Definition at line 184 of file MaterialAccountingGroup.h.
Referenced by inside(), and MaterialAccountingGroup().