CMS 3D CMS Logo

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

#include <TrackingMaterialAnalyser.h>

Inheritance diagram for TrackingMaterialAnalyser:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 TrackingMaterialAnalyser (const edm::ParameterSet &)
 
virtual ~TrackingMaterialAnalyser ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

enum  SplitMode { NEAREST_LAYER, INNER_LAYER, OUTER_LAYER, UNDEFINED }
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
 
void beginJob ()
 
void endJob ()
 
int findLayer (const MaterialAccountingDetector &detector)
 
void saveLayerPlots ()
 
void saveParameters (const char *name)
 
void saveXml (const char *name)
 
void split (MaterialAccountingTrack &track)
 

Private Attributes

std::vector< std::string > m_groupNames
 
std::vector
< MaterialAccountingGroup * > 
m_groups
 
edm::InputTag m_material
 
TrackingMaterialPlotterm_plotter
 
bool m_saveDetailedPlots
 
bool m_saveParameters
 
bool m_saveSummaryPlot
 
bool m_saveXml
 
bool m_skipAfterLastDetector
 
bool m_skipBeforeFirstDetector
 
SplitMode m_splitMode
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 16 of file TrackingMaterialAnalyser.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

TrackingMaterialAnalyser::TrackingMaterialAnalyser ( const edm::ParameterSet iPSet)

Definition at line 30 of file TrackingMaterialAnalyser.cc.

References edm::hlt::Exception, edm::ParameterSet::getParameter(), INNER_LAYER, edm::errors::LogicError, m_groupNames, m_material, m_plotter, m_saveDetailedPlots, m_saveParameters, m_saveSummaryPlot, m_saveXml, m_skipAfterLastDetector, m_skipBeforeFirstDetector, m_splitMode, NEAREST_LAYER, NULL, OUTER_LAYER, AlCaHLTBitMon_QueryRunRegistry::string, and UNDEFINED.

31 {
32  m_material = iPSet.getParameter<edm::InputTag>("MaterialAccounting");
33  m_groupNames = iPSet.getParameter<std::vector<std::string> >("Groups");
34  const std::string & splitmode = iPSet.getParameter<std::string>("SplitMode");
35  if (strcasecmp(splitmode.c_str(), "NearestLayer") == 0) {
37  } else if (strcasecmp(splitmode.c_str(), "InnerLayer") == 0) {
39  } else if (strcasecmp(splitmode.c_str(), "OuterLayer") == 0) {
41  } else {
43  throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode \"" << splitmode << "\". Acceptable values are \"NearestLayer\", \"InnerLayer\", \"OuterLayer\".";
44  }
45  m_skipAfterLastDetector = iPSet.getParameter<bool>("SkipAfterLastDetector");
46  m_skipBeforeFirstDetector = iPSet.getParameter<bool>("SkipBeforeFirstDetector");
47  m_saveSummaryPlot = iPSet.getParameter<bool>("SaveSummaryPlot");
48  m_saveDetailedPlots = iPSet.getParameter<bool>("SaveDetailedPlots");
49  m_saveParameters = iPSet.getParameter<bool>("SaveParameters");
50  m_saveXml = iPSet.getParameter<bool>("SaveXML");
51  if (m_saveSummaryPlot)
52  m_plotter = new TrackingMaterialPlotter( 300., 120., 10 ); // 10x10 points per cm2
53  else
54  m_plotter = NULL;
55 }
T getParameter(std::string const &) const
#define NULL
Definition: scimark2.h:8
std::vector< std::string > m_groupNames
TrackingMaterialPlotter * m_plotter
TrackingMaterialAnalyser::~TrackingMaterialAnalyser ( void  )
virtual

Definition at line 58 of file TrackingMaterialAnalyser.cc.

References m_plotter.

59 {
60  if (m_plotter)
61  delete m_plotter;
62 }
TrackingMaterialPlotter * m_plotter

Member Function Documentation

void TrackingMaterialAnalyser::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
privatevirtual

Implements edm::EDAnalyzer.

Definition at line 136 of file TrackingMaterialAnalyser.cc.

References gather_cfg::cout, end, edm::EventSetup::get(), i, info, m_groupNames, m_groups, m_material, split(), and lumiQTWidget::t.

137 {
139  setup.get<IdealGeometryRecord>().get( hDDD );
140 
141  m_groups.reserve( m_groupNames.size() );
142  for (unsigned int i = 0; i < m_groupNames.size(); ++i)
143  m_groups.push_back( new MaterialAccountingGroup( m_groupNames[i], * hDDD) );
144 
145  // INFO
146  std::cout << "TrackingMaterialAnalyser: List of the tracker groups: " << std::endl;
147  for (unsigned int i = 0; i < m_groups.size(); ++i)
148  std::cout << '\t' << m_groups[i]->info() << std::endl;
149  std::cout << std::endl;
150 
152  event.getByLabel(m_material, h_tracks);
153 
154  for (std::vector<MaterialAccountingTrack>::const_iterator t = h_tracks->begin(), end = h_tracks->end(); t != end; ++t) {
155  MaterialAccountingTrack track(*t);
156  split( track );
157  }
158 }
int i
Definition: DBlmapReader.cc:9
void split(MaterialAccountingTrack &track)
#define end
Definition: vmac.h:38
std::vector< std::string > m_groupNames
std::vector< MaterialAccountingGroup * > m_groups
const T & get() const
Definition: EventSetup.h:55
tuple cout
Definition: gather_cfg.py:121
void TrackingMaterialAnalyser::beginJob ( void  )
inlineprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 31 of file TrackingMaterialAnalyser.h.

31 {}
void TrackingMaterialAnalyser::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 116 of file TrackingMaterialAnalyser.cc.

References m_plotter, m_saveDetailedPlots, m_saveParameters, m_saveSummaryPlot, m_saveXml, saveLayerPlots(), saveParameters(), and saveXml().

117 {
118  if (m_saveParameters)
119  saveParameters("parameters");
120 
121  if (m_saveXml)
122  saveXml("parameters.xml");
123 
125  saveLayerPlots();
126 
127  if (m_saveSummaryPlot and m_plotter) {
128  m_plotter->normalize();
129  m_plotter->draw();
130  }
131 }
void saveXml(const char *name)
void saveParameters(const char *name)
TrackingMaterialPlotter * m_plotter
int TrackingMaterialAnalyser::findLayer ( const MaterialAccountingDetector detector)
private

Definition at line 343 of file TrackingMaterialAnalyser.cc.

References dtNoiseDBValidation_cfg::cerr, i, getHLTprescales::index, m_groups, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), MaterialAccountingDetector::position(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by split().

344 {
345  int index = 0;
346  size_t inside = 0;
347  for (size_t i = 0; i < m_groups.size(); ++i)
348  if (m_groups[i]->inside(detector)) {
349  ++inside;
350  index = i+1;
351  }
352  if (inside == 0) {
353  index = 0;
354  std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector does not belong to any DetLayer" << std::endl;
355  std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
356  << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
357  << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
358  << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")"
359  << std::endl;
360  }
361  if (inside > 1) {
362  index = 0;
363  std::cerr << "TrackingMaterialAnalyser::findLayer(...): ERROR: detector belongs to " << inside << "DetLayers" << std::endl;
364  std::cerr << "TrackingMaterialAnalyser::findLayer(...): detector position: " << std::fixed
365  << " (r: " << std::setprecision(1) << std::setw(5) << detector.position().perp()
366  << ", z: " << std::setprecision(1) << std::setw(6) << detector.position().z()
367  << ", phi: " << std::setprecision(3) << std::setw(6) << detector.position().phi() << ")"
368  << std::endl;
369  }
370 
371  return index;
372 }
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:72
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T z() const
Definition: PV3DBase.h:64
std::vector< MaterialAccountingGroup * > m_groups
const GlobalPoint & position() const
void TrackingMaterialAnalyser::saveLayerPlots ( )
private

Definition at line 107 of file TrackingMaterialAnalyser.cc.

References i, m_groups, and MaterialAccountingGroup::savePlots().

Referenced by endJob().

108 {
109  for (unsigned int i = 0; i < m_groups.size(); ++i) {
110  MaterialAccountingGroup & layer = *(m_groups[i]);
111  layer.savePlots();
112  }
113 }
int i
Definition: DBlmapReader.cc:9
std::vector< MaterialAccountingGroup * > m_groups
void savePlots(void)
save the plots
void TrackingMaterialAnalyser::saveParameters ( const char *  name)
private

Definition at line 65 of file TrackingMaterialAnalyser.cc.

References MaterialAccountingGroup::averageEnergyLoss(), MaterialAccountingGroup::averageLength(), MaterialAccountingGroup::averageRadiationLengths(), gather_cfg::cout, cfg-viewer::format(), i, m_groups, MaterialAccountingGroup::name(), Parameters::parameters, MaterialAccountingGroup::sigmaEnergyLoss(), MaterialAccountingGroup::sigmaLength(), MaterialAccountingGroup::sigmaRadiationLengths(), and MaterialAccountingGroup::tracks().

Referenced by endJob().

66 {
67  std::ofstream parameters(name);
68  std::cout << std::endl;
69  for (unsigned int i = 0; i < m_groups.size(); ++i) {
70  MaterialAccountingGroup & layer = *(m_groups[i]);
71  std::cout << layer.name() << std::endl;
72  std::cout << boost::format("\tnumber of hits: %9d") % layer.tracks() << std::endl;
73  std::cout << boost::format("\tnormalized segment length: %9.1f ± %9.1f cm") % layer.averageLength() % layer.sigmaLength() << std::endl;
74  std::cout << boost::format("\tnormalized radiation lengths: %9.3f ± %9.3f") % layer.averageRadiationLengths() % layer.sigmaRadiationLengths() << std::endl;
75  std::cout << boost::format("\tnormalized energy loss: %9.3f ± %9.3f MeV") % layer.averageEnergyLoss() % layer.sigmaEnergyLoss() << std::endl;
76  parameters << boost::format("%-20s\t%7d\t%5.1f ± %5.1f cm\t%6.4f ± %6.4f \t%6.4fe-03 ± %6.4fe-03 GeV")
77  % layer.name()
78  % layer.tracks()
79  % layer.averageLength() % layer.sigmaLength()
81  % layer.averageEnergyLoss() % layer.sigmaEnergyLoss()
82  << std::endl;
83  }
84  std::cout << std::endl;
85 
86  parameters.close();
87 }
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
double averageRadiationLengths(void) const
return the average normalized number of radiation lengths
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
double sigmaEnergyLoss(void) const
return the sigma of the normalized energy loss density factor for Bethe-Bloch
unsigned int tracks(void) const
return the number of tracks that hit this layer
double averageLength(void) const
return the average normalized layer thickness
double averageEnergyLoss(void) const
return the average normalized energy loss density factor for Bethe-Bloch
const std::string & name(void) const
get the layer name
std::vector< MaterialAccountingGroup * > m_groups
tuple cout
Definition: gather_cfg.py:121
void TrackingMaterialAnalyser::saveXml ( const char *  name)
private

Definition at line 90 of file TrackingMaterialAnalyser.cc.

References MaterialAccountingGroup::averageEnergyLoss(), MaterialAccountingGroup::averageRadiationLengths(), i, m_groups, and MaterialAccountingGroup::name().

Referenced by endJob().

91 {
92  std::ofstream xml(name);
93  xml << "<?xml version=\"1.0\" encoding=\"utf-8\"?>" << std::endl;
94  xml << "<Groups>" << std::endl;
95  for (unsigned int i = 0; i < m_groups.size(); ++i) {
96  MaterialAccountingGroup & layer = *(m_groups[i]);
97  xml << " <Group name=\"" << layer.name() << "\">\n"
98  << " <Parameter name=\"TrackerRadLength\" value=\"" << layer.averageRadiationLengths() << "\"/>\n"
99  << " <Parameter name=\"TrackerXi\" value=\"" << layer.averageEnergyLoss() << "\"/>\n"
100  << " </Group>\n"
101  << std::endl;
102  }
103  xml << "</Groups>" << std::endl;
104 }
int i
Definition: DBlmapReader.cc:9
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
const std::string & name(void) const
get the layer name
std::vector< MaterialAccountingGroup * > m_groups
void TrackingMaterialAnalyser::split ( MaterialAccountingTrack track)
private

Definition at line 167 of file TrackingMaterialAnalyser.cc.

References begin, dtNoiseDBValidation_cfg::cerr, MuonGeometrySanityCheck_cfi::detectors(), end, findLayer(), i, getHLTprescales::index, INNER_LAYER, MaterialAccountingStep::length(), edm::errors::LogicError, MaterialAccountingTrack::m_detectors, m_groups, m_plotter, m_skipAfterLastDetector, m_skipBeforeFirstDetector, m_splitMode, MaterialAccountingTrack::m_steps, MaterialAccountingTrack::m_total, NEAREST_LAYER, or, OUTER_LAYER, CfgNavigationSchool_cfi::parts, TrackingMaterialPlotter::plotSegmentInLayer(), TrackingMaterialPlotter::plotSegmentUnassigned(), MaterialAccountingStep::split(), relval_parameters_module::step, and UNDEFINED.

Referenced by analyze().

168 {
169  // group sensitive detectors by their DetLayer
170  std::vector<int> group( track.m_detectors.size() );
171  for (unsigned int i = 0; i < track.m_detectors.size(); ++i)
172  group[i] = findLayer( track.m_detectors[i] );
173 
174  unsigned int detectors = track.m_detectors.size();
175  if (detectors == 0) {
176  // the track doesn't cross any active detector:
177  // keep al material as unassigned
178  if (m_plotter)
179  for (unsigned int i = 1; i < track.m_steps.size(); ++i)
181  } else {
182  const double TOLERANCE = 0.0001; // 1 um tolerance
183  std::vector<double> limits(detectors + 2);
184 
185  // define the trivial limits
187  limits[0] = track.m_detectors[0].m_curvilinearIn - TOLERANCE;
188  else
189  limits[0] = - TOLERANCE;
191  limits[detectors] = track.m_detectors[detectors-1].m_curvilinearOut + TOLERANCE;
192  else
193  limits[detectors] = track.m_total.length() + TOLERANCE;
194  limits[detectors+1] = INFINITY; // this is probably no more needed, but doesn't harm...
195 
196  // pick the algorithm to define the non-trivial limits
197  switch (m_splitMode) {
198  // assign each segment to the the nearest layer
199  // e.g. the material between pixel barrel 3 and TIB 1 will be split among the two
200  case NEAREST_LAYER:
201  for (unsigned int i = 1; i < detectors; ++i)
202  limits[i] = (track.m_detectors[i-1].m_curvilinearOut + track.m_detectors[i].m_curvilinearIn) / 2.;
203  break;
204 
205  // assign each segment to the the inner layer
206  // e.g. all material between pixel barrel 3 and TIB 1 will go into the pixel barrel
207  case INNER_LAYER:
208  for (unsigned int i = 1; i < detectors; ++i)
209  limits[i] = track.m_detectors[i].m_curvilinearIn - TOLERANCE;
210  break;
211 
212  // assign each segment to the the outer layer
213  // e.g. all material between pixel barrel 3 and TIB 1 will go into the TIB
214  case OUTER_LAYER:
215  for (unsigned int i = 1; i < detectors; ++i)
216  limits[i] = track.m_detectors[i-1].m_curvilinearOut + TOLERANCE;
217  break;
218 
219  case UNDEFINED:
220  default:
221  // throw something
222  throw edm::Exception(edm::errors::LogicError) << "Invalid SplitMode";
223  }
224 
225  //for (unsigned int i = 0; i < detectors; ++i)
226  // std::cout << "MaterialAccountingTrack::split(): detector region boundaries: [" << limits[i] << ", " << limits[i+1] << "] along track" << std::endl;
227 
228  double begin = 0.; // begginning of step, along the track
229  double end = 0.; // end of step, along the track
230  unsigned int i = 1; // step conter
231 
232  // skip the material before the first layer
233  //std::cout << "before first layer, skipping" << std::endl;
234  while (end < limits[0]) {
235  const MaterialAccountingStep & step = track.m_steps[i++];
236  end = begin + step.length();
237 
238  // do not account material before the first layer
239  if (m_plotter)
241 
242  begin = end;
243  //std::cout << '.';
244  }
245  //std::cout << std::endl;
246 
247  // optionally split a step across the first layer boundary
248  //std::cout << "first layer (0): " << limits[0] << ".." << limits[1] << std::endl;
249  if (begin < limits[0] and end > limits[0]) {
250  const MaterialAccountingStep & step = track.m_steps[i++];
251  end = begin + step.length();
252 
253  double fraction = (limits[0] - begin) / (end - begin);
254  std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
255 
256  //std::cout << '!' << std::endl;
257  track.m_detectors[0].account( parts.second, limits[1], end );
258 
259  if (m_plotter) {
260  // step partially before first layer, keep first part as unassocated
261  m_plotter->plotSegmentUnassigned( parts.first );
262 
263  // associate second part to first layer
264  m_plotter->plotSegmentInLayer( parts.second, group[0] );
265  }
266  begin = end;
267  }
268 
269  unsigned int index = 0; // which detector
270  while (i < track.m_steps.size()) {
271  const MaterialAccountingStep & step = track.m_steps[i++];
272 
273  end = begin + step.length();
274 
275  if (begin > limits[detectors]) {
276  // segment after last layer and skipping requested in configuation
277  if (m_plotter)
279  begin = end;
280  continue;
281  }
282 
283  // from here onwards we should be in the accountable region, either completely in a single layer:
284  // limits[index] <= begin < end <= limits[index+1]
285  // or possibly split between 2 layers
286  // limits[index] < begin < limits[index+1] < end < limits[index+2]
287  if (begin < limits[index] or end > limits[index+2]) {
288  // sanity check
289  std::cerr << "MaterialAccountingTrack::split(): ERROR: internal logic error, expected " << limits[index] << " < " << begin << " < " << limits[index+1] << std::endl;
290  break;
291  }
292 
293  //std::cout << '.';
294  if (limits[index] <= begin and end <= limits[index+1]) {
295  // step completely inside current detector range
296  track.m_detectors[index].account( step, begin, end );
297  if (m_plotter)
298  m_plotter->plotSegmentInLayer( step, group[index] );
299  } else {
300  // step shared beteewn two detectors, transition at limits[index+1]
301  double fraction = (limits[index+1] - begin) / (end - begin);
302  std::pair<MaterialAccountingStep, MaterialAccountingStep> parts = step.split(fraction);
303 
304  if (m_plotter) {
305  if (index > 0)
306  m_plotter->plotSegmentInLayer( parts.first, group[index] );
307  else
308  // track outside acceptance, keep as unassocated
309  m_plotter->plotSegmentUnassigned( parts.first );
310 
311  if (index+1 < detectors)
312  m_plotter->plotSegmentInLayer( parts.second, group[index+1] );
313  else
314  // track outside acceptance, keep as unassocated
315  m_plotter->plotSegmentUnassigned( parts.second );
316  }
317 
318  track.m_detectors[index].account( parts.first, begin, limits[index+1] );
319  ++index; // next layer
320  //std::cout << '!' << std::endl;
321  //std::cout << "next layer (" << index << "): " << limits[index] << ".." << limits[index+1] << std::endl;
322  if (index < detectors)
323  track.m_detectors[index].account( parts.second, limits[index+1], end );
324  }
325  begin = end;
326  }
327 
328  }
329  //std::cout << std::endl;
330 
331  // add the material from each detector to its layer (if there is one and only one)
332  for (unsigned int i = 0; i < track.m_detectors.size(); ++i)
333  if (group[i] != 0)
334  m_groups[group[i]-1]->addDetector( track.m_detectors[i] );
335 
336  // end of track: commit internal buffers and reset the m_groups internal state for a new track
337  for (unsigned int i = 0; i < m_groups.size(); ++i)
338  m_groups[i]->endOfTrack();
339 }
void plotSegmentUnassigned(const MaterialAccountingStep &step)
int i
Definition: DBlmapReader.cc:9
std::vector< MaterialAccountingDetector > m_detectors
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::vector< MaterialAccountingStep > m_steps
double length(void) const
MaterialAccountingStep m_total
int findLayer(const MaterialAccountingDetector &detector)
#define end
Definition: vmac.h:38
TrackingMaterialPlotter * m_plotter
std::vector< MaterialAccountingGroup * > m_groups
std::pair< MaterialAccountingStep, MaterialAccountingStep > split(double fraction) const
split the step (0..1) in (0..f) + (f..1) using linear interpolation
#define begin
Definition: vmac.h:31
void plotSegmentInLayer(const MaterialAccountingStep &step, int layer)

Member Data Documentation

std::vector<std::string> TrackingMaterialAnalyser::m_groupNames
private

Definition at line 50 of file TrackingMaterialAnalyser.h.

Referenced by analyze(), and TrackingMaterialAnalyser().

std::vector<MaterialAccountingGroup *> TrackingMaterialAnalyser::m_groups
private
edm::InputTag TrackingMaterialAnalyser::m_material
private

Definition at line 41 of file TrackingMaterialAnalyser.h.

Referenced by analyze(), and TrackingMaterialAnalyser().

TrackingMaterialPlotter* TrackingMaterialAnalyser::m_plotter
private
bool TrackingMaterialAnalyser::m_saveDetailedPlots
private

Definition at line 46 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_saveParameters
private

Definition at line 47 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_saveSummaryPlot
private

Definition at line 45 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_saveXml
private

Definition at line 48 of file TrackingMaterialAnalyser.h.

Referenced by endJob(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_skipAfterLastDetector
private

Definition at line 43 of file TrackingMaterialAnalyser.h.

Referenced by split(), and TrackingMaterialAnalyser().

bool TrackingMaterialAnalyser::m_skipBeforeFirstDetector
private

Definition at line 44 of file TrackingMaterialAnalyser.h.

Referenced by split(), and TrackingMaterialAnalyser().

SplitMode TrackingMaterialAnalyser::m_splitMode
private

Definition at line 42 of file TrackingMaterialAnalyser.h.

Referenced by split(), and TrackingMaterialAnalyser().